xref: /petsc/src/mat/tests/ex104.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown static char help[] = "Test MatMatMult(), MatTranspose(), MatTransposeMatMult() for Dense and Elemental matrices.\n\n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown  Example:
4c4762a1bSJed Brown    mpiexec -n <np> ./ex104 -mat_type elemental
5c4762a1bSJed Brown */
6c4762a1bSJed Brown 
7c4762a1bSJed Brown #include <petscmat.h>
8c4762a1bSJed Brown 
9c4762a1bSJed Brown int main(int argc,char **argv)
10c4762a1bSJed Brown {
11c4762a1bSJed Brown   Mat            A,B,C,D;
12c4762a1bSJed Brown   PetscInt       i,M=10,N=5,j,nrows,ncols,am,an,rstart,rend;
13c4762a1bSJed Brown   PetscRandom    r;
14c20d7725SJed Brown   PetscBool      equal,Aiselemental;
15c4762a1bSJed Brown   PetscReal      fill = 1.0;
16c4762a1bSJed Brown   IS             isrows,iscols;
17c4762a1bSJed Brown   const PetscInt *rows,*cols;
18c4762a1bSJed Brown   PetscScalar    *v,rval;
19c4762a1bSJed Brown #if defined(PETSC_HAVE_ELEMENTAL)
20c4762a1bSJed Brown   PetscBool      Test_MatMatMult=PETSC_TRUE;
21c4762a1bSJed Brown #else
22c4762a1bSJed Brown   PetscBool      Test_MatMatMult=PETSC_FALSE;
23c4762a1bSJed Brown #endif
24c4762a1bSJed Brown   PetscMPIInt    size;
25c4762a1bSJed Brown 
26*327415f7SBarry Smith   PetscFunctionBeginUser;
279566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
29c4762a1bSJed Brown 
309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL));
329566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
339566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
349566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATDENSE));
359566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
369566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
379566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&r));
389566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(r));
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   /* Set local matrix entries */
419566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipIS(A,&isrows,&iscols));
429566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isrows,&nrows));
439566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrows,&rows));
449566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(iscols,&ncols));
459566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(iscols,&cols));
469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrows*ncols,&v));
47c4762a1bSJed Brown   for (i=0; i<nrows; i++) {
48c4762a1bSJed Brown     for (j=0; j<ncols; j++) {
499566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(r,&rval));
50c4762a1bSJed Brown       v[i*ncols+j] = rval;
51c4762a1bSJed Brown     }
52c4762a1bSJed Brown   }
539566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES));
549566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
559566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
569566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrows,&rows));
579566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscols,&cols));
589566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrows));
599566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscols));
609566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&r));
61c4762a1bSJed Brown 
629566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A,MATELEMENTAL,&Aiselemental));
63c20d7725SJed Brown 
64c20d7725SJed Brown   /* Test MatCreateTranspose() and MatTranspose() */
659566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(A,&C));
669566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX,&B)); /* B = A^T */
679566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(C,B,10,&equal));
6828b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"A^T*x != (x^T*A)^T");
699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
70c20d7725SJed Brown 
719566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&B));
72c20d7725SJed Brown   if (!Aiselemental) {
739566063dSJacob Faibussowitsch     PetscCall(MatTranspose(B,MAT_INPLACE_MATRIX,&B));
749566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(C,B,10,&equal));
7528b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"C*x != B*x");
76c20d7725SJed Brown   }
779566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
78c20d7725SJed Brown 
79c20d7725SJed Brown   /* Test B = C*A for matrix type transpose and seqdense */
80c20d7725SJed Brown   if (size == 1 && !Aiselemental) {
819566063dSJacob Faibussowitsch     PetscCall(MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&B));
829566063dSJacob Faibussowitsch     PetscCall(MatMatMultEqual(C,A,B,10,&equal));
8328b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B != C*A for matrix type transpose and seqdense");
849566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
85c20d7725SJed Brown   }
869566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /* Test MatMatMult() */
89c4762a1bSJed Brown   if (Test_MatMatMult) {
90c4762a1bSJed Brown #if !defined(PETSC_HAVE_ELEMENTAL)
91be096a46SBarry Smith     PetscCheck(size == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"This test requires ELEMENTAL");
92c4762a1bSJed Brown #endif
939566063dSJacob Faibussowitsch     PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX,&B)); /* B = A^T */
949566063dSJacob Faibussowitsch     PetscCall(MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C)); /* C = B*A = A^T*A */
959566063dSJacob Faibussowitsch     PetscCall(MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C));
969566063dSJacob Faibussowitsch     PetscCall(MatMatMultEqual(B,A,C,10,&equal));
9728b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B*A*x != C*x");
98c4762a1bSJed Brown 
99c4762a1bSJed Brown     /* Test MatDuplicate for matrix product */
1009566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(C,MAT_COPY_VALUES,&D));
101c20d7725SJed Brown 
1029566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&D));
1039566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
1049566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
105c4762a1bSJed Brown   }
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /* Test MatTransposeMatMult() */
108c20d7725SJed Brown   if (!Aiselemental) {
1099566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&D)); /* D = A^T*A */
1109566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMult(A,A,MAT_REUSE_MATRIX,fill,&D));
1119566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMultEqual(A,A,D,10,&equal));
11228b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*A*x");
113c4762a1bSJed Brown 
114c4762a1bSJed Brown     /* Test MatDuplicate for matrix product */
1159566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(D,MAT_COPY_VALUES,&C));
1169566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
1179566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&D));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown     /* Test D*x = A^T*C*A*x, where C is in AIJ format */
1209566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A,&am,&an));
1219566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
122c4762a1bSJed Brown     if (size == 1) {
1239566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,am,am));
124c4762a1bSJed Brown     } else {
1259566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(C,am,am,PETSC_DECIDE,PETSC_DECIDE));
126c4762a1bSJed Brown     }
1279566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(C));
1289566063dSJacob Faibussowitsch     PetscCall(MatSetUp(C));
1299566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(C,&rstart,&rend));
130c4762a1bSJed Brown     v[0] = 1.0;
131c4762a1bSJed Brown     for (i=rstart; i<rend; i++) {
1329566063dSJacob Faibussowitsch       PetscCall(MatSetValues(C,1,&i,1,&i,v,INSERT_VALUES));
133c4762a1bSJed Brown     }
1349566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
1359566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown     /* B = C*A, D = A^T*B */
1389566063dSJacob Faibussowitsch     PetscCall(MatMatMult(C,A,MAT_INITIAL_MATRIX,1.0,&B));
1399566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMult(A,B,MAT_INITIAL_MATRIX,fill,&D));
1409566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMultEqual(A,B,D,10,&equal));
14128b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*B*x");
142c4762a1bSJed Brown 
1439566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&D));
1449566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
1459566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
146c4762a1bSJed Brown   }
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* Test MatMatTransposeMult() */
149c20d7725SJed Brown   if (!Aiselemental) {
150c4762a1bSJed Brown     PetscReal diff, scale;
151c4762a1bSJed Brown     PetscInt  am, an, aM, aN;
152c4762a1bSJed Brown 
1539566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A, &am, &an));
1549566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A, &aM, &aN));
1559566063dSJacob Faibussowitsch     PetscCall(MatCreateDense(PetscObjectComm((PetscObject)A),PETSC_DECIDE, an, aM + 10, aN, NULL, &B));
1569566063dSJacob Faibussowitsch     PetscCall(MatSetRandom(B, NULL));
1579566063dSJacob Faibussowitsch     PetscCall(MatMatTransposeMult(A,B,MAT_INITIAL_MATRIX,fill,&D)); /* D = A*A^T */
158c4762a1bSJed Brown 
159c4762a1bSJed Brown     /* Test MatDuplicate for matrix product */
1609566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(D,MAT_COPY_VALUES,&C));
161c4762a1bSJed Brown 
1629566063dSJacob Faibussowitsch     PetscCall(MatMatTransposeMult(A,B,MAT_REUSE_MATRIX,fill,&D));
1639566063dSJacob Faibussowitsch     PetscCall(MatAXPY(C, -1., D, SAME_NONZERO_PATTERN));
164c4762a1bSJed Brown 
1659566063dSJacob Faibussowitsch     PetscCall(MatNorm(C, NORM_FROBENIUS, &diff));
1669566063dSJacob Faibussowitsch     PetscCall(MatNorm(D, NORM_FROBENIUS, &scale));
167e00437b9SBarry Smith     PetscCheck(diff <= PETSC_SMALL * scale,PetscObjectComm((PetscObject)D), PETSC_ERR_PLIB, "MatMatTransposeMult() differs between MAT_INITIAL_MATRIX and MAT_REUSE_MATRIX");
1689566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
169c4762a1bSJed Brown 
1709566063dSJacob Faibussowitsch     PetscCall(MatMatTransposeMultEqual(A,B,D,10,&equal));
17128b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*A*x");
1729566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&D));
1739566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   }
176c4762a1bSJed Brown 
1779566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1789566063dSJacob Faibussowitsch   PetscCall(PetscFree(v));
1799566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
180b122ec5aSJacob Faibussowitsch   return 0;
181c4762a1bSJed Brown }
182c4762a1bSJed Brown 
183c4762a1bSJed Brown /*TEST
184c4762a1bSJed Brown 
185c4762a1bSJed Brown     test:
186c4762a1bSJed Brown       output_file: output/ex104.out
187c4762a1bSJed Brown 
188c4762a1bSJed Brown     test:
189c4762a1bSJed Brown       suffix: 2
190c4762a1bSJed Brown       nsize: 2
191c4762a1bSJed Brown       output_file: output/ex104.out
192c4762a1bSJed Brown 
193c4762a1bSJed Brown     test:
194c4762a1bSJed Brown       suffix: 3
195c4762a1bSJed Brown       nsize: 4
196c4762a1bSJed Brown       output_file: output/ex104.out
197c4762a1bSJed Brown       args: -M 23 -N 31
198c4762a1bSJed Brown 
199c4762a1bSJed Brown     test:
200c4762a1bSJed Brown       suffix: 4
201c4762a1bSJed Brown       nsize: 4
202c4762a1bSJed Brown       output_file: output/ex104.out
203c4762a1bSJed Brown       args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via cyclic
204c4762a1bSJed Brown 
205c4762a1bSJed Brown     test:
206c4762a1bSJed Brown       suffix: 5
207c4762a1bSJed Brown       nsize: 4
208c4762a1bSJed Brown       output_file: output/ex104.out
209c4762a1bSJed Brown       args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via allgatherv
210c4762a1bSJed Brown 
211c20d7725SJed Brown     test:
212c20d7725SJed Brown       suffix: 6
213c20d7725SJed Brown       args: -mat_type elemental
214c20d7725SJed Brown       requires: elemental
215c20d7725SJed Brown       output_file: output/ex104.out
216c20d7725SJed Brown 
217c20d7725SJed Brown     test:
218c20d7725SJed Brown       suffix: 7
219c20d7725SJed Brown       nsize: 2
220c20d7725SJed Brown       args: -mat_type elemental
221c20d7725SJed Brown       requires: elemental
222c20d7725SJed Brown       output_file: output/ex104.out
223c20d7725SJed Brown 
224c4762a1bSJed Brown TEST*/
225