xref: /petsc/src/mat/tests/ex104.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
275f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
28c4762a1bSJed Brown 
295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATDENSE));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&r));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(r));
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   /* Set local matrix entries */
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipIS(A,&isrows,&iscols));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(isrows,&nrows));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(isrows,&rows));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(iscols,&ncols));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(iscols,&cols));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrows*ncols,&v));
46c4762a1bSJed Brown   for (i=0; i<nrows; i++) {
47c4762a1bSJed Brown     for (j=0; j<ncols; j++) {
485f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(r,&rval));
49c4762a1bSJed Brown       v[i*ncols+j] = rval;
50c4762a1bSJed Brown     }
51c4762a1bSJed Brown   }
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(isrows,&rows));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(iscols,&cols));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isrows));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iscols));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&r));
60c4762a1bSJed Brown 
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATELEMENTAL,&Aiselemental));
62c20d7725SJed Brown 
63c20d7725SJed Brown   /* Test MatCreateTranspose() and MatTranspose() */
645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateTranspose(A,&C));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX,&B)); /* B = A^T */
665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultEqual(C,B,10,&equal));
6728b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"A^T*x != (x^T*A)^T");
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
69c20d7725SJed Brown 
705f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&B));
71c20d7725SJed Brown   if (!Aiselemental) {
725f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTranspose(B,MAT_INPLACE_MATRIX,&B));
735f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultEqual(C,B,10,&equal));
7428b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"C*x != B*x");
75c20d7725SJed Brown   }
765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
77c20d7725SJed Brown 
78c20d7725SJed Brown   /* Test B = C*A for matrix type transpose and seqdense */
79c20d7725SJed Brown   if (size == 1 && !Aiselemental) {
805f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&B));
815f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMultEqual(C,A,B,10,&equal));
8228b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B != C*A for matrix type transpose and seqdense");
835f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&B));
84c20d7725SJed Brown   }
855f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /* Test MatMatMult() */
88c4762a1bSJed Brown   if (Test_MatMatMult) {
89c4762a1bSJed Brown #if !defined(PETSC_HAVE_ELEMENTAL)
902c71b3e2SJacob Faibussowitsch     PetscCheckFalse(size > 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"This test requires ELEMENTAL");
91c4762a1bSJed Brown #endif
925f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX,&B)); /* B = A^T */
935f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C)); /* C = B*A = A^T*A */
945f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C));
955f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMultEqual(B,A,C,10,&equal));
9628b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B*A*x != C*x");
97c4762a1bSJed Brown 
98c4762a1bSJed Brown     /* Test MatDuplicate for matrix product */
995f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&D));
100c20d7725SJed Brown 
1015f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&D));
1025f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&C));
1035f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&B));
104c4762a1bSJed Brown   }
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* Test MatTransposeMatMult() */
107c20d7725SJed Brown   if (!Aiselemental) {
1085f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&D)); /* D = A^T*A */
1095f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTransposeMatMult(A,A,MAT_REUSE_MATRIX,fill,&D));
1105f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTransposeMatMultEqual(A,A,D,10,&equal));
11128b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*A*x");
112c4762a1bSJed Brown 
113c4762a1bSJed Brown     /* Test MatDuplicate for matrix product */
1145f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDuplicate(D,MAT_COPY_VALUES,&C));
1155f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&C));
1165f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&D));
117c4762a1bSJed Brown 
118c4762a1bSJed Brown     /* Test D*x = A^T*C*A*x, where C is in AIJ format */
1195f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetLocalSize(A,&am,&an));
1205f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
121c4762a1bSJed Brown     if (size == 1) {
1225f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,am,am));
123c4762a1bSJed Brown     } else {
1245f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetSizes(C,am,am,PETSC_DECIDE,PETSC_DECIDE));
125c4762a1bSJed Brown     }
1265f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetFromOptions(C));
1275f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetUp(C));
1285f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetOwnershipRange(C,&rstart,&rend));
129c4762a1bSJed Brown     v[0] = 1.0;
130c4762a1bSJed Brown     for (i=rstart; i<rend; i++) {
1315f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(C,1,&i,1,&i,v,INSERT_VALUES));
132c4762a1bSJed Brown     }
1335f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
1345f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown     /* B = C*A, D = A^T*B */
1375f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(C,A,MAT_INITIAL_MATRIX,1.0,&B));
1385f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTransposeMatMult(A,B,MAT_INITIAL_MATRIX,fill,&D));
1395f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTransposeMatMultEqual(A,B,D,10,&equal));
14028b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*B*x");
141c4762a1bSJed Brown 
1425f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&D));
1435f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&C));
1445f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&B));
145c4762a1bSJed Brown   }
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /* Test MatMatTransposeMult() */
148c20d7725SJed Brown   if (!Aiselemental) {
149c4762a1bSJed Brown     PetscReal diff, scale;
150c4762a1bSJed Brown     PetscInt  am, an, aM, aN;
151c4762a1bSJed Brown 
1525f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetLocalSize(A, &am, &an));
1535f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetSize(A, &aM, &aN));
1545f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateDense(PetscObjectComm((PetscObject)A),PETSC_DECIDE, an, aM + 10, aN, NULL, &B));
1555f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetRandom(B, NULL));
1565f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1575f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
1585f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatTransposeMult(A,B,MAT_INITIAL_MATRIX,fill,&D)); /* D = A*A^T */
159c4762a1bSJed Brown 
160c4762a1bSJed Brown     /* Test MatDuplicate for matrix product */
1615f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDuplicate(D,MAT_COPY_VALUES,&C));
162c4762a1bSJed Brown 
1635f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatTransposeMult(A,B,MAT_REUSE_MATRIX,fill,&D));
1645f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAXPY(C, -1., D, SAME_NONZERO_PATTERN));
165c4762a1bSJed Brown 
1665f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNorm(C, NORM_FROBENIUS, &diff));
1675f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNorm(D, NORM_FROBENIUS, &scale));
1682c71b3e2SJacob Faibussowitsch     PetscCheckFalse(diff > PETSC_SMALL * scale,PetscObjectComm((PetscObject)D), PETSC_ERR_PLIB, "MatMatTransposeMult() differs between MAT_INITIAL_MATRIX and MAT_REUSE_MATRIX");
1695f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&C));
170c4762a1bSJed Brown 
1715f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatTransposeMultEqual(A,B,D,10,&equal));
17228b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*A*x");
1735f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&D));
1745f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&B));
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   }
177c4762a1bSJed Brown 
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(v));
180*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
181*b122ec5aSJacob Faibussowitsch   return 0;
182c4762a1bSJed Brown }
183c4762a1bSJed Brown 
184c4762a1bSJed Brown /*TEST
185c4762a1bSJed Brown 
186c4762a1bSJed Brown     test:
187c4762a1bSJed Brown       output_file: output/ex104.out
188c4762a1bSJed Brown 
189c4762a1bSJed Brown     test:
190c4762a1bSJed Brown       suffix: 2
191c4762a1bSJed Brown       nsize: 2
192c4762a1bSJed Brown       output_file: output/ex104.out
193c4762a1bSJed Brown 
194c4762a1bSJed Brown     test:
195c4762a1bSJed Brown       suffix: 3
196c4762a1bSJed Brown       nsize: 4
197c4762a1bSJed Brown       output_file: output/ex104.out
198c4762a1bSJed Brown       args: -M 23 -N 31
199c4762a1bSJed Brown 
200c4762a1bSJed Brown     test:
201c4762a1bSJed Brown       suffix: 4
202c4762a1bSJed Brown       nsize: 4
203c4762a1bSJed Brown       output_file: output/ex104.out
204c4762a1bSJed Brown       args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via cyclic
205c4762a1bSJed Brown 
206c4762a1bSJed Brown     test:
207c4762a1bSJed Brown       suffix: 5
208c4762a1bSJed Brown       nsize: 4
209c4762a1bSJed Brown       output_file: output/ex104.out
210c4762a1bSJed Brown       args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via allgatherv
211c4762a1bSJed Brown 
212c20d7725SJed Brown     test:
213c20d7725SJed Brown       suffix: 6
214c20d7725SJed Brown       args: -mat_type elemental
215c20d7725SJed Brown       requires: elemental
216c20d7725SJed Brown       output_file: output/ex104.out
217c20d7725SJed Brown 
218c20d7725SJed Brown     test:
219c20d7725SJed Brown       suffix: 7
220c20d7725SJed Brown       nsize: 2
221c20d7725SJed Brown       args: -mat_type elemental
222c20d7725SJed Brown       requires: elemental
223c20d7725SJed Brown       output_file: output/ex104.out
224c20d7725SJed Brown 
225c4762a1bSJed Brown TEST*/
226