xref: /petsc/src/mat/tests/ex109.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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,D,AT;
8c4762a1bSJed Brown   PetscInt       i,M,N,Istart,Iend,n=7,j,J,Ii,m=8,am,an;
9c4762a1bSJed Brown   PetscScalar    v;
10c4762a1bSJed Brown   PetscRandom    r;
113fbca975SStefano Zampini   PetscBool      equal=PETSC_FALSE,flg;
12c20d7725SJed Brown   PetscReal      fill = 1.0,norm;
13c4762a1bSJed Brown   PetscMPIInt    size;
14c4762a1bSJed Brown 
15*327415f7SBarry Smith   PetscFunctionBeginUser;
169566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL));
20c4762a1bSJed Brown 
219566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&r));
229566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(r));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   /* Create a aij matrix A */
25c4762a1bSJed Brown   M    = N = m*n;
269566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
279566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
289566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATAIJ));
299566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
309566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL));
319566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(A,5,NULL));
32c4762a1bSJed Brown 
339566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
34c4762a1bSJed Brown   am   = Iend - Istart;
35c4762a1bSJed Brown   for (Ii=Istart; Ii<Iend; Ii++) {
36c4762a1bSJed Brown     v = -1.0; i = Ii/n; j = Ii - i*n;
379566063dSJacob Faibussowitsch     if (i>0)   {J = Ii - n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
389566063dSJacob Faibussowitsch     if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
399566063dSJacob Faibussowitsch     if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
409566063dSJacob Faibussowitsch     if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
419566063dSJacob Faibussowitsch     v = 4.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
42c4762a1bSJed Brown   }
439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   /* Create a dense matrix B */
479566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A,&am,&an));
489566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
499566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B,an,PETSC_DECIDE,PETSC_DECIDE,M));
509566063dSJacob Faibussowitsch   PetscCall(MatSetType(B,MATDENSE));
519566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(B,NULL));
529566063dSJacob Faibussowitsch   PetscCall(MatMPIDenseSetPreallocation(B,NULL));
539566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(B));
549566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(B,r));
559566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&r));
56c4762a1bSJed Brown 
57c20d7725SJed Brown   /* Test reuse of user-provided dense C (unassembled) -- not recommended usage */
589566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
599566063dSJacob Faibussowitsch   PetscCall(MatSetType(C,MATDENSE));
609566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C,am,PETSC_DECIDE,PETSC_DECIDE,N));
619566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
629566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
639566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(C));
649566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C));
659566063dSJacob Faibussowitsch   PetscCall(MatNorm(C,NORM_INFINITY,&norm));
669566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
67c20d7725SJed Brown 
68c4762a1bSJed Brown   /* Test C = A*B (aij*dense) */
699566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A,B,MAT_INITIAL_MATRIX,fill,&C));
709566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C));
71c4762a1bSJed Brown 
72c20d7725SJed Brown   /* Test developer API */
739566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(A,B,NULL,&D));
749566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(D,MATPRODUCT_AB));
759566063dSJacob Faibussowitsch   PetscCall(MatProductSetAlgorithm(D,"default"));
769566063dSJacob Faibussowitsch   PetscCall(MatProductSetFill(D,fill));
779566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(D));
789566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(D));
79c4762a1bSJed Brown   for (i=0; i<2; i++) {
809566063dSJacob Faibussowitsch     PetscCall(MatProductNumeric(D));
81c4762a1bSJed Brown   }
829566063dSJacob Faibussowitsch   PetscCall(MatEqual(C,D,&equal));
8328b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != D");
849566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
85c4762a1bSJed Brown 
86c20d7725SJed Brown   /* Test D = AT*B (transpose(aij)*dense) */
879566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(A,&AT));
889566063dSJacob Faibussowitsch   PetscCall(MatMatMult(AT,B,MAT_INITIAL_MATRIX,fill,&D));
899566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(AT,B,D,10,&equal));
9028b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != AT*B (transpose(aij)*dense)");
919566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
929566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&AT));
93c20d7725SJed Brown 
94c4762a1bSJed Brown   /* Test D = C*A (dense*aij) */
959566063dSJacob Faibussowitsch   PetscCall(MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&D));
969566063dSJacob Faibussowitsch   PetscCall(MatMatMult(C,A,MAT_REUSE_MATRIX,fill,&D));
979566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(C,A,D,10,&equal));
9828b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != C*A (dense*aij)");
999566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /* Test D = A*C (aij*dense) */
1029566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A,C,MAT_INITIAL_MATRIX,fill,&D));
1039566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A,C,MAT_REUSE_MATRIX,fill,&D));
1049566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(A,C,D,10,&equal));
10528b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != A*C (aij*dense)");
1069566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   /* Test D = B*C (dense*dense) */
1099566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
110c4762a1bSJed Brown   if (size == 1) {
1119566063dSJacob Faibussowitsch     PetscCall(MatMatMult(B,C,MAT_INITIAL_MATRIX,fill,&D));
1129566063dSJacob Faibussowitsch     PetscCall(MatMatMult(B,C,MAT_REUSE_MATRIX,fill,&D));
1139566063dSJacob Faibussowitsch     PetscCall(MatMatMultEqual(B,C,D,10,&equal));
11428b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != B*C (dense*dense)");
1159566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&D));
116c4762a1bSJed Brown   }
117c4762a1bSJed Brown 
118c20d7725SJed Brown   /* Test D = B*C^T (dense*dense) */
1199566063dSJacob Faibussowitsch   PetscCall(MatMatTransposeMult(B,C,MAT_INITIAL_MATRIX,fill,&D));
1209566063dSJacob Faibussowitsch   PetscCall(MatMatTransposeMult(B,C,MAT_REUSE_MATRIX,fill,&D));
1219566063dSJacob Faibussowitsch   PetscCall(MatMatTransposeMultEqual(B,C,D,10,&equal));
12228b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != B*C^T (dense*dense)");
1239566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
124c20d7725SJed Brown 
1253fbca975SStefano Zampini   /* Test MatProductCreateWithMat() and reuse C and B for B = A*C */
1263fbca975SStefano Zampini   flg = PETSC_FALSE;
1279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-test_userAPI",&flg));
1283fbca975SStefano Zampini   if (flg) {
1293fbca975SStefano Zampini     /* user driver */
1309566063dSJacob Faibussowitsch     PetscCall(MatMatMult(A,C,MAT_REUSE_MATRIX,fill,&B));
1313fbca975SStefano Zampini   } else {
1323fbca975SStefano Zampini     /* clear internal data structures related with previous products to avoid circular references */
1339566063dSJacob Faibussowitsch     PetscCall(MatProductClear(A));
1349566063dSJacob Faibussowitsch     PetscCall(MatProductClear(B));
1359566063dSJacob Faibussowitsch     PetscCall(MatProductClear(C));
1369566063dSJacob Faibussowitsch     PetscCall(MatProductCreateWithMat(A,C,NULL,B));
1379566063dSJacob Faibussowitsch     PetscCall(MatProductSetType(B,MATPRODUCT_AB));
1389566063dSJacob Faibussowitsch     PetscCall(MatProductSetFromOptions(B));
1399566063dSJacob Faibussowitsch     PetscCall(MatProductSymbolic(B));
1409566063dSJacob Faibussowitsch     PetscCall(MatProductNumeric(B));
1413fbca975SStefano Zampini   }
1423fbca975SStefano Zampini 
1439566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
1449566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1459566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1469566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
147b122ec5aSJacob Faibussowitsch   return 0;
148c4762a1bSJed Brown }
149c4762a1bSJed Brown 
150c4762a1bSJed Brown /*TEST
151c4762a1bSJed Brown 
152c4762a1bSJed Brown    test:
153c4762a1bSJed Brown       args: -M 10 -N 10
154c4762a1bSJed Brown       output_file: output/ex109.out
155c4762a1bSJed Brown 
156c4762a1bSJed Brown    test:
157c4762a1bSJed Brown       suffix: 2
158c4762a1bSJed Brown       nsize: 3
159c4762a1bSJed Brown       output_file: output/ex109.out
160c4762a1bSJed Brown 
161c20d7725SJed Brown    test:
162c20d7725SJed Brown       suffix: 3
163c20d7725SJed Brown       nsize: 2
164c20d7725SJed Brown       args: -matmattransmult_mpidense_mpidense_via cyclic
165c20d7725SJed Brown       output_file: output/ex109.out
166c20d7725SJed Brown 
1673fbca975SStefano Zampini    test:
1683fbca975SStefano Zampini       suffix: 4
1693fbca975SStefano Zampini       args: -test_userAPI
1703fbca975SStefano Zampini       output_file: output/ex109.out
1713fbca975SStefano Zampini 
1723fbca975SStefano Zampini    test:
1733fbca975SStefano Zampini       suffix: 5
1743fbca975SStefano Zampini       nsize: 3
1753fbca975SStefano Zampini       args: -test_userAPI
1763fbca975SStefano Zampini       output_file: output/ex109.out
1773fbca975SStefano Zampini 
178c4762a1bSJed Brown TEST*/
179