xref: /petsc/src/mat/tests/ex109.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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   PetscErrorCode ierr;
11c4762a1bSJed Brown   PetscRandom    r;
123fbca975SStefano Zampini   PetscBool      equal=PETSC_FALSE,flg;
13c20d7725SJed Brown   PetscReal      fill = 1.0,norm;
14c4762a1bSJed Brown   PetscMPIInt    size;
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL));
20c4762a1bSJed Brown 
215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&r));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(r));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   /* Create a aij matrix A */
25c4762a1bSJed Brown   M    = N = m*n;
265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATAIJ));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(A,5,NULL));
32c4762a1bSJed Brown 
335f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
375f80ce2aSJacob Faibussowitsch     if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
385f80ce2aSJacob Faibussowitsch     if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
395f80ce2aSJacob Faibussowitsch     if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
405f80ce2aSJacob Faibussowitsch     if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
415f80ce2aSJacob Faibussowitsch     v = 4.0; CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
42c4762a1bSJed Brown   }
435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   /* Create a dense matrix B */
475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(A,&am,&an));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(B,an,PETSC_DECIDE,PETSC_DECIDE,M));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(B,MATDENSE));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqDenseSetPreallocation(B,NULL));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIDenseSetPreallocation(B,NULL));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(B));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(B,r));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&r));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
58c4762a1bSJed Brown 
59c20d7725SJed Brown   /* Test reuse of user-provided dense C (unassembled) -- not recommended usage */
605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(C,MATDENSE));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,am,PETSC_DECIDE,PETSC_DECIDE,N));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(C));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroEntries(C));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(C,NORM_INFINITY,&norm));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
69c20d7725SJed Brown 
70c4762a1bSJed Brown   /* Test C = A*B (aij*dense) */
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,B,MAT_INITIAL_MATRIX,fill,&C));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C));
73c4762a1bSJed Brown 
74c20d7725SJed Brown   /* Test developer API */
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatProductCreate(A,B,NULL,&D));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatProductSetType(D,MATPRODUCT_AB));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatProductSetAlgorithm(D,"default"));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatProductSetFill(D,fill));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatProductSetFromOptions(D));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(MatProductSymbolic(D));
81c4762a1bSJed Brown   for (i=0; i<2; i++) {
825f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductNumeric(D));
83c4762a1bSJed Brown   }
845f80ce2aSJacob Faibussowitsch   CHKERRQ(MatEqual(C,D,&equal));
85*28b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != D");
865f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&D));
87c4762a1bSJed Brown 
88c20d7725SJed Brown   /* Test D = AT*B (transpose(aij)*dense) */
895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateTranspose(A,&AT));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(AT,B,MAT_INITIAL_MATRIX,fill,&D));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual(AT,B,D,10,&equal));
92*28b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != AT*B (transpose(aij)*dense)");
935f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&D));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&AT));
95c20d7725SJed Brown 
96c4762a1bSJed Brown   /* Test D = C*A (dense*aij) */
975f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&D));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(C,A,MAT_REUSE_MATRIX,fill,&D));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual(C,A,D,10,&equal));
100*28b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != C*A (dense*aij)");
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&D));
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   /* Test D = A*C (aij*dense) */
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,C,MAT_INITIAL_MATRIX,fill,&D));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,C,MAT_REUSE_MATRIX,fill,&D));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual(A,C,D,10,&equal));
107*28b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != A*C (aij*dense)");
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&D));
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   /* Test D = B*C (dense*dense) */
1115f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
112c4762a1bSJed Brown   if (size == 1) {
1135f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(B,C,MAT_INITIAL_MATRIX,fill,&D));
1145f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(B,C,MAT_REUSE_MATRIX,fill,&D));
1155f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMultEqual(B,C,D,10,&equal));
116*28b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != B*C (dense*dense)");
1175f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&D));
118c4762a1bSJed Brown   }
119c4762a1bSJed Brown 
120c20d7725SJed Brown   /* Test D = B*C^T (dense*dense) */
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatTransposeMult(B,C,MAT_INITIAL_MATRIX,fill,&D));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatTransposeMult(B,C,MAT_REUSE_MATRIX,fill,&D));
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatTransposeMultEqual(B,C,D,10,&equal));
124*28b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != B*C^T (dense*dense)");
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&D));
126c20d7725SJed Brown 
1273fbca975SStefano Zampini   /* Test MatProductCreateWithMat() and reuse C and B for B = A*C */
1283fbca975SStefano Zampini   flg = PETSC_FALSE;
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-test_userAPI",&flg));
1303fbca975SStefano Zampini   if (flg) {
1313fbca975SStefano Zampini     /* user driver */
1325f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(A,C,MAT_REUSE_MATRIX,fill,&B));
1333fbca975SStefano Zampini   } else {
1343fbca975SStefano Zampini     /* clear internal data structures related with previous products to avoid circular references */
1355f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductClear(A));
1365f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductClear(B));
1375f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductClear(C));
1385f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductCreateWithMat(A,C,NULL,B));
1395f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductSetType(B,MATPRODUCT_AB));
1405f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductSetFromOptions(B));
1415f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductSymbolic(B));
1425f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductNumeric(B));
1433fbca975SStefano Zampini   }
1443fbca975SStefano Zampini 
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
148c4762a1bSJed Brown   ierr = PetscFinalize();
149c4762a1bSJed Brown   return ierr;
150c4762a1bSJed Brown }
151c4762a1bSJed Brown 
152c4762a1bSJed Brown /*TEST
153c4762a1bSJed Brown 
154c4762a1bSJed Brown    test:
155c4762a1bSJed Brown       args: -M 10 -N 10
156c4762a1bSJed Brown       output_file: output/ex109.out
157c4762a1bSJed Brown 
158c4762a1bSJed Brown    test:
159c4762a1bSJed Brown       suffix: 2
160c4762a1bSJed Brown       nsize: 3
161c4762a1bSJed Brown       output_file: output/ex109.out
162c4762a1bSJed Brown 
163c20d7725SJed Brown    test:
164c20d7725SJed Brown       suffix: 3
165c20d7725SJed Brown       nsize: 2
166c20d7725SJed Brown       args: -matmattransmult_mpidense_mpidense_via cyclic
167c20d7725SJed Brown       output_file: output/ex109.out
168c20d7725SJed Brown 
1693fbca975SStefano Zampini    test:
1703fbca975SStefano Zampini       suffix: 4
1713fbca975SStefano Zampini       args: -test_userAPI
1723fbca975SStefano Zampini       output_file: output/ex109.out
1733fbca975SStefano Zampini 
1743fbca975SStefano Zampini    test:
1753fbca975SStefano Zampini       suffix: 5
1763fbca975SStefano Zampini       nsize: 3
1773fbca975SStefano Zampini       args: -test_userAPI
1783fbca975SStefano Zampini       output_file: output/ex109.out
1793fbca975SStefano Zampini 
180c4762a1bSJed Brown TEST*/
181