xref: /petsc/src/mat/tests/ex104.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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   PetscErrorCode ierr;
14c4762a1bSJed Brown   PetscRandom    r;
15c20d7725SJed Brown   PetscBool      equal,Aiselemental;
16c4762a1bSJed Brown   PetscReal      fill = 1.0;
17c4762a1bSJed Brown   IS             isrows,iscols;
18c4762a1bSJed Brown   const PetscInt *rows,*cols;
19c4762a1bSJed Brown   PetscScalar    *v,rval;
20c4762a1bSJed Brown #if defined(PETSC_HAVE_ELEMENTAL)
21c4762a1bSJed Brown   PetscBool      Test_MatMatMult=PETSC_TRUE;
22c4762a1bSJed Brown #else
23c4762a1bSJed Brown   PetscBool      Test_MatMatMult=PETSC_FALSE;
24c4762a1bSJed Brown #endif
25c4762a1bSJed Brown   PetscMPIInt    size;
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
285f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
29c4762a1bSJed Brown 
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATDENSE));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&r));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(r));
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   /* Set local matrix entries */
415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipIS(A,&isrows,&iscols));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(isrows,&nrows));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(isrows,&rows));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(iscols,&ncols));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(iscols,&cols));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrows*ncols,&v));
47c4762a1bSJed Brown   for (i=0; i<nrows; i++) {
48c4762a1bSJed Brown     for (j=0; j<ncols; j++) {
495f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(r,&rval));
50c4762a1bSJed Brown       v[i*ncols+j] = rval;
51c4762a1bSJed Brown     }
52c4762a1bSJed Brown   }
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(isrows,&rows));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(iscols,&cols));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isrows));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iscols));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&r));
61c4762a1bSJed Brown 
625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATELEMENTAL,&Aiselemental));
63c20d7725SJed Brown 
64c20d7725SJed Brown   /* Test MatCreateTranspose() and MatTranspose() */
655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateTranspose(A,&C));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX,&B)); /* B = A^T */
675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultEqual(C,B,10,&equal));
68*28b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"A^T*x != (x^T*A)^T");
695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
70c20d7725SJed Brown 
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&B));
72c20d7725SJed Brown   if (!Aiselemental) {
735f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTranspose(B,MAT_INPLACE_MATRIX,&B));
745f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultEqual(C,B,10,&equal));
75*28b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"C*x != B*x");
76c20d7725SJed Brown   }
775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
78c20d7725SJed Brown 
79c20d7725SJed Brown   /* Test B = C*A for matrix type transpose and seqdense */
80c20d7725SJed Brown   if (size == 1 && !Aiselemental) {
815f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&B));
825f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMultEqual(C,A,B,10,&equal));
83*28b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B != C*A for matrix type transpose and seqdense");
845f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&B));
85c20d7725SJed Brown   }
865f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /* Test MatMatMult() */
89c4762a1bSJed Brown   if (Test_MatMatMult) {
90c4762a1bSJed Brown #if !defined(PETSC_HAVE_ELEMENTAL)
912c71b3e2SJacob Faibussowitsch     PetscCheckFalse(size > 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"This test requires ELEMENTAL");
92c4762a1bSJed Brown #endif
935f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX,&B)); /* B = A^T */
945f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C)); /* C = B*A = A^T*A */
955f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C));
965f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMultEqual(B,A,C,10,&equal));
97*28b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B*A*x != C*x");
98c4762a1bSJed Brown 
99c4762a1bSJed Brown     /* Test MatDuplicate for matrix product */
1005f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&D));
101c20d7725SJed Brown 
1025f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&D));
1035f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&C));
1045f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&B));
105c4762a1bSJed Brown   }
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /* Test MatTransposeMatMult() */
108c20d7725SJed Brown   if (!Aiselemental) {
1095f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&D)); /* D = A^T*A */
1105f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTransposeMatMult(A,A,MAT_REUSE_MATRIX,fill,&D));
1115f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTransposeMatMultEqual(A,A,D,10,&equal));
112*28b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*A*x");
113c4762a1bSJed Brown 
114c4762a1bSJed Brown     /* Test MatDuplicate for matrix product */
1155f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDuplicate(D,MAT_COPY_VALUES,&C));
1165f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&C));
1175f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&D));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown     /* Test D*x = A^T*C*A*x, where C is in AIJ format */
1205f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetLocalSize(A,&am,&an));
1215f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
122c4762a1bSJed Brown     if (size == 1) {
1235f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,am,am));
124c4762a1bSJed Brown     } else {
1255f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetSizes(C,am,am,PETSC_DECIDE,PETSC_DECIDE));
126c4762a1bSJed Brown     }
1275f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetFromOptions(C));
1285f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetUp(C));
1295f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetOwnershipRange(C,&rstart,&rend));
130c4762a1bSJed Brown     v[0] = 1.0;
131c4762a1bSJed Brown     for (i=rstart; i<rend; i++) {
1325f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(C,1,&i,1,&i,v,INSERT_VALUES));
133c4762a1bSJed Brown     }
1345f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
1355f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown     /* B = C*A, D = A^T*B */
1385f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(C,A,MAT_INITIAL_MATRIX,1.0,&B));
1395f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTransposeMatMult(A,B,MAT_INITIAL_MATRIX,fill,&D));
1405f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTransposeMatMultEqual(A,B,D,10,&equal));
141*28b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*B*x");
142c4762a1bSJed Brown 
1435f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&D));
1445f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&C));
1455f80ce2aSJacob Faibussowitsch     CHKERRQ(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 
1535f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetLocalSize(A, &am, &an));
1545f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetSize(A, &aM, &aN));
1555f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateDense(PetscObjectComm((PetscObject)A),PETSC_DECIDE, an, aM + 10, aN, NULL, &B));
1565f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetRandom(B, NULL));
1575f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1585f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
1595f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatTransposeMult(A,B,MAT_INITIAL_MATRIX,fill,&D)); /* D = A*A^T */
160c4762a1bSJed Brown 
161c4762a1bSJed Brown     /* Test MatDuplicate for matrix product */
1625f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDuplicate(D,MAT_COPY_VALUES,&C));
163c4762a1bSJed Brown 
1645f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatTransposeMult(A,B,MAT_REUSE_MATRIX,fill,&D));
1655f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAXPY(C, -1., D, SAME_NONZERO_PATTERN));
166c4762a1bSJed Brown 
1675f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNorm(C, NORM_FROBENIUS, &diff));
1685f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNorm(D, NORM_FROBENIUS, &scale));
1692c71b3e2SJacob Faibussowitsch     PetscCheckFalse(diff > PETSC_SMALL * scale,PetscObjectComm((PetscObject)D), PETSC_ERR_PLIB, "MatMatTransposeMult() differs between MAT_INITIAL_MATRIX and MAT_REUSE_MATRIX");
1705f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&C));
171c4762a1bSJed Brown 
1725f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatTransposeMultEqual(A,B,D,10,&equal));
173*28b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*A*x");
1745f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&D));
1755f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&B));
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   }
178c4762a1bSJed Brown 
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(v));
181c4762a1bSJed Brown   ierr = PetscFinalize();
182c4762a1bSJed Brown   return ierr;
183c4762a1bSJed Brown }
184c4762a1bSJed Brown 
185c4762a1bSJed Brown /*TEST
186c4762a1bSJed Brown 
187c4762a1bSJed Brown     test:
188c4762a1bSJed Brown       output_file: output/ex104.out
189c4762a1bSJed Brown 
190c4762a1bSJed Brown     test:
191c4762a1bSJed Brown       suffix: 2
192c4762a1bSJed Brown       nsize: 2
193c4762a1bSJed Brown       output_file: output/ex104.out
194c4762a1bSJed Brown 
195c4762a1bSJed Brown     test:
196c4762a1bSJed Brown       suffix: 3
197c4762a1bSJed Brown       nsize: 4
198c4762a1bSJed Brown       output_file: output/ex104.out
199c4762a1bSJed Brown       args: -M 23 -N 31
200c4762a1bSJed Brown 
201c4762a1bSJed Brown     test:
202c4762a1bSJed Brown       suffix: 4
203c4762a1bSJed Brown       nsize: 4
204c4762a1bSJed Brown       output_file: output/ex104.out
205c4762a1bSJed Brown       args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via cyclic
206c4762a1bSJed Brown 
207c4762a1bSJed Brown     test:
208c4762a1bSJed Brown       suffix: 5
209c4762a1bSJed Brown       nsize: 4
210c4762a1bSJed Brown       output_file: output/ex104.out
211c4762a1bSJed Brown       args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via allgatherv
212c4762a1bSJed Brown 
213c20d7725SJed Brown     test:
214c20d7725SJed Brown       suffix: 6
215c20d7725SJed Brown       args: -mat_type elemental
216c20d7725SJed Brown       requires: elemental
217c20d7725SJed Brown       output_file: output/ex104.out
218c20d7725SJed Brown 
219c20d7725SJed Brown     test:
220c20d7725SJed Brown       suffix: 7
221c20d7725SJed Brown       nsize: 2
222c20d7725SJed Brown       args: -mat_type elemental
223c20d7725SJed Brown       requires: elemental
224c20d7725SJed Brown       output_file: output/ex104.out
225c20d7725SJed Brown 
226c4762a1bSJed Brown TEST*/
227