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