xref: /petsc/src/mat/tests/ex23.c (revision a50ef18c006c8158869f1d5a0ae3686c1cb9c29f)
1d0dbe9f7SStefano Zampini static char help[] = "Tests the use of interface functions for MATIS matrices and conversion routines.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown PetscErrorCode TestMatZeroRows(Mat, Mat, PetscBool, IS, PetscScalar);
6c4762a1bSJed Brown PetscErrorCode CheckMat(Mat, Mat, PetscBool, const char *);
7c4762a1bSJed Brown 
8d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
9d71ae5a4SJacob Faibussowitsch {
10c4762a1bSJed Brown   Mat                    A, B, A2, B2, T;
11c4762a1bSJed Brown   Mat                    Aee, Aeo, Aoe, Aoo;
12d0dbe9f7SStefano Zampini   Mat                   *mats, *Asub, *Bsub;
13c4762a1bSJed Brown   Vec                    x, y;
14c4762a1bSJed Brown   MatInfo                info;
15c4762a1bSJed Brown   ISLocalToGlobalMapping cmap, rmap;
16c4762a1bSJed Brown   IS                     is, is2, reven, rodd, ceven, codd;
17c4762a1bSJed Brown   IS                    *rows, *cols;
18d0dbe9f7SStefano Zampini   IS                     irow[2], icol[2];
19d0dbe9f7SStefano Zampini   PetscLayout            rlayout, clayout;
20d0dbe9f7SStefano Zampini   const PetscInt        *rrange, *crange;
21c4762a1bSJed Brown   MatType                lmtype;
22c4762a1bSJed Brown   PetscScalar            diag = 2.;
23e432b41dSStefano Zampini   PetscInt               n, m, i, lm, ln;
24*a50ef18cSStefano Zampini   PetscInt               rst, ren, cst, cen, nr, nc, rbs = 1, cbs = 1;
25d0dbe9f7SStefano Zampini   PetscMPIInt            rank, size, lrank, rrank;
26c4762a1bSJed Brown   PetscBool              testT, squaretest, isaij;
274f58015eSStefano Zampini   PetscBool              permute = PETSC_FALSE, negmap = PETSC_FALSE, repmap = PETSC_FALSE, allow_repeated = PETSC_TRUE;
285042aa92SStefano Zampini   PetscBool              diffmap = PETSC_TRUE, symmetric = PETSC_FALSE, issymmetric, test_matlab = PETSC_FALSE;
29c4762a1bSJed Brown 
30327415f7SBarry Smith   PetscFunctionBeginUser;
319566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
339566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
34c4762a1bSJed Brown   m = n = 2 * size;
359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-symmetric", &symmetric, NULL));
369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-negmap", &negmap, NULL));
399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-repmap", &repmap, NULL));
409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-permmap", &permute, NULL));
419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-diffmap", &diffmap, NULL));
424f58015eSStefano Zampini   PetscCall(PetscOptionsGetBool(NULL, NULL, "-allow_repeated", &allow_repeated, NULL));
435042aa92SStefano Zampini   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_matlab", &test_matlab, NULL));
44*a50ef18cSStefano Zampini   PetscCall(PetscOptionsGetInt(NULL, NULL, "-rbs", &rbs, NULL));
45*a50ef18cSStefano Zampini   PetscCall(PetscOptionsGetInt(NULL, NULL, "-cbs", &cbs, NULL));
46e00437b9SBarry Smith   PetscCheck(size == 1 || m >= 4, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Number of rows should be larger or equal 4 for parallel runs");
47e00437b9SBarry Smith   PetscCheck(size != 1 || m >= 2, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Number of rows should be larger or equal 2 for uniprocessor runs");
48e00437b9SBarry Smith   PetscCheck(n >= 2, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Number of cols should be larger or equal 2");
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* create a MATIS matrix */
519566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
529566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m, n));
539566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATIS));
549566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
55e432b41dSStefano Zampini   if (!negmap && !repmap) {
56fc989267SStefano Zampini     /* This is not the proper setting for MATIS for finite elements, it is just used to test the routines
57fc989267SStefano Zampini        Here we use a one-to-one correspondence between local row/column spaces and global row/column spaces
58fc989267SStefano Zampini        Equivalent to passing NULL for the mapping */
599566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_WORLD, n, 0, 1, &is));
60e432b41dSStefano Zampini   } else if (negmap && !repmap) { /* non repeated but with negative indices */
619566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_WORLD, n + 2, -2, 1, &is));
62e432b41dSStefano Zampini   } else if (!negmap && repmap) { /* non negative but repeated indices */
63e432b41dSStefano Zampini     IS isl[2];
64e432b41dSStefano Zampini 
659566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_WORLD, n, 0, 1, &isl[0]));
669566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_WORLD, n, n - 1, -1, &isl[1]));
679566063dSJacob Faibussowitsch     PetscCall(ISConcatenate(PETSC_COMM_WORLD, 2, isl, &is));
689566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isl[0]));
699566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isl[1]));
70e432b41dSStefano Zampini   } else { /* negative and repeated indices */
71e432b41dSStefano Zampini     IS isl[2];
72e432b41dSStefano Zampini 
739566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_WORLD, n + 1, -1, 1, &isl[0]));
749566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_WORLD, n + 1, n - 1, -1, &isl[1]));
759566063dSJacob Faibussowitsch     PetscCall(ISConcatenate(PETSC_COMM_WORLD, 2, isl, &is));
769566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isl[0]));
779566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isl[1]));
78e432b41dSStefano Zampini   }
799566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreateIS(is, &cmap));
809566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
81c4762a1bSJed Brown 
82e432b41dSStefano Zampini   if (m != n || diffmap) {
839566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_WORLD, m, permute ? m - 1 : 0, permute ? -1 : 1, &is));
849566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreateIS(is, &rmap));
859566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is));
86c4762a1bSJed Brown   } else {
879566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)cmap));
88c4762a1bSJed Brown     rmap = cmap;
89c4762a1bSJed Brown   }
90e432b41dSStefano Zampini 
914f58015eSStefano Zampini   PetscCall(MatISSetAllowRepeated(A, allow_repeated));
929566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(A, rmap, cmap));
93*a50ef18cSStefano Zampini   PetscCall(MatSetBlockSizes(A, rbs, cbs));
949566063dSJacob Faibussowitsch   PetscCall(MatISStoreL2L(A, PETSC_FALSE));
959566063dSJacob Faibussowitsch   PetscCall(MatISSetPreallocation(A, 3, NULL, 3, NULL));
969566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, (PetscBool) !(repmap || negmap))); /* I do not want to precompute the pattern */
979566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetSize(rmap, &lm));
989566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetSize(cmap, &ln));
99e432b41dSStefano Zampini   for (i = 0; i < lm; i++) {
100c4762a1bSJed Brown     PetscScalar v[3];
101c4762a1bSJed Brown     PetscInt    cols[3];
102c4762a1bSJed Brown 
103c4762a1bSJed Brown     cols[0] = (i - 1 + n) % n;
104c4762a1bSJed Brown     cols[1] = i % n;
105c4762a1bSJed Brown     cols[2] = (i + 1) % n;
106c4762a1bSJed Brown     v[0]    = -1. * (symmetric ? PetscMin(i + 1, cols[0] + 1) : i + 1);
107c4762a1bSJed Brown     v[1]    = 2. * (symmetric ? PetscMin(i + 1, cols[1] + 1) : i + 1);
108c4762a1bSJed Brown     v[2]    = -1. * (symmetric ? PetscMin(i + 1, cols[2] + 1) : i + 1);
1099566063dSJacob Faibussowitsch     PetscCall(ISGlobalToLocalMappingApply(cmap, IS_GTOLM_MASK, 3, cols, NULL, cols));
1109566063dSJacob Faibussowitsch     PetscCall(MatSetValuesLocal(A, 1, &i, 3, cols, v, ADD_VALUES));
111c4762a1bSJed Brown   }
1129566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1139566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
114c4762a1bSJed Brown 
115e432b41dSStefano Zampini   /* activate tests for square matrices with same maps only */
1169566063dSJacob Faibussowitsch   PetscCall(MatHasCongruentLayouts(A, &squaretest));
117e432b41dSStefano Zampini   if (squaretest && rmap != cmap) {
118e432b41dSStefano Zampini     PetscInt nr, nc;
119e432b41dSStefano Zampini 
1209566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(rmap, &nr));
1219566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(cmap, &nc));
122e432b41dSStefano Zampini     if (nr != nc) squaretest = PETSC_FALSE;
123e432b41dSStefano Zampini     else {
124e432b41dSStefano Zampini       const PetscInt *idxs1, *idxs2;
125e432b41dSStefano Zampini 
1269566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingGetIndices(rmap, &idxs1));
1279566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingGetIndices(cmap, &idxs2));
1289566063dSJacob Faibussowitsch       PetscCall(PetscArraycmp(idxs1, idxs2, nr, &squaretest));
1299566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingRestoreIndices(rmap, &idxs1));
1309566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingRestoreIndices(cmap, &idxs2));
131e432b41dSStefano Zampini     }
1321c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &squaretest, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A)));
133e432b41dSStefano Zampini   }
1344f58015eSStefano Zampini   if (negmap && repmap) squaretest = PETSC_FALSE;
135e432b41dSStefano Zampini 
136e432b41dSStefano Zampini   /* test MatISGetLocalMat */
1379566063dSJacob Faibussowitsch   PetscCall(MatISGetLocalMat(A, &B));
1389566063dSJacob Faibussowitsch   PetscCall(MatGetType(B, &lmtype));
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* test MatGetInfo */
1419566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatGetInfo\n"));
1429566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(A, MAT_LOCAL, &info));
1439566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
1449371c9d4SSatish Balay   PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "Process  %2d: %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", PetscGlobalRank, (PetscInt)info.nz_used, (PetscInt)info.nz_allocated,
1459371c9d4SSatish Balay                                                (PetscInt)info.nz_unneeded, (PetscInt)info.assemblies, (PetscInt)info.mallocs));
1469566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
1479566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(A, MAT_GLOBAL_MAX, &info));
1489371c9d4SSatish Balay   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "GlobalMax  : %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", (PetscInt)info.nz_used, (PetscInt)info.nz_allocated, (PetscInt)info.nz_unneeded,
1499371c9d4SSatish Balay                                    (PetscInt)info.assemblies, (PetscInt)info.mallocs));
1509566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(A, MAT_GLOBAL_SUM, &info));
1519371c9d4SSatish Balay   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "GlobalSum  : %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", (PetscInt)info.nz_used, (PetscInt)info.nz_allocated, (PetscInt)info.nz_unneeded,
1529371c9d4SSatish Balay                                    (PetscInt)info.assemblies, (PetscInt)info.mallocs));
153c4762a1bSJed Brown 
154e432b41dSStefano Zampini   /* test MatIsSymmetric */
1559566063dSJacob Faibussowitsch   PetscCall(MatIsSymmetric(A, 0.0, &issymmetric));
1569566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatIsSymmetric: %d\n", issymmetric));
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   /* Create a MPIAIJ matrix, same as A */
1599566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
1609566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n));
161*a50ef18cSStefano Zampini   PetscCall(MatSetBlockSizes(B, rbs, cbs));
1629566063dSJacob Faibussowitsch   PetscCall(MatSetType(B, MATAIJ));
1639566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(B));
1649566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(B, rmap, cmap));
1659566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(B, 3, NULL, 3, NULL));
1669566063dSJacob Faibussowitsch   PetscCall(MatMPIBAIJSetPreallocation(B, 1, 3, NULL, 3, NULL));
167c4762a1bSJed Brown #if defined(PETSC_HAVE_HYPRE)
1689566063dSJacob Faibussowitsch   PetscCall(MatHYPRESetPreallocation(B, 3, NULL, 3, NULL));
169c4762a1bSJed Brown #endif
1709566063dSJacob Faibussowitsch   PetscCall(MatISSetPreallocation(B, 3, NULL, 3, NULL));
1719566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, (PetscBool) !(repmap || negmap))); /* I do not want to precompute the pattern */
172e432b41dSStefano Zampini   for (i = 0; i < lm; i++) {
173c4762a1bSJed Brown     PetscScalar v[3];
174c4762a1bSJed Brown     PetscInt    cols[3];
175c4762a1bSJed Brown 
176c4762a1bSJed Brown     cols[0] = (i - 1 + n) % n;
177c4762a1bSJed Brown     cols[1] = i % n;
178c4762a1bSJed Brown     cols[2] = (i + 1) % n;
179c4762a1bSJed Brown     v[0]    = -1. * (symmetric ? PetscMin(i + 1, cols[0] + 1) : i + 1);
180c4762a1bSJed Brown     v[1]    = 2. * (symmetric ? PetscMin(i + 1, cols[1] + 1) : i + 1);
181c4762a1bSJed Brown     v[2]    = -1. * (symmetric ? PetscMin(i + 1, cols[2] + 1) : i + 1);
1829566063dSJacob Faibussowitsch     PetscCall(ISGlobalToLocalMappingApply(cmap, IS_GTOLM_MASK, 3, cols, NULL, cols));
1839566063dSJacob Faibussowitsch     PetscCall(MatSetValuesLocal(B, 1, &i, 3, cols, v, ADD_VALUES));
184c4762a1bSJed Brown   }
1859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1869566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
187e432b41dSStefano Zampini 
188e432b41dSStefano Zampini   /* test MatView */
1899566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatView\n"));
1909566063dSJacob Faibussowitsch   PetscCall(MatView(A, NULL));
1919566063dSJacob Faibussowitsch   PetscCall(MatView(B, NULL));
192c4762a1bSJed Brown 
1935042aa92SStefano Zampini   /* test MATLAB ASCII view */
1945042aa92SStefano Zampini   if (test_matlab) { /* output is different when using real or complex numbers */
1955042aa92SStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatView ASCII MATLAB\n"));
1965042aa92SStefano Zampini     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_MATLAB));
1975042aa92SStefano Zampini     PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
1985042aa92SStefano Zampini     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
1995042aa92SStefano Zampini   }
2005042aa92SStefano Zampini 
201c4762a1bSJed Brown   /* test CheckMat */
2029566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test CheckMat\n"));
2039566063dSJacob Faibussowitsch   PetscCall(CheckMat(A, B, PETSC_FALSE, "CheckMat"));
204c4762a1bSJed Brown 
2055042aa92SStefano Zampini   /* test binary MatView/MatLoad */
2065042aa92SStefano Zampini   {
2075042aa92SStefano Zampini     PetscMPIInt color = rank % 2;
2085042aa92SStefano Zampini     MPI_Comm    comm;
2095042aa92SStefano Zampini     char        name[PETSC_MAX_PATH_LEN];
2105042aa92SStefano Zampini     PetscViewer wview, cview, sview, view;
2115042aa92SStefano Zampini     Mat         A2;
2125042aa92SStefano Zampini 
2135042aa92SStefano Zampini     PetscCallMPI(MPI_Comm_split(PETSC_COMM_WORLD, color, rank, &comm));
2145042aa92SStefano Zampini 
2155042aa92SStefano Zampini     PetscCall(PetscSNPrintf(name, PETSC_STATIC_ARRAY_LENGTH(name), "world_is"));
2165042aa92SStefano Zampini     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, name, FILE_MODE_WRITE, &wview));
2175042aa92SStefano Zampini     PetscCall(PetscSNPrintf(name, PETSC_STATIC_ARRAY_LENGTH(name), "seq_is_%d", rank));
2185042aa92SStefano Zampini     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, name, FILE_MODE_WRITE, &sview));
2195042aa92SStefano Zampini     PetscCall(PetscSNPrintf(name, PETSC_STATIC_ARRAY_LENGTH(name), "color_is_%d", color));
2205042aa92SStefano Zampini     PetscCall(PetscViewerBinaryOpen(comm, name, FILE_MODE_WRITE, &cview));
2215042aa92SStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatView on binary world\n"));
2225042aa92SStefano Zampini     PetscCall(MatView(A, wview));
2235042aa92SStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatView on binary self\n"));
2245042aa92SStefano Zampini     PetscCall(MatView(A, sview));
2255042aa92SStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatView on binary subcomm\n"));
2265042aa92SStefano Zampini     PetscCall(MatView(A, cview));
2275042aa92SStefano Zampini     PetscCall(PetscViewerDestroy(&wview));
2285042aa92SStefano Zampini     PetscCall(PetscViewerDestroy(&cview));
2295042aa92SStefano Zampini     PetscCall(PetscViewerDestroy(&sview));
2305042aa92SStefano Zampini 
2315042aa92SStefano Zampini     /* Load a world matrix */
2325042aa92SStefano Zampini     PetscCall(MatCreate(PETSC_COMM_WORLD, &A2));
2335042aa92SStefano Zampini     PetscCall(MatSetType(A2, MATIS));
2345042aa92SStefano Zampini     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL));
2355042aa92SStefano Zampini 
2365042aa92SStefano Zampini     /* Read back the same matrix and check */
2375042aa92SStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatLoad from world\n"));
2385042aa92SStefano Zampini     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "world_is", FILE_MODE_READ, &view));
2395042aa92SStefano Zampini     PetscCall(MatLoad(A2, view));
2405042aa92SStefano Zampini     PetscCall(CheckMat(A, A2, PETSC_TRUE, "Load"));
2415042aa92SStefano Zampini     PetscCall(MatView(A2, PETSC_VIEWER_STDOUT_WORLD));
2425042aa92SStefano Zampini     PetscCall(PetscViewerDestroy(&view));
2435042aa92SStefano Zampini 
2445042aa92SStefano Zampini     /* Read the matrix from rank 0 only */
2455042aa92SStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatLoad from self\n"));
2465042aa92SStefano Zampini     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "seq_is_0", FILE_MODE_READ, &view));
2475042aa92SStefano Zampini     PetscCall(MatLoad(A2, view));
2485042aa92SStefano Zampini     PetscCall(MatView(A2, PETSC_VIEWER_STDOUT_WORLD));
2495042aa92SStefano Zampini     PetscCall(PetscViewerDestroy(&view));
2505042aa92SStefano Zampini 
2515042aa92SStefano Zampini     /* Read the matrix from subcomm */
2525042aa92SStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatLoad from subcomm\n"));
2535042aa92SStefano Zampini     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "color_is_0", FILE_MODE_READ, &view));
2545042aa92SStefano Zampini     PetscCall(MatLoad(A2, view));
2555042aa92SStefano Zampini     PetscCall(MatView(A2, PETSC_VIEWER_STDOUT_WORLD));
2565042aa92SStefano Zampini     PetscCall(PetscViewerDestroy(&view));
2575042aa92SStefano Zampini 
2585042aa92SStefano Zampini     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
2595042aa92SStefano Zampini     PetscCall(MatDestroy(&A2));
2605042aa92SStefano Zampini 
2615042aa92SStefano Zampini     /* now load the original matrix from color 0 only processes */
2625042aa92SStefano Zampini     if (!color) {
2635042aa92SStefano Zampini       PetscCall(PetscPrintf(comm, "Test subcomm MatLoad from world\n"));
2645042aa92SStefano Zampini       PetscCall(MatCreate(comm, &A2));
2655042aa92SStefano Zampini       PetscCall(MatSetType(A2, MATIS));
2665042aa92SStefano Zampini       PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(comm), PETSC_VIEWER_ASCII_INFO_DETAIL));
2675042aa92SStefano Zampini       PetscCall(PetscViewerBinaryOpen(comm, "world_is", FILE_MODE_READ, &view));
2685042aa92SStefano Zampini       PetscCall(MatLoad(A2, view));
2695042aa92SStefano Zampini       PetscCall(MatView(A2, PETSC_VIEWER_STDOUT_(comm)));
2705042aa92SStefano Zampini       PetscCall(PetscViewerDestroy(&view));
2715042aa92SStefano Zampini       PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(comm)));
2725042aa92SStefano Zampini       PetscCall(MatDestroy(&A2));
2735042aa92SStefano Zampini     }
2745042aa92SStefano Zampini 
2755042aa92SStefano Zampini     PetscCallMPI(MPI_Comm_free(&comm));
2765042aa92SStefano Zampini   }
2775042aa92SStefano Zampini 
278c4762a1bSJed Brown   /* test MatDuplicate and MatAXPY */
2799566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatDuplicate and MatAXPY\n"));
2809566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2));
2819566063dSJacob Faibussowitsch   PetscCall(CheckMat(A, A2, PETSC_FALSE, "MatDuplicate and MatAXPY"));
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   /* test MatConvert */
2849566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatConvert_IS_XAIJ\n"));
2859566063dSJacob Faibussowitsch   PetscCall(MatConvert(A2, MATAIJ, MAT_INITIAL_MATRIX, &B2));
2869566063dSJacob Faibussowitsch   PetscCall(CheckMat(B, B2, PETSC_TRUE, "MatConvert_IS_XAIJ MAT_INITIAL_MATRIX"));
2879566063dSJacob Faibussowitsch   PetscCall(MatConvert(A2, MATAIJ, MAT_REUSE_MATRIX, &B2));
2889566063dSJacob Faibussowitsch   PetscCall(CheckMat(B, B2, PETSC_TRUE, "MatConvert_IS_XAIJ MAT_REUSE_MATRIX"));
2899566063dSJacob Faibussowitsch   PetscCall(MatConvert(A2, MATAIJ, MAT_INPLACE_MATRIX, &A2));
2909566063dSJacob Faibussowitsch   PetscCall(CheckMat(B, A2, PETSC_TRUE, "MatConvert_IS_XAIJ MAT_INPLACE_MATRIX"));
2919566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A2));
2929566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B2));
2939566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatConvert_XAIJ_IS\n"));
2949566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &B2));
2959566063dSJacob Faibussowitsch   PetscCall(MatConvert(B2, MATIS, MAT_INITIAL_MATRIX, &A2));
2969566063dSJacob Faibussowitsch   PetscCall(CheckMat(A, A2, PETSC_TRUE, "MatConvert_XAIJ_IS MAT_INITIAL_MATRIX"));
2979566063dSJacob Faibussowitsch   PetscCall(MatConvert(B2, MATIS, MAT_REUSE_MATRIX, &A2));
2989566063dSJacob Faibussowitsch   PetscCall(CheckMat(A, A2, PETSC_TRUE, "MatConvert_XAIJ_IS MAT_REUSE_MATRIX"));
2999566063dSJacob Faibussowitsch   PetscCall(MatConvert(B2, MATIS, MAT_INPLACE_MATRIX, &B2));
3009566063dSJacob Faibussowitsch   PetscCall(CheckMat(A, B2, PETSC_TRUE, "MatConvert_XAIJ_IS MAT_INPLACE_MATRIX"));
3019566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A2));
3029566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B2));
3039566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(lmtype, MATSEQAIJ, &isaij));
304c4762a1bSJed Brown   if (size == 1 && isaij) { /* tests special code paths in MatConvert_IS_XAIJ */
305c4762a1bSJed Brown     PetscInt               ri, ci, rr[3] = {0, 1, 0}, cr[4] = {1, 2, 0, 1}, rk[3] = {0, 2, 1}, ck[4] = {1, 0, 3, 2};
306e432b41dSStefano Zampini     ISLocalToGlobalMapping tcmap, trmap;
307c4762a1bSJed Brown 
308c4762a1bSJed Brown     for (ri = 0; ri < 2; ri++) {
309c4762a1bSJed Brown       PetscInt *r;
310c4762a1bSJed Brown 
311c4762a1bSJed Brown       r = (PetscInt *)(ri == 0 ? rr : rk);
312c4762a1bSJed Brown       for (ci = 0; ci < 2; ci++) {
313c4762a1bSJed Brown         PetscInt *c, rb, cb;
314c4762a1bSJed Brown 
315c4762a1bSJed Brown         c = (PetscInt *)(ci == 0 ? cr : ck);
316c4762a1bSJed Brown         for (rb = 1; rb < 4; rb++) {
3179566063dSJacob Faibussowitsch           PetscCall(ISCreateBlock(PETSC_COMM_SELF, rb, 3, r, PETSC_COPY_VALUES, &is));
3189566063dSJacob Faibussowitsch           PetscCall(ISLocalToGlobalMappingCreateIS(is, &trmap));
3199566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&is));
320c4762a1bSJed Brown           for (cb = 1; cb < 4; cb++) {
321c4762a1bSJed Brown             Mat  T, lT, T2;
322c4762a1bSJed Brown             char testname[256];
323c4762a1bSJed Brown 
3249566063dSJacob Faibussowitsch             PetscCall(PetscSNPrintf(testname, sizeof(testname), "MatConvert_IS_XAIJ special case (%" PetscInt_FMT " %" PetscInt_FMT ", bs %" PetscInt_FMT " %" PetscInt_FMT ")", ri, ci, rb, cb));
3259566063dSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test %s\n", testname));
326c4762a1bSJed Brown 
3279566063dSJacob Faibussowitsch             PetscCall(ISCreateBlock(PETSC_COMM_SELF, cb, 4, c, PETSC_COPY_VALUES, &is));
3289566063dSJacob Faibussowitsch             PetscCall(ISLocalToGlobalMappingCreateIS(is, &tcmap));
3299566063dSJacob Faibussowitsch             PetscCall(ISDestroy(&is));
330c4762a1bSJed Brown 
3319566063dSJacob Faibussowitsch             PetscCall(MatCreate(PETSC_COMM_SELF, &T));
3329566063dSJacob Faibussowitsch             PetscCall(MatSetSizes(T, PETSC_DECIDE, PETSC_DECIDE, rb * 3, cb * 4));
3339566063dSJacob Faibussowitsch             PetscCall(MatSetType(T, MATIS));
3349566063dSJacob Faibussowitsch             PetscCall(MatSetLocalToGlobalMapping(T, trmap, tcmap));
3359566063dSJacob Faibussowitsch             PetscCall(ISLocalToGlobalMappingDestroy(&tcmap));
3369566063dSJacob Faibussowitsch             PetscCall(MatISGetLocalMat(T, &lT));
3379566063dSJacob Faibussowitsch             PetscCall(MatSetType(lT, MATSEQAIJ));
3389566063dSJacob Faibussowitsch             PetscCall(MatSeqAIJSetPreallocation(lT, cb * 4, NULL));
3399566063dSJacob Faibussowitsch             PetscCall(MatSetRandom(lT, NULL));
3409566063dSJacob Faibussowitsch             PetscCall(MatConvert(lT, lmtype, MAT_INPLACE_MATRIX, &lT));
3419566063dSJacob Faibussowitsch             PetscCall(MatISRestoreLocalMat(T, &lT));
3429566063dSJacob Faibussowitsch             PetscCall(MatAssemblyBegin(T, MAT_FINAL_ASSEMBLY));
3439566063dSJacob Faibussowitsch             PetscCall(MatAssemblyEnd(T, MAT_FINAL_ASSEMBLY));
344c4762a1bSJed Brown 
3459566063dSJacob Faibussowitsch             PetscCall(MatConvert(T, MATAIJ, MAT_INITIAL_MATRIX, &T2));
3469566063dSJacob Faibussowitsch             PetscCall(CheckMat(T, T2, PETSC_TRUE, "MAT_INITIAL_MATRIX"));
3479566063dSJacob Faibussowitsch             PetscCall(MatConvert(T, MATAIJ, MAT_REUSE_MATRIX, &T2));
3489566063dSJacob Faibussowitsch             PetscCall(CheckMat(T, T2, PETSC_TRUE, "MAT_REUSE_MATRIX"));
3499566063dSJacob Faibussowitsch             PetscCall(MatDestroy(&T2));
3509566063dSJacob Faibussowitsch             PetscCall(MatDuplicate(T, MAT_COPY_VALUES, &T2));
3519566063dSJacob Faibussowitsch             PetscCall(MatConvert(T2, MATAIJ, MAT_INPLACE_MATRIX, &T2));
3529566063dSJacob Faibussowitsch             PetscCall(CheckMat(T, T2, PETSC_TRUE, "MAT_INPLACE_MATRIX"));
3539566063dSJacob Faibussowitsch             PetscCall(MatDestroy(&T));
3549566063dSJacob Faibussowitsch             PetscCall(MatDestroy(&T2));
355c4762a1bSJed Brown           }
3569566063dSJacob Faibussowitsch           PetscCall(ISLocalToGlobalMappingDestroy(&trmap));
357c4762a1bSJed Brown         }
358c4762a1bSJed Brown       }
359c4762a1bSJed Brown     }
360c4762a1bSJed Brown   }
361c4762a1bSJed Brown 
362c4762a1bSJed Brown   /* test MatDiagonalScale */
3639566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatDiagonalScale\n"));
3649566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2));
3659566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &B2));
3669566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &x, &y));
3679566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(x, NULL));
368e432b41dSStefano Zampini   if (issymmetric) {
3699566063dSJacob Faibussowitsch     PetscCall(VecCopy(x, y));
370c4762a1bSJed Brown   } else {
3719566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(y, NULL));
3729566063dSJacob Faibussowitsch     PetscCall(VecScale(y, 8.));
373c4762a1bSJed Brown   }
3749566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(A2, y, x));
3759566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(B2, y, x));
3769566063dSJacob Faibussowitsch   PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatDiagonalScale"));
3779566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A2));
3789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B2));
3799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
3809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
381c4762a1bSJed Brown 
382c4762a1bSJed Brown   /* test MatPtAP (A IS and B AIJ) */
383c4762a1bSJed Brown   if (isaij && m == n) {
3849566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatPtAP\n"));
3854f58015eSStefano Zampini     /* There's a bug in MatCreateSubMatrices_MPIAIJ I cannot figure out */
3864f58015eSStefano Zampini     if (!allow_repeated || !repmap || size == 1) {
3879566063dSJacob Faibussowitsch       PetscCall(MatISStoreL2L(A, PETSC_TRUE));
3889566063dSJacob Faibussowitsch       PetscCall(MatPtAP(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &A2));
3899566063dSJacob Faibussowitsch       PetscCall(MatPtAP(B, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &B2));
3909566063dSJacob Faibussowitsch       PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatPtAP MAT_INITIAL_MATRIX"));
3919566063dSJacob Faibussowitsch       PetscCall(MatPtAP(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &A2));
3929566063dSJacob Faibussowitsch       PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatPtAP MAT_REUSE_MATRIX"));
3939566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A2));
3949566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&B2));
395c4762a1bSJed Brown     }
3964f58015eSStefano Zampini   }
397c4762a1bSJed Brown 
398c4762a1bSJed Brown   /* test MatGetLocalSubMatrix */
3999566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatGetLocalSubMatrix\n"));
4009566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A, MAT_DO_NOT_COPY_VALUES, &A2));
4019566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, lm / 2 + lm % 2, 0, 2, &reven));
4029566063dSJacob Faibussowitsch   PetscCall(ISComplement(reven, 0, lm, &rodd));
4039566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, ln / 2 + ln % 2, 0, 2, &ceven));
4049566063dSJacob Faibussowitsch   PetscCall(ISComplement(ceven, 0, ln, &codd));
4059566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSubMatrix(A2, reven, ceven, &Aee));
4069566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSubMatrix(A2, reven, codd, &Aeo));
4079566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSubMatrix(A2, rodd, ceven, &Aoe));
4089566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSubMatrix(A2, rodd, codd, &Aoo));
409e432b41dSStefano Zampini   for (i = 0; i < lm; i++) {
410c4762a1bSJed Brown     PetscInt    j, je, jo, colse[3], colso[3];
411c4762a1bSJed Brown     PetscScalar ve[3], vo[3];
412c4762a1bSJed Brown     PetscScalar v[3];
413c4762a1bSJed Brown     PetscInt    cols[3];
414e432b41dSStefano Zampini     PetscInt    row = i / 2;
415c4762a1bSJed Brown 
416c4762a1bSJed Brown     cols[0] = (i - 1 + n) % n;
417c4762a1bSJed Brown     cols[1] = i % n;
418c4762a1bSJed Brown     cols[2] = (i + 1) % n;
419c4762a1bSJed Brown     v[0]    = -1. * (symmetric ? PetscMin(i + 1, cols[0] + 1) : i + 1);
420c4762a1bSJed Brown     v[1]    = 2. * (symmetric ? PetscMin(i + 1, cols[1] + 1) : i + 1);
421c4762a1bSJed Brown     v[2]    = -1. * (symmetric ? PetscMin(i + 1, cols[2] + 1) : i + 1);
4229566063dSJacob Faibussowitsch     PetscCall(ISGlobalToLocalMappingApply(cmap, IS_GTOLM_MASK, 3, cols, NULL, cols));
423c4762a1bSJed Brown     for (j = 0, je = 0, jo = 0; j < 3; j++) {
424c4762a1bSJed Brown       if (cols[j] % 2) {
425c4762a1bSJed Brown         vo[jo]      = v[j];
426c4762a1bSJed Brown         colso[jo++] = cols[j] / 2;
427c4762a1bSJed Brown       } else {
428c4762a1bSJed Brown         ve[je]      = v[j];
429c4762a1bSJed Brown         colse[je++] = cols[j] / 2;
430c4762a1bSJed Brown       }
431c4762a1bSJed Brown     }
432c4762a1bSJed Brown     if (i % 2) {
4339566063dSJacob Faibussowitsch       PetscCall(MatSetValuesLocal(Aoe, 1, &row, je, colse, ve, ADD_VALUES));
4349566063dSJacob Faibussowitsch       PetscCall(MatSetValuesLocal(Aoo, 1, &row, jo, colso, vo, ADD_VALUES));
435c4762a1bSJed Brown     } else {
4369566063dSJacob Faibussowitsch       PetscCall(MatSetValuesLocal(Aee, 1, &row, je, colse, ve, ADD_VALUES));
4379566063dSJacob Faibussowitsch       PetscCall(MatSetValuesLocal(Aeo, 1, &row, jo, colso, vo, ADD_VALUES));
438c4762a1bSJed Brown     }
439c4762a1bSJed Brown   }
4409566063dSJacob Faibussowitsch   PetscCall(MatRestoreLocalSubMatrix(A2, reven, ceven, &Aee));
4419566063dSJacob Faibussowitsch   PetscCall(MatRestoreLocalSubMatrix(A2, reven, codd, &Aeo));
4429566063dSJacob Faibussowitsch   PetscCall(MatRestoreLocalSubMatrix(A2, rodd, ceven, &Aoe));
4439566063dSJacob Faibussowitsch   PetscCall(MatRestoreLocalSubMatrix(A2, rodd, codd, &Aoo));
4449566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&reven));
4459566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ceven));
4469566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&rodd));
4479566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&codd));
4489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A2, MAT_FINAL_ASSEMBLY));
4499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A2, MAT_FINAL_ASSEMBLY));
4509566063dSJacob Faibussowitsch   PetscCall(MatAXPY(A2, -1., A, SAME_NONZERO_PATTERN));
4519566063dSJacob Faibussowitsch   PetscCall(CheckMat(A2, NULL, PETSC_FALSE, "MatGetLocalSubMatrix"));
4529566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A2));
453c4762a1bSJed Brown 
454c4762a1bSJed Brown   /* test MatConvert_Nest_IS */
455c4762a1bSJed Brown   testT = PETSC_FALSE;
4569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_trans", &testT, NULL));
457c4762a1bSJed Brown 
4589566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatConvert_Nest_IS\n"));
459c4762a1bSJed Brown   nr = 2;
460c4762a1bSJed Brown   nc = 2;
4619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nr", &nr, NULL));
4629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nc", &nc, NULL));
463c4762a1bSJed Brown   if (testT) {
4649566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(A, &cst, &cen));
4659566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRangeColumn(A, &rst, &ren));
466c4762a1bSJed Brown   } else {
4679566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(A, &rst, &ren));
4689566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRangeColumn(A, &cst, &cen));
469c4762a1bSJed Brown   }
4709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(nr, &rows, nc, &cols, 2 * nr * nc, &mats));
471c4762a1bSJed Brown   for (i = 0; i < nr * nc; i++) {
472c4762a1bSJed Brown     if (testT) {
4739566063dSJacob Faibussowitsch       PetscCall(MatCreateTranspose(A, &mats[i]));
4749566063dSJacob Faibussowitsch       PetscCall(MatTranspose(B, MAT_INITIAL_MATRIX, &mats[i + nr * nc]));
475c4762a1bSJed Brown     } else {
4769566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &mats[i]));
4779566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &mats[i + nr * nc]));
478c4762a1bSJed Brown     }
479c4762a1bSJed Brown   }
48048a46eb9SPierre Jolivet   for (i = 0; i < nr; i++) PetscCall(ISCreateStride(PETSC_COMM_WORLD, ren - rst, i + rst, nr, &rows[i]));
48148a46eb9SPierre Jolivet   for (i = 0; i < nc; i++) PetscCall(ISCreateStride(PETSC_COMM_WORLD, cen - cst, i + cst, nc, &cols[i]));
4829566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(PETSC_COMM_WORLD, nr, rows, nc, cols, mats, &A2));
4839566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(PETSC_COMM_WORLD, nr, rows, nc, cols, mats + nr * nc, &B2));
48448a46eb9SPierre Jolivet   for (i = 0; i < nr; i++) PetscCall(ISDestroy(&rows[i]));
48548a46eb9SPierre Jolivet   for (i = 0; i < nc; i++) PetscCall(ISDestroy(&cols[i]));
48648a46eb9SPierre Jolivet   for (i = 0; i < 2 * nr * nc; i++) PetscCall(MatDestroy(&mats[i]));
4879566063dSJacob Faibussowitsch   PetscCall(PetscFree3(rows, cols, mats));
4889566063dSJacob Faibussowitsch   PetscCall(MatConvert(B2, MATAIJ, MAT_INITIAL_MATRIX, &T));
4899566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B2));
4909566063dSJacob Faibussowitsch   PetscCall(MatConvert(A2, MATIS, MAT_INITIAL_MATRIX, &B2));
4919566063dSJacob Faibussowitsch   PetscCall(CheckMat(B2, T, PETSC_TRUE, "MatConvert_Nest_IS MAT_INITIAL_MATRIX"));
4929566063dSJacob Faibussowitsch   PetscCall(MatConvert(A2, MATIS, MAT_REUSE_MATRIX, &B2));
4939566063dSJacob Faibussowitsch   PetscCall(CheckMat(B2, T, PETSC_TRUE, "MatConvert_Nest_IS MAT_REUSE_MATRIX"));
4949566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B2));
4959566063dSJacob Faibussowitsch   PetscCall(MatConvert(A2, MATIS, MAT_INPLACE_MATRIX, &A2));
4969566063dSJacob Faibussowitsch   PetscCall(CheckMat(A2, T, PETSC_TRUE, "MatConvert_Nest_IS MAT_INPLACE_MATRIX"));
4979566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&T));
4989566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A2));
499c4762a1bSJed Brown 
500c4762a1bSJed Brown   /* test MatCreateSubMatrix */
5019566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatCreateSubMatrix\n"));
502dd400576SPatrick Sanan   if (rank == 0) {
5039566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_WORLD, 1, 1, 1, &is));
5049566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_WORLD, 2, 0, 1, &is2));
505c4762a1bSJed Brown   } else if (rank == 1) {
5069566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_WORLD, 1, 0, 1, &is));
507c4762a1bSJed Brown     if (n > 3) {
5089566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_WORLD, 1, 3, 1, &is2));
509c4762a1bSJed Brown     } else {
5109566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_WORLD, 0, 0, 1, &is2));
511c4762a1bSJed Brown     }
512c4762a1bSJed Brown   } else if (rank == 2 && n > 4) {
5139566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_WORLD, 0, 0, 1, &is));
5149566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_WORLD, n - 4, 4, 1, &is2));
515c4762a1bSJed Brown   } else {
5169566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_WORLD, 0, 0, 1, &is));
5179566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_WORLD, 0, 0, 1, &is2));
518c4762a1bSJed Brown   }
5199566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A, is, is, MAT_INITIAL_MATRIX, &A2));
5209566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(B, is, is, MAT_INITIAL_MATRIX, &B2));
5219566063dSJacob Faibussowitsch   PetscCall(CheckMat(A2, B2, PETSC_TRUE, "first MatCreateSubMatrix"));
522c4762a1bSJed Brown 
5239566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(A, is, is, MAT_REUSE_MATRIX, &A2));
5249566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(B, is, is, MAT_REUSE_MATRIX, &B2));
5259566063dSJacob Faibussowitsch   PetscCall(CheckMat(A2, B2, PETSC_FALSE, "reuse MatCreateSubMatrix"));
5269566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A2));
5279566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B2));
528c4762a1bSJed Brown 
529e432b41dSStefano Zampini   if (!issymmetric) {
5309566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(A, is, is2, MAT_INITIAL_MATRIX, &A2));
5319566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(B, is, is2, MAT_INITIAL_MATRIX, &B2));
5329566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(A, is, is2, MAT_REUSE_MATRIX, &A2));
5339566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(B, is, is2, MAT_REUSE_MATRIX, &B2));
5349566063dSJacob Faibussowitsch     PetscCall(CheckMat(A2, B2, PETSC_FALSE, "second MatCreateSubMatrix"));
535c4762a1bSJed Brown   }
536c4762a1bSJed Brown 
5379566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A2));
5389566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B2));
5399566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
5409566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is2));
541c4762a1bSJed Brown 
542d0dbe9f7SStefano Zampini   /* test MatCreateSubMatrices */
543d0dbe9f7SStefano Zampini   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatCreateSubMatrices\n"));
544d0dbe9f7SStefano Zampini   PetscCall(MatGetLayouts(A, &rlayout, &clayout));
545d0dbe9f7SStefano Zampini   PetscCall(PetscLayoutGetRanges(rlayout, &rrange));
546d0dbe9f7SStefano Zampini   PetscCall(PetscLayoutGetRanges(clayout, &crange));
547d0dbe9f7SStefano Zampini   lrank = (size + rank - 1) % size;
548d0dbe9f7SStefano Zampini   rrank = (rank + 1) % size;
549d0dbe9f7SStefano Zampini   PetscCall(ISCreateStride(PETSC_COMM_SELF, (rrange[lrank + 1] - rrange[lrank]), rrange[lrank], 1, &irow[0]));
550d0dbe9f7SStefano Zampini   PetscCall(ISCreateStride(PETSC_COMM_SELF, (crange[rrank + 1] - crange[rrank]), crange[rrank], 1, &icol[0]));
551d0dbe9f7SStefano Zampini   PetscCall(ISCreateStride(PETSC_COMM_SELF, (rrange[rrank + 1] - rrange[rrank]), rrange[rrank], 1, &irow[1]));
552d0dbe9f7SStefano Zampini   PetscCall(ISCreateStride(PETSC_COMM_SELF, (crange[lrank + 1] - crange[lrank]), crange[lrank], 1, &icol[1]));
553d0dbe9f7SStefano Zampini   PetscCall(MatCreateSubMatrices(A, 2, irow, icol, MAT_INITIAL_MATRIX, &Asub));
554d0dbe9f7SStefano Zampini   PetscCall(MatCreateSubMatrices(B, 2, irow, icol, MAT_INITIAL_MATRIX, &Bsub));
555d0dbe9f7SStefano Zampini   PetscCall(CheckMat(Asub[0], Bsub[0], PETSC_FALSE, "MatCreateSubMatrices[0]"));
556d0dbe9f7SStefano Zampini   PetscCall(CheckMat(Asub[1], Bsub[1], PETSC_FALSE, "MatCreateSubMatrices[1]"));
557d0dbe9f7SStefano Zampini   PetscCall(MatCreateSubMatrices(A, 2, irow, icol, MAT_REUSE_MATRIX, &Asub));
558d0dbe9f7SStefano Zampini   PetscCall(MatCreateSubMatrices(B, 2, irow, icol, MAT_REUSE_MATRIX, &Bsub));
559d0dbe9f7SStefano Zampini   PetscCall(CheckMat(Asub[0], Bsub[0], PETSC_FALSE, "MatCreateSubMatrices[0]"));
560d0dbe9f7SStefano Zampini   PetscCall(CheckMat(Asub[1], Bsub[1], PETSC_FALSE, "MatCreateSubMatrices[1]"));
561d0dbe9f7SStefano Zampini   PetscCall(MatDestroySubMatrices(2, &Asub));
562d0dbe9f7SStefano Zampini   PetscCall(MatDestroySubMatrices(2, &Bsub));
563d0dbe9f7SStefano Zampini   PetscCall(ISDestroy(&irow[0]));
564d0dbe9f7SStefano Zampini   PetscCall(ISDestroy(&irow[1]));
565d0dbe9f7SStefano Zampini   PetscCall(ISDestroy(&icol[0]));
566d0dbe9f7SStefano Zampini   PetscCall(ISDestroy(&icol[1]));
567d0dbe9f7SStefano Zampini 
568c4762a1bSJed Brown   /* Create an IS required by MatZeroRows(): just rank zero provides the rows to be eliminated */
569c4762a1bSJed Brown   if (size > 1) {
570dd400576SPatrick Sanan     if (rank == 0) {
571c4762a1bSJed Brown       PetscInt st, len;
572c4762a1bSJed Brown 
573c4762a1bSJed Brown       st  = (m + 1) / 2;
574c4762a1bSJed Brown       len = PetscMin(m / 2, PetscMax(m - (m + 1) / 2 - 1, 0));
5759566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_WORLD, len, st, 1, &is));
576c4762a1bSJed Brown     } else {
5779566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_WORLD, 0, 0, 1, &is));
578c4762a1bSJed Brown     }
579c4762a1bSJed Brown   } else {
5809566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_WORLD, 1, 0, 1, &is));
581c4762a1bSJed Brown   }
582c4762a1bSJed Brown 
583c4762a1bSJed Brown   if (squaretest) { /* tests for square matrices only, with same maps for rows and columns */
584d0dbe9f7SStefano Zampini     PetscInt *idx0, *idx1, n0, n1;
585d0dbe9f7SStefano Zampini     IS        Ais[2], Bis[2];
586d0dbe9f7SStefano Zampini 
587c4762a1bSJed Brown     /* test MatDiagonalSet */
5889566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatDiagonalSet\n"));
5899566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2));
5909566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &B2));
5919566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(A, NULL, &x));
5929566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x, NULL));
5934f58015eSStefano Zampini     PetscCall(MatDiagonalSet(A2, x, allow_repeated ? ADD_VALUES : INSERT_VALUES));
5944f58015eSStefano Zampini     PetscCall(MatDiagonalSet(B2, x, allow_repeated ? ADD_VALUES : INSERT_VALUES));
5959566063dSJacob Faibussowitsch     PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatDiagonalSet"));
5969566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&x));
5979566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A2));
5989566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B2));
599c4762a1bSJed Brown 
600c4762a1bSJed Brown     /* test MatShift (MatShift_IS internally uses MatDiagonalSet_IS with ADD_VALUES) */
6019566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatShift\n"));
6029566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2));
6039566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &B2));
6049566063dSJacob Faibussowitsch     PetscCall(MatShift(A2, 2.0));
6059566063dSJacob Faibussowitsch     PetscCall(MatShift(B2, 2.0));
6069566063dSJacob Faibussowitsch     PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatShift"));
6079566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A2));
6089566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B2));
609c4762a1bSJed Brown 
610c4762a1bSJed Brown     /* nonzero diag value is supported for square matrices only */
6119566063dSJacob Faibussowitsch     PetscCall(TestMatZeroRows(A, B, PETSC_TRUE, is, diag));
612d0dbe9f7SStefano Zampini 
613d0dbe9f7SStefano Zampini     /* test MatIncreaseOverlap */
614d0dbe9f7SStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatIncreaseOverlap\n"));
615d0dbe9f7SStefano Zampini     PetscCall(MatGetOwnershipRange(A, &rst, &ren));
616d0dbe9f7SStefano Zampini     n0 = (ren - rst) / 2;
617d0dbe9f7SStefano Zampini     n1 = (ren - rst) / 3;
618d0dbe9f7SStefano Zampini     PetscCall(PetscMalloc1(n0, &idx0));
619d0dbe9f7SStefano Zampini     PetscCall(PetscMalloc1(n1, &idx1));
620d0dbe9f7SStefano Zampini     for (i = 0; i < n0; i++) idx0[i] = ren - i - 1;
621d0dbe9f7SStefano Zampini     for (i = 0; i < n1; i++) idx1[i] = rst + i;
622d0dbe9f7SStefano Zampini     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n0, idx0, PETSC_OWN_POINTER, &Ais[0]));
623d0dbe9f7SStefano Zampini     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n1, idx1, PETSC_OWN_POINTER, &Ais[1]));
624d0dbe9f7SStefano Zampini     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n0, idx0, PETSC_COPY_VALUES, &Bis[0]));
625d0dbe9f7SStefano Zampini     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n1, idx1, PETSC_COPY_VALUES, &Bis[1]));
626d0dbe9f7SStefano Zampini     PetscCall(MatIncreaseOverlap(A, 2, Ais, 3));
627d0dbe9f7SStefano Zampini     PetscCall(MatIncreaseOverlap(B, 2, Bis, 3));
628d0dbe9f7SStefano Zampini     /* Non deterministic output! */
629d0dbe9f7SStefano Zampini     PetscCall(ISSort(Ais[0]));
630d0dbe9f7SStefano Zampini     PetscCall(ISSort(Ais[1]));
631d0dbe9f7SStefano Zampini     PetscCall(ISSort(Bis[0]));
632d0dbe9f7SStefano Zampini     PetscCall(ISSort(Bis[1]));
633d0dbe9f7SStefano Zampini     PetscCall(ISView(Ais[0], NULL));
634d0dbe9f7SStefano Zampini     PetscCall(ISView(Bis[0], NULL));
635d0dbe9f7SStefano Zampini     PetscCall(ISView(Ais[1], NULL));
636d0dbe9f7SStefano Zampini     PetscCall(ISView(Bis[1], NULL));
637d0dbe9f7SStefano Zampini     PetscCall(MatCreateSubMatrices(A, 2, Ais, Ais, MAT_INITIAL_MATRIX, &Asub));
638d0dbe9f7SStefano Zampini     PetscCall(MatCreateSubMatrices(B, 2, Bis, Ais, MAT_INITIAL_MATRIX, &Bsub));
639d0dbe9f7SStefano Zampini     PetscCall(CheckMat(Asub[0], Bsub[0], PETSC_FALSE, "MatIncreaseOverlap[0]"));
640d0dbe9f7SStefano Zampini     PetscCall(CheckMat(Asub[1], Bsub[1], PETSC_FALSE, "MatIncreaseOverlap[1]"));
641d0dbe9f7SStefano Zampini     PetscCall(MatDestroySubMatrices(2, &Asub));
642d0dbe9f7SStefano Zampini     PetscCall(MatDestroySubMatrices(2, &Bsub));
643d0dbe9f7SStefano Zampini     PetscCall(ISDestroy(&Ais[0]));
644d0dbe9f7SStefano Zampini     PetscCall(ISDestroy(&Ais[1]));
645d0dbe9f7SStefano Zampini     PetscCall(ISDestroy(&Bis[0]));
646d0dbe9f7SStefano Zampini     PetscCall(ISDestroy(&Bis[1]));
647c4762a1bSJed Brown   }
6489566063dSJacob Faibussowitsch   PetscCall(TestMatZeroRows(A, B, squaretest, is, 0.0));
6499566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
650c4762a1bSJed Brown 
651c4762a1bSJed Brown   /* test MatTranspose */
6529566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatTranspose\n"));
6539566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &A2));
6549566063dSJacob Faibussowitsch   PetscCall(MatTranspose(B, MAT_INITIAL_MATRIX, &B2));
6559566063dSJacob Faibussowitsch   PetscCall(CheckMat(A2, B2, PETSC_FALSE, "initial matrix MatTranspose"));
656c4762a1bSJed Brown 
6579566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A, MAT_REUSE_MATRIX, &A2));
6589566063dSJacob Faibussowitsch   PetscCall(CheckMat(A2, B2, PETSC_FALSE, "reuse matrix (not in place) MatTranspose"));
6599566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A2));
660c4762a1bSJed Brown 
6619566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2));
6629566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A2, MAT_INPLACE_MATRIX, &A2));
6639566063dSJacob Faibussowitsch   PetscCall(CheckMat(A2, B2, PETSC_FALSE, "reuse matrix (in place) MatTranspose"));
6649566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A2));
665c4762a1bSJed Brown 
6669566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &A2));
6679566063dSJacob Faibussowitsch   PetscCall(CheckMat(A2, B2, PETSC_FALSE, "reuse matrix (different type) MatTranspose"));
6689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A2));
6699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B2));
670c4762a1bSJed Brown 
671c4762a1bSJed Brown   /* test MatISFixLocalEmpty */
672c4762a1bSJed Brown   if (isaij) {
673c4762a1bSJed Brown     PetscInt r[2];
674c4762a1bSJed Brown 
675c4762a1bSJed Brown     r[0] = 0;
676c4762a1bSJed Brown     r[1] = PetscMin(m, n) - 1;
6779566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatISFixLocalEmpty\n"));
6789566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2));
679e432b41dSStefano Zampini 
6809566063dSJacob Faibussowitsch     PetscCall(MatISFixLocalEmpty(A2, PETSC_TRUE));
6819566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A2, MAT_FINAL_ASSEMBLY));
6829566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A2, MAT_FINAL_ASSEMBLY));
6839566063dSJacob Faibussowitsch     PetscCall(CheckMat(A2, B, PETSC_FALSE, "MatISFixLocalEmpty (null)"));
684c4762a1bSJed Brown 
6859566063dSJacob Faibussowitsch     PetscCall(MatZeroRows(A2, 2, r, 0.0, NULL, NULL));
6869566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(A2, NULL, "-fixempty_view"));
6879566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &B2));
6889566063dSJacob Faibussowitsch     PetscCall(MatZeroRows(B2, 2, r, 0.0, NULL, NULL));
6899566063dSJacob Faibussowitsch     PetscCall(MatISFixLocalEmpty(A2, PETSC_TRUE));
6909566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A2, MAT_FINAL_ASSEMBLY));
6919566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A2, MAT_FINAL_ASSEMBLY));
6929566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(A2, NULL, "-fixempty_view"));
6939566063dSJacob Faibussowitsch     PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatISFixLocalEmpty (rows)"));
6949566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A2));
695c4762a1bSJed Brown 
6969566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2));
6979566063dSJacob Faibussowitsch     PetscCall(MatZeroRows(A2, 2, r, 0.0, NULL, NULL));
6989566063dSJacob Faibussowitsch     PetscCall(MatTranspose(A2, MAT_INPLACE_MATRIX, &A2));
6999566063dSJacob Faibussowitsch     PetscCall(MatTranspose(B2, MAT_INPLACE_MATRIX, &B2));
7009566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(A2, NULL, "-fixempty_view"));
7019566063dSJacob Faibussowitsch     PetscCall(MatISFixLocalEmpty(A2, PETSC_TRUE));
7029566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A2, MAT_FINAL_ASSEMBLY));
7039566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A2, MAT_FINAL_ASSEMBLY));
7049566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(A2, NULL, "-fixempty_view"));
7059566063dSJacob Faibussowitsch     PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatISFixLocalEmpty (cols)"));
706c4762a1bSJed Brown 
7079566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A2));
7089566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B2));
709c4762a1bSJed Brown 
710c4762a1bSJed Brown     if (squaretest) {
7119566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2));
7129566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &B2));
7139566063dSJacob Faibussowitsch       PetscCall(MatZeroRowsColumns(A2, 2, r, 0.0, NULL, NULL));
7149566063dSJacob Faibussowitsch       PetscCall(MatZeroRowsColumns(B2, 2, r, 0.0, NULL, NULL));
7159566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(A2, NULL, "-fixempty_view"));
7169566063dSJacob Faibussowitsch       PetscCall(MatISFixLocalEmpty(A2, PETSC_TRUE));
7179566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(A2, MAT_FINAL_ASSEMBLY));
7189566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(A2, MAT_FINAL_ASSEMBLY));
7199566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(A2, NULL, "-fixempty_view"));
7209566063dSJacob Faibussowitsch       PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatISFixLocalEmpty (rows+cols)"));
7219566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A2));
7229566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&B2));
723c4762a1bSJed Brown     }
724c4762a1bSJed Brown   }
725c4762a1bSJed Brown 
726c4762a1bSJed Brown   /* test MatInvertBlockDiagonal
727c4762a1bSJed Brown        special cases for block-diagonal matrices */
728c4762a1bSJed Brown   if (m == n) {
729c4762a1bSJed Brown     ISLocalToGlobalMapping map;
730c4762a1bSJed Brown     Mat                    Abd, Bbd;
731c4762a1bSJed Brown     IS                     is, bis;
732c4762a1bSJed Brown     const PetscScalar     *isbd, *aijbd;
733c4762a1bSJed Brown     PetscScalar           *vals;
734c4762a1bSJed Brown     const PetscInt        *sts, *idxs;
735c4762a1bSJed Brown     PetscInt              *idxs2, diff, perm, nl, bs, st, en, in;
736c4762a1bSJed Brown     PetscBool              ok;
737c4762a1bSJed Brown 
738c4762a1bSJed Brown     for (diff = 0; diff < 3; diff++) {
739c4762a1bSJed Brown       for (perm = 0; perm < 3; perm++) {
740c4762a1bSJed Brown         for (bs = 1; bs < 4; bs++) {
7419566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatInvertBlockDiagonal blockdiag %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", n, diff, perm, bs));
7429566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(bs * bs, &vals));
7439566063dSJacob Faibussowitsch           PetscCall(MatGetOwnershipRanges(A, &sts));
744c4762a1bSJed Brown           switch (diff) {
745c4762a1bSJed Brown           case 1: /* inverted layout by processes */
746c4762a1bSJed Brown             in = 1;
747c4762a1bSJed Brown             st = sts[size - rank - 1];
748c4762a1bSJed Brown             en = sts[size - rank];
749c4762a1bSJed Brown             nl = en - st;
750c4762a1bSJed Brown             break;
751c4762a1bSJed Brown           case 2: /* round-robin layout */
752c4762a1bSJed Brown             in = size;
753c4762a1bSJed Brown             st = rank;
754c4762a1bSJed Brown             nl = n / size;
755c4762a1bSJed Brown             if (rank < n % size) nl++;
756c4762a1bSJed Brown             break;
757c4762a1bSJed Brown           default: /* same layout */
758c4762a1bSJed Brown             in = 1;
759c4762a1bSJed Brown             st = sts[rank];
760c4762a1bSJed Brown             en = sts[rank + 1];
761c4762a1bSJed Brown             nl = en - st;
762c4762a1bSJed Brown             break;
763c4762a1bSJed Brown           }
7649566063dSJacob Faibussowitsch           PetscCall(ISCreateStride(PETSC_COMM_WORLD, nl, st, in, &is));
7659566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(is, &nl));
7669566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(is, &idxs));
7679566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(nl, &idxs2));
768c4762a1bSJed Brown           for (i = 0; i < nl; i++) {
769c4762a1bSJed Brown             switch (perm) { /* invert some of the indices */
770d71ae5a4SJacob Faibussowitsch             case 2:
771d71ae5a4SJacob Faibussowitsch               idxs2[i] = rank % 2 ? idxs[i] : idxs[nl - i - 1];
772d71ae5a4SJacob Faibussowitsch               break;
773d71ae5a4SJacob Faibussowitsch             case 1:
774d71ae5a4SJacob Faibussowitsch               idxs2[i] = rank % 2 ? idxs[nl - i - 1] : idxs[i];
775d71ae5a4SJacob Faibussowitsch               break;
776d71ae5a4SJacob Faibussowitsch             default:
777d71ae5a4SJacob Faibussowitsch               idxs2[i] = idxs[i];
778d71ae5a4SJacob Faibussowitsch               break;
779c4762a1bSJed Brown             }
780c4762a1bSJed Brown           }
7819566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(is, &idxs));
7829566063dSJacob Faibussowitsch           PetscCall(ISCreateBlock(PETSC_COMM_WORLD, bs, nl, idxs2, PETSC_OWN_POINTER, &bis));
7839566063dSJacob Faibussowitsch           PetscCall(ISLocalToGlobalMappingCreateIS(bis, &map));
7849566063dSJacob Faibussowitsch           PetscCall(MatCreateIS(PETSC_COMM_WORLD, bs, PETSC_DECIDE, PETSC_DECIDE, bs * n, bs * n, map, map, &Abd));
7859566063dSJacob Faibussowitsch           PetscCall(ISLocalToGlobalMappingDestroy(&map));
7869566063dSJacob Faibussowitsch           PetscCall(MatISSetPreallocation(Abd, bs, NULL, 0, NULL));
787c4762a1bSJed Brown           for (i = 0; i < nl; i++) {
788c4762a1bSJed Brown             PetscInt b1, b2;
789c4762a1bSJed Brown 
790c4762a1bSJed Brown             for (b1 = 0; b1 < bs; b1++) {
791ad540459SPierre Jolivet               for (b2 = 0; b2 < bs; b2++) vals[b1 * bs + b2] = i * bs * bs + b1 * bs + b2 + 1 + (b1 == b2 ? 1.0 : 0);
792c4762a1bSJed Brown             }
7939566063dSJacob Faibussowitsch             PetscCall(MatSetValuesBlockedLocal(Abd, 1, &i, 1, &i, vals, INSERT_VALUES));
794c4762a1bSJed Brown           }
7959566063dSJacob Faibussowitsch           PetscCall(MatAssemblyBegin(Abd, MAT_FINAL_ASSEMBLY));
7969566063dSJacob Faibussowitsch           PetscCall(MatAssemblyEnd(Abd, MAT_FINAL_ASSEMBLY));
7979566063dSJacob Faibussowitsch           PetscCall(MatConvert(Abd, MATAIJ, MAT_INITIAL_MATRIX, &Bbd));
7989566063dSJacob Faibussowitsch           PetscCall(MatInvertBlockDiagonal(Abd, &isbd));
7999566063dSJacob Faibussowitsch           PetscCall(MatInvertBlockDiagonal(Bbd, &aijbd));
8009566063dSJacob Faibussowitsch           PetscCall(MatGetLocalSize(Bbd, &nl, NULL));
801c4762a1bSJed Brown           ok = PETSC_TRUE;
802c4762a1bSJed Brown           for (i = 0; i < nl / bs; i++) {
803c4762a1bSJed Brown             PetscInt b1, b2;
804c4762a1bSJed Brown 
805c4762a1bSJed Brown             for (b1 = 0; b1 < bs; b1++) {
806c4762a1bSJed Brown               for (b2 = 0; b2 < bs; b2++) {
807c4762a1bSJed Brown                 if (PetscAbsScalar(isbd[i * bs * bs + b1 * bs + b2] - aijbd[i * bs * bs + b1 * bs + b2]) > PETSC_SMALL) ok = PETSC_FALSE;
808c4762a1bSJed Brown                 if (!ok) {
8099566063dSJacob Faibussowitsch                   PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] ERROR block %" PetscInt_FMT ", entry %" PetscInt_FMT " %" PetscInt_FMT ": %g %g\n", rank, i, b1, b2, (double)PetscAbsScalar(isbd[i * bs * bs + b1 * bs + b2]), (double)PetscAbsScalar(aijbd[i * bs * bs + b1 * bs + b2])));
810c4762a1bSJed Brown                   break;
811c4762a1bSJed Brown                 }
812c4762a1bSJed Brown               }
813c4762a1bSJed Brown               if (!ok) break;
814c4762a1bSJed Brown             }
815c4762a1bSJed Brown             if (!ok) break;
816c4762a1bSJed Brown           }
8179566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&Abd));
8189566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&Bbd));
8199566063dSJacob Faibussowitsch           PetscCall(PetscFree(vals));
8209566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&is));
8219566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&bis));
822c4762a1bSJed Brown         }
823c4762a1bSJed Brown       }
824c4762a1bSJed Brown     }
825c4762a1bSJed Brown   }
826d0dbe9f7SStefano Zampini 
827d0dbe9f7SStefano Zampini   /* test MatGetDiagonalBlock */
828d0dbe9f7SStefano Zampini   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatGetDiagonalBlock\n"));
829d0dbe9f7SStefano Zampini   PetscCall(MatGetDiagonalBlock(A, &A2));
830d0dbe9f7SStefano Zampini   PetscCall(MatGetDiagonalBlock(B, &B2));
831d0dbe9f7SStefano Zampini   PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatGetDiagonalBlock"));
832d0dbe9f7SStefano Zampini   PetscCall(MatScale(A, 2.0));
833d0dbe9f7SStefano Zampini   PetscCall(MatScale(B, 2.0));
834d0dbe9f7SStefano Zampini   PetscCall(MatGetDiagonalBlock(A, &A2));
835d0dbe9f7SStefano Zampini   PetscCall(MatGetDiagonalBlock(B, &B2));
836d0dbe9f7SStefano Zampini   PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatGetDiagonalBlock"));
837d0dbe9f7SStefano Zampini 
8384f58015eSStefano Zampini   /* test MatISSetAllowRepeated on a MATIS */
8394f58015eSStefano Zampini   PetscCall(MatISSetAllowRepeated(A, allow_repeated));
8404f58015eSStefano Zampini   if (allow_repeated) { /* original MATIS maybe with repeated entries, test assembling of local matrices */
8414f58015eSStefano Zampini     Mat lA, lA2;
8424f58015eSStefano Zampini 
8434f58015eSStefano Zampini     for (PetscInt i = 0; i < 1; i++) { /* TODO: make MatScatter inherit from MATSHELL and support MatProducts */
8444f58015eSStefano Zampini       PetscBool usemult = PETSC_FALSE;
8454f58015eSStefano Zampini 
8464f58015eSStefano Zampini       PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2));
8474f58015eSStefano Zampini       if (i) {
8484f58015eSStefano Zampini         Mat tA;
8494f58015eSStefano Zampini 
8504f58015eSStefano Zampini         usemult = PETSC_TRUE;
8514f58015eSStefano Zampini         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatISSetAllowRepeated(false) with possibly repeated entries and local shell matrices\n"));
8524f58015eSStefano Zampini         PetscCall(MatISGetLocalMat(A2, &lA2));
8534f58015eSStefano Zampini         PetscCall(MatConvert(lA2, MATSHELL, MAT_INITIAL_MATRIX, &tA));
8544f58015eSStefano Zampini         PetscCall(MatISRestoreLocalMat(A2, &lA2));
8554f58015eSStefano Zampini         PetscCall(MatISSetLocalMat(A2, tA));
8564f58015eSStefano Zampini         PetscCall(MatDestroy(&tA));
8574f58015eSStefano Zampini       } else {
8584f58015eSStefano Zampini         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatISSetAllowRepeated(false) with possibly repeated entries\n"));
8594f58015eSStefano Zampini       }
8604f58015eSStefano Zampini       PetscCall(MatISSetAllowRepeated(A2, PETSC_FALSE));
8614f58015eSStefano Zampini       PetscCall(MatISGetLocalMat(A, &lA));
8624f58015eSStefano Zampini       PetscCall(MatISGetLocalMat(A2, &lA2));
8634f58015eSStefano Zampini       if (!repmap) PetscCall(CheckMat(lA, lA2, usemult, "MatISSetAllowRepeated(false) with non-repeated entries"));
8644f58015eSStefano Zampini       PetscCall(MatISRestoreLocalMat(A, &lA));
8654f58015eSStefano Zampini       PetscCall(MatISRestoreLocalMat(A2, &lA2));
8664f58015eSStefano Zampini       if (repmap) PetscCall(CheckMat(A2, B, usemult, "MatISSetAllowRepeated(false) with repeated entries"));
8674f58015eSStefano Zampini       else PetscCall(CheckMat(A2, B, usemult, "MatISSetAllowRepeated(false) with non-repeated entries"));
8684f58015eSStefano Zampini       PetscCall(MatDestroy(&A2));
8694f58015eSStefano Zampini     }
8704f58015eSStefano Zampini   } else { /* original matis with non-repeated entries, this should only recreate the local matrices */
8714f58015eSStefano Zampini     Mat       lA;
8724f58015eSStefano Zampini     PetscBool flg;
8734f58015eSStefano Zampini 
8744f58015eSStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatISSetAllowRepeated(true) with non repeated entries\n"));
8754f58015eSStefano Zampini     PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2));
8764f58015eSStefano Zampini     PetscCall(MatISSetAllowRepeated(A2, PETSC_TRUE));
8774f58015eSStefano Zampini     PetscCall(MatISGetLocalMat(A2, &lA));
8784f58015eSStefano Zampini     PetscCall(MatAssembled(lA, &flg));
8794f58015eSStefano Zampini     PetscCheck(!flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Local mat should be unassembled");
8804f58015eSStefano Zampini     PetscCall(MatISRestoreLocalMat(A2, &lA));
8814f58015eSStefano Zampini     PetscCall(MatDestroy(&A2));
8824f58015eSStefano Zampini   }
8834f58015eSStefano Zampini 
884c4762a1bSJed Brown   /* free testing matrices */
8859566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&cmap));
8869566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&rmap));
8879566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
8889566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
8899566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
890b122ec5aSJacob Faibussowitsch   return 0;
891c4762a1bSJed Brown }
892c4762a1bSJed Brown 
893d71ae5a4SJacob Faibussowitsch PetscErrorCode CheckMat(Mat A, Mat B, PetscBool usemult, const char *func)
894d71ae5a4SJacob Faibussowitsch {
895c4762a1bSJed Brown   Mat       Bcheck;
896c4762a1bSJed Brown   PetscReal error;
897c4762a1bSJed Brown 
898c4762a1bSJed Brown   PetscFunctionBeginUser;
899c4762a1bSJed Brown   if (!usemult && B) {
900c4762a1bSJed Brown     PetscBool hasnorm;
901c4762a1bSJed Brown 
9029566063dSJacob Faibussowitsch     PetscCall(MatHasOperation(B, MATOP_NORM, &hasnorm));
903c4762a1bSJed Brown     if (!hasnorm) usemult = PETSC_TRUE;
904c4762a1bSJed Brown   }
905c4762a1bSJed Brown   if (!usemult) {
906c4762a1bSJed Brown     if (B) {
907c4762a1bSJed Brown       MatType Btype;
908c4762a1bSJed Brown 
9099566063dSJacob Faibussowitsch       PetscCall(MatGetType(B, &Btype));
9109566063dSJacob Faibussowitsch       PetscCall(MatConvert(A, Btype, MAT_INITIAL_MATRIX, &Bcheck));
911c4762a1bSJed Brown     } else {
9129566063dSJacob Faibussowitsch       PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &Bcheck));
913c4762a1bSJed Brown     }
914c4762a1bSJed Brown     if (B) { /* if B is present, subtract it */
9159566063dSJacob Faibussowitsch       PetscCall(MatAXPY(Bcheck, -1., B, DIFFERENT_NONZERO_PATTERN));
916c4762a1bSJed Brown     }
9179566063dSJacob Faibussowitsch     PetscCall(MatNorm(Bcheck, NORM_INFINITY, &error));
918c4762a1bSJed Brown     if (error > PETSC_SQRT_MACHINE_EPSILON) {
919c4762a1bSJed Brown       ISLocalToGlobalMapping rl2g, cl2g;
920c4762a1bSJed Brown 
921d0dbe9f7SStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)Bcheck, "Bcheck"));
9229566063dSJacob Faibussowitsch       PetscCall(MatView(Bcheck, NULL));
923c4762a1bSJed Brown       if (B) {
924d0dbe9f7SStefano Zampini         PetscCall(PetscObjectSetName((PetscObject)B, "B"));
9259566063dSJacob Faibussowitsch         PetscCall(MatView(B, NULL));
9269566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&Bcheck));
9279566063dSJacob Faibussowitsch         PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &Bcheck));
928d0dbe9f7SStefano Zampini         PetscCall(PetscObjectSetName((PetscObject)Bcheck, "Assembled A"));
9299566063dSJacob Faibussowitsch         PetscCall(MatView(Bcheck, NULL));
930c4762a1bSJed Brown       }
9319566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Bcheck));
932d0dbe9f7SStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)A, "A"));
9339566063dSJacob Faibussowitsch       PetscCall(MatView(A, NULL));
9349566063dSJacob Faibussowitsch       PetscCall(MatGetLocalToGlobalMapping(A, &rl2g, &cl2g));
935d0dbe9f7SStefano Zampini       if (rl2g) PetscCall(ISLocalToGlobalMappingView(rl2g, NULL));
936d0dbe9f7SStefano Zampini       if (cl2g) PetscCall(ISLocalToGlobalMappingView(cl2g, NULL));
937d0dbe9f7SStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "ERROR ON %s: %g", func, (double)error);
938c4762a1bSJed Brown     }
9399566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Bcheck));
940c4762a1bSJed Brown   } else {
941c4762a1bSJed Brown     PetscBool ok, okt;
942c4762a1bSJed Brown 
9439566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(A, B, 3, &ok));
9449566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeEqual(A, B, 3, &okt));
945e00437b9SBarry Smith     PetscCheck(ok && okt, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "ERROR ON %s: mult ok ?  %d, multtranspose ok ? %d", func, ok, okt);
946c4762a1bSJed Brown   }
9473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
948c4762a1bSJed Brown }
949c4762a1bSJed Brown 
950d71ae5a4SJacob Faibussowitsch PetscErrorCode TestMatZeroRows(Mat A, Mat Afull, PetscBool squaretest, IS is, PetscScalar diag)
951d71ae5a4SJacob Faibussowitsch {
952c4762a1bSJed Brown   Mat                    B, Bcheck, B2 = NULL, lB;
953c4762a1bSJed Brown   Vec                    x = NULL, b = NULL, b2 = NULL;
954c4762a1bSJed Brown   ISLocalToGlobalMapping l2gr, l2gc;
955c4762a1bSJed Brown   PetscReal              error;
956c4762a1bSJed Brown   char                   diagstr[16];
957c4762a1bSJed Brown   const PetscInt        *idxs;
958c4762a1bSJed Brown   PetscInt               rst, ren, i, n, N, d;
959c4762a1bSJed Brown   PetscMPIInt            rank;
960c4762a1bSJed Brown   PetscBool              miss, haszerorows;
961c4762a1bSJed Brown 
962c4762a1bSJed Brown   PetscFunctionBeginUser;
963c4762a1bSJed Brown   if (diag == 0.) {
964c6a7a370SJeremy L Thompson     PetscCall(PetscStrncpy(diagstr, "zero", sizeof(diagstr)));
965c4762a1bSJed Brown   } else {
966c6a7a370SJeremy L Thompson     PetscCall(PetscStrncpy(diagstr, "nonzero", sizeof(diagstr)));
967c4762a1bSJed Brown   }
9689566063dSJacob Faibussowitsch   PetscCall(ISView(is, NULL));
9699566063dSJacob Faibussowitsch   PetscCall(MatGetLocalToGlobalMapping(A, &l2gr, &l2gc));
970c4762a1bSJed Brown   /* tests MatDuplicate and MatCopy */
971c4762a1bSJed Brown   if (diag == 0.) {
9729566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B));
973c4762a1bSJed Brown   } else {
9749566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A, MAT_DO_NOT_COPY_VALUES, &B));
9759566063dSJacob Faibussowitsch     PetscCall(MatCopy(A, B, SAME_NONZERO_PATTERN));
976c4762a1bSJed Brown   }
9779566063dSJacob Faibussowitsch   PetscCall(MatISGetLocalMat(B, &lB));
9789566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(lB, MATOP_ZERO_ROWS, &haszerorows));
979c4762a1bSJed Brown   if (squaretest && haszerorows) {
9809566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(B, &x, &b));
9819566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &B2));
9829566063dSJacob Faibussowitsch     PetscCall(VecSetLocalToGlobalMapping(b, l2gr));
9839566063dSJacob Faibussowitsch     PetscCall(VecSetLocalToGlobalMapping(x, l2gc));
9849566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x, NULL));
9859566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(b, NULL));
986c4762a1bSJed Brown     /* mimic b[is] = x[is] */
9879566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(b, &b2));
9889566063dSJacob Faibussowitsch     PetscCall(VecSetLocalToGlobalMapping(b2, l2gr));
9899566063dSJacob Faibussowitsch     PetscCall(VecCopy(b, b2));
9909566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(is, &n));
9919566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is, &idxs));
9929566063dSJacob Faibussowitsch     PetscCall(VecGetSize(x, &N));
993c4762a1bSJed Brown     for (i = 0; i < n; i++) {
994c4762a1bSJed Brown       if (0 <= idxs[i] && idxs[i] < N) {
9959566063dSJacob Faibussowitsch         PetscCall(VecSetValue(b2, idxs[i], diag, INSERT_VALUES));
9969566063dSJacob Faibussowitsch         PetscCall(VecSetValue(x, idxs[i], 1., INSERT_VALUES));
997c4762a1bSJed Brown       }
998c4762a1bSJed Brown     }
9999566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(b2));
10009566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(b2));
10019566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(x));
10029566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(x));
10039566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is, &idxs));
1004c4762a1bSJed Brown     /*  test ZeroRows on MATIS */
10059566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatZeroRows (diag %s)\n", diagstr));
10069566063dSJacob Faibussowitsch     PetscCall(MatZeroRowsIS(B, is, diag, x, b));
10079566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatZeroRowsColumns (diag %s)\n", diagstr));
10089566063dSJacob Faibussowitsch     PetscCall(MatZeroRowsColumnsIS(B2, is, diag, NULL, NULL));
1009c4762a1bSJed Brown   } else if (haszerorows) {
1010c4762a1bSJed Brown     /*  test ZeroRows on MATIS */
10119566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatZeroRows (diag %s)\n", diagstr));
10129566063dSJacob Faibussowitsch     PetscCall(MatZeroRowsIS(B, is, diag, NULL, NULL));
1013c4762a1bSJed Brown     b = b2 = x = NULL;
1014c4762a1bSJed Brown   } else {
10159566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Skipping MatZeroRows (diag %s)\n", diagstr));
1016c4762a1bSJed Brown     b = b2 = x = NULL;
1017c4762a1bSJed Brown   }
1018c4762a1bSJed Brown 
1019c4762a1bSJed Brown   if (squaretest && haszerorows) {
10209566063dSJacob Faibussowitsch     PetscCall(VecAXPY(b2, -1., b));
10219566063dSJacob Faibussowitsch     PetscCall(VecNorm(b2, NORM_INFINITY, &error));
1022e00437b9SBarry Smith     PetscCheck(error <= PETSC_SQRT_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "ERROR IN ZEROROWS ON B %g (diag %s)", (double)error, diagstr);
1023c4762a1bSJed Brown   }
1024c4762a1bSJed Brown 
1025c4762a1bSJed Brown   /* test MatMissingDiagonal */
10269566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatMissingDiagonal\n"));
10279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
10289566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(B, &miss, &d));
10299566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(B, &rst, &ren));
10309566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
10319566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "[%d] [%" PetscInt_FMT ",%" PetscInt_FMT ") Missing %d, row %" PetscInt_FMT " (diag %s)\n", rank, rst, ren, (int)miss, d, diagstr));
10329566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
10339566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
1034c4762a1bSJed Brown 
10359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
10369566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
10379566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b2));
1038c4762a1bSJed Brown 
1039c4762a1bSJed Brown   /* check the result of ZeroRows with that from MPIAIJ routines
1040c4762a1bSJed Brown      assuming that MatConvert_IS_XAIJ and MatZeroRows_MPIAIJ work fine */
1041c4762a1bSJed Brown   if (haszerorows) {
10429566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(Afull, MAT_COPY_VALUES, &Bcheck));
10439566063dSJacob Faibussowitsch     PetscCall(MatSetOption(Bcheck, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
10449566063dSJacob Faibussowitsch     PetscCall(MatZeroRowsIS(Bcheck, is, diag, NULL, NULL));
10459566063dSJacob Faibussowitsch     PetscCall(CheckMat(B, Bcheck, PETSC_FALSE, "Zerorows"));
10469566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Bcheck));
1047c4762a1bSJed Brown   }
10489566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1049c4762a1bSJed Brown 
1050c4762a1bSJed Brown   if (B2) { /* test MatZeroRowsColumns */
10519566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(Afull, MAT_COPY_VALUES, &B));
10529566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
10539566063dSJacob Faibussowitsch     PetscCall(MatZeroRowsColumnsIS(B, is, diag, NULL, NULL));
10549566063dSJacob Faibussowitsch     PetscCall(CheckMat(B2, B, PETSC_FALSE, "MatZeroRowsColumns"));
10559566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
10569566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B2));
1057c4762a1bSJed Brown   }
10583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1059c4762a1bSJed Brown }
1060c4762a1bSJed Brown 
1061c4762a1bSJed Brown /*TEST
1062c4762a1bSJed Brown 
1063c4762a1bSJed Brown    test:
10645042aa92SStefano Zampini       requires: !complex
10655042aa92SStefano Zampini       args: -test_matlab -test_trans
1066c4762a1bSJed Brown 
1067c4762a1bSJed Brown    test:
1068c4762a1bSJed Brown       suffix: 2
1069c4762a1bSJed Brown       nsize: 4
10704f58015eSStefano Zampini       args: -mat_is_convert_local_nest -nr 3 -nc 4
1071c4762a1bSJed Brown 
1072c4762a1bSJed Brown    test:
1073c4762a1bSJed Brown       suffix: 3
1074c4762a1bSJed Brown       nsize: 5
1075*a50ef18cSStefano Zampini       args: -m 11 -n 10 -mat_is_convert_local_nest -nr 2 -nc 1 -cbs 2
1076c4762a1bSJed Brown 
1077c4762a1bSJed Brown    test:
1078c4762a1bSJed Brown       suffix: 4
1079c4762a1bSJed Brown       nsize: 6
1080c4762a1bSJed Brown       args: -m 9 -n 12 -test_trans -nr 2 -nc 7
1081c4762a1bSJed Brown 
1082c4762a1bSJed Brown    test:
1083c4762a1bSJed Brown       suffix: 5
1084c4762a1bSJed Brown       nsize: 6
1085*a50ef18cSStefano Zampini       args: -m 12 -n 12 -test_trans -nr 3 -nc 1 -rbs 2
1086c4762a1bSJed Brown 
1087c4762a1bSJed Brown    test:
1088c4762a1bSJed Brown       suffix: 6
1089*a50ef18cSStefano Zampini       args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap -rbs 6 -cbs 3
1090c4762a1bSJed Brown 
1091c4762a1bSJed Brown    test:
1092c4762a1bSJed Brown       suffix: 7
1093c4762a1bSJed Brown       args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap -permmap
1094c4762a1bSJed Brown 
1095c4762a1bSJed Brown    test:
1096c4762a1bSJed Brown       suffix: 8
1097c4762a1bSJed Brown       args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -permmap
1098c4762a1bSJed Brown 
1099c4762a1bSJed Brown    test:
1100c4762a1bSJed Brown       suffix: 9
1101c4762a1bSJed Brown       nsize: 5
1102c4762a1bSJed Brown       args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap
1103c4762a1bSJed Brown 
1104c4762a1bSJed Brown    test:
1105c4762a1bSJed Brown       suffix: 10
1106c4762a1bSJed Brown       nsize: 5
1107c4762a1bSJed Brown       args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap -permmap
1108c4762a1bSJed Brown 
1109c20d7725SJed Brown    test:
1110c20d7725SJed Brown       suffix: vscat_default
1111c4762a1bSJed Brown       nsize: 5
1112c4762a1bSJed Brown       args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -permmap
1113c4762a1bSJed Brown       output_file: output/ex23_11.out
1114c4762a1bSJed Brown 
1115c4762a1bSJed Brown    test:
1116c4762a1bSJed Brown       suffix: 12
1117c4762a1bSJed Brown       nsize: 3
11184f58015eSStefano Zampini       args: -m 12 -n 12 -symmetric -mat_is_localmat_type sbaij -test_trans -nr 2 -nc 3
1119c4762a1bSJed Brown 
1120c4762a1bSJed Brown    testset:
1121c4762a1bSJed Brown       output_file: output/ex23_13.out
1122c4762a1bSJed Brown       nsize: 3
1123c4762a1bSJed Brown       args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -diffmap -permmap
11245042aa92SStefano Zampini       filter: grep -v "type:" | grep -v "block size is 1" | grep -v "not using I-node routines"
1125c4762a1bSJed Brown       test:
1126c4762a1bSJed Brown         suffix: baij
11274f58015eSStefano Zampini         args: -mat_is_localmat_type baij
1128c4762a1bSJed Brown       test:
1129c4762a1bSJed Brown         requires: viennacl
1130c4762a1bSJed Brown         suffix: viennacl
11314f58015eSStefano Zampini         args: -mat_is_localmat_type aijviennacl
1132c4762a1bSJed Brown       test:
1133c4762a1bSJed Brown         requires: cuda
1134c4762a1bSJed Brown         suffix: cusparse
11354f58015eSStefano Zampini         args: -mat_is_localmat_type aijcusparse
11365042aa92SStefano Zampini       test:
11375042aa92SStefano Zampini         requires: kokkos_kernels
11385042aa92SStefano Zampini         suffix: kokkos
11395042aa92SStefano Zampini         args: -mat_is_localmat_type aijkokkos
1140c4762a1bSJed Brown 
1141e432b41dSStefano Zampini    test:
1142e432b41dSStefano Zampini       suffix: negrep
1143e432b41dSStefano Zampini       nsize: {{1 3}separate output}
11444f58015eSStefano Zampini       args: -m {{5 7}separate output} -n {{5 7}separate output} -test_trans -nr 2 -nc 3 -negmap {{0 1}separate output} -repmap {{0 1}separate output} -permmap -diffmap {{0 1}separate output} -allow_repeated 0
11454f58015eSStefano Zampini 
11464f58015eSStefano Zampini    test:
11474f58015eSStefano Zampini       suffix: negrep_allowrep
11484f58015eSStefano Zampini       nsize: {{1 3}separate output}
11494f58015eSStefano Zampini       args: -m {{5 7}separate output} -n {{5 7}separate output} -test_trans -nr 2 -nc 3 -negmap {{0 1}separate output} -repmap {{0 1}separate output} -permmap -diffmap {{0 1}separate output} -allow_repeated
1150e432b41dSStefano Zampini 
1151c4762a1bSJed Brown TEST*/
1152