xref: /petsc/src/mat/tests/ex104.c (revision 6d073425f519d7d5901a193c565dd78299bee587)
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 
main(int argc,char ** argv)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;
15d016bddeSToby 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;
21d016bddeSToby 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));
26d016bddeSToby 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));
30d016bddeSToby Isaac   PetscCall(PetscOptionsGetBool(NULL, NULL, "-columns_on_one_rank", &columns_on_one_rank, NULL));
319566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
32d016bddeSToby Isaac   if (!columns_on_one_rank) {
339566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N));
34d016bddeSToby Isaac   } else {
35d016bddeSToby Isaac     PetscCall(MatSetSizes(A, PETSC_DECIDE, rank == 0 ? N : 0, M, N));
36d016bddeSToby 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) {
84*93c18bbdSPierre Jolivet     PetscCall(MatScale(C, -1.0));
859566063dSJacob Faibussowitsch     PetscCall(MatMatMult(C, A, MAT_INITIAL_MATRIX, fill, &B));
86*93c18bbdSPierre Jolivet     PetscCall(MatScale(C, -1.0));
87*93c18bbdSPierre Jolivet     PetscCall(MatScale(B, -1.0));
889566063dSJacob Faibussowitsch     PetscCall(MatMatMultEqual(C, A, B, 10, &equal));
8928b400f6SJacob Faibussowitsch     PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B != C*A for matrix type transpose and seqdense");
909566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
91c20d7725SJed Brown   }
929566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   /* Test MatMatMult() */
95c4762a1bSJed Brown   if (Test_MatMatMult) {
969566063dSJacob Faibussowitsch     PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &B));        /* B = A^T */
979566063dSJacob Faibussowitsch     PetscCall(MatMatMult(B, A, MAT_INITIAL_MATRIX, fill, &C)); /* C = B*A = A^T*A */
989566063dSJacob Faibussowitsch     PetscCall(MatMatMult(B, A, MAT_REUSE_MATRIX, fill, &C));
999566063dSJacob Faibussowitsch     PetscCall(MatMatMultEqual(B, A, C, 10, &equal));
10028b400f6SJacob Faibussowitsch     PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B*A*x != C*x");
101c4762a1bSJed Brown 
102c4762a1bSJed Brown     /* Test MatDuplicate for matrix product */
1039566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &D));
104c20d7725SJed Brown 
1059566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&D));
1069566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
1079566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
108c4762a1bSJed Brown   }
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   /* Test MatTransposeMatMult() */
111c20d7725SJed Brown   if (!Aiselemental) {
112d016bddeSToby Isaac     Mat E;
113d016bddeSToby Isaac 
1149566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMult(A, A, MAT_INITIAL_MATRIX, fill, &D)); /* D = A^T*A */
1159566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMult(A, A, MAT_REUSE_MATRIX, fill, &D));
1169566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMultEqual(A, A, D, 10, &equal));
11728b400f6SJacob Faibussowitsch     PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "D*x != A^T*A*x");
118c4762a1bSJed Brown 
119c4762a1bSJed Brown     /* Test MatDuplicate for matrix product */
1209566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(D, MAT_COPY_VALUES, &C));
1219566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
122d016bddeSToby Isaac 
123d016bddeSToby Isaac     /* Test A*D for fast path when D is on one process */
124d016bddeSToby Isaac     PetscCall(MatSetRandom(D, NULL));
125d016bddeSToby Isaac     PetscCall(MatMatMult(A, D, MAT_INITIAL_MATRIX, fill, &E));
126d016bddeSToby Isaac     PetscCall(MatMatMult(A, D, MAT_REUSE_MATRIX, fill, &E));
127d016bddeSToby Isaac     PetscCall(MatMatMultEqual(A, D, E, 10, &equal));
128d016bddeSToby Isaac     PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "E*x != A*D*x");
129d016bddeSToby Isaac     PetscCall(MatDestroy(&E));
130d016bddeSToby Isaac 
1319566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&D));
132c4762a1bSJed Brown 
133c4762a1bSJed Brown     /* Test D*x = A^T*C*A*x, where C is in AIJ format */
1349566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A, &am, &an));
1359566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
136c4762a1bSJed Brown     if (size == 1) {
1379566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, am, am));
138c4762a1bSJed Brown     } else {
1399566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(C, am, am, PETSC_DECIDE, PETSC_DECIDE));
140c4762a1bSJed Brown     }
1419566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(C));
1429566063dSJacob Faibussowitsch     PetscCall(MatSetUp(C));
1439566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(C, &rstart, &rend));
144c4762a1bSJed Brown     v[0] = 1.0;
14548a46eb9SPierre Jolivet     for (i = rstart; i < rend; i++) PetscCall(MatSetValues(C, 1, &i, 1, &i, v, INSERT_VALUES));
1469566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
1479566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
148c4762a1bSJed Brown 
149c4762a1bSJed Brown     /* B = C*A, D = A^T*B */
1509566063dSJacob Faibussowitsch     PetscCall(MatMatMult(C, A, MAT_INITIAL_MATRIX, 1.0, &B));
1519566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMult(A, B, MAT_INITIAL_MATRIX, fill, &D));
1529566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMultEqual(A, B, D, 10, &equal));
15328b400f6SJacob Faibussowitsch     PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "D*x != A^T*B*x");
154c4762a1bSJed Brown 
1559566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&D));
1569566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
1579566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
158c4762a1bSJed Brown   }
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   /* Test MatMatTransposeMult() */
161c20d7725SJed Brown   if (!Aiselemental) {
162c4762a1bSJed Brown     PetscReal diff, scale;
163c4762a1bSJed Brown     PetscInt  am, an, aM, aN;
164c4762a1bSJed Brown 
1659566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A, &am, &an));
1669566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A, &aM, &aN));
1679566063dSJacob Faibussowitsch     PetscCall(MatCreateDense(PetscObjectComm((PetscObject)A), PETSC_DECIDE, an, aM + 10, aN, NULL, &B));
1689566063dSJacob Faibussowitsch     PetscCall(MatSetRandom(B, NULL));
1699566063dSJacob Faibussowitsch     PetscCall(MatMatTransposeMult(A, B, MAT_INITIAL_MATRIX, fill, &D)); /* D = A*A^T */
170c4762a1bSJed Brown 
171c4762a1bSJed Brown     /* Test MatDuplicate for matrix product */
1729566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(D, MAT_COPY_VALUES, &C));
173c4762a1bSJed Brown 
1749566063dSJacob Faibussowitsch     PetscCall(MatMatTransposeMult(A, B, MAT_REUSE_MATRIX, fill, &D));
1759566063dSJacob Faibussowitsch     PetscCall(MatAXPY(C, -1., D, SAME_NONZERO_PATTERN));
176c4762a1bSJed Brown 
1779566063dSJacob Faibussowitsch     PetscCall(MatNorm(C, NORM_FROBENIUS, &diff));
1789566063dSJacob Faibussowitsch     PetscCall(MatNorm(D, NORM_FROBENIUS, &scale));
179e00437b9SBarry Smith     PetscCheck(diff <= PETSC_SMALL * scale, PetscObjectComm((PetscObject)D), PETSC_ERR_PLIB, "MatMatTransposeMult() differs between MAT_INITIAL_MATRIX and MAT_REUSE_MATRIX");
1809566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
181c4762a1bSJed Brown 
1829566063dSJacob Faibussowitsch     PetscCall(MatMatTransposeMultEqual(A, B, D, 10, &equal));
18328b400f6SJacob Faibussowitsch     PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "D*x != A^T*A*x");
1849566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&D));
1859566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
186c4762a1bSJed Brown   }
187c4762a1bSJed Brown 
1889566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1899566063dSJacob Faibussowitsch   PetscCall(PetscFree(v));
1909566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
191b122ec5aSJacob Faibussowitsch   return 0;
192c4762a1bSJed Brown }
193c4762a1bSJed Brown 
194c4762a1bSJed Brown /*TEST
195c4762a1bSJed Brown 
196c4762a1bSJed Brown     test:
1973886731fSPierre Jolivet       output_file: output/empty.out
198c4762a1bSJed Brown 
199c4762a1bSJed Brown     test:
200c4762a1bSJed Brown       suffix: 2
201c4762a1bSJed Brown       nsize: 2
2023886731fSPierre Jolivet       output_file: output/empty.out
203c4762a1bSJed Brown 
204c4762a1bSJed Brown     test:
205c4762a1bSJed Brown       suffix: 3
206c4762a1bSJed Brown       nsize: 4
2073886731fSPierre Jolivet       output_file: output/empty.out
208c4762a1bSJed Brown       args: -M 23 -N 31
209c4762a1bSJed Brown 
210c4762a1bSJed Brown     test:
211c4762a1bSJed Brown       suffix: 4
212c4762a1bSJed Brown       nsize: 4
2133886731fSPierre Jolivet       output_file: output/empty.out
214c4762a1bSJed Brown       args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via cyclic
215c4762a1bSJed Brown 
216c4762a1bSJed Brown     test:
217c4762a1bSJed Brown       suffix: 5
218c4762a1bSJed Brown       nsize: 4
2193886731fSPierre Jolivet       output_file: output/empty.out
220c4762a1bSJed Brown       args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via allgatherv
221c4762a1bSJed Brown 
222c20d7725SJed Brown     test:
223c20d7725SJed Brown       suffix: 6
224c20d7725SJed Brown       args: -mat_type elemental
225c20d7725SJed Brown       requires: elemental
2263886731fSPierre Jolivet       output_file: output/empty.out
227c20d7725SJed Brown 
2285d8c7819SPierre Jolivet     testset:
229c20d7725SJed Brown       nsize: 2
2303886731fSPierre Jolivet       output_file: output/empty.out
2315d8c7819SPierre Jolivet       requires: elemental
2325d8c7819SPierre Jolivet       test:
2335d8c7819SPierre Jolivet         suffix: 7_dense
2345d8c7819SPierre Jolivet         args: -mat_type dense -mat_product_algorithm elemental
2355d8c7819SPierre Jolivet       test:
2365d8c7819SPierre Jolivet         suffix: 7_elemental
2375d8c7819SPierre Jolivet         args: -mat_type elemental
238c20d7725SJed Brown 
239d016bddeSToby Isaac     test:
240d016bddeSToby Isaac       suffix: 8
241d016bddeSToby Isaac       nsize: 4
242d016bddeSToby Isaac       args: -columns_on_one_rank
2433886731fSPierre Jolivet       output_file: output/empty.out
244d016bddeSToby Isaac 
245d016bddeSToby Isaac     test:
246d016bddeSToby Isaac       suffix: 9
247d016bddeSToby Isaac       nsize: 4
248d016bddeSToby Isaac       requires: cuda
249d016bddeSToby Isaac       args: -columns_on_one_rank -mat_type densecuda
2503886731fSPierre Jolivet       output_file: output/empty.out
251d016bddeSToby Isaac 
252c4762a1bSJed Brown TEST*/
253