xref: /petsc/src/mat/tests/ex221.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
1c4762a1bSJed Brown static char help[] = "Tests various routines for MATSHELL\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown typedef struct _n_User *User;
6c4762a1bSJed Brown struct _n_User {
7c4762a1bSJed Brown   Mat B;
8c4762a1bSJed Brown };
9c4762a1bSJed Brown 
10d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_User(Mat A, Vec X)
11d71ae5a4SJacob Faibussowitsch {
12c4762a1bSJed Brown   User user;
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   PetscFunctionBegin;
159566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &user));
169566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(user->B, X));
17*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18c4762a1bSJed Brown }
19c4762a1bSJed Brown 
20d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_User(Mat A, Vec X, Vec Y)
21d71ae5a4SJacob Faibussowitsch {
22c4762a1bSJed Brown   User user;
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   PetscFunctionBegin;
259566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &user));
269566063dSJacob Faibussowitsch   PetscCall(MatMult(user->B, X, Y));
27*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28c4762a1bSJed Brown }
29c4762a1bSJed Brown 
30d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_User(Mat A, Vec X, Vec Y)
31d71ae5a4SJacob Faibussowitsch {
32c4762a1bSJed Brown   User user;
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   PetscFunctionBegin;
359566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &user));
369566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(user->B, X, Y));
37*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38c4762a1bSJed Brown }
39c4762a1bSJed Brown 
40d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCopy_User(Mat A, Mat X, MatStructure str)
41d71ae5a4SJacob Faibussowitsch {
42c4762a1bSJed Brown   User user, userX;
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   PetscFunctionBegin;
459566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &user));
469566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(X, &userX));
4708401ef6SPierre Jolivet   PetscCheck(user == userX, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "This should not happen");
489566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)user->B));
49*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
50c4762a1bSJed Brown }
51c4762a1bSJed Brown 
52d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_User(Mat A)
53d71ae5a4SJacob Faibussowitsch {
54c4762a1bSJed Brown   User user;
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   PetscFunctionBegin;
579566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &user));
589566063dSJacob Faibussowitsch   PetscCall(PetscObjectDereference((PetscObject)user->B));
59*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
60c4762a1bSJed Brown }
61c4762a1bSJed Brown 
62d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
63d71ae5a4SJacob Faibussowitsch {
64c4762a1bSJed Brown   User         user;
65c4762a1bSJed Brown   Mat          A, S;
66c4762a1bSJed Brown   PetscScalar *data, diag = 1.3;
67c4762a1bSJed Brown   PetscReal    tol = PETSC_SMALL;
68c4762a1bSJed Brown   PetscInt     i, j, m = PETSC_DECIDE, n = PETSC_DECIDE, M = 17, N = 15, s1, s2;
69c4762a1bSJed Brown   PetscInt     test, ntest = 2;
70c4762a1bSJed Brown   PetscMPIInt  rank, size;
71c4762a1bSJed Brown   PetscBool    nc        = PETSC_FALSE, cong, flg;
72c4762a1bSJed Brown   PetscBool    ronl      = PETSC_TRUE;
73c4762a1bSJed Brown   PetscBool    randomize = PETSC_FALSE, submat = PETSC_FALSE;
74c4762a1bSJed Brown   PetscBool    keep         = PETSC_FALSE;
75c4762a1bSJed Brown   PetscBool    testzerorows = PETSC_TRUE, testdiagscale = PETSC_TRUE, testgetdiag = PETSC_TRUE, testsubmat = PETSC_TRUE;
76c4762a1bSJed Brown   PetscBool    testshift = PETSC_TRUE, testscale = PETSC_TRUE, testdup = PETSC_TRUE, testreset = PETSC_TRUE;
77c4762a1bSJed Brown   PetscBool    testaxpy = PETSC_TRUE, testaxpyd = PETSC_TRUE, testaxpyerr = PETSC_FALSE;
78c4762a1bSJed Brown 
79327415f7SBarry Smith   PetscFunctionBeginUser;
809566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
819566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
829566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL));
849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL));
859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ml", &m, NULL));
869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nl", &n, NULL));
879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-square_nc", &nc, NULL));
889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-rows_only", &ronl, NULL));
899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-randomize", &randomize, NULL));
909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-submat", &submat, NULL));
919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_zerorows", &testzerorows, NULL));
929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_diagscale", &testdiagscale, NULL));
939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_getdiag", &testgetdiag, NULL));
949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_shift", &testshift, NULL));
959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_scale", &testscale, NULL));
969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_dup", &testdup, NULL));
979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_reset", &testreset, NULL));
989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_submat", &testsubmat, NULL));
999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_axpy", &testaxpy, NULL));
1009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_axpy_different", &testaxpyd, NULL));
1019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_axpy_error", &testaxpyerr, NULL));
1029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-loop", &ntest, NULL));
1039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-tol", &tol, NULL));
1049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-diag", &diag, NULL));
1059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-keep", &keep, NULL));
106c4762a1bSJed Brown   /* This tests square matrices with different row/col layout */
107c4762a1bSJed Brown   if (nc && size > 1) {
108c4762a1bSJed Brown     M = PetscMax(PetscMax(N, M), 1);
109c4762a1bSJed Brown     N = M;
110c4762a1bSJed Brown     m = n = 0;
1119371c9d4SSatish Balay     if (rank == 0) {
1129371c9d4SSatish Balay       m = M - 1;
1139371c9d4SSatish Balay       n = 1;
1149371c9d4SSatish Balay     } else if (rank == 1) {
1159371c9d4SSatish Balay       m = 1;
1169371c9d4SSatish Balay       n = N - 1;
1179371c9d4SSatish Balay     }
118c4762a1bSJed Brown   }
1199566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD, m, n, M, N, NULL, &A));
1209566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
1219566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
1229566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &s1, NULL));
123c4762a1bSJed Brown   s2 = 1;
124c4762a1bSJed Brown   while (s2 < M) s2 *= 10;
1259566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(A, &data));
126c4762a1bSJed Brown   for (j = 0; j < N; j++) {
127ad540459SPierre Jolivet     for (i = 0; i < m; i++) data[j * m + i] = s2 * j + i + s1 + 1;
128c4762a1bSJed Brown   }
1299566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   if (submat) {
133c4762a1bSJed Brown     Mat      A2;
134c4762a1bSJed Brown     IS       r, c;
135c4762a1bSJed Brown     PetscInt rst, ren, cst, cen;
136c4762a1bSJed Brown 
1379566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(A, &rst, &ren));
1389566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRangeColumn(A, &cst, &cen));
1399566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)A), (ren - rst) / 2, rst, 1, &r));
1409566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)A), (cen - cst) / 2, cst, 1, &c));
1419566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(A, r, c, MAT_INITIAL_MATRIX, &A2));
1429566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&r));
1439566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&c));
1449566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
145c4762a1bSJed Brown     A = A2;
146c4762a1bSJed Brown   }
147c4762a1bSJed Brown 
1489566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
1499566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
1509566063dSJacob Faibussowitsch   PetscCall(MatHasCongruentLayouts(A, &cong));
151c4762a1bSJed Brown 
1529566063dSJacob Faibussowitsch   PetscCall(MatConvert(A, MATAIJ, MAT_INPLACE_MATRIX, &A));
1539566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_KEEP_NONZERO_PATTERN, keep));
1549566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)A, "initial"));
1559566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(A, NULL, "-view_mat"));
156c4762a1bSJed Brown 
1579566063dSJacob Faibussowitsch   PetscCall(PetscNew(&user));
1589566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD, m, n, M, N, user, &S));
1599566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(S, MATOP_MULT, (void (*)(void))MatMult_User));
1609566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(S, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_User));
1611baa6e33SBarry Smith   if (cong) PetscCall(MatShellSetOperation(S, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_User));
1629566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(S, MATOP_COPY, (void (*)(void))MatCopy_User));
1639566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(S, MATOP_DESTROY, (void (*)(void))MatDestroy_User));
1649566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &user->B));
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   /* Square and rows only scaling */
167c4762a1bSJed Brown   ronl = cong ? ronl : PETSC_TRUE;
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   for (test = 0; test < ntest; test++) {
170c4762a1bSJed Brown     PetscReal err;
171c4762a1bSJed Brown 
1729566063dSJacob Faibussowitsch     PetscCall(MatMultAddEqual(A, S, 10, &flg));
17348a46eb9SPierre Jolivet     if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[test %" PetscInt_FMT "] Error mult add\n", test));
1749566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAddEqual(A, S, 10, &flg));
17548a46eb9SPierre Jolivet     if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[test %" PetscInt_FMT "] Error mult add (T)\n", test));
176c4762a1bSJed Brown     if (testzerorows) {
177c4762a1bSJed Brown       Mat       ST, B, C, BT, BTT;
178c4762a1bSJed Brown       IS        zr;
179c4762a1bSJed Brown       Vec       x = NULL, b1 = NULL, b2 = NULL;
180c4762a1bSJed Brown       PetscInt *idxs = NULL, nr = 0;
181c4762a1bSJed Brown 
182c4762a1bSJed Brown       if (rank == (test % size)) {
183c4762a1bSJed Brown         nr = 1;
1849566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(nr, &idxs));
185c4762a1bSJed Brown         if (test % 2) {
186c4762a1bSJed Brown           idxs[0] = (2 * M - 1 - test / 2) % M;
187c4762a1bSJed Brown         } else {
188c4762a1bSJed Brown           idxs[0] = (test / 2) % M;
189c4762a1bSJed Brown         }
190c4762a1bSJed Brown         idxs[0] = PetscMax(idxs[0], 0);
191c4762a1bSJed Brown       }
1929566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, nr, idxs, PETSC_OWN_POINTER, &zr));
1939566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)zr, "ZR"));
1949566063dSJacob Faibussowitsch       PetscCall(ISViewFromOptions(zr, NULL, "-view_is"));
1959566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(A, &x, &b1));
196c4762a1bSJed Brown       if (randomize) {
1979566063dSJacob Faibussowitsch         PetscCall(VecSetRandom(x, NULL));
1989566063dSJacob Faibussowitsch         PetscCall(VecSetRandom(b1, NULL));
199c4762a1bSJed Brown       } else {
2009566063dSJacob Faibussowitsch         PetscCall(VecSet(x, 11.4));
2019566063dSJacob Faibussowitsch         PetscCall(VecSet(b1, -14.2));
202c4762a1bSJed Brown       }
2039566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(b1, &b2));
2049566063dSJacob Faibussowitsch       PetscCall(VecCopy(b1, b2));
2059566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)b1, "A_B1"));
2069566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)b2, "A_B2"));
207c4762a1bSJed Brown       if (size > 1 && !cong) { /* MATMPIAIJ ZeroRows and ZeroRowsColumns are buggy in this case */
2089566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&b1));
209c4762a1bSJed Brown       }
210c4762a1bSJed Brown       if (ronl) {
2119566063dSJacob Faibussowitsch         PetscCall(MatZeroRowsIS(A, zr, diag, x, b1));
2129566063dSJacob Faibussowitsch         PetscCall(MatZeroRowsIS(S, zr, diag, x, b2));
213c4762a1bSJed Brown       } else {
2149566063dSJacob Faibussowitsch         PetscCall(MatZeroRowsColumnsIS(A, zr, diag, x, b1));
2159566063dSJacob Faibussowitsch         PetscCall(MatZeroRowsColumnsIS(S, zr, diag, x, b2));
2169566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&zr));
217c4762a1bSJed Brown         /* Mix zerorows and zerorowscols */
218c4762a1bSJed Brown         nr   = 0;
219c4762a1bSJed Brown         idxs = NULL;
220dd400576SPatrick Sanan         if (rank == 0) {
221c4762a1bSJed Brown           nr = 1;
2229566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(nr, &idxs));
223c4762a1bSJed Brown           if (test % 2) {
224c4762a1bSJed Brown             idxs[0] = (3 * M - 2 - test / 2) % M;
225c4762a1bSJed Brown           } else {
226c4762a1bSJed Brown             idxs[0] = (test / 2 + 1) % M;
227c4762a1bSJed Brown           }
228c4762a1bSJed Brown           idxs[0] = PetscMax(idxs[0], 0);
229c4762a1bSJed Brown         }
2309566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, nr, idxs, PETSC_OWN_POINTER, &zr));
2319566063dSJacob Faibussowitsch         PetscCall(PetscObjectSetName((PetscObject)zr, "ZR2"));
2329566063dSJacob Faibussowitsch         PetscCall(ISViewFromOptions(zr, NULL, "-view_is"));
2339566063dSJacob Faibussowitsch         PetscCall(MatZeroRowsIS(A, zr, diag * 2.0 + PETSC_SMALL, NULL, NULL));
2349566063dSJacob Faibussowitsch         PetscCall(MatZeroRowsIS(S, zr, diag * 2.0 + PETSC_SMALL, NULL, NULL));
235c4762a1bSJed Brown       }
2369566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&zr));
237c4762a1bSJed Brown 
238c4762a1bSJed Brown       if (b1) {
239c4762a1bSJed Brown         Vec b;
240c4762a1bSJed Brown 
2419566063dSJacob Faibussowitsch         PetscCall(VecViewFromOptions(b1, NULL, "-view_b"));
2429566063dSJacob Faibussowitsch         PetscCall(VecViewFromOptions(b2, NULL, "-view_b"));
2439566063dSJacob Faibussowitsch         PetscCall(VecDuplicate(b1, &b));
2449566063dSJacob Faibussowitsch         PetscCall(VecCopy(b1, b));
2459566063dSJacob Faibussowitsch         PetscCall(VecAXPY(b, -1.0, b2));
2469566063dSJacob Faibussowitsch         PetscCall(VecNorm(b, NORM_INFINITY, &err));
24748a46eb9SPierre Jolivet         if (err >= tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[test %" PetscInt_FMT "] Error b %g\n", test, (double)err));
2489566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&b));
249c4762a1bSJed Brown       }
2509566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&b1));
2519566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&b2));
2529566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&x));
2539566063dSJacob Faibussowitsch       PetscCall(MatConvert(S, MATDENSE, MAT_INITIAL_MATRIX, &B));
254c4762a1bSJed Brown 
2559566063dSJacob Faibussowitsch       PetscCall(MatCreateTranspose(S, &ST));
2569566063dSJacob Faibussowitsch       PetscCall(MatComputeOperator(ST, MATDENSE, &BT));
2579566063dSJacob Faibussowitsch       PetscCall(MatTranspose(BT, MAT_INITIAL_MATRIX, &BTT));
2589566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)B, "S"));
2599566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)BTT, "STT"));
2609566063dSJacob Faibussowitsch       PetscCall(MatConvert(A, MATDENSE, MAT_INITIAL_MATRIX, &C));
2619566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)C, "A"));
262c4762a1bSJed Brown 
2639566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(C, NULL, "-view_mat"));
2649566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(B, NULL, "-view_mat"));
2659566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(BTT, NULL, "-view_mat"));
266c4762a1bSJed Brown 
2679566063dSJacob Faibussowitsch       PetscCall(MatAXPY(C, -1.0, B, SAME_NONZERO_PATTERN));
2689566063dSJacob Faibussowitsch       PetscCall(MatNorm(C, NORM_FROBENIUS, &err));
26948a46eb9SPierre Jolivet       if (err >= tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[test %" PetscInt_FMT "] Error mat mult after %s %g\n", test, ronl ? "MatZeroRows" : "MatZeroRowsColumns", (double)err));
270c4762a1bSJed Brown 
2719566063dSJacob Faibussowitsch       PetscCall(MatConvert(A, MATDENSE, MAT_REUSE_MATRIX, &C));
2729566063dSJacob Faibussowitsch       PetscCall(MatAXPY(C, -1.0, BTT, SAME_NONZERO_PATTERN));
2739566063dSJacob Faibussowitsch       PetscCall(MatNorm(C, NORM_FROBENIUS, &err));
27448a46eb9SPierre Jolivet       if (err >= tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[test %" PetscInt_FMT "] Error mat mult transpose after %s %g\n", test, ronl ? "MatZeroRows" : "MatZeroRowsColumns", (double)err));
275c4762a1bSJed Brown 
2769566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&ST));
2779566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&BTT));
2789566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&BT));
2799566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&B));
2809566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&C));
281c4762a1bSJed Brown     }
282c4762a1bSJed Brown     if (testdiagscale) { /* MatDiagonalScale() */
283c4762a1bSJed Brown       Vec vr, vl;
284c4762a1bSJed Brown 
2859566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(A, &vr, &vl));
286c4762a1bSJed Brown       if (randomize) {
2879566063dSJacob Faibussowitsch         PetscCall(VecSetRandom(vr, NULL));
2889566063dSJacob Faibussowitsch         PetscCall(VecSetRandom(vl, NULL));
289c4762a1bSJed Brown       } else {
2909566063dSJacob Faibussowitsch         PetscCall(VecSet(vr, test % 2 ? 0.15 : 1.0 / 0.15));
2919566063dSJacob Faibussowitsch         PetscCall(VecSet(vl, test % 2 ? -1.2 : 1.0 / -1.2));
292c4762a1bSJed Brown       }
2939566063dSJacob Faibussowitsch       PetscCall(MatDiagonalScale(A, vl, vr));
2949566063dSJacob Faibussowitsch       PetscCall(MatDiagonalScale(S, vl, vr));
2959566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&vr));
2969566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&vl));
297c4762a1bSJed Brown     }
298c4762a1bSJed Brown 
299c4762a1bSJed Brown     if (testscale) { /* MatScale() */
3009566063dSJacob Faibussowitsch       PetscCall(MatScale(A, test % 2 ? 1.4 : 1.0 / 1.4));
3019566063dSJacob Faibussowitsch       PetscCall(MatScale(S, test % 2 ? 1.4 : 1.0 / 1.4));
302c4762a1bSJed Brown     }
303c4762a1bSJed Brown 
304c4762a1bSJed Brown     if (testshift && cong) { /* MatShift() : MATSHELL shift is broken when row/cols layout are not congruent and left/right scaling have been applied */
3059566063dSJacob Faibussowitsch       PetscCall(MatShift(A, test % 2 ? -77.5 : 77.5));
3069566063dSJacob Faibussowitsch       PetscCall(MatShift(S, test % 2 ? -77.5 : 77.5));
307c4762a1bSJed Brown     }
308c4762a1bSJed Brown 
309c4762a1bSJed Brown     if (testgetdiag && cong) { /* MatGetDiagonal() */
310c4762a1bSJed Brown       Vec dA, dS;
311c4762a1bSJed Brown 
3129566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(A, &dA, NULL));
3139566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(S, &dS, NULL));
3149566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(A, dA));
3159566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(S, dS));
3169566063dSJacob Faibussowitsch       PetscCall(VecAXPY(dA, -1.0, dS));
3179566063dSJacob Faibussowitsch       PetscCall(VecNorm(dA, NORM_INFINITY, &err));
31848a46eb9SPierre Jolivet       if (err >= tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[test %" PetscInt_FMT "] Error diag %g\n", test, (double)err));
3199566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&dA));
3209566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&dS));
321c4762a1bSJed Brown     }
322c4762a1bSJed Brown 
323c4762a1bSJed Brown     if (testdup && !test) {
324c4762a1bSJed Brown       Mat A2, S2;
325c4762a1bSJed Brown 
3269566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2));
3279566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(S, MAT_COPY_VALUES, &S2));
3289566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A));
3299566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&S));
330c4762a1bSJed Brown       A = A2;
331c4762a1bSJed Brown       S = S2;
332c4762a1bSJed Brown     }
333c4762a1bSJed Brown 
334c4762a1bSJed Brown     if (testsubmat) {
335c4762a1bSJed Brown       Mat      sA, sS, dA, dS, At, St;
336c4762a1bSJed Brown       IS       r, c;
337c4762a1bSJed Brown       PetscInt rst, ren, cst, cen;
338c4762a1bSJed Brown 
3399566063dSJacob Faibussowitsch       PetscCall(MatGetOwnershipRange(A, &rst, &ren));
3409566063dSJacob Faibussowitsch       PetscCall(MatGetOwnershipRangeColumn(A, &cst, &cen));
3419566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PetscObjectComm((PetscObject)A), (ren - rst) / 2, rst, 1, &r));
3429566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PetscObjectComm((PetscObject)A), (cen - cst) / 2, cst, 1, &c));
3439566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(A, r, c, MAT_INITIAL_MATRIX, &sA));
3449566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(S, r, c, MAT_INITIAL_MATRIX, &sS));
3459566063dSJacob Faibussowitsch       PetscCall(MatMultAddEqual(sA, sS, 10, &flg));
34648a46eb9SPierre Jolivet       if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[test %" PetscInt_FMT "] Error submatrix mult add\n", test));
3479566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAddEqual(sA, sS, 10, &flg));
34848a46eb9SPierre Jolivet       if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[test %" PetscInt_FMT "] Error submatrix mult add (T)\n", test));
3499566063dSJacob Faibussowitsch       PetscCall(MatConvert(sA, MATDENSE, MAT_INITIAL_MATRIX, &dA));
3509566063dSJacob Faibussowitsch       PetscCall(MatConvert(sS, MATDENSE, MAT_INITIAL_MATRIX, &dS));
3519566063dSJacob Faibussowitsch       PetscCall(MatAXPY(dA, -1.0, dS, SAME_NONZERO_PATTERN));
3529566063dSJacob Faibussowitsch       PetscCall(MatNorm(dA, NORM_FROBENIUS, &err));
35348a46eb9SPierre Jolivet       if (err >= tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[test %" PetscInt_FMT "] Error mat submatrix %g\n", test, (double)err));
3549566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&sA));
3559566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&sS));
3569566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&dA));
3579566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&dS));
3589566063dSJacob Faibussowitsch       PetscCall(MatCreateTranspose(A, &At));
3599566063dSJacob Faibussowitsch       PetscCall(MatCreateTranspose(S, &St));
3609566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(At, c, r, MAT_INITIAL_MATRIX, &sA));
3619566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(St, c, r, MAT_INITIAL_MATRIX, &sS));
3629566063dSJacob Faibussowitsch       PetscCall(MatMultAddEqual(sA, sS, 10, &flg));
36348a46eb9SPierre Jolivet       if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[test %" PetscInt_FMT "] Error submatrix (T) mult add\n", test));
3649566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAddEqual(sA, sS, 10, &flg));
36548a46eb9SPierre Jolivet       if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[test %" PetscInt_FMT "] Error submatrix (T) mult add (T)\n", test));
3669566063dSJacob Faibussowitsch       PetscCall(MatConvert(sA, MATDENSE, MAT_INITIAL_MATRIX, &dA));
3679566063dSJacob Faibussowitsch       PetscCall(MatConvert(sS, MATDENSE, MAT_INITIAL_MATRIX, &dS));
3689566063dSJacob Faibussowitsch       PetscCall(MatAXPY(dA, -1.0, dS, SAME_NONZERO_PATTERN));
3699566063dSJacob Faibussowitsch       PetscCall(MatNorm(dA, NORM_FROBENIUS, &err));
37048a46eb9SPierre Jolivet       if (err >= tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[test %" PetscInt_FMT "] Error mat submatrix (T) %g\n", test, (double)err));
3719566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&sA));
3729566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&sS));
3739566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&dA));
3749566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&dS));
3759566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&At));
3769566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&St));
3779566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&r));
3789566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&c));
379c4762a1bSJed Brown     }
380c4762a1bSJed Brown 
381c4762a1bSJed Brown     if (testaxpy) {
382c4762a1bSJed Brown       Mat          tA, tS, dA, dS;
383c4762a1bSJed Brown       MatStructure str[3] = {SAME_NONZERO_PATTERN, SUBSET_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN};
384c4762a1bSJed Brown 
3859566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &tA));
386c4762a1bSJed Brown       if (testaxpyd && !(test % 2)) {
3879566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)tA));
388c4762a1bSJed Brown         tS = tA;
389c4762a1bSJed Brown       } else {
3909566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)S));
391c4762a1bSJed Brown         tS = S;
392c4762a1bSJed Brown       }
3939566063dSJacob Faibussowitsch       PetscCall(MatAXPY(A, 0.5, tA, str[test % 3]));
3949566063dSJacob Faibussowitsch       PetscCall(MatAXPY(S, 0.5, tS, str[test % 3]));
395c4762a1bSJed Brown       /* this will trigger an error the next MatMult or MatMultTranspose call for S */
3969566063dSJacob Faibussowitsch       if (testaxpyerr) PetscCall(MatScale(tA, 0));
3979566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&tA));
3989566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&tS));
3999566063dSJacob Faibussowitsch       PetscCall(MatMultAddEqual(A, S, 10, &flg));
40048a46eb9SPierre Jolivet       if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[test %" PetscInt_FMT "] Error axpy mult add\n", test));
4019566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAddEqual(A, S, 10, &flg));
40248a46eb9SPierre Jolivet       if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[test %" PetscInt_FMT "] Error axpy mult add (T)\n", test));
4039566063dSJacob Faibussowitsch       PetscCall(MatConvert(A, MATDENSE, MAT_INITIAL_MATRIX, &dA));
4049566063dSJacob Faibussowitsch       PetscCall(MatConvert(S, MATDENSE, MAT_INITIAL_MATRIX, &dS));
4059566063dSJacob Faibussowitsch       PetscCall(MatAXPY(dA, -1.0, dS, SAME_NONZERO_PATTERN));
4069566063dSJacob Faibussowitsch       PetscCall(MatNorm(dA, NORM_FROBENIUS, &err));
40748a46eb9SPierre Jolivet       if (err >= tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[test %" PetscInt_FMT "] Error mat submatrix %g\n", test, (double)err));
4089566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&dA));
4099566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&dS));
410c4762a1bSJed Brown     }
411c4762a1bSJed Brown 
412c4762a1bSJed Brown     if (testreset && (ntest == 1 || test == ntest - 2)) {
413c4762a1bSJed Brown       /* reset MATSHELL */
4149566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY));
4159566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY));
416c4762a1bSJed Brown       /* reset A */
4179566063dSJacob Faibussowitsch       PetscCall(MatCopy(user->B, A, DIFFERENT_NONZERO_PATTERN));
418c4762a1bSJed Brown     }
419c4762a1bSJed Brown   }
420c4762a1bSJed Brown 
4219566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
4229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&S));
4239566063dSJacob Faibussowitsch   PetscCall(PetscFree(user));
4249566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
425b122ec5aSJacob Faibussowitsch   return 0;
426c4762a1bSJed Brown }
427c4762a1bSJed Brown 
428c4762a1bSJed Brown /*TEST
429c4762a1bSJed Brown 
430c4762a1bSJed Brown    testset:
431c4762a1bSJed Brown      suffix: rect
432c4762a1bSJed Brown      requires: !single
433c4762a1bSJed Brown      output_file: output/ex221_1.out
434c4762a1bSJed Brown      nsize: {{1 3}}
435c4762a1bSJed Brown      args: -loop 3 -keep {{0 1}} -M {{12 19}} -N {{19 12}} -submat {{0 1}} -test_axpy_different {{0 1}}
436c4762a1bSJed Brown 
437c4762a1bSJed Brown    testset:
438c4762a1bSJed Brown      suffix: square
439c4762a1bSJed Brown      requires: !single
440c4762a1bSJed Brown      output_file: output/ex221_1.out
441c4762a1bSJed Brown      nsize: {{1 3}}
442c4762a1bSJed Brown      args: -M 21 -N 21 -loop 4 -rows_only {{0 1}} -keep {{0 1}} -submat {{0 1}} -test_axpy_different {{0 1}}
443c4762a1bSJed Brown TEST*/
444