1c4762a1bSJed Brown static char help[] = "Tests sequential and parallel MatMatMult() and MatPtAP(), MatTransposeMatMult(), sequential MatMatTransposeMult(), MatRARt()\n\
2c4762a1bSJed Brown Input arguments are:\n\
3c4762a1bSJed Brown -f0 <input_file> -f1 <input_file> -f2 <input_file> -f3 <input_file> : file to load\n\n";
4c4762a1bSJed Brown /* Example of usage:
5c4762a1bSJed Brown ./ex94 -f0 <A_binary> -f1 <B_binary> -matmatmult_mat_view ascii::ascii_info -matmatmulttr_mat_view
6c4762a1bSJed Brown mpiexec -n 3 ./ex94 -f0 medium -f1 medium -f2 arco1 -f3 arco1 -matmatmult_mat_view
7c4762a1bSJed Brown */
8c4762a1bSJed Brown
9c4762a1bSJed Brown #include <petscmat.h>
10c4762a1bSJed Brown
11c4762a1bSJed Brown /*
12c4762a1bSJed Brown B = A - B
13c4762a1bSJed Brown norm = norm(B)
14c4762a1bSJed Brown */
MatNormDifference(Mat A,Mat B,PetscReal * norm)15d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNormDifference(Mat A, Mat B, PetscReal *norm)
16d71ae5a4SJacob Faibussowitsch {
17c4762a1bSJed Brown PetscFunctionBegin;
189566063dSJacob Faibussowitsch PetscCall(MatAXPY(B, -1.0, A, DIFFERENT_NONZERO_PATTERN));
199566063dSJacob Faibussowitsch PetscCall(MatNorm(B, NORM_FROBENIUS, norm));
203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
21c4762a1bSJed Brown }
22c4762a1bSJed Brown
main(int argc,char ** args)23d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
24d71ae5a4SJacob Faibussowitsch {
25c4762a1bSJed Brown Mat A, A_save, B, AT, ATT, BT, BTT, P, R, C, C1;
26c4762a1bSJed Brown Vec x, v1, v2, v3, v4;
27c4762a1bSJed Brown PetscViewer viewer;
28c4762a1bSJed Brown PetscMPIInt size, rank;
29c4762a1bSJed Brown PetscInt i, m, n, j, *idxn, M, N, nzp, rstart, rend;
30c4762a1bSJed Brown PetscReal norm, norm_abs, norm_tmp, fill = 4.0;
31c4762a1bSJed Brown PetscRandom rdm;
32c4762a1bSJed Brown char file[4][128];
33c4762a1bSJed Brown PetscBool flg, preload = PETSC_TRUE;
34c4762a1bSJed Brown PetscScalar *a, rval, alpha, none = -1.0;
35c4762a1bSJed Brown PetscBool Test_MatMatMult = PETSC_TRUE, Test_MatMatTr = PETSC_TRUE, Test_MatPtAP = PETSC_TRUE, Test_MatRARt = PETSC_TRUE, Test_MatMatMatMult = PETSC_TRUE;
36c4762a1bSJed Brown PetscBool Test_MatAXPY = PETSC_FALSE, view = PETSC_FALSE;
37c4762a1bSJed Brown PetscInt pm, pn, pM, pN;
38c4762a1bSJed Brown MatInfo info;
39c4762a1bSJed Brown PetscBool seqaij;
40c4762a1bSJed Brown MatType mattype;
41c4762a1bSJed Brown Mat Cdensetest, Pdense, Cdense, Adense;
42c4762a1bSJed Brown
43327415f7SBarry Smith PetscFunctionBeginUser;
44c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
47c4762a1bSJed Brown
489566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-fill", &fill, NULL));
499566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-matops_view", &view, NULL));
509566063dSJacob Faibussowitsch if (view) PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO));
51c4762a1bSJed Brown
52c4762a1bSJed Brown /* Load the matrices A_save and B */
539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f0", file[0], sizeof(file[0]), &flg));
5428b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate a file name for small matrix A with the -f0 option.");
559566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f1", file[1], sizeof(file[1]), &flg));
5628b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate a file name for small matrix B with the -f1 option.");
579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f2", file[2], sizeof(file[2]), &flg));
58c4762a1bSJed Brown if (!flg) {
59c4762a1bSJed Brown preload = PETSC_FALSE;
60c4762a1bSJed Brown } else {
619566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f3", file[3], sizeof(file[3]), &flg));
6228b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate a file name for test matrix B with the -f3 option.");
63c4762a1bSJed Brown }
64c4762a1bSJed Brown
65c4762a1bSJed Brown PetscPreLoadBegin(preload, "Load system");
669566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file[2 * PetscPreLoadIt], FILE_MODE_READ, &viewer));
679566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A_save));
689566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A_save));
699566063dSJacob Faibussowitsch PetscCall(MatLoad(A_save, viewer));
709566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer));
71c4762a1bSJed Brown
729566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file[2 * PetscPreLoadIt + 1], FILE_MODE_READ, &viewer));
739566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
749566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B));
759566063dSJacob Faibussowitsch PetscCall(MatLoad(B, viewer));
769566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer));
77c4762a1bSJed Brown
789566063dSJacob Faibussowitsch PetscCall(MatGetType(B, &mattype));
79c4762a1bSJed Brown
809566063dSJacob Faibussowitsch PetscCall(MatGetSize(B, &M, &N));
81c4762a1bSJed Brown nzp = PetscMax((PetscInt)(0.1 * M), 5);
82f1e99121SPierre Jolivet PetscCall(PetscMalloc2(nzp + 1, &idxn, nzp + 1, &a));
83c4762a1bSJed Brown
84c4762a1bSJed Brown /* Create vectors v1 and v2 that are compatible with A_save */
859566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &v1));
869566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A_save, &m, NULL));
879566063dSJacob Faibussowitsch PetscCall(VecSetSizes(v1, m, PETSC_DECIDE));
889566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(v1));
899566063dSJacob Faibussowitsch PetscCall(VecDuplicate(v1, &v2));
90c4762a1bSJed Brown
919566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
929566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm));
939566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-fill", &fill, NULL));
94c4762a1bSJed Brown
95c4762a1bSJed Brown /* Test MatAXPY() */
96c4762a1bSJed Brown /*-------------------*/
979566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-test_MatAXPY", &Test_MatAXPY));
98c4762a1bSJed Brown if (Test_MatAXPY) {
99c4762a1bSJed Brown Mat Btmp;
1009566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A_save, MAT_COPY_VALUES, &A));
1019566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &Btmp));
1029566063dSJacob Faibussowitsch PetscCall(MatAXPY(A, -1.0, B, DIFFERENT_NONZERO_PATTERN)); /* A = -B + A_save */
103c4762a1bSJed Brown
1049566063dSJacob Faibussowitsch PetscCall(MatScale(A, -1.0)); /* A = -A = B - A_save */
1059566063dSJacob Faibussowitsch PetscCall(MatAXPY(Btmp, -1.0, A, DIFFERENT_NONZERO_PATTERN)); /* Btmp = -A + B = A_save */
1069566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A_save, Btmp, 10, &flg));
10728b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatAXPY() is incorrect");
1089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
1099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Btmp));
110c4762a1bSJed Brown
111c4762a1bSJed Brown Test_MatMatMult = PETSC_FALSE;
112c4762a1bSJed Brown Test_MatMatTr = PETSC_FALSE;
113c4762a1bSJed Brown Test_MatPtAP = PETSC_FALSE;
114c4762a1bSJed Brown Test_MatRARt = PETSC_FALSE;
115c4762a1bSJed Brown Test_MatMatMatMult = PETSC_FALSE;
116c4762a1bSJed Brown }
117c4762a1bSJed Brown
118c4762a1bSJed Brown /* 1) Test MatMatMult() */
119c4762a1bSJed Brown /* ---------------------*/
120c4762a1bSJed Brown if (Test_MatMatMult) {
1219566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A_save, MAT_COPY_VALUES, &A));
1229566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(A, &AT));
1239566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(AT, &ATT));
1249566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(B, &BT));
1259566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(BT, &BTT));
126c20d7725SJed Brown
1279566063dSJacob Faibussowitsch PetscCall(MatMatMult(AT, B, MAT_INITIAL_MATRIX, fill, &C));
1289566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(AT, B, C, 10, &flg));
12928b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error: MatMatMult() for C=AT*B");
1309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
131c20d7725SJed Brown
1329566063dSJacob Faibussowitsch PetscCall(MatMatMult(ATT, B, MAT_INITIAL_MATRIX, fill, &C));
1339566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(ATT, B, C, 10, &flg));
13428b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error: MatMatMult() for C=ATT*B");
1359566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
136c20d7725SJed Brown
1379566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, B, MAT_INITIAL_MATRIX, fill, &C));
1389566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A, B, C, 10, &flg));
13928b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error: MatMatMult() for reuse C=A*B");
140c20d7725SJed Brown /* ATT has different matrix type as A (although they have same internal data structure),
141c20d7725SJed Brown we cannot call MatProductReplaceMats(ATT,NULL,NULL,C) and MatMatMult(ATT,B,MAT_REUSE_MATRIX,fill,&C) */
1429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
143c20d7725SJed Brown
1449566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, BTT, MAT_INITIAL_MATRIX, fill, &C));
1459566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A, BTT, C, 10, &flg));
14628b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error: MatMatMult() for C=A*BTT");
1479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
148c20d7725SJed Brown
1499566063dSJacob Faibussowitsch PetscCall(MatMatMult(ATT, BTT, MAT_INITIAL_MATRIX, fill, &C));
1509566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A, B, C, 10, &flg));
15128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error: MatMatMult()");
1529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
153c20d7725SJed Brown
1549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BTT));
1559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BT));
1569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ATT));
1579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AT));
158c20d7725SJed Brown
1599566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, B, MAT_INITIAL_MATRIX, fill, &C));
1609566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(C, "matmatmult_")); /* enable option '-matmatmult_' for matrix C */
1619566063dSJacob Faibussowitsch PetscCall(MatGetInfo(C, MAT_GLOBAL_SUM, &info));
162c4762a1bSJed Brown
163c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */
164c4762a1bSJed Brown alpha = 1.0;
165c4762a1bSJed Brown for (i = 0; i < 2; i++) {
166c4762a1bSJed Brown alpha -= 0.1;
1679566063dSJacob Faibussowitsch PetscCall(MatScale(A, alpha));
1689566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, fill, &C));
169c4762a1bSJed Brown }
1709566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A, B, C, 10, &flg));
17128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error: MatMatMult()");
1729566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
173c4762a1bSJed Brown
174c4762a1bSJed Brown /* Test MatDuplicate() of C=A*B */
1759566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &C1));
1769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C1));
1779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
178c4762a1bSJed Brown } /* if (Test_MatMatMult) */
179c4762a1bSJed Brown
180c4762a1bSJed Brown /* 2) Test MatTransposeMatMult() and MatMatTransposeMult() */
181c4762a1bSJed Brown /* ------------------------------------------------------- */
182c4762a1bSJed Brown if (Test_MatMatTr) {
183c4762a1bSJed Brown /* Create P */
184c4762a1bSJed Brown PetscInt PN, rstart, rend;
185c4762a1bSJed Brown PN = M / 2;
186c4762a1bSJed Brown nzp = 5; /* num of nonzeros in each row of P */
1879566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &P));
1889566063dSJacob Faibussowitsch PetscCall(MatSetSizes(P, PETSC_DECIDE, PETSC_DECIDE, M, PN));
1899566063dSJacob Faibussowitsch PetscCall(MatSetType(P, mattype));
1909566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(P, nzp, NULL));
1919566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(P, nzp, NULL, nzp, NULL));
1929566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(P, &rstart, &rend));
19348a46eb9SPierre Jolivet for (i = 0; i < nzp; i++) PetscCall(PetscRandomGetValue(rdm, &a[i]));
194c4762a1bSJed Brown for (i = rstart; i < rend; i++) {
195c4762a1bSJed Brown for (j = 0; j < nzp; j++) {
1969566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval));
197c4762a1bSJed Brown idxn[j] = (PetscInt)(PetscRealPart(rval) * PN);
198c4762a1bSJed Brown }
1999566063dSJacob Faibussowitsch PetscCall(MatSetValues(P, 1, &i, nzp, idxn, a, ADD_VALUES));
200c4762a1bSJed Brown }
2019566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY));
2029566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY));
203c4762a1bSJed Brown
204c4762a1bSJed Brown /* Create R = P^T */
2059566063dSJacob Faibussowitsch PetscCall(MatTranspose(P, MAT_INITIAL_MATRIX, &R));
206c4762a1bSJed Brown
207c4762a1bSJed Brown { /* Test R = P^T, C1 = R*B */
2089566063dSJacob Faibussowitsch PetscCall(MatMatMult(R, B, MAT_INITIAL_MATRIX, fill, &C1));
2099566063dSJacob Faibussowitsch PetscCall(MatTranspose(P, MAT_REUSE_MATRIX, &R));
2109566063dSJacob Faibussowitsch PetscCall(MatMatMult(R, B, MAT_REUSE_MATRIX, fill, &C1));
2119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C1));
212c4762a1bSJed Brown }
213c4762a1bSJed Brown
214c4762a1bSJed Brown /* C = P^T*B */
2159566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(P, B, MAT_INITIAL_MATRIX, fill, &C));
2169566063dSJacob Faibussowitsch PetscCall(MatGetInfo(C, MAT_GLOBAL_SUM, &info));
217c4762a1bSJed Brown
218c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */
2199566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(P, B, MAT_REUSE_MATRIX, fill, &C));
220c4762a1bSJed Brown if (view) {
2219566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "C = P^T * B:\n"));
2229566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
223c4762a1bSJed Brown }
2249566063dSJacob Faibussowitsch PetscCall(MatProductClear(C));
225c4762a1bSJed Brown if (view) {
2269566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nC = P^T * B after MatProductClear():\n"));
2279566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
228c4762a1bSJed Brown }
229c4762a1bSJed Brown
230c4762a1bSJed Brown /* Compare P^T*B and R*B */
2319566063dSJacob Faibussowitsch PetscCall(MatMatMult(R, B, MAT_INITIAL_MATRIX, fill, &C1));
2329566063dSJacob Faibussowitsch PetscCall(MatNormDifference(C, C1, &norm));
23308401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in MatTransposeMatMult(): %g", (double)norm);
2349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C1));
235c4762a1bSJed Brown
236c4762a1bSJed Brown /* Test MatDuplicate() of C=P^T*B */
2379566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &C1));
2389566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C1));
2399566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
240c4762a1bSJed Brown
241c4762a1bSJed Brown /* C = B*R^T */
2429566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)B, MATSEQAIJ, &seqaij));
243c4762a1bSJed Brown if (size == 1 && seqaij) {
2449566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(B, R, MAT_INITIAL_MATRIX, fill, &C));
2459566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(C, "matmatmulttr_")); /* enable '-matmatmulttr_' for matrix C */
2469566063dSJacob Faibussowitsch PetscCall(MatGetInfo(C, MAT_GLOBAL_SUM, &info));
247c4762a1bSJed Brown
248c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */
2499566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(B, R, MAT_REUSE_MATRIX, fill, &C));
250c4762a1bSJed Brown
251c4762a1bSJed Brown /* Check */
2529566063dSJacob Faibussowitsch PetscCall(MatMatMult(B, P, MAT_INITIAL_MATRIX, fill, &C1));
2539566063dSJacob Faibussowitsch PetscCall(MatNormDifference(C, C1, &norm));
25408401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in MatMatTransposeMult() %g", (double)norm);
2559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C1));
2569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
257c4762a1bSJed Brown }
2589566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P));
2599566063dSJacob Faibussowitsch PetscCall(MatDestroy(&R));
260c4762a1bSJed Brown }
261c4762a1bSJed Brown
262c4762a1bSJed Brown /* 3) Test MatPtAP() */
263c4762a1bSJed Brown /*-------------------*/
264c4762a1bSJed Brown if (Test_MatPtAP) {
265c4762a1bSJed Brown PetscInt PN;
266c4762a1bSJed Brown Mat Cdup;
267c4762a1bSJed Brown
2689566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A_save, MAT_COPY_VALUES, &A));
2699566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N));
2709566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n));
271c4762a1bSJed Brown
272c4762a1bSJed Brown PN = M / 2;
273bbea24aaSStefano Zampini nzp = (PetscInt)(0.1 * PN + 1); /* num of nonzeros in each row of P */
2749566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &P));
2759566063dSJacob Faibussowitsch PetscCall(MatSetSizes(P, PETSC_DECIDE, PETSC_DECIDE, N, PN));
2769566063dSJacob Faibussowitsch PetscCall(MatSetType(P, mattype));
2779566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(P, nzp, NULL));
2789566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(P, nzp, NULL, nzp, NULL));
27948a46eb9SPierre Jolivet for (i = 0; i < nzp; i++) PetscCall(PetscRandomGetValue(rdm, &a[i]));
2809566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(P, &rstart, &rend));
281c4762a1bSJed Brown for (i = rstart; i < rend; i++) {
282c4762a1bSJed Brown for (j = 0; j < nzp; j++) {
2839566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval));
284c4762a1bSJed Brown idxn[j] = (PetscInt)(PetscRealPart(rval) * PN);
285c4762a1bSJed Brown }
2869566063dSJacob Faibussowitsch PetscCall(MatSetValues(P, 1, &i, nzp, idxn, a, ADD_VALUES));
287c4762a1bSJed Brown }
2889566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY));
2899566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY));
290c4762a1bSJed Brown
2919566063dSJacob Faibussowitsch /* PetscCall(MatView(P,PETSC_VIEWER_STDOUT_WORLD)); */
2929566063dSJacob Faibussowitsch PetscCall(MatGetSize(P, &pM, &pN));
2939566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(P, &pm, &pn));
2949566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, P, MAT_INITIAL_MATRIX, fill, &C));
295c4762a1bSJed Brown
296c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */
297c4762a1bSJed Brown alpha = 1.0;
298c4762a1bSJed Brown for (i = 0; i < 2; i++) {
299c4762a1bSJed Brown alpha -= 0.1;
3009566063dSJacob Faibussowitsch PetscCall(MatScale(A, alpha));
3019566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, P, MAT_REUSE_MATRIX, fill, &C));
302c4762a1bSJed Brown }
303c4762a1bSJed Brown
304c4762a1bSJed Brown /* Test PtAP ops with P Dense and A either AIJ or SeqDense (it assumes MatPtAP_XAIJ_XAIJ is fine) */
3059566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &seqaij));
306c4762a1bSJed Brown if (seqaij) {
3079566063dSJacob Faibussowitsch PetscCall(MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cdensetest));
3089566063dSJacob Faibussowitsch PetscCall(MatConvert(P, MATSEQDENSE, MAT_INITIAL_MATRIX, &Pdense));
309c4762a1bSJed Brown } else {
3109566063dSJacob Faibussowitsch PetscCall(MatConvert(C, MATMPIDENSE, MAT_INITIAL_MATRIX, &Cdensetest));
3119566063dSJacob Faibussowitsch PetscCall(MatConvert(P, MATMPIDENSE, MAT_INITIAL_MATRIX, &Pdense));
312c4762a1bSJed Brown }
313c4762a1bSJed Brown
314c20d7725SJed Brown /* test with A(AIJ), Pdense -- call MatPtAP_Basic() when np>1 */
3159566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, Pdense, MAT_INITIAL_MATRIX, fill, &Cdense));
3169566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, Pdense, MAT_REUSE_MATRIX, fill, &Cdense));
3179566063dSJacob Faibussowitsch PetscCall(MatPtAPMultEqual(A, Pdense, Cdense, 10, &flg));
31828b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in MatPtAP with A AIJ and P Dense");
3199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cdense));
320c4762a1bSJed Brown
321c4762a1bSJed Brown /* test with A SeqDense */
322c4762a1bSJed Brown if (seqaij) {
3239566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSEQDENSE, MAT_INITIAL_MATRIX, &Adense));
3249566063dSJacob Faibussowitsch PetscCall(MatPtAP(Adense, Pdense, MAT_INITIAL_MATRIX, fill, &Cdense));
3259566063dSJacob Faibussowitsch PetscCall(MatPtAP(Adense, Pdense, MAT_REUSE_MATRIX, fill, &Cdense));
3269566063dSJacob Faibussowitsch PetscCall(MatPtAPMultEqual(Adense, Pdense, Cdense, 10, &flg));
32728b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in MatPtAP with A SeqDense and P SeqDense");
3289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cdense));
3299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Adense));
330c4762a1bSJed Brown }
3319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cdensetest));
3329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Pdense));
333c4762a1bSJed Brown
334c4762a1bSJed Brown /* Test MatDuplicate() of C=PtAP and MatView(Cdup,...) */
3359566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &Cdup));
336c4762a1bSJed Brown if (view) {
3379566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nC = P^T * A * P:\n"));
3389566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
339c4762a1bSJed Brown
3409566063dSJacob Faibussowitsch PetscCall(MatProductClear(C));
3419566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nC = P^T * A * P after MatProductClear():\n"));
3429566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
343c4762a1bSJed Brown
3449566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nCdup:\n"));
3459566063dSJacob Faibussowitsch PetscCall(MatView(Cdup, PETSC_VIEWER_STDOUT_WORLD));
346c4762a1bSJed Brown }
3479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cdup));
348c4762a1bSJed Brown
349c4762a1bSJed Brown if (size > 1 || !seqaij) Test_MatRARt = PETSC_FALSE;
350c4762a1bSJed Brown /* 4) Test MatRARt() */
351c4762a1bSJed Brown /* ----------------- */
352c4762a1bSJed Brown if (Test_MatRARt) {
353c20d7725SJed Brown Mat R, RARt, Rdense, RARtdense;
3549566063dSJacob Faibussowitsch PetscCall(MatTranspose(P, MAT_INITIAL_MATRIX, &R));
355c20d7725SJed Brown
356c20d7725SJed Brown /* Test MatRARt_Basic(), MatMatMatMult_Basic() */
3579566063dSJacob Faibussowitsch PetscCall(MatConvert(R, MATDENSE, MAT_INITIAL_MATRIX, &Rdense));
3589566063dSJacob Faibussowitsch PetscCall(MatRARt(A, Rdense, MAT_INITIAL_MATRIX, 2.0, &RARtdense));
3599566063dSJacob Faibussowitsch PetscCall(MatRARt(A, Rdense, MAT_REUSE_MATRIX, 2.0, &RARtdense));
360c20d7725SJed Brown
3619566063dSJacob Faibussowitsch PetscCall(MatConvert(RARtdense, MATAIJ, MAT_INITIAL_MATRIX, &RARt));
3629566063dSJacob Faibussowitsch PetscCall(MatNormDifference(C, RARt, &norm));
36308401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "|PtAP - RARtdense| = %g", (double)norm);
3649566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Rdense));
3659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RARtdense));
3669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RARt));
367c20d7725SJed Brown
368c20d7725SJed Brown /* Test MatRARt() for aij matrices */
3699566063dSJacob Faibussowitsch PetscCall(MatRARt(A, R, MAT_INITIAL_MATRIX, 2.0, &RARt));
3709566063dSJacob Faibussowitsch PetscCall(MatRARt(A, R, MAT_REUSE_MATRIX, 2.0, &RARt));
3719566063dSJacob Faibussowitsch PetscCall(MatNormDifference(C, RARt, &norm));
37208401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "|PtAP - RARt| = %g", (double)norm);
3739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&R));
3749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RARt));
375c4762a1bSJed Brown }
376c4762a1bSJed Brown
377c4762a1bSJed Brown if (Test_MatMatMatMult && size == 1) {
378c4762a1bSJed Brown Mat R, RAP;
3799566063dSJacob Faibussowitsch PetscCall(MatTranspose(P, MAT_INITIAL_MATRIX, &R));
3809566063dSJacob Faibussowitsch PetscCall(MatMatMatMult(R, A, P, MAT_INITIAL_MATRIX, 2.0, &RAP));
3819566063dSJacob Faibussowitsch PetscCall(MatMatMatMult(R, A, P, MAT_REUSE_MATRIX, 2.0, &RAP));
3829566063dSJacob Faibussowitsch PetscCall(MatNormDifference(C, RAP, &norm));
38308401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PtAP != RAP %g", (double)norm);
3849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&R));
3859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RAP));
386c4762a1bSJed Brown }
387c4762a1bSJed Brown
388c4762a1bSJed Brown /* Create vector x that is compatible with P */
3899566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
3909566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(P, &m, &n));
3919566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, n, PETSC_DECIDE));
3929566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x));
393c4762a1bSJed Brown
3949566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &v3));
3959566063dSJacob Faibussowitsch PetscCall(VecSetSizes(v3, n, PETSC_DECIDE));
3969566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(v3));
3979566063dSJacob Faibussowitsch PetscCall(VecDuplicate(v3, &v4));
398c4762a1bSJed Brown
399c4762a1bSJed Brown norm = 0.0;
400c4762a1bSJed Brown for (i = 0; i < 10; i++) {
4019566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm));
4029566063dSJacob Faibussowitsch PetscCall(MatMult(P, x, v1));
4039566063dSJacob Faibussowitsch PetscCall(MatMult(A, v1, v2)); /* v2 = A*P*x */
404c4762a1bSJed Brown
4059566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(P, v2, v3)); /* v3 = Pt*A*P*x */
4069566063dSJacob Faibussowitsch PetscCall(MatMult(C, x, v4)); /* v3 = C*x */
4079566063dSJacob Faibussowitsch PetscCall(VecNorm(v4, NORM_2, &norm_abs));
4089566063dSJacob Faibussowitsch PetscCall(VecAXPY(v4, none, v3));
4099566063dSJacob Faibussowitsch PetscCall(VecNorm(v4, NORM_2, &norm_tmp));
410c4762a1bSJed Brown
411c4762a1bSJed Brown norm_tmp /= norm_abs;
412c4762a1bSJed Brown if (norm_tmp > norm) norm = norm_tmp;
413c4762a1bSJed Brown }
41408401ef6SPierre Jolivet PetscCheck(norm < PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error: MatPtAP(), |v1 - v2|: %g", (double)norm);
415c4762a1bSJed Brown
4169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
4179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P));
4189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
4199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v3));
4209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v4));
4219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
422c4762a1bSJed Brown }
423c4762a1bSJed Brown
424c4762a1bSJed Brown /* Destroy objects */
4259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v1));
4269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v2));
4279566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm));
428f1e99121SPierre Jolivet PetscCall(PetscFree2(idxn, a));
429c4762a1bSJed Brown
4309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A_save));
4319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B));
432c4762a1bSJed Brown
433c4762a1bSJed Brown PetscPreLoadEnd();
4349566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
435b122ec5aSJacob Faibussowitsch return 0;
436c4762a1bSJed Brown }
437c4762a1bSJed Brown
438c4762a1bSJed Brown /*TEST
439c4762a1bSJed Brown
440c4762a1bSJed Brown test:
441c4762a1bSJed Brown suffix: 2_mattransposematmult_matmatmult
442c4762a1bSJed Brown nsize: 3
443dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
444c20d7725SJed Brown args: -f0 ${DATAFILESPATH}/matrices/medium -f1 ${DATAFILESPATH}/matrices/medium -mattransposematmult_via at*b> ex94_2.tmp 2>&1
445*3886731fSPierre Jolivet output_file: output/empty.out
446c4762a1bSJed Brown
447c4762a1bSJed Brown test:
448c4762a1bSJed Brown suffix: 2_mattransposematmult_scalable
449c4762a1bSJed Brown nsize: 3
450dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
451c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/medium -f1 ${DATAFILESPATH}/matrices/medium -mattransposematmult_via scalable> ex94_2.tmp 2>&1
452*3886731fSPierre Jolivet output_file: output/empty.out
453c4762a1bSJed Brown
454c4762a1bSJed Brown test:
455c4762a1bSJed Brown suffix: axpy_mpiaij
456dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
457c4762a1bSJed Brown nsize: 8
458c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/poisson_2d5p -f1 ${DATAFILESPATH}/matrices/poisson_2d13p -test_MatAXPY
459*3886731fSPierre Jolivet output_file: output/empty.out
460c4762a1bSJed Brown
461c4762a1bSJed Brown test:
462c4762a1bSJed Brown suffix: axpy_mpibaij
463dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
464c4762a1bSJed Brown nsize: 8
465c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/poisson_2d5p -f1 ${DATAFILESPATH}/matrices/poisson_2d13p -test_MatAXPY -mat_type baij
466*3886731fSPierre Jolivet output_file: output/empty.out
467c4762a1bSJed Brown
468c4762a1bSJed Brown test:
469c4762a1bSJed Brown suffix: axpy_mpisbaij
470dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
471c4762a1bSJed Brown nsize: 8
472c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/poisson_2d5p -f1 ${DATAFILESPATH}/matrices/poisson_2d13p -test_MatAXPY -mat_type sbaij
473*3886731fSPierre Jolivet output_file: output/empty.out
474c4762a1bSJed Brown
475c4762a1bSJed Brown test:
476c4762a1bSJed Brown suffix: matmatmult
477dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
478c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -viewer_binary_skip_info
479*3886731fSPierre Jolivet output_file: output/empty.out
480c4762a1bSJed Brown
481c4762a1bSJed Brown test:
482c4762a1bSJed Brown suffix: matmatmult_2
483dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
484c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -mat_type mpiaij -viewer_binary_skip_info
485*3886731fSPierre Jolivet output_file: output/empty.out
486c4762a1bSJed Brown
487c4762a1bSJed Brown test:
488c4762a1bSJed Brown suffix: matmatmult_scalable
489c4762a1bSJed Brown nsize: 4
490dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
491c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -matmatmult_via scalable
492*3886731fSPierre Jolivet output_file: output/empty.out
493c4762a1bSJed Brown
494c4762a1bSJed Brown test:
495c4762a1bSJed Brown suffix: ptap
496c4762a1bSJed Brown nsize: 3
497dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
498c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/medium -f1 ${DATAFILESPATH}/matrices/medium -matptap_via scalable
499*3886731fSPierre Jolivet output_file: output/empty.out
500c4762a1bSJed Brown
501c4762a1bSJed Brown test:
502c4762a1bSJed Brown suffix: rap
503c4762a1bSJed Brown nsize: 3
504dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
505c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/medium -f1 ${DATAFILESPATH}/matrices/medium
506*3886731fSPierre Jolivet output_file: output/empty.out
507c4762a1bSJed Brown
508c4762a1bSJed Brown test:
509c4762a1bSJed Brown suffix: scalable0
510dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
511c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -viewer_binary_skip_info
512*3886731fSPierre Jolivet output_file: output/empty.out
513c4762a1bSJed Brown
514c4762a1bSJed Brown test:
515c4762a1bSJed Brown suffix: scalable1
516dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
517c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -viewer_binary_skip_info -matptap_via scalable
518*3886731fSPierre Jolivet output_file: output/empty.out
519c4762a1bSJed Brown
520c4762a1bSJed Brown test:
521c4762a1bSJed Brown suffix: view
522c4762a1bSJed Brown nsize: 2
523dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
524c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/tiny -f1 ${DATAFILESPATH}/matrices/tiny -viewer_binary_skip_info -matops_view
525c4762a1bSJed Brown output_file: output/ex94_2.out
526c4762a1bSJed Brown
527c4762a1bSJed Brown TEST*/
528