1c4762a1bSJed Brown static char help[] = "Test MatMatMult() for AIJ and Dense matrices.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown int main(int argc,char **argv) 6c4762a1bSJed Brown { 7c20d7725SJed Brown Mat A,B,C; 8c20d7725SJed Brown PetscInt M=10,N=5; 9c4762a1bSJed Brown PetscErrorCode ierr; 10c4762a1bSJed Brown PetscRandom r; 11c4762a1bSJed Brown PetscBool equal=PETSC_FALSE; 12c4762a1bSJed Brown PetscReal fill = 1.0; 13c4762a1bSJed Brown PetscInt nza,am,an; 14c4762a1bSJed Brown 15c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 16c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); 17c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr); 18c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL);CHKERRQ(ierr); 19c4762a1bSJed Brown 20c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&r);CHKERRQ(ierr); 21c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr); 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* create a aij matrix A */ 24c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 25c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,M);CHKERRQ(ierr); 26c4762a1bSJed Brown ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 27c4762a1bSJed Brown nza = (PetscInt)(.3*M); /* num of nozeros in each row of A */ 28c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(A,nza,NULL);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(A,nza,NULL,nza,NULL);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = MatSetRandom(A,r);CHKERRQ(ierr); 31c4762a1bSJed Brown 32c4762a1bSJed Brown /* create a dense matrix B */ 33c4762a1bSJed Brown ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 34c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr); 35c4762a1bSJed Brown ierr = MatSetSizes(B,PETSC_DECIDE,am,N,PETSC_DECIDE);CHKERRQ(ierr); 36c4762a1bSJed Brown ierr = MatSetType(B,MATDENSE);CHKERRQ(ierr); 37c4762a1bSJed Brown ierr = MatSeqDenseSetPreallocation(B,NULL);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = MatMPIDenseSetPreallocation(B,NULL);CHKERRQ(ierr); 39c4762a1bSJed Brown ierr = MatSetRandom(B,r);CHKERRQ(ierr); 40c4762a1bSJed Brown ierr = PetscRandomDestroy(&r);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 42c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* Test MatMatMult() */ 45c4762a1bSJed Brown ierr = MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); 46c20d7725SJed Brown ierr = MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr); 47c20d7725SJed Brown ierr = MatMatMultEqual(B,A,C,10,&equal);CHKERRQ(ierr); 48*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != B*A"); 49c4762a1bSJed Brown 50c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = PetscFinalize(); 54c4762a1bSJed Brown return ierr; 55c4762a1bSJed Brown } 56c4762a1bSJed Brown 57c4762a1bSJed Brown /*TEST 58c4762a1bSJed Brown 59c4762a1bSJed Brown test: 60c4762a1bSJed Brown output_file: output/ex122.out 61c4762a1bSJed Brown 62c4762a1bSJed Brown TEST*/ 63