xref: /petsc/src/mat/tests/ex104.c (revision d016bdde269de9549a736fe23cc3868ea52c341b)
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 
9d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
10d71ae5a4SJacob Faibussowitsch {
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;
15*d016bddeSToby Isaac   PetscBool       columns_on_one_rank = PETSC_FALSE;
16c4762a1bSJed Brown   PetscReal       fill                = 1.0;
17c4762a1bSJed Brown   IS              isrows, iscols;
18c4762a1bSJed Brown   const PetscInt *rows, *cols;
19c4762a1bSJed Brown   PetscScalar    *v, rval;
20c4762a1bSJed Brown   PetscBool       Test_MatMatMult = PETSC_TRUE;
21*d016bddeSToby Isaac   PetscMPIInt     size, rank;
22c4762a1bSJed Brown 
23327415f7SBarry Smith   PetscFunctionBeginUser;
24c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
26*d016bddeSToby Isaac   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
27c4762a1bSJed Brown 
289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL));
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL));
30*d016bddeSToby Isaac   PetscCall(PetscOptionsGetBool(NULL, NULL, "-columns_on_one_rank", &columns_on_one_rank, NULL));
319566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
32*d016bddeSToby Isaac   if (!columns_on_one_rank) {
339566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N));
34*d016bddeSToby Isaac   } else {
35*d016bddeSToby Isaac     PetscCall(MatSetSizes(A, PETSC_DECIDE, rank == 0 ? N : 0, M, N));
36*d016bddeSToby Isaac   }
379566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATDENSE));
389566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
399566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
409566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &r));
419566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(r));
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   /* Set local matrix entries */
449566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipIS(A, &isrows, &iscols));
459566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isrows, &nrows));
469566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrows, &rows));
479566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(iscols, &ncols));
489566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(iscols, &cols));
499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrows * ncols, &v));
50c4762a1bSJed Brown   for (i = 0; i < nrows; i++) {
51c4762a1bSJed Brown     for (j = 0; j < ncols; j++) {
529566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(r, &rval));
53c4762a1bSJed Brown       v[i * ncols + j] = rval;
54c4762a1bSJed Brown     }
55c4762a1bSJed Brown   }
569566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, nrows, rows, ncols, cols, v, INSERT_VALUES));
579566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
599566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrows, &rows));
609566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscols, &cols));
619566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrows));
629566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscols));
639566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&r));
64c4762a1bSJed Brown 
659566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATELEMENTAL, &Aiselemental));
66c20d7725SJed Brown 
67c20d7725SJed Brown   /* Test MatCreateTranspose() and MatTranspose() */
689566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(A, &C));
699566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &B)); /* B = A^T */
709566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(C, B, 10, &equal));
7128b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "A^T*x != (x^T*A)^T");
729566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
73c20d7725SJed Brown 
749566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B));
75c20d7725SJed Brown   if (!Aiselemental) {
769566063dSJacob Faibussowitsch     PetscCall(MatTranspose(B, MAT_INPLACE_MATRIX, &B));
779566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(C, B, 10, &equal));
7828b400f6SJacob Faibussowitsch     PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "C*x != B*x");
79c20d7725SJed Brown   }
809566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
81c20d7725SJed Brown 
82c20d7725SJed Brown   /* Test B = C*A for matrix type transpose and seqdense */
83c20d7725SJed Brown   if (size == 1 && !Aiselemental) {
849566063dSJacob Faibussowitsch     PetscCall(MatMatMult(C, A, MAT_INITIAL_MATRIX, fill, &B));
859566063dSJacob Faibussowitsch     PetscCall(MatMatMultEqual(C, A, B, 10, &equal));
8628b400f6SJacob Faibussowitsch     PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B != C*A for matrix type transpose and seqdense");
879566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
88c20d7725SJed Brown   }
899566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /* Test MatMatMult() */
92c4762a1bSJed Brown   if (Test_MatMatMult) {
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) {
109*d016bddeSToby Isaac     Mat E;
110*d016bddeSToby Isaac 
1119566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMult(A, A, MAT_INITIAL_MATRIX, fill, &D)); /* D = A^T*A */
1129566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMult(A, A, MAT_REUSE_MATRIX, fill, &D));
1139566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMultEqual(A, A, D, 10, &equal));
11428b400f6SJacob Faibussowitsch     PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "D*x != A^T*A*x");
115c4762a1bSJed Brown 
116c4762a1bSJed Brown     /* Test MatDuplicate for matrix product */
1179566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(D, MAT_COPY_VALUES, &C));
1189566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
119*d016bddeSToby Isaac 
120*d016bddeSToby Isaac     /* Test A*D for fast path when D is on one process */
121*d016bddeSToby Isaac     PetscCall(MatSetRandom(D, NULL));
122*d016bddeSToby Isaac     PetscCall(MatMatMult(A, D, MAT_INITIAL_MATRIX, fill, &E));
123*d016bddeSToby Isaac     PetscCall(MatMatMult(A, D, MAT_REUSE_MATRIX, fill, &E));
124*d016bddeSToby Isaac     PetscCall(MatMatMultEqual(A, D, E, 10, &equal));
125*d016bddeSToby Isaac     PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "E*x != A*D*x");
126*d016bddeSToby Isaac     PetscCall(MatDestroy(&E));
127*d016bddeSToby Isaac 
1289566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&D));
129c4762a1bSJed Brown 
130c4762a1bSJed Brown     /* Test D*x = A^T*C*A*x, where C is in AIJ format */
1319566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A, &am, &an));
1329566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
133c4762a1bSJed Brown     if (size == 1) {
1349566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, am, am));
135c4762a1bSJed Brown     } else {
1369566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(C, am, am, PETSC_DECIDE, PETSC_DECIDE));
137c4762a1bSJed Brown     }
1389566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(C));
1399566063dSJacob Faibussowitsch     PetscCall(MatSetUp(C));
1409566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(C, &rstart, &rend));
141c4762a1bSJed Brown     v[0] = 1.0;
14248a46eb9SPierre Jolivet     for (i = rstart; i < rend; i++) PetscCall(MatSetValues(C, 1, &i, 1, &i, v, INSERT_VALUES));
1439566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
1449566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
145c4762a1bSJed Brown 
146c4762a1bSJed Brown     /* B = C*A, D = A^T*B */
1479566063dSJacob Faibussowitsch     PetscCall(MatMatMult(C, A, MAT_INITIAL_MATRIX, 1.0, &B));
1489566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMult(A, B, MAT_INITIAL_MATRIX, fill, &D));
1499566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMultEqual(A, B, D, 10, &equal));
15028b400f6SJacob Faibussowitsch     PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "D*x != A^T*B*x");
151c4762a1bSJed Brown 
1529566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&D));
1539566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
1549566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
155c4762a1bSJed Brown   }
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   /* Test MatMatTransposeMult() */
158c20d7725SJed Brown   if (!Aiselemental) {
159c4762a1bSJed Brown     PetscReal diff, scale;
160c4762a1bSJed Brown     PetscInt  am, an, aM, aN;
161c4762a1bSJed Brown 
1629566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A, &am, &an));
1639566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A, &aM, &aN));
1649566063dSJacob Faibussowitsch     PetscCall(MatCreateDense(PetscObjectComm((PetscObject)A), PETSC_DECIDE, an, aM + 10, aN, NULL, &B));
1659566063dSJacob Faibussowitsch     PetscCall(MatSetRandom(B, NULL));
1669566063dSJacob Faibussowitsch     PetscCall(MatMatTransposeMult(A, B, MAT_INITIAL_MATRIX, fill, &D)); /* D = A*A^T */
167c4762a1bSJed Brown 
168c4762a1bSJed Brown     /* Test MatDuplicate for matrix product */
1699566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(D, MAT_COPY_VALUES, &C));
170c4762a1bSJed Brown 
1719566063dSJacob Faibussowitsch     PetscCall(MatMatTransposeMult(A, B, MAT_REUSE_MATRIX, fill, &D));
1729566063dSJacob Faibussowitsch     PetscCall(MatAXPY(C, -1., D, SAME_NONZERO_PATTERN));
173c4762a1bSJed Brown 
1749566063dSJacob Faibussowitsch     PetscCall(MatNorm(C, NORM_FROBENIUS, &diff));
1759566063dSJacob Faibussowitsch     PetscCall(MatNorm(D, NORM_FROBENIUS, &scale));
176e00437b9SBarry Smith     PetscCheck(diff <= PETSC_SMALL * scale, PetscObjectComm((PetscObject)D), PETSC_ERR_PLIB, "MatMatTransposeMult() differs between MAT_INITIAL_MATRIX and MAT_REUSE_MATRIX");
1779566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
178c4762a1bSJed Brown 
1799566063dSJacob Faibussowitsch     PetscCall(MatMatTransposeMultEqual(A, B, D, 10, &equal));
18028b400f6SJacob Faibussowitsch     PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "D*x != A^T*A*x");
1819566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&D));
1829566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
183c4762a1bSJed Brown   }
184c4762a1bSJed Brown 
1859566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1869566063dSJacob Faibussowitsch   PetscCall(PetscFree(v));
1879566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
188b122ec5aSJacob Faibussowitsch   return 0;
189c4762a1bSJed Brown }
190c4762a1bSJed Brown 
191c4762a1bSJed Brown /*TEST
192c4762a1bSJed Brown 
193c4762a1bSJed Brown     test:
194c4762a1bSJed Brown       output_file: output/ex104.out
195c4762a1bSJed Brown 
196c4762a1bSJed Brown     test:
197c4762a1bSJed Brown       suffix: 2
198c4762a1bSJed Brown       nsize: 2
199c4762a1bSJed Brown       output_file: output/ex104.out
200c4762a1bSJed Brown 
201c4762a1bSJed Brown     test:
202c4762a1bSJed Brown       suffix: 3
203c4762a1bSJed Brown       nsize: 4
204c4762a1bSJed Brown       output_file: output/ex104.out
205c4762a1bSJed Brown       args: -M 23 -N 31
206c4762a1bSJed Brown 
207c4762a1bSJed Brown     test:
208c4762a1bSJed Brown       suffix: 4
209c4762a1bSJed Brown       nsize: 4
210c4762a1bSJed Brown       output_file: output/ex104.out
211c4762a1bSJed Brown       args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via cyclic
212c4762a1bSJed Brown 
213c4762a1bSJed Brown     test:
214c4762a1bSJed Brown       suffix: 5
215c4762a1bSJed Brown       nsize: 4
216c4762a1bSJed Brown       output_file: output/ex104.out
217c4762a1bSJed Brown       args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via allgatherv
218c4762a1bSJed Brown 
219c20d7725SJed Brown     test:
220c20d7725SJed Brown       suffix: 6
221c20d7725SJed Brown       args: -mat_type elemental
222c20d7725SJed Brown       requires: elemental
223c20d7725SJed Brown       output_file: output/ex104.out
224c20d7725SJed Brown 
2255d8c7819SPierre Jolivet     testset:
226c20d7725SJed Brown       nsize: 2
227c20d7725SJed Brown       output_file: output/ex104.out
2285d8c7819SPierre Jolivet       requires: elemental
2295d8c7819SPierre Jolivet       test:
2305d8c7819SPierre Jolivet         suffix: 7_dense
2315d8c7819SPierre Jolivet         args: -mat_type dense -mat_product_algorithm elemental
2325d8c7819SPierre Jolivet       test:
2335d8c7819SPierre Jolivet         suffix: 7_elemental
2345d8c7819SPierre Jolivet         args: -mat_type elemental
235c20d7725SJed Brown 
236*d016bddeSToby Isaac     test:
237*d016bddeSToby Isaac       suffix: 8
238*d016bddeSToby Isaac       nsize: 4
239*d016bddeSToby Isaac       args: -columns_on_one_rank
240*d016bddeSToby Isaac       output_file: output/ex104.out
241*d016bddeSToby Isaac 
242*d016bddeSToby Isaac     test:
243*d016bddeSToby Isaac       suffix: 9
244*d016bddeSToby Isaac       nsize: 4
245*d016bddeSToby Isaac       requires: cuda
246*d016bddeSToby Isaac       args: -columns_on_one_rank -mat_type densecuda
247*d016bddeSToby Isaac       output_file: output/ex104.out
248*d016bddeSToby Isaac 
249c4762a1bSJed Brown TEST*/
250