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 PetscRandom r; 10c4762a1bSJed Brown PetscBool equal=PETSC_FALSE; 11c4762a1bSJed Brown PetscReal fill = 1.0; 12c4762a1bSJed Brown PetscInt nza,am,an; 13c4762a1bSJed Brown 14*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 15*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 16*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); 17*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL)); 18c4762a1bSJed Brown 19*9566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&r)); 20*9566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(r)); 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* create a aij matrix A */ 23*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 24*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,M)); 25*9566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATAIJ)); 26c4762a1bSJed Brown nza = (PetscInt)(.3*M); /* num of nozeros in each row of A */ 27*9566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,nza,NULL)); 28*9566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A,nza,NULL,nza,NULL)); 29*9566063dSJacob Faibussowitsch PetscCall(MatSetRandom(A,r)); 30c4762a1bSJed Brown 31c4762a1bSJed Brown /* create a dense matrix B */ 32*9566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&am,&an)); 33*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); 34*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,PETSC_DECIDE,am,N,PETSC_DECIDE)); 35*9566063dSJacob Faibussowitsch PetscCall(MatSetType(B,MATDENSE)); 36*9566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(B,NULL)); 37*9566063dSJacob Faibussowitsch PetscCall(MatMPIDenseSetPreallocation(B,NULL)); 38*9566063dSJacob Faibussowitsch PetscCall(MatSetRandom(B,r)); 39*9566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&r)); 40*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 41*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Test MatMatMult() */ 44*9566063dSJacob Faibussowitsch PetscCall(MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C)); 45*9566063dSJacob Faibussowitsch PetscCall(MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C)); 46*9566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(B,A,C,10,&equal)); 4728b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != B*A"); 48c4762a1bSJed Brown 49*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 50*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 51*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 52*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 53b122ec5aSJacob Faibussowitsch return 0; 54c4762a1bSJed Brown } 55c4762a1bSJed Brown 56c4762a1bSJed Brown /*TEST 57c4762a1bSJed Brown 58c4762a1bSJed Brown test: 59c4762a1bSJed Brown output_file: output/ex122.out 60c4762a1bSJed Brown 61c4762a1bSJed Brown TEST*/ 62