xref: /petsc/src/mat/tests/ex109.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL));
19c4762a1bSJed Brown 
205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&r));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(r));
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   /* Create a aij matrix A */
24c4762a1bSJed Brown   M    = N = m*n;
255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATAIJ));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(A,5,NULL));
31c4762a1bSJed Brown 
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(A,&Istart,&Iend));
33c4762a1bSJed Brown   am   = Iend - Istart;
34c4762a1bSJed Brown   for (Ii=Istart; Ii<Iend; Ii++) {
35c4762a1bSJed Brown     v = -1.0; i = Ii/n; j = Ii - i*n;
365f80ce2aSJacob Faibussowitsch     if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
375f80ce2aSJacob Faibussowitsch     if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
385f80ce2aSJacob Faibussowitsch     if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
395f80ce2aSJacob Faibussowitsch     if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
405f80ce2aSJacob Faibussowitsch     v = 4.0; CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
41c4762a1bSJed Brown   }
425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* Create a dense matrix B */
465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(A,&am,&an));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(B,an,PETSC_DECIDE,PETSC_DECIDE,M));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(B,MATDENSE));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqDenseSetPreallocation(B,NULL));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIDenseSetPreallocation(B,NULL));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(B));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(B,r));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&r));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
57c4762a1bSJed Brown 
58c20d7725SJed Brown   /* Test reuse of user-provided dense C (unassembled) -- not recommended usage */
595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(C,MATDENSE));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,am,PETSC_DECIDE,PETSC_DECIDE,N));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(C));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroEntries(C));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(C,NORM_INFINITY,&norm));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
68c20d7725SJed Brown 
69c4762a1bSJed Brown   /* Test C = A*B (aij*dense) */
705f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,B,MAT_INITIAL_MATRIX,fill,&C));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C));
72c4762a1bSJed Brown 
73c20d7725SJed Brown   /* Test developer API */
745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatProductCreate(A,B,NULL,&D));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatProductSetType(D,MATPRODUCT_AB));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatProductSetAlgorithm(D,"default"));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatProductSetFill(D,fill));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatProductSetFromOptions(D));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatProductSymbolic(D));
80c4762a1bSJed Brown   for (i=0; i<2; i++) {
815f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductNumeric(D));
82c4762a1bSJed Brown   }
835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatEqual(C,D,&equal));
8428b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != D");
855f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&D));
86c4762a1bSJed Brown 
87c20d7725SJed Brown   /* Test D = AT*B (transpose(aij)*dense) */
885f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateTranspose(A,&AT));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(AT,B,MAT_INITIAL_MATRIX,fill,&D));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual(AT,B,D,10,&equal));
9128b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != AT*B (transpose(aij)*dense)");
925f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&D));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&AT));
94c20d7725SJed Brown 
95c4762a1bSJed Brown   /* Test D = C*A (dense*aij) */
965f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&D));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(C,A,MAT_REUSE_MATRIX,fill,&D));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual(C,A,D,10,&equal));
9928b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != C*A (dense*aij)");
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&D));
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   /* Test D = A*C (aij*dense) */
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,C,MAT_INITIAL_MATRIX,fill,&D));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,C,MAT_REUSE_MATRIX,fill,&D));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual(A,C,D,10,&equal));
10628b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != A*C (aij*dense)");
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&D));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /* Test D = B*C (dense*dense) */
1105f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
111c4762a1bSJed Brown   if (size == 1) {
1125f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(B,C,MAT_INITIAL_MATRIX,fill,&D));
1135f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(B,C,MAT_REUSE_MATRIX,fill,&D));
1145f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMultEqual(B,C,D,10,&equal));
11528b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != B*C (dense*dense)");
1165f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&D));
117c4762a1bSJed Brown   }
118c4762a1bSJed Brown 
119c20d7725SJed Brown   /* Test D = B*C^T (dense*dense) */
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatTransposeMult(B,C,MAT_INITIAL_MATRIX,fill,&D));
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatTransposeMult(B,C,MAT_REUSE_MATRIX,fill,&D));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatTransposeMultEqual(B,C,D,10,&equal));
12328b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != B*C^T (dense*dense)");
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&D));
125c20d7725SJed Brown 
1263fbca975SStefano Zampini   /* Test MatProductCreateWithMat() and reuse C and B for B = A*C */
1273fbca975SStefano Zampini   flg = PETSC_FALSE;
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-test_userAPI",&flg));
1293fbca975SStefano Zampini   if (flg) {
1303fbca975SStefano Zampini     /* user driver */
1315f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(A,C,MAT_REUSE_MATRIX,fill,&B));
1323fbca975SStefano Zampini   } else {
1333fbca975SStefano Zampini     /* clear internal data structures related with previous products to avoid circular references */
1345f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductClear(A));
1355f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductClear(B));
1365f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductClear(C));
1375f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductCreateWithMat(A,C,NULL,B));
1385f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductSetType(B,MATPRODUCT_AB));
1395f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductSetFromOptions(B));
1405f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductSymbolic(B));
1415f80ce2aSJacob Faibussowitsch     CHKERRQ(MatProductNumeric(B));
1423fbca975SStefano Zampini   }
1433fbca975SStefano Zampini 
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
147*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
148*b122ec5aSJacob Faibussowitsch   return 0;
149c4762a1bSJed Brown }
150c4762a1bSJed Brown 
151c4762a1bSJed Brown /*TEST
152c4762a1bSJed Brown 
153c4762a1bSJed Brown    test:
154c4762a1bSJed Brown       args: -M 10 -N 10
155c4762a1bSJed Brown       output_file: output/ex109.out
156c4762a1bSJed Brown 
157c4762a1bSJed Brown    test:
158c4762a1bSJed Brown       suffix: 2
159c4762a1bSJed Brown       nsize: 3
160c4762a1bSJed Brown       output_file: output/ex109.out
161c4762a1bSJed Brown 
162c20d7725SJed Brown    test:
163c20d7725SJed Brown       suffix: 3
164c20d7725SJed Brown       nsize: 2
165c20d7725SJed Brown       args: -matmattransmult_mpidense_mpidense_via cyclic
166c20d7725SJed Brown       output_file: output/ex109.out
167c20d7725SJed Brown 
1683fbca975SStefano Zampini    test:
1693fbca975SStefano Zampini       suffix: 4
1703fbca975SStefano Zampini       args: -test_userAPI
1713fbca975SStefano Zampini       output_file: output/ex109.out
1723fbca975SStefano Zampini 
1733fbca975SStefano Zampini    test:
1743fbca975SStefano Zampini       suffix: 5
1753fbca975SStefano Zampini       nsize: 3
1763fbca975SStefano Zampini       args: -test_userAPI
1773fbca975SStefano Zampini       output_file: output/ex109.out
1783fbca975SStefano Zampini 
179c4762a1bSJed Brown TEST*/
180