1c4762a1bSJed Brown 2d0dbe9f7SStefano Zampini static char help[] = "Tests the use of interface functions for MATIS matrices and conversion routines.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown PetscErrorCode TestMatZeroRows(Mat, Mat, PetscBool, IS, PetscScalar); 7c4762a1bSJed Brown PetscErrorCode CheckMat(Mat, Mat, PetscBool, const char *); 8c4762a1bSJed Brown 9d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 10d71ae5a4SJacob Faibussowitsch { 11c4762a1bSJed Brown Mat A, B, A2, B2, T; 12c4762a1bSJed Brown Mat Aee, Aeo, Aoe, Aoo; 13d0dbe9f7SStefano Zampini Mat *mats, *Asub, *Bsub; 14c4762a1bSJed Brown Vec x, y; 15c4762a1bSJed Brown MatInfo info; 16c4762a1bSJed Brown ISLocalToGlobalMapping cmap, rmap; 17c4762a1bSJed Brown IS is, is2, reven, rodd, ceven, codd; 18c4762a1bSJed Brown IS *rows, *cols; 19d0dbe9f7SStefano Zampini IS irow[2], icol[2]; 20d0dbe9f7SStefano Zampini PetscLayout rlayout, clayout; 21d0dbe9f7SStefano Zampini const PetscInt *rrange, *crange; 22c4762a1bSJed Brown MatType lmtype; 23c4762a1bSJed Brown PetscScalar diag = 2.; 24e432b41dSStefano Zampini PetscInt n, m, i, lm, ln; 25c4762a1bSJed Brown PetscInt rst, ren, cst, cen, nr, nc; 26d0dbe9f7SStefano Zampini PetscMPIInt rank, size, lrank, rrank; 27c4762a1bSJed Brown PetscBool testT, squaretest, isaij; 28e432b41dSStefano Zampini PetscBool permute = PETSC_FALSE, negmap = PETSC_FALSE, repmap = PETSC_FALSE; 29e432b41dSStefano Zampini PetscBool diffmap = PETSC_TRUE, symmetric = PETSC_FALSE, issymmetric; 30c4762a1bSJed Brown 31327415f7SBarry Smith PetscFunctionBeginUser; 329566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 35c4762a1bSJed Brown m = n = 2 * size; 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-symmetric", &symmetric, NULL)); 379566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 389566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-negmap", &negmap, NULL)); 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-repmap", &repmap, NULL)); 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-permmap", &permute, NULL)); 429566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-diffmap", &diffmap, NULL)); 43e00437b9SBarry 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"); 44e00437b9SBarry 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"); 45e00437b9SBarry Smith PetscCheck(n >= 2, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Number of cols should be larger or equal 2"); 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* create a MATIS matrix */ 489566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 499566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m, n)); 509566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATIS)); 519566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 52e432b41dSStefano Zampini if (!negmap && !repmap) { 53fc989267SStefano Zampini /* This is not the proper setting for MATIS for finite elements, it is just used to test the routines 54fc989267SStefano Zampini Here we use a one-to-one correspondence between local row/column spaces and global row/column spaces 55fc989267SStefano Zampini Equivalent to passing NULL for the mapping */ 569566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, n, 0, 1, &is)); 57e432b41dSStefano Zampini } else if (negmap && !repmap) { /* non repeated but with negative indices */ 589566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, n + 2, -2, 1, &is)); 59e432b41dSStefano Zampini } else if (!negmap && repmap) { /* non negative but repeated indices */ 60e432b41dSStefano Zampini IS isl[2]; 61e432b41dSStefano Zampini 629566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, n, 0, 1, &isl[0])); 639566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, n, n - 1, -1, &isl[1])); 649566063dSJacob Faibussowitsch PetscCall(ISConcatenate(PETSC_COMM_WORLD, 2, isl, &is)); 659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl[0])); 669566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl[1])); 67e432b41dSStefano Zampini } else { /* negative and repeated indices */ 68e432b41dSStefano Zampini IS isl[2]; 69e432b41dSStefano Zampini 709566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, n + 1, -1, 1, &isl[0])); 719566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, n + 1, n - 1, -1, &isl[1])); 729566063dSJacob Faibussowitsch PetscCall(ISConcatenate(PETSC_COMM_WORLD, 2, isl, &is)); 739566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl[0])); 749566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl[1])); 75e432b41dSStefano Zampini } 769566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(is, &cmap)); 779566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 78c4762a1bSJed Brown 79e432b41dSStefano Zampini if (m != n || diffmap) { 809566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, m, permute ? m - 1 : 0, permute ? -1 : 1, &is)); 819566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(is, &rmap)); 829566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 83c4762a1bSJed Brown } else { 849566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)cmap)); 85c4762a1bSJed Brown rmap = cmap; 86c4762a1bSJed Brown } 87e432b41dSStefano Zampini 889566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(A, rmap, cmap)); 899566063dSJacob Faibussowitsch PetscCall(MatISStoreL2L(A, PETSC_FALSE)); 909566063dSJacob Faibussowitsch PetscCall(MatISSetPreallocation(A, 3, NULL, 3, NULL)); 919566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, (PetscBool) !(repmap || negmap))); /* I do not want to precompute the pattern */ 929566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(rmap, &lm)); 939566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(cmap, &ln)); 94e432b41dSStefano Zampini for (i = 0; i < lm; i++) { 95c4762a1bSJed Brown PetscScalar v[3]; 96c4762a1bSJed Brown PetscInt cols[3]; 97c4762a1bSJed Brown 98c4762a1bSJed Brown cols[0] = (i - 1 + n) % n; 99c4762a1bSJed Brown cols[1] = i % n; 100c4762a1bSJed Brown cols[2] = (i + 1) % n; 101c4762a1bSJed Brown v[0] = -1. * (symmetric ? PetscMin(i + 1, cols[0] + 1) : i + 1); 102c4762a1bSJed Brown v[1] = 2. * (symmetric ? PetscMin(i + 1, cols[1] + 1) : i + 1); 103c4762a1bSJed Brown v[2] = -1. * (symmetric ? PetscMin(i + 1, cols[2] + 1) : i + 1); 1049566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(cmap, IS_GTOLM_MASK, 3, cols, NULL, cols)); 1059566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A, 1, &i, 3, cols, v, ADD_VALUES)); 106c4762a1bSJed Brown } 1079566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1089566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 109c4762a1bSJed Brown 110e432b41dSStefano Zampini /* activate tests for square matrices with same maps only */ 1119566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(A, &squaretest)); 112e432b41dSStefano Zampini if (squaretest && rmap != cmap) { 113e432b41dSStefano Zampini PetscInt nr, nc; 114e432b41dSStefano Zampini 1159566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(rmap, &nr)); 1169566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(cmap, &nc)); 117e432b41dSStefano Zampini if (nr != nc) squaretest = PETSC_FALSE; 118e432b41dSStefano Zampini else { 119e432b41dSStefano Zampini const PetscInt *idxs1, *idxs2; 120e432b41dSStefano Zampini 1219566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(rmap, &idxs1)); 1229566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(cmap, &idxs2)); 1239566063dSJacob Faibussowitsch PetscCall(PetscArraycmp(idxs1, idxs2, nr, &squaretest)); 1249566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(rmap, &idxs1)); 1259566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(cmap, &idxs2)); 126e432b41dSStefano Zampini } 1271c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &squaretest, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A))); 128e432b41dSStefano Zampini } 129e432b41dSStefano Zampini 130e432b41dSStefano Zampini /* test MatISGetLocalMat */ 1319566063dSJacob Faibussowitsch PetscCall(MatISGetLocalMat(A, &B)); 1329566063dSJacob Faibussowitsch PetscCall(MatGetType(B, &lmtype)); 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* test MatGetInfo */ 1359566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatGetInfo\n")); 1369566063dSJacob Faibussowitsch PetscCall(MatGetInfo(A, MAT_LOCAL, &info)); 1379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 1389371c9d4SSatish 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, 1399371c9d4SSatish Balay (PetscInt)info.nz_unneeded, (PetscInt)info.assemblies, (PetscInt)info.mallocs)); 1409566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 1419566063dSJacob Faibussowitsch PetscCall(MatGetInfo(A, MAT_GLOBAL_MAX, &info)); 1429371c9d4SSatish 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, 1439371c9d4SSatish Balay (PetscInt)info.assemblies, (PetscInt)info.mallocs)); 1449566063dSJacob Faibussowitsch PetscCall(MatGetInfo(A, MAT_GLOBAL_SUM, &info)); 1459371c9d4SSatish 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, 1469371c9d4SSatish Balay (PetscInt)info.assemblies, (PetscInt)info.mallocs)); 147c4762a1bSJed Brown 148e432b41dSStefano Zampini /* test MatIsSymmetric */ 1499566063dSJacob Faibussowitsch PetscCall(MatIsSymmetric(A, 0.0, &issymmetric)); 1509566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatIsSymmetric: %d\n", issymmetric)); 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* Create a MPIAIJ matrix, same as A */ 1539566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 1549566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n)); 1559566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATAIJ)); 1569566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 1579566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 1589566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(B, rmap, cmap)); 1599566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(B, 3, NULL, 3, NULL)); 1609566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(B, 1, 3, NULL, 3, NULL)); 161c4762a1bSJed Brown #if defined(PETSC_HAVE_HYPRE) 1629566063dSJacob Faibussowitsch PetscCall(MatHYPRESetPreallocation(B, 3, NULL, 3, NULL)); 163c4762a1bSJed Brown #endif 1649566063dSJacob Faibussowitsch PetscCall(MatISSetPreallocation(B, 3, NULL, 3, NULL)); 1659566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, (PetscBool) !(repmap || negmap))); /* I do not want to precompute the pattern */ 166e432b41dSStefano Zampini for (i = 0; i < lm; i++) { 167c4762a1bSJed Brown PetscScalar v[3]; 168c4762a1bSJed Brown PetscInt cols[3]; 169c4762a1bSJed Brown 170c4762a1bSJed Brown cols[0] = (i - 1 + n) % n; 171c4762a1bSJed Brown cols[1] = i % n; 172c4762a1bSJed Brown cols[2] = (i + 1) % n; 173c4762a1bSJed Brown v[0] = -1. * (symmetric ? PetscMin(i + 1, cols[0] + 1) : i + 1); 174c4762a1bSJed Brown v[1] = 2. * (symmetric ? PetscMin(i + 1, cols[1] + 1) : i + 1); 175c4762a1bSJed Brown v[2] = -1. * (symmetric ? PetscMin(i + 1, cols[2] + 1) : i + 1); 1769566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(cmap, IS_GTOLM_MASK, 3, cols, NULL, cols)); 1779566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(B, 1, &i, 3, cols, v, ADD_VALUES)); 178c4762a1bSJed Brown } 1799566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 1809566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 181e432b41dSStefano Zampini 182e432b41dSStefano Zampini /* test MatView */ 1839566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatView\n")); 1849566063dSJacob Faibussowitsch PetscCall(MatView(A, NULL)); 1859566063dSJacob Faibussowitsch PetscCall(MatView(B, NULL)); 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* test CheckMat */ 1889566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test CheckMat\n")); 1899566063dSJacob Faibussowitsch PetscCall(CheckMat(A, B, PETSC_FALSE, "CheckMat")); 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* test MatDuplicate and MatAXPY */ 1929566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatDuplicate and MatAXPY\n")); 1939566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2)); 1949566063dSJacob Faibussowitsch PetscCall(CheckMat(A, A2, PETSC_FALSE, "MatDuplicate and MatAXPY")); 195c4762a1bSJed Brown 196c4762a1bSJed Brown /* test MatConvert */ 1979566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatConvert_IS_XAIJ\n")); 1989566063dSJacob Faibussowitsch PetscCall(MatConvert(A2, MATAIJ, MAT_INITIAL_MATRIX, &B2)); 1999566063dSJacob Faibussowitsch PetscCall(CheckMat(B, B2, PETSC_TRUE, "MatConvert_IS_XAIJ MAT_INITIAL_MATRIX")); 2009566063dSJacob Faibussowitsch PetscCall(MatConvert(A2, MATAIJ, MAT_REUSE_MATRIX, &B2)); 2019566063dSJacob Faibussowitsch PetscCall(CheckMat(B, B2, PETSC_TRUE, "MatConvert_IS_XAIJ MAT_REUSE_MATRIX")); 2029566063dSJacob Faibussowitsch PetscCall(MatConvert(A2, MATAIJ, MAT_INPLACE_MATRIX, &A2)); 2039566063dSJacob Faibussowitsch PetscCall(CheckMat(B, A2, PETSC_TRUE, "MatConvert_IS_XAIJ MAT_INPLACE_MATRIX")); 2049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 2059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 2069566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatConvert_XAIJ_IS\n")); 2079566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &B2)); 2089566063dSJacob Faibussowitsch PetscCall(MatConvert(B2, MATIS, MAT_INITIAL_MATRIX, &A2)); 2099566063dSJacob Faibussowitsch PetscCall(CheckMat(A, A2, PETSC_TRUE, "MatConvert_XAIJ_IS MAT_INITIAL_MATRIX")); 2109566063dSJacob Faibussowitsch PetscCall(MatConvert(B2, MATIS, MAT_REUSE_MATRIX, &A2)); 2119566063dSJacob Faibussowitsch PetscCall(CheckMat(A, A2, PETSC_TRUE, "MatConvert_XAIJ_IS MAT_REUSE_MATRIX")); 2129566063dSJacob Faibussowitsch PetscCall(MatConvert(B2, MATIS, MAT_INPLACE_MATRIX, &B2)); 2139566063dSJacob Faibussowitsch PetscCall(CheckMat(A, B2, PETSC_TRUE, "MatConvert_XAIJ_IS MAT_INPLACE_MATRIX")); 2149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 2159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 2169566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(lmtype, MATSEQAIJ, &isaij)); 217c4762a1bSJed Brown if (size == 1 && isaij) { /* tests special code paths in MatConvert_IS_XAIJ */ 218c4762a1bSJed 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}; 219e432b41dSStefano Zampini ISLocalToGlobalMapping tcmap, trmap; 220c4762a1bSJed Brown 221c4762a1bSJed Brown for (ri = 0; ri < 2; ri++) { 222c4762a1bSJed Brown PetscInt *r; 223c4762a1bSJed Brown 224c4762a1bSJed Brown r = (PetscInt *)(ri == 0 ? rr : rk); 225c4762a1bSJed Brown for (ci = 0; ci < 2; ci++) { 226c4762a1bSJed Brown PetscInt *c, rb, cb; 227c4762a1bSJed Brown 228c4762a1bSJed Brown c = (PetscInt *)(ci == 0 ? cr : ck); 229c4762a1bSJed Brown for (rb = 1; rb < 4; rb++) { 2309566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, rb, 3, r, PETSC_COPY_VALUES, &is)); 2319566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(is, &trmap)); 2329566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 233c4762a1bSJed Brown for (cb = 1; cb < 4; cb++) { 234c4762a1bSJed Brown Mat T, lT, T2; 235c4762a1bSJed Brown char testname[256]; 236c4762a1bSJed Brown 2379566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(testname, sizeof(testname), "MatConvert_IS_XAIJ special case (%" PetscInt_FMT " %" PetscInt_FMT ", bs %" PetscInt_FMT " %" PetscInt_FMT ")", ri, ci, rb, cb)); 2389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test %s\n", testname)); 239c4762a1bSJed Brown 2409566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, cb, 4, c, PETSC_COPY_VALUES, &is)); 2419566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(is, &tcmap)); 2429566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 243c4762a1bSJed Brown 2449566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &T)); 2459566063dSJacob Faibussowitsch PetscCall(MatSetSizes(T, PETSC_DECIDE, PETSC_DECIDE, rb * 3, cb * 4)); 2469566063dSJacob Faibussowitsch PetscCall(MatSetType(T, MATIS)); 2479566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(T, trmap, tcmap)); 2489566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&tcmap)); 2499566063dSJacob Faibussowitsch PetscCall(MatISGetLocalMat(T, &lT)); 2509566063dSJacob Faibussowitsch PetscCall(MatSetType(lT, MATSEQAIJ)); 2519566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(lT, cb * 4, NULL)); 2529566063dSJacob Faibussowitsch PetscCall(MatSetRandom(lT, NULL)); 2539566063dSJacob Faibussowitsch PetscCall(MatConvert(lT, lmtype, MAT_INPLACE_MATRIX, &lT)); 2549566063dSJacob Faibussowitsch PetscCall(MatISRestoreLocalMat(T, &lT)); 2559566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(T, MAT_FINAL_ASSEMBLY)); 2569566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(T, MAT_FINAL_ASSEMBLY)); 257c4762a1bSJed Brown 2589566063dSJacob Faibussowitsch PetscCall(MatConvert(T, MATAIJ, MAT_INITIAL_MATRIX, &T2)); 2599566063dSJacob Faibussowitsch PetscCall(CheckMat(T, T2, PETSC_TRUE, "MAT_INITIAL_MATRIX")); 2609566063dSJacob Faibussowitsch PetscCall(MatConvert(T, MATAIJ, MAT_REUSE_MATRIX, &T2)); 2619566063dSJacob Faibussowitsch PetscCall(CheckMat(T, T2, PETSC_TRUE, "MAT_REUSE_MATRIX")); 2629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T2)); 2639566063dSJacob Faibussowitsch PetscCall(MatDuplicate(T, MAT_COPY_VALUES, &T2)); 2649566063dSJacob Faibussowitsch PetscCall(MatConvert(T2, MATAIJ, MAT_INPLACE_MATRIX, &T2)); 2659566063dSJacob Faibussowitsch PetscCall(CheckMat(T, T2, PETSC_TRUE, "MAT_INPLACE_MATRIX")); 2669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T)); 2679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T2)); 268c4762a1bSJed Brown } 2699566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&trmap)); 270c4762a1bSJed Brown } 271c4762a1bSJed Brown } 272c4762a1bSJed Brown } 273c4762a1bSJed Brown } 274c4762a1bSJed Brown 275c4762a1bSJed Brown /* test MatDiagonalScale */ 2769566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatDiagonalScale\n")); 2779566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2)); 2789566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &B2)); 2799566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, &y)); 2809566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, NULL)); 281e432b41dSStefano Zampini if (issymmetric) { 2829566063dSJacob Faibussowitsch PetscCall(VecCopy(x, y)); 283c4762a1bSJed Brown } else { 2849566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y, NULL)); 2859566063dSJacob Faibussowitsch PetscCall(VecScale(y, 8.)); 286c4762a1bSJed Brown } 2879566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A2, y, x)); 2889566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(B2, y, x)); 2899566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatDiagonalScale")); 2909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 2919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 2929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 294c4762a1bSJed Brown 295c4762a1bSJed Brown /* test MatPtAP (A IS and B AIJ) */ 296c4762a1bSJed Brown if (isaij && m == n) { 2979566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatPtAP\n")); 2989566063dSJacob Faibussowitsch PetscCall(MatISStoreL2L(A, PETSC_TRUE)); 2999566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &A2)); 3009566063dSJacob Faibussowitsch PetscCall(MatPtAP(B, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &B2)); 3019566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatPtAP MAT_INITIAL_MATRIX")); 3029566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &A2)); 3039566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatPtAP MAT_REUSE_MATRIX")); 3049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 3059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 306c4762a1bSJed Brown } 307c4762a1bSJed Brown 308c4762a1bSJed Brown /* test MatGetLocalSubMatrix */ 3099566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatGetLocalSubMatrix\n")); 3109566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_DO_NOT_COPY_VALUES, &A2)); 3119566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, lm / 2 + lm % 2, 0, 2, &reven)); 3129566063dSJacob Faibussowitsch PetscCall(ISComplement(reven, 0, lm, &rodd)); 3139566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, ln / 2 + ln % 2, 0, 2, &ceven)); 3149566063dSJacob Faibussowitsch PetscCall(ISComplement(ceven, 0, ln, &codd)); 3159566063dSJacob Faibussowitsch PetscCall(MatGetLocalSubMatrix(A2, reven, ceven, &Aee)); 3169566063dSJacob Faibussowitsch PetscCall(MatGetLocalSubMatrix(A2, reven, codd, &Aeo)); 3179566063dSJacob Faibussowitsch PetscCall(MatGetLocalSubMatrix(A2, rodd, ceven, &Aoe)); 3189566063dSJacob Faibussowitsch PetscCall(MatGetLocalSubMatrix(A2, rodd, codd, &Aoo)); 319e432b41dSStefano Zampini for (i = 0; i < lm; i++) { 320c4762a1bSJed Brown PetscInt j, je, jo, colse[3], colso[3]; 321c4762a1bSJed Brown PetscScalar ve[3], vo[3]; 322c4762a1bSJed Brown PetscScalar v[3]; 323c4762a1bSJed Brown PetscInt cols[3]; 324e432b41dSStefano Zampini PetscInt row = i / 2; 325c4762a1bSJed Brown 326c4762a1bSJed Brown cols[0] = (i - 1 + n) % n; 327c4762a1bSJed Brown cols[1] = i % n; 328c4762a1bSJed Brown cols[2] = (i + 1) % n; 329c4762a1bSJed Brown v[0] = -1. * (symmetric ? PetscMin(i + 1, cols[0] + 1) : i + 1); 330c4762a1bSJed Brown v[1] = 2. * (symmetric ? PetscMin(i + 1, cols[1] + 1) : i + 1); 331c4762a1bSJed Brown v[2] = -1. * (symmetric ? PetscMin(i + 1, cols[2] + 1) : i + 1); 3329566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(cmap, IS_GTOLM_MASK, 3, cols, NULL, cols)); 333c4762a1bSJed Brown for (j = 0, je = 0, jo = 0; j < 3; j++) { 334c4762a1bSJed Brown if (cols[j] % 2) { 335c4762a1bSJed Brown vo[jo] = v[j]; 336c4762a1bSJed Brown colso[jo++] = cols[j] / 2; 337c4762a1bSJed Brown } else { 338c4762a1bSJed Brown ve[je] = v[j]; 339c4762a1bSJed Brown colse[je++] = cols[j] / 2; 340c4762a1bSJed Brown } 341c4762a1bSJed Brown } 342c4762a1bSJed Brown if (i % 2) { 3439566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(Aoe, 1, &row, je, colse, ve, ADD_VALUES)); 3449566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(Aoo, 1, &row, jo, colso, vo, ADD_VALUES)); 345c4762a1bSJed Brown } else { 3469566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(Aee, 1, &row, je, colse, ve, ADD_VALUES)); 3479566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(Aeo, 1, &row, jo, colso, vo, ADD_VALUES)); 348c4762a1bSJed Brown } 349c4762a1bSJed Brown } 3509566063dSJacob Faibussowitsch PetscCall(MatRestoreLocalSubMatrix(A2, reven, ceven, &Aee)); 3519566063dSJacob Faibussowitsch PetscCall(MatRestoreLocalSubMatrix(A2, reven, codd, &Aeo)); 3529566063dSJacob Faibussowitsch PetscCall(MatRestoreLocalSubMatrix(A2, rodd, ceven, &Aoe)); 3539566063dSJacob Faibussowitsch PetscCall(MatRestoreLocalSubMatrix(A2, rodd, codd, &Aoo)); 3549566063dSJacob Faibussowitsch PetscCall(ISDestroy(&reven)); 3559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ceven)); 3569566063dSJacob Faibussowitsch PetscCall(ISDestroy(&rodd)); 3579566063dSJacob Faibussowitsch PetscCall(ISDestroy(&codd)); 3589566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A2, MAT_FINAL_ASSEMBLY)); 3599566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A2, MAT_FINAL_ASSEMBLY)); 3609566063dSJacob Faibussowitsch PetscCall(MatAXPY(A2, -1., A, SAME_NONZERO_PATTERN)); 3619566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, NULL, PETSC_FALSE, "MatGetLocalSubMatrix")); 3629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 363c4762a1bSJed Brown 364c4762a1bSJed Brown /* test MatConvert_Nest_IS */ 365c4762a1bSJed Brown testT = PETSC_FALSE; 3669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_trans", &testT, NULL)); 367c4762a1bSJed Brown 3689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatConvert_Nest_IS\n")); 369c4762a1bSJed Brown nr = 2; 370c4762a1bSJed Brown nc = 2; 3719566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nr", &nr, NULL)); 3729566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nc", &nc, NULL)); 373c4762a1bSJed Brown if (testT) { 3749566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &cst, &cen)); 3759566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(A, &rst, &ren)); 376c4762a1bSJed Brown } else { 3779566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rst, &ren)); 3789566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(A, &cst, &cen)); 379c4762a1bSJed Brown } 3809566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(nr, &rows, nc, &cols, 2 * nr * nc, &mats)); 381c4762a1bSJed Brown for (i = 0; i < nr * nc; i++) { 382c4762a1bSJed Brown if (testT) { 3839566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(A, &mats[i])); 3849566063dSJacob Faibussowitsch PetscCall(MatTranspose(B, MAT_INITIAL_MATRIX, &mats[i + nr * nc])); 385c4762a1bSJed Brown } else { 3869566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &mats[i])); 3879566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &mats[i + nr * nc])); 388c4762a1bSJed Brown } 389c4762a1bSJed Brown } 39048a46eb9SPierre Jolivet for (i = 0; i < nr; i++) PetscCall(ISCreateStride(PETSC_COMM_WORLD, ren - rst, i + rst, nr, &rows[i])); 39148a46eb9SPierre Jolivet for (i = 0; i < nc; i++) PetscCall(ISCreateStride(PETSC_COMM_WORLD, cen - cst, i + cst, nc, &cols[i])); 3929566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PETSC_COMM_WORLD, nr, rows, nc, cols, mats, &A2)); 3939566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PETSC_COMM_WORLD, nr, rows, nc, cols, mats + nr * nc, &B2)); 39448a46eb9SPierre Jolivet for (i = 0; i < nr; i++) PetscCall(ISDestroy(&rows[i])); 39548a46eb9SPierre Jolivet for (i = 0; i < nc; i++) PetscCall(ISDestroy(&cols[i])); 39648a46eb9SPierre Jolivet for (i = 0; i < 2 * nr * nc; i++) PetscCall(MatDestroy(&mats[i])); 3979566063dSJacob Faibussowitsch PetscCall(PetscFree3(rows, cols, mats)); 3989566063dSJacob Faibussowitsch PetscCall(MatConvert(B2, MATAIJ, MAT_INITIAL_MATRIX, &T)); 3999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 4009566063dSJacob Faibussowitsch PetscCall(MatConvert(A2, MATIS, MAT_INITIAL_MATRIX, &B2)); 4019566063dSJacob Faibussowitsch PetscCall(CheckMat(B2, T, PETSC_TRUE, "MatConvert_Nest_IS MAT_INITIAL_MATRIX")); 4029566063dSJacob Faibussowitsch PetscCall(MatConvert(A2, MATIS, MAT_REUSE_MATRIX, &B2)); 4039566063dSJacob Faibussowitsch PetscCall(CheckMat(B2, T, PETSC_TRUE, "MatConvert_Nest_IS MAT_REUSE_MATRIX")); 4049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 4059566063dSJacob Faibussowitsch PetscCall(MatConvert(A2, MATIS, MAT_INPLACE_MATRIX, &A2)); 4069566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, T, PETSC_TRUE, "MatConvert_Nest_IS MAT_INPLACE_MATRIX")); 4079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T)); 4089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 409c4762a1bSJed Brown 410c4762a1bSJed Brown /* test MatCreateSubMatrix */ 4119566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatCreateSubMatrix\n")); 412dd400576SPatrick Sanan if (rank == 0) { 4139566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, 1, 1, 1, &is)); 4149566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, 2, 0, 1, &is2)); 415c4762a1bSJed Brown } else if (rank == 1) { 4169566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, 1, 0, 1, &is)); 417c4762a1bSJed Brown if (n > 3) { 4189566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, 1, 3, 1, &is2)); 419c4762a1bSJed Brown } else { 4209566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, 0, 0, 1, &is2)); 421c4762a1bSJed Brown } 422c4762a1bSJed Brown } else if (rank == 2 && n > 4) { 4239566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, 0, 0, 1, &is)); 4249566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, n - 4, 4, 1, &is2)); 425c4762a1bSJed Brown } else { 4269566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, 0, 0, 1, &is)); 4279566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, 0, 0, 1, &is2)); 428c4762a1bSJed Brown } 4299566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, is, is, MAT_INITIAL_MATRIX, &A2)); 4309566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(B, is, is, MAT_INITIAL_MATRIX, &B2)); 4319566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, B2, PETSC_TRUE, "first MatCreateSubMatrix")); 432c4762a1bSJed Brown 4339566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, is, is, MAT_REUSE_MATRIX, &A2)); 4349566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(B, is, is, MAT_REUSE_MATRIX, &B2)); 4359566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, B2, PETSC_FALSE, "reuse MatCreateSubMatrix")); 4369566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 4379566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 438c4762a1bSJed Brown 439e432b41dSStefano Zampini if (!issymmetric) { 4409566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, is, is2, MAT_INITIAL_MATRIX, &A2)); 4419566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(B, is, is2, MAT_INITIAL_MATRIX, &B2)); 4429566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, is, is2, MAT_REUSE_MATRIX, &A2)); 4439566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(B, is, is2, MAT_REUSE_MATRIX, &B2)); 4449566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, B2, PETSC_FALSE, "second MatCreateSubMatrix")); 445c4762a1bSJed Brown } 446c4762a1bSJed Brown 4479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 4489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 4499566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 4509566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2)); 451c4762a1bSJed Brown 452d0dbe9f7SStefano Zampini /* test MatCreateSubMatrices */ 453d0dbe9f7SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatCreateSubMatrices\n")); 454d0dbe9f7SStefano Zampini PetscCall(MatGetLayouts(A, &rlayout, &clayout)); 455d0dbe9f7SStefano Zampini PetscCall(PetscLayoutGetRanges(rlayout, &rrange)); 456d0dbe9f7SStefano Zampini PetscCall(PetscLayoutGetRanges(clayout, &crange)); 457d0dbe9f7SStefano Zampini lrank = (size + rank - 1) % size; 458d0dbe9f7SStefano Zampini rrank = (rank + 1) % size; 459d0dbe9f7SStefano Zampini PetscCall(ISCreateStride(PETSC_COMM_SELF, (rrange[lrank + 1] - rrange[lrank]), rrange[lrank], 1, &irow[0])); 460d0dbe9f7SStefano Zampini PetscCall(ISCreateStride(PETSC_COMM_SELF, (crange[rrank + 1] - crange[rrank]), crange[rrank], 1, &icol[0])); 461d0dbe9f7SStefano Zampini PetscCall(ISCreateStride(PETSC_COMM_SELF, (rrange[rrank + 1] - rrange[rrank]), rrange[rrank], 1, &irow[1])); 462d0dbe9f7SStefano Zampini PetscCall(ISCreateStride(PETSC_COMM_SELF, (crange[lrank + 1] - crange[lrank]), crange[lrank], 1, &icol[1])); 463d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(A, 2, irow, icol, MAT_INITIAL_MATRIX, &Asub)); 464d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(B, 2, irow, icol, MAT_INITIAL_MATRIX, &Bsub)); 465d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[0], Bsub[0], PETSC_FALSE, "MatCreateSubMatrices[0]")); 466d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[1], Bsub[1], PETSC_FALSE, "MatCreateSubMatrices[1]")); 467d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(A, 2, irow, icol, MAT_REUSE_MATRIX, &Asub)); 468d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(B, 2, irow, icol, MAT_REUSE_MATRIX, &Bsub)); 469d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[0], Bsub[0], PETSC_FALSE, "MatCreateSubMatrices[0]")); 470d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[1], Bsub[1], PETSC_FALSE, "MatCreateSubMatrices[1]")); 471d0dbe9f7SStefano Zampini PetscCall(MatDestroySubMatrices(2, &Asub)); 472d0dbe9f7SStefano Zampini PetscCall(MatDestroySubMatrices(2, &Bsub)); 473d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&irow[0])); 474d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&irow[1])); 475d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&icol[0])); 476d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&icol[1])); 477d0dbe9f7SStefano Zampini 478c4762a1bSJed Brown /* Create an IS required by MatZeroRows(): just rank zero provides the rows to be eliminated */ 479c4762a1bSJed Brown if (size > 1) { 480dd400576SPatrick Sanan if (rank == 0) { 481c4762a1bSJed Brown PetscInt st, len; 482c4762a1bSJed Brown 483c4762a1bSJed Brown st = (m + 1) / 2; 484c4762a1bSJed Brown len = PetscMin(m / 2, PetscMax(m - (m + 1) / 2 - 1, 0)); 4859566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, len, st, 1, &is)); 486c4762a1bSJed Brown } else { 4879566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, 0, 0, 1, &is)); 488c4762a1bSJed Brown } 489c4762a1bSJed Brown } else { 4909566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, 1, 0, 1, &is)); 491c4762a1bSJed Brown } 492c4762a1bSJed Brown 493c4762a1bSJed Brown if (squaretest) { /* tests for square matrices only, with same maps for rows and columns */ 494d0dbe9f7SStefano Zampini PetscInt *idx0, *idx1, n0, n1; 495d0dbe9f7SStefano Zampini IS Ais[2], Bis[2]; 496d0dbe9f7SStefano Zampini 497c4762a1bSJed Brown /* test MatDiagonalSet */ 4989566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatDiagonalSet\n")); 4999566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2)); 5009566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &B2)); 5019566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, NULL, &x)); 5029566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, NULL)); 5039566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(A2, x, INSERT_VALUES)); 5049566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(B2, x, INSERT_VALUES)); 5059566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatDiagonalSet")); 5069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 5079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 5089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 509c4762a1bSJed Brown 510c4762a1bSJed Brown /* test MatShift (MatShift_IS internally uses MatDiagonalSet_IS with ADD_VALUES) */ 5119566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatShift\n")); 5129566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2)); 5139566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &B2)); 5149566063dSJacob Faibussowitsch PetscCall(MatShift(A2, 2.0)); 5159566063dSJacob Faibussowitsch PetscCall(MatShift(B2, 2.0)); 5169566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatShift")); 5179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 5189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 519c4762a1bSJed Brown 520c4762a1bSJed Brown /* nonzero diag value is supported for square matrices only */ 5219566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows(A, B, PETSC_TRUE, is, diag)); 522d0dbe9f7SStefano Zampini 523d0dbe9f7SStefano Zampini /* test MatIncreaseOverlap */ 524d0dbe9f7SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatIncreaseOverlap\n")); 525d0dbe9f7SStefano Zampini PetscCall(MatGetOwnershipRange(A, &rst, &ren)); 526d0dbe9f7SStefano Zampini n0 = (ren - rst) / 2; 527d0dbe9f7SStefano Zampini n1 = (ren - rst) / 3; 528d0dbe9f7SStefano Zampini PetscCall(PetscMalloc1(n0, &idx0)); 529d0dbe9f7SStefano Zampini PetscCall(PetscMalloc1(n1, &idx1)); 530d0dbe9f7SStefano Zampini for (i = 0; i < n0; i++) idx0[i] = ren - i - 1; 531d0dbe9f7SStefano Zampini for (i = 0; i < n1; i++) idx1[i] = rst + i; 532d0dbe9f7SStefano Zampini PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n0, idx0, PETSC_OWN_POINTER, &Ais[0])); 533d0dbe9f7SStefano Zampini PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n1, idx1, PETSC_OWN_POINTER, &Ais[1])); 534d0dbe9f7SStefano Zampini PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n0, idx0, PETSC_COPY_VALUES, &Bis[0])); 535d0dbe9f7SStefano Zampini PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n1, idx1, PETSC_COPY_VALUES, &Bis[1])); 536d0dbe9f7SStefano Zampini PetscCall(MatIncreaseOverlap(A, 2, Ais, 3)); 537d0dbe9f7SStefano Zampini PetscCall(MatIncreaseOverlap(B, 2, Bis, 3)); 538d0dbe9f7SStefano Zampini /* Non deterministic output! */ 539d0dbe9f7SStefano Zampini PetscCall(ISSort(Ais[0])); 540d0dbe9f7SStefano Zampini PetscCall(ISSort(Ais[1])); 541d0dbe9f7SStefano Zampini PetscCall(ISSort(Bis[0])); 542d0dbe9f7SStefano Zampini PetscCall(ISSort(Bis[1])); 543d0dbe9f7SStefano Zampini PetscCall(ISView(Ais[0], NULL)); 544d0dbe9f7SStefano Zampini PetscCall(ISView(Bis[0], NULL)); 545d0dbe9f7SStefano Zampini PetscCall(ISView(Ais[1], NULL)); 546d0dbe9f7SStefano Zampini PetscCall(ISView(Bis[1], NULL)); 547d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(A, 2, Ais, Ais, MAT_INITIAL_MATRIX, &Asub)); 548d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(B, 2, Bis, Ais, MAT_INITIAL_MATRIX, &Bsub)); 549d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[0], Bsub[0], PETSC_FALSE, "MatIncreaseOverlap[0]")); 550d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[1], Bsub[1], PETSC_FALSE, "MatIncreaseOverlap[1]")); 551d0dbe9f7SStefano Zampini PetscCall(MatDestroySubMatrices(2, &Asub)); 552d0dbe9f7SStefano Zampini PetscCall(MatDestroySubMatrices(2, &Bsub)); 553d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&Ais[0])); 554d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&Ais[1])); 555d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&Bis[0])); 556d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&Bis[1])); 557c4762a1bSJed Brown } 5589566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows(A, B, squaretest, is, 0.0)); 5599566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 560c4762a1bSJed Brown 561c4762a1bSJed Brown /* test MatTranspose */ 5629566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatTranspose\n")); 5639566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &A2)); 5649566063dSJacob Faibussowitsch PetscCall(MatTranspose(B, MAT_INITIAL_MATRIX, &B2)); 5659566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, B2, PETSC_FALSE, "initial matrix MatTranspose")); 566c4762a1bSJed Brown 5679566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_REUSE_MATRIX, &A2)); 5689566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, B2, PETSC_FALSE, "reuse matrix (not in place) MatTranspose")); 5699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 570c4762a1bSJed Brown 5719566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2)); 5729566063dSJacob Faibussowitsch PetscCall(MatTranspose(A2, MAT_INPLACE_MATRIX, &A2)); 5739566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, B2, PETSC_FALSE, "reuse matrix (in place) MatTranspose")); 5749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 575c4762a1bSJed Brown 5769566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &A2)); 5779566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, B2, PETSC_FALSE, "reuse matrix (different type) MatTranspose")); 5789566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 5799566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 580c4762a1bSJed Brown 581c4762a1bSJed Brown /* test MatISFixLocalEmpty */ 582c4762a1bSJed Brown if (isaij) { 583c4762a1bSJed Brown PetscInt r[2]; 584c4762a1bSJed Brown 585c4762a1bSJed Brown r[0] = 0; 586c4762a1bSJed Brown r[1] = PetscMin(m, n) - 1; 5879566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatISFixLocalEmpty\n")); 5889566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2)); 589e432b41dSStefano Zampini 5909566063dSJacob Faibussowitsch PetscCall(MatISFixLocalEmpty(A2, PETSC_TRUE)); 5919566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A2, MAT_FINAL_ASSEMBLY)); 5929566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A2, MAT_FINAL_ASSEMBLY)); 5939566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, B, PETSC_FALSE, "MatISFixLocalEmpty (null)")); 594c4762a1bSJed Brown 5959566063dSJacob Faibussowitsch PetscCall(MatZeroRows(A2, 2, r, 0.0, NULL, NULL)); 5969566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2, NULL, "-fixempty_view")); 5979566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &B2)); 5989566063dSJacob Faibussowitsch PetscCall(MatZeroRows(B2, 2, r, 0.0, NULL, NULL)); 5999566063dSJacob Faibussowitsch PetscCall(MatISFixLocalEmpty(A2, PETSC_TRUE)); 6009566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A2, MAT_FINAL_ASSEMBLY)); 6019566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A2, MAT_FINAL_ASSEMBLY)); 6029566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2, NULL, "-fixempty_view")); 6039566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatISFixLocalEmpty (rows)")); 6049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 605c4762a1bSJed Brown 6069566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2)); 6079566063dSJacob Faibussowitsch PetscCall(MatZeroRows(A2, 2, r, 0.0, NULL, NULL)); 6089566063dSJacob Faibussowitsch PetscCall(MatTranspose(A2, MAT_INPLACE_MATRIX, &A2)); 6099566063dSJacob Faibussowitsch PetscCall(MatTranspose(B2, MAT_INPLACE_MATRIX, &B2)); 6109566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2, NULL, "-fixempty_view")); 6119566063dSJacob Faibussowitsch PetscCall(MatISFixLocalEmpty(A2, PETSC_TRUE)); 6129566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A2, MAT_FINAL_ASSEMBLY)); 6139566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A2, MAT_FINAL_ASSEMBLY)); 6149566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2, NULL, "-fixempty_view")); 6159566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatISFixLocalEmpty (cols)")); 616c4762a1bSJed Brown 6179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 6189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 619c4762a1bSJed Brown 620c4762a1bSJed Brown if (squaretest) { 6219566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2)); 6229566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &B2)); 6239566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumns(A2, 2, r, 0.0, NULL, NULL)); 6249566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumns(B2, 2, r, 0.0, NULL, NULL)); 6259566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2, NULL, "-fixempty_view")); 6269566063dSJacob Faibussowitsch PetscCall(MatISFixLocalEmpty(A2, PETSC_TRUE)); 6279566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A2, MAT_FINAL_ASSEMBLY)); 6289566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A2, MAT_FINAL_ASSEMBLY)); 6299566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2, NULL, "-fixempty_view")); 6309566063dSJacob Faibussowitsch PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatISFixLocalEmpty (rows+cols)")); 6319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 6329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 633c4762a1bSJed Brown } 634c4762a1bSJed Brown } 635c4762a1bSJed Brown 636c4762a1bSJed Brown /* test MatInvertBlockDiagonal 637c4762a1bSJed Brown special cases for block-diagonal matrices */ 638c4762a1bSJed Brown if (m == n) { 639c4762a1bSJed Brown ISLocalToGlobalMapping map; 640c4762a1bSJed Brown Mat Abd, Bbd; 641c4762a1bSJed Brown IS is, bis; 642c4762a1bSJed Brown const PetscScalar *isbd, *aijbd; 643c4762a1bSJed Brown PetscScalar *vals; 644c4762a1bSJed Brown const PetscInt *sts, *idxs; 645c4762a1bSJed Brown PetscInt *idxs2, diff, perm, nl, bs, st, en, in; 646c4762a1bSJed Brown PetscBool ok; 647c4762a1bSJed Brown 648c4762a1bSJed Brown for (diff = 0; diff < 3; diff++) { 649c4762a1bSJed Brown for (perm = 0; perm < 3; perm++) { 650c4762a1bSJed Brown for (bs = 1; bs < 4; bs++) { 6519566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatInvertBlockDiagonal blockdiag %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", n, diff, perm, bs)); 6529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * bs, &vals)); 6539566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(A, &sts)); 654c4762a1bSJed Brown switch (diff) { 655c4762a1bSJed Brown case 1: /* inverted layout by processes */ 656c4762a1bSJed Brown in = 1; 657c4762a1bSJed Brown st = sts[size - rank - 1]; 658c4762a1bSJed Brown en = sts[size - rank]; 659c4762a1bSJed Brown nl = en - st; 660c4762a1bSJed Brown break; 661c4762a1bSJed Brown case 2: /* round-robin layout */ 662c4762a1bSJed Brown in = size; 663c4762a1bSJed Brown st = rank; 664c4762a1bSJed Brown nl = n / size; 665c4762a1bSJed Brown if (rank < n % size) nl++; 666c4762a1bSJed Brown break; 667c4762a1bSJed Brown default: /* same layout */ 668c4762a1bSJed Brown in = 1; 669c4762a1bSJed Brown st = sts[rank]; 670c4762a1bSJed Brown en = sts[rank + 1]; 671c4762a1bSJed Brown nl = en - st; 672c4762a1bSJed Brown break; 673c4762a1bSJed Brown } 6749566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, nl, st, in, &is)); 6759566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &nl)); 6769566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &idxs)); 6779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nl, &idxs2)); 678c4762a1bSJed Brown for (i = 0; i < nl; i++) { 679c4762a1bSJed Brown switch (perm) { /* invert some of the indices */ 680d71ae5a4SJacob Faibussowitsch case 2: 681d71ae5a4SJacob Faibussowitsch idxs2[i] = rank % 2 ? idxs[i] : idxs[nl - i - 1]; 682d71ae5a4SJacob Faibussowitsch break; 683d71ae5a4SJacob Faibussowitsch case 1: 684d71ae5a4SJacob Faibussowitsch idxs2[i] = rank % 2 ? idxs[nl - i - 1] : idxs[i]; 685d71ae5a4SJacob Faibussowitsch break; 686d71ae5a4SJacob Faibussowitsch default: 687d71ae5a4SJacob Faibussowitsch idxs2[i] = idxs[i]; 688d71ae5a4SJacob Faibussowitsch break; 689c4762a1bSJed Brown } 690c4762a1bSJed Brown } 6919566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &idxs)); 6929566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_WORLD, bs, nl, idxs2, PETSC_OWN_POINTER, &bis)); 6939566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(bis, &map)); 6949566063dSJacob Faibussowitsch PetscCall(MatCreateIS(PETSC_COMM_WORLD, bs, PETSC_DECIDE, PETSC_DECIDE, bs * n, bs * n, map, map, &Abd)); 6959566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&map)); 6969566063dSJacob Faibussowitsch PetscCall(MatISSetPreallocation(Abd, bs, NULL, 0, NULL)); 697c4762a1bSJed Brown for (i = 0; i < nl; i++) { 698c4762a1bSJed Brown PetscInt b1, b2; 699c4762a1bSJed Brown 700c4762a1bSJed Brown for (b1 = 0; b1 < bs; b1++) { 701ad540459SPierre Jolivet for (b2 = 0; b2 < bs; b2++) vals[b1 * bs + b2] = i * bs * bs + b1 * bs + b2 + 1 + (b1 == b2 ? 1.0 : 0); 702c4762a1bSJed Brown } 7039566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedLocal(Abd, 1, &i, 1, &i, vals, INSERT_VALUES)); 704c4762a1bSJed Brown } 7059566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Abd, MAT_FINAL_ASSEMBLY)); 7069566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Abd, MAT_FINAL_ASSEMBLY)); 7079566063dSJacob Faibussowitsch PetscCall(MatConvert(Abd, MATAIJ, MAT_INITIAL_MATRIX, &Bbd)); 7089566063dSJacob Faibussowitsch PetscCall(MatInvertBlockDiagonal(Abd, &isbd)); 7099566063dSJacob Faibussowitsch PetscCall(MatInvertBlockDiagonal(Bbd, &aijbd)); 7109566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Bbd, &nl, NULL)); 711c4762a1bSJed Brown ok = PETSC_TRUE; 712c4762a1bSJed Brown for (i = 0; i < nl / bs; i++) { 713c4762a1bSJed Brown PetscInt b1, b2; 714c4762a1bSJed Brown 715c4762a1bSJed Brown for (b1 = 0; b1 < bs; b1++) { 716c4762a1bSJed Brown for (b2 = 0; b2 < bs; b2++) { 717c4762a1bSJed Brown if (PetscAbsScalar(isbd[i * bs * bs + b1 * bs + b2] - aijbd[i * bs * bs + b1 * bs + b2]) > PETSC_SMALL) ok = PETSC_FALSE; 718c4762a1bSJed Brown if (!ok) { 7199566063dSJacob 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]))); 720c4762a1bSJed Brown break; 721c4762a1bSJed Brown } 722c4762a1bSJed Brown } 723c4762a1bSJed Brown if (!ok) break; 724c4762a1bSJed Brown } 725c4762a1bSJed Brown if (!ok) break; 726c4762a1bSJed Brown } 7279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Abd)); 7289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bbd)); 7299566063dSJacob Faibussowitsch PetscCall(PetscFree(vals)); 7309566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 7319566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bis)); 732c4762a1bSJed Brown } 733c4762a1bSJed Brown } 734c4762a1bSJed Brown } 735c4762a1bSJed Brown } 736d0dbe9f7SStefano Zampini 737d0dbe9f7SStefano Zampini /* test MatGetDiagonalBlock */ 738d0dbe9f7SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatGetDiagonalBlock\n")); 739d0dbe9f7SStefano Zampini PetscCall(MatGetDiagonalBlock(A, &A2)); 740d0dbe9f7SStefano Zampini PetscCall(MatGetDiagonalBlock(B, &B2)); 741d0dbe9f7SStefano Zampini PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatGetDiagonalBlock")); 742d0dbe9f7SStefano Zampini PetscCall(MatScale(A, 2.0)); 743d0dbe9f7SStefano Zampini PetscCall(MatScale(B, 2.0)); 744d0dbe9f7SStefano Zampini PetscCall(MatGetDiagonalBlock(A, &A2)); 745d0dbe9f7SStefano Zampini PetscCall(MatGetDiagonalBlock(B, &B2)); 746d0dbe9f7SStefano Zampini PetscCall(CheckMat(A2, B2, PETSC_FALSE, "MatGetDiagonalBlock")); 747d0dbe9f7SStefano Zampini 748c4762a1bSJed Brown /* free testing matrices */ 7499566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&cmap)); 7509566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&rmap)); 7519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 7529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 7539566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 754b122ec5aSJacob Faibussowitsch return 0; 755c4762a1bSJed Brown } 756c4762a1bSJed Brown 757d71ae5a4SJacob Faibussowitsch PetscErrorCode CheckMat(Mat A, Mat B, PetscBool usemult, const char *func) 758d71ae5a4SJacob Faibussowitsch { 759c4762a1bSJed Brown Mat Bcheck; 760c4762a1bSJed Brown PetscReal error; 761c4762a1bSJed Brown 762c4762a1bSJed Brown PetscFunctionBeginUser; 763c4762a1bSJed Brown if (!usemult && B) { 764c4762a1bSJed Brown PetscBool hasnorm; 765c4762a1bSJed Brown 7669566063dSJacob Faibussowitsch PetscCall(MatHasOperation(B, MATOP_NORM, &hasnorm)); 767c4762a1bSJed Brown if (!hasnorm) usemult = PETSC_TRUE; 768c4762a1bSJed Brown } 769c4762a1bSJed Brown if (!usemult) { 770c4762a1bSJed Brown if (B) { 771c4762a1bSJed Brown MatType Btype; 772c4762a1bSJed Brown 7739566063dSJacob Faibussowitsch PetscCall(MatGetType(B, &Btype)); 7749566063dSJacob Faibussowitsch PetscCall(MatConvert(A, Btype, MAT_INITIAL_MATRIX, &Bcheck)); 775c4762a1bSJed Brown } else { 7769566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &Bcheck)); 777c4762a1bSJed Brown } 778c4762a1bSJed Brown if (B) { /* if B is present, subtract it */ 7799566063dSJacob Faibussowitsch PetscCall(MatAXPY(Bcheck, -1., B, DIFFERENT_NONZERO_PATTERN)); 780c4762a1bSJed Brown } 7819566063dSJacob Faibussowitsch PetscCall(MatNorm(Bcheck, NORM_INFINITY, &error)); 782c4762a1bSJed Brown if (error > PETSC_SQRT_MACHINE_EPSILON) { 783c4762a1bSJed Brown ISLocalToGlobalMapping rl2g, cl2g; 784c4762a1bSJed Brown 785d0dbe9f7SStefano Zampini PetscCall(PetscObjectSetName((PetscObject)Bcheck, "Bcheck")); 7869566063dSJacob Faibussowitsch PetscCall(MatView(Bcheck, NULL)); 787c4762a1bSJed Brown if (B) { 788d0dbe9f7SStefano Zampini PetscCall(PetscObjectSetName((PetscObject)B, "B")); 7899566063dSJacob Faibussowitsch PetscCall(MatView(B, NULL)); 7909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bcheck)); 7919566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &Bcheck)); 792d0dbe9f7SStefano Zampini PetscCall(PetscObjectSetName((PetscObject)Bcheck, "Assembled A")); 7939566063dSJacob Faibussowitsch PetscCall(MatView(Bcheck, NULL)); 794c4762a1bSJed Brown } 7959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bcheck)); 796d0dbe9f7SStefano Zampini PetscCall(PetscObjectSetName((PetscObject)A, "A")); 7979566063dSJacob Faibussowitsch PetscCall(MatView(A, NULL)); 7989566063dSJacob Faibussowitsch PetscCall(MatGetLocalToGlobalMapping(A, &rl2g, &cl2g)); 799d0dbe9f7SStefano Zampini if (rl2g) PetscCall(ISLocalToGlobalMappingView(rl2g, NULL)); 800d0dbe9f7SStefano Zampini if (cl2g) PetscCall(ISLocalToGlobalMappingView(cl2g, NULL)); 801d0dbe9f7SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "ERROR ON %s: %g", func, (double)error); 802c4762a1bSJed Brown } 8039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bcheck)); 804c4762a1bSJed Brown } else { 805c4762a1bSJed Brown PetscBool ok, okt; 806c4762a1bSJed Brown 8079566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A, B, 3, &ok)); 8089566063dSJacob Faibussowitsch PetscCall(MatMultTransposeEqual(A, B, 3, &okt)); 809e00437b9SBarry Smith PetscCheck(ok && okt, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "ERROR ON %s: mult ok ? %d, multtranspose ok ? %d", func, ok, okt); 810c4762a1bSJed Brown } 811*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 812c4762a1bSJed Brown } 813c4762a1bSJed Brown 814d71ae5a4SJacob Faibussowitsch PetscErrorCode TestMatZeroRows(Mat A, Mat Afull, PetscBool squaretest, IS is, PetscScalar diag) 815d71ae5a4SJacob Faibussowitsch { 816c4762a1bSJed Brown Mat B, Bcheck, B2 = NULL, lB; 817c4762a1bSJed Brown Vec x = NULL, b = NULL, b2 = NULL; 818c4762a1bSJed Brown ISLocalToGlobalMapping l2gr, l2gc; 819c4762a1bSJed Brown PetscReal error; 820c4762a1bSJed Brown char diagstr[16]; 821c4762a1bSJed Brown const PetscInt *idxs; 822c4762a1bSJed Brown PetscInt rst, ren, i, n, N, d; 823c4762a1bSJed Brown PetscMPIInt rank; 824c4762a1bSJed Brown PetscBool miss, haszerorows; 825c4762a1bSJed Brown 826c4762a1bSJed Brown PetscFunctionBeginUser; 827c4762a1bSJed Brown if (diag == 0.) { 8289566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(diagstr, "zero")); 829c4762a1bSJed Brown } else { 8309566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(diagstr, "nonzero")); 831c4762a1bSJed Brown } 8329566063dSJacob Faibussowitsch PetscCall(ISView(is, NULL)); 8339566063dSJacob Faibussowitsch PetscCall(MatGetLocalToGlobalMapping(A, &l2gr, &l2gc)); 834c4762a1bSJed Brown /* tests MatDuplicate and MatCopy */ 835c4762a1bSJed Brown if (diag == 0.) { 8369566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B)); 837c4762a1bSJed Brown } else { 8389566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_DO_NOT_COPY_VALUES, &B)); 8399566063dSJacob Faibussowitsch PetscCall(MatCopy(A, B, SAME_NONZERO_PATTERN)); 840c4762a1bSJed Brown } 8419566063dSJacob Faibussowitsch PetscCall(MatISGetLocalMat(B, &lB)); 8429566063dSJacob Faibussowitsch PetscCall(MatHasOperation(lB, MATOP_ZERO_ROWS, &haszerorows)); 843c4762a1bSJed Brown if (squaretest && haszerorows) { 8449566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &x, &b)); 8459566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &B2)); 8469566063dSJacob Faibussowitsch PetscCall(VecSetLocalToGlobalMapping(b, l2gr)); 8479566063dSJacob Faibussowitsch PetscCall(VecSetLocalToGlobalMapping(x, l2gc)); 8489566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, NULL)); 8499566063dSJacob Faibussowitsch PetscCall(VecSetRandom(b, NULL)); 850c4762a1bSJed Brown /* mimic b[is] = x[is] */ 8519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(b, &b2)); 8529566063dSJacob Faibussowitsch PetscCall(VecSetLocalToGlobalMapping(b2, l2gr)); 8539566063dSJacob Faibussowitsch PetscCall(VecCopy(b, b2)); 8549566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &n)); 8559566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &idxs)); 8569566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &N)); 857c4762a1bSJed Brown for (i = 0; i < n; i++) { 858c4762a1bSJed Brown if (0 <= idxs[i] && idxs[i] < N) { 8599566063dSJacob Faibussowitsch PetscCall(VecSetValue(b2, idxs[i], diag, INSERT_VALUES)); 8609566063dSJacob Faibussowitsch PetscCall(VecSetValue(x, idxs[i], 1., INSERT_VALUES)); 861c4762a1bSJed Brown } 862c4762a1bSJed Brown } 8639566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(b2)); 8649566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(b2)); 8659566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x)); 8669566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x)); 8679566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &idxs)); 868c4762a1bSJed Brown /* test ZeroRows on MATIS */ 8699566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatZeroRows (diag %s)\n", diagstr)); 8709566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(B, is, diag, x, b)); 8719566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatZeroRowsColumns (diag %s)\n", diagstr)); 8729566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumnsIS(B2, is, diag, NULL, NULL)); 873c4762a1bSJed Brown } else if (haszerorows) { 874c4762a1bSJed Brown /* test ZeroRows on MATIS */ 8759566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatZeroRows (diag %s)\n", diagstr)); 8769566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(B, is, diag, NULL, NULL)); 877c4762a1bSJed Brown b = b2 = x = NULL; 878c4762a1bSJed Brown } else { 8799566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Skipping MatZeroRows (diag %s)\n", diagstr)); 880c4762a1bSJed Brown b = b2 = x = NULL; 881c4762a1bSJed Brown } 882c4762a1bSJed Brown 883c4762a1bSJed Brown if (squaretest && haszerorows) { 8849566063dSJacob Faibussowitsch PetscCall(VecAXPY(b2, -1., b)); 8859566063dSJacob Faibussowitsch PetscCall(VecNorm(b2, NORM_INFINITY, &error)); 886e00437b9SBarry Smith PetscCheck(error <= PETSC_SQRT_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "ERROR IN ZEROROWS ON B %g (diag %s)", (double)error, diagstr); 887c4762a1bSJed Brown } 888c4762a1bSJed Brown 889c4762a1bSJed Brown /* test MatMissingDiagonal */ 8909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test MatMissingDiagonal\n")); 8919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 8929566063dSJacob Faibussowitsch PetscCall(MatMissingDiagonal(B, &miss, &d)); 8939566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(B, &rst, &ren)); 8949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 8959566063dSJacob 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)); 8969566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 8979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 898c4762a1bSJed Brown 8999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 9009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 9019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b2)); 902c4762a1bSJed Brown 903c4762a1bSJed Brown /* check the result of ZeroRows with that from MPIAIJ routines 904c4762a1bSJed Brown assuming that MatConvert_IS_XAIJ and MatZeroRows_MPIAIJ work fine */ 905c4762a1bSJed Brown if (haszerorows) { 9069566063dSJacob Faibussowitsch PetscCall(MatDuplicate(Afull, MAT_COPY_VALUES, &Bcheck)); 9079566063dSJacob Faibussowitsch PetscCall(MatSetOption(Bcheck, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 9089566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(Bcheck, is, diag, NULL, NULL)); 9099566063dSJacob Faibussowitsch PetscCall(CheckMat(B, Bcheck, PETSC_FALSE, "Zerorows")); 9109566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bcheck)); 911c4762a1bSJed Brown } 9129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 913c4762a1bSJed Brown 914c4762a1bSJed Brown if (B2) { /* test MatZeroRowsColumns */ 9159566063dSJacob Faibussowitsch PetscCall(MatDuplicate(Afull, MAT_COPY_VALUES, &B)); 9169566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 9179566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumnsIS(B, is, diag, NULL, NULL)); 9189566063dSJacob Faibussowitsch PetscCall(CheckMat(B2, B, PETSC_FALSE, "MatZeroRowsColumns")); 9199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 9209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 921c4762a1bSJed Brown } 922*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 923c4762a1bSJed Brown } 924c4762a1bSJed Brown 925c4762a1bSJed Brown /*TEST 926c4762a1bSJed Brown 927c4762a1bSJed Brown test: 928c4762a1bSJed Brown args: -test_trans 929c4762a1bSJed Brown 930c4762a1bSJed Brown test: 931c4762a1bSJed Brown suffix: 2 932c4762a1bSJed Brown nsize: 4 933c4762a1bSJed Brown args: -matis_convert_local_nest -nr 3 -nc 4 934c4762a1bSJed Brown 935c4762a1bSJed Brown test: 936c4762a1bSJed Brown suffix: 3 937c4762a1bSJed Brown nsize: 5 938c4762a1bSJed Brown args: -m 11 -n 10 -matis_convert_local_nest -nr 2 -nc 1 939c4762a1bSJed Brown 940c4762a1bSJed Brown test: 941c4762a1bSJed Brown suffix: 4 942c4762a1bSJed Brown nsize: 6 943c4762a1bSJed Brown args: -m 9 -n 12 -test_trans -nr 2 -nc 7 944c4762a1bSJed Brown 945c4762a1bSJed Brown test: 946c4762a1bSJed Brown suffix: 5 947c4762a1bSJed Brown nsize: 6 948c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 3 -nc 1 949c4762a1bSJed Brown 950c4762a1bSJed Brown test: 951c4762a1bSJed Brown suffix: 6 952c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap 953c4762a1bSJed Brown 954c4762a1bSJed Brown test: 955c4762a1bSJed Brown suffix: 7 956c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap -permmap 957c4762a1bSJed Brown 958c4762a1bSJed Brown test: 959c4762a1bSJed Brown suffix: 8 960c4762a1bSJed Brown args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -permmap 961c4762a1bSJed Brown 962c4762a1bSJed Brown test: 963c4762a1bSJed Brown suffix: 9 964c4762a1bSJed Brown nsize: 5 965c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap 966c4762a1bSJed Brown 967c4762a1bSJed Brown test: 968c4762a1bSJed Brown suffix: 10 969c4762a1bSJed Brown nsize: 5 970c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap -permmap 971c4762a1bSJed Brown 972c20d7725SJed Brown test: 973c20d7725SJed Brown suffix: vscat_default 974c4762a1bSJed Brown nsize: 5 975c4762a1bSJed Brown args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -permmap 976c4762a1bSJed Brown output_file: output/ex23_11.out 977c4762a1bSJed Brown 978c4762a1bSJed Brown test: 979c4762a1bSJed Brown suffix: 12 980c4762a1bSJed Brown nsize: 3 981c4762a1bSJed Brown args: -m 12 -n 12 -symmetric -matis_localmat_type sbaij -test_trans -nr 2 -nc 3 982c4762a1bSJed Brown 983c4762a1bSJed Brown testset: 984c4762a1bSJed Brown output_file: output/ex23_13.out 985c4762a1bSJed Brown nsize: 3 986c4762a1bSJed Brown args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -diffmap -permmap 987c4762a1bSJed Brown filter: grep -v "type:" 988c4762a1bSJed Brown test: 989c4762a1bSJed Brown suffix: baij 990c4762a1bSJed Brown args: -matis_localmat_type baij 991c4762a1bSJed Brown test: 992c4762a1bSJed Brown requires: viennacl 993c4762a1bSJed Brown suffix: viennacl 994c4762a1bSJed Brown args: -matis_localmat_type aijviennacl 995c4762a1bSJed Brown test: 996c4762a1bSJed Brown requires: cuda 997c4762a1bSJed Brown suffix: cusparse 998c4762a1bSJed Brown args: -matis_localmat_type aijcusparse 999c4762a1bSJed Brown 1000e432b41dSStefano Zampini test: 1001e432b41dSStefano Zampini suffix: negrep 1002e432b41dSStefano Zampini nsize: {{1 3}separate output} 1003e432b41dSStefano 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} 1004e432b41dSStefano Zampini 1005c4762a1bSJed Brown TEST*/ 1006