1c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for parallel MatSBAIJ format.\n"; 2c4762a1bSJed Brown /* Example of usage: 3c4762a1bSJed Brown mpiexec -n 2 ./ex92 -nd 2 -ov 3 -mat_block_size 2 -view_id 0 -test_overlap -test_submat 4c4762a1bSJed Brown */ 5c4762a1bSJed Brown #include <petscmat.h> 6c4762a1bSJed Brown 7d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 8d71ae5a4SJacob Faibussowitsch { 9c4762a1bSJed Brown Mat A, Atrans, sA, *submatA, *submatsA; 10c4762a1bSJed Brown PetscMPIInt size, rank; 11c4762a1bSJed Brown PetscInt bs = 1, mbs = 10, ov = 1, i, j, k, *rows, *cols, nd = 2, *idx, rstart, rend, sz, M, N, Mbs; 12c4762a1bSJed Brown PetscScalar *vals, rval, one = 1.0; 13c4762a1bSJed Brown IS *is1, *is2; 14c4762a1bSJed Brown PetscRandom rand; 15c4762a1bSJed Brown PetscBool flg, TestOverlap, TestSubMat, TestAllcols, test_sorted = PETSC_FALSE; 16c4762a1bSJed Brown PetscInt vid = -1; 17c4762a1bSJed Brown PetscLogStage stages[2]; 18c4762a1bSJed Brown 19327415f7SBarry Smith PetscFunctionBeginUser; 20c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 23c4762a1bSJed Brown 249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL)); 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_mbs", &mbs, NULL)); 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL)); 279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL)); 289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-view_id", &vid, NULL)); 299566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-test_overlap", &TestOverlap)); 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-test_submat", &TestSubMat)); 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-test_allcols", &TestAllcols)); 329566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_sorted", &test_sorted, NULL)); 33c4762a1bSJed Brown 349566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 359566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, mbs * bs, mbs * bs, PETSC_DECIDE, PETSC_DECIDE)); 369566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATBAIJ)); 379566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(A, bs, PETSC_DEFAULT, NULL)); 389566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(A, bs, PETSC_DEFAULT, NULL, PETSC_DEFAULT, NULL)); 39c4762a1bSJed Brown 409566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand)); 419566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 42c4762a1bSJed Brown 439566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 449566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 45c4762a1bSJed Brown Mbs = M / bs; 46c4762a1bSJed Brown 479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs, &rows)); 489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs, &cols)); 499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * bs, &vals)); 509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M, &idx)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Now set blocks of values */ 53c4762a1bSJed Brown for (j = 0; j < bs * bs; j++) vals[j] = 0.0; 54c4762a1bSJed Brown for (i = 0; i < Mbs; i++) { 559371c9d4SSatish Balay cols[0] = i * bs; 569371c9d4SSatish Balay rows[0] = i * bs; 57c4762a1bSJed Brown for (j = 1; j < bs; j++) { 58c4762a1bSJed Brown rows[j] = rows[j - 1] + 1; 59c4762a1bSJed Brown cols[j] = cols[j - 1] + 1; 60c4762a1bSJed Brown } 619566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, bs, rows, bs, cols, vals, ADD_VALUES)); 62c4762a1bSJed Brown } 63c4762a1bSJed Brown /* second, add random blocks */ 64c4762a1bSJed Brown for (i = 0; i < 20 * bs; i++) { 659566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &rval)); 66c4762a1bSJed Brown cols[0] = bs * (PetscInt)(PetscRealPart(rval) * Mbs); 679566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &rval)); 68c4762a1bSJed Brown rows[0] = rstart + bs * (PetscInt)(PetscRealPart(rval) * mbs); 69c4762a1bSJed Brown for (j = 1; j < bs; j++) { 70c4762a1bSJed Brown rows[j] = rows[j - 1] + 1; 71c4762a1bSJed Brown cols[j] = cols[j - 1] + 1; 72c4762a1bSJed Brown } 73c4762a1bSJed Brown 74c4762a1bSJed Brown for (j = 0; j < bs * bs; j++) { 759566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &rval)); 76c4762a1bSJed Brown vals[j] = rval; 77c4762a1bSJed Brown } 789566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, bs, rows, bs, cols, vals, ADD_VALUES)); 79c4762a1bSJed Brown } 80c4762a1bSJed Brown 819566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 829566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* make A a symmetric matrix: A <- A^T + A */ 859566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &Atrans)); 869566063dSJacob Faibussowitsch PetscCall(MatAXPY(A, one, Atrans, DIFFERENT_NONZERO_PATTERN)); 879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Atrans)); 889566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &Atrans)); 899566063dSJacob Faibussowitsch PetscCall(MatEqual(A, Atrans, &flg)); 90*966bd95aSPierre Jolivet PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "A+A^T is non-symmetric"); 919566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE)); 929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Atrans)); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* create a SeqSBAIJ matrix sA (= A) */ 959566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSBAIJ, MAT_INITIAL_MATRIX, &sA)); 96c4762a1bSJed Brown if (vid >= 0 && vid < size) { 979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "A:\n")); 989566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 999566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "sA:\n")); 1009566063dSJacob Faibussowitsch PetscCall(MatView(sA, PETSC_VIEWER_STDOUT_WORLD)); 101c4762a1bSJed Brown } 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* Test sA==A through MatMult() */ 1049566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A, sA, 10, &flg)); 10528b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Error in MatConvert(): A != sA"); 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* Test MatIncreaseOverlap() */ 1089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is1)); 1099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is2)); 110c4762a1bSJed Brown 111c4762a1bSJed Brown for (i = 0; i < nd; i++) { 112c4762a1bSJed Brown if (!TestAllcols) { 1139566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &rval)); 114c4762a1bSJed Brown sz = (PetscInt)((0.5 + 0.2 * PetscRealPart(rval)) * mbs); /* 0.5*mbs < sz < 0.7*mbs */ 115c4762a1bSJed Brown 116c4762a1bSJed Brown for (j = 0; j < sz; j++) { 1179566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &rval)); 118c4762a1bSJed Brown idx[j * bs] = bs * (PetscInt)(PetscRealPart(rval) * Mbs); 119c4762a1bSJed Brown for (k = 1; k < bs; k++) idx[j * bs + k] = idx[j * bs] + k; 120c4762a1bSJed Brown } 1219566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, sz * bs, idx, PETSC_COPY_VALUES, is1 + i)); 1229566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, sz * bs, idx, PETSC_COPY_VALUES, is2 + i)); 123c4762a1bSJed Brown if (rank == vid) { 1249566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " [%d] IS sz[%" PetscInt_FMT "]: %" PetscInt_FMT "\n", rank, i, sz)); 1259566063dSJacob Faibussowitsch PetscCall(ISView(is2[i], PETSC_VIEWER_STDOUT_SELF)); 126c4762a1bSJed Brown } 127a5b23f4aSJose E. Roman } else { /* Test all rows and columns */ 128c4762a1bSJed Brown sz = M; 1299566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, sz, 0, 1, is1 + i)); 1309566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, sz, 0, 1, is2 + i)); 131c4762a1bSJed Brown 132c4762a1bSJed Brown if (rank == vid) { 133c4762a1bSJed Brown PetscBool colflag; 1349566063dSJacob Faibussowitsch PetscCall(ISIdentity(is2[i], &colflag)); 1359566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] is2[%" PetscInt_FMT "], colflag %d\n", rank, i, colflag)); 1369566063dSJacob Faibussowitsch PetscCall(ISView(is2[i], PETSC_VIEWER_STDOUT_SELF)); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown } 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 1419566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("MatOv_SBAIJ", &stages[0])); 1429566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("MatOv_BAIJ", &stages[1])); 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* Test MatIncreaseOverlap */ 145c4762a1bSJed Brown if (TestOverlap) { 1469566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stages[0])); 1479566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(sA, nd, is2, ov)); 1489566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 149c4762a1bSJed Brown 1509566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stages[1])); 1519566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A, nd, is1, ov)); 1529566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 153c4762a1bSJed Brown 154c4762a1bSJed Brown if (rank == vid) { 1559566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n[%d] IS from BAIJ:\n", rank)); 1569566063dSJacob Faibussowitsch PetscCall(ISView(is1[0], PETSC_VIEWER_STDOUT_SELF)); 1579566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n[%d] IS from SBAIJ:\n", rank)); 1589566063dSJacob Faibussowitsch PetscCall(ISView(is2[0], PETSC_VIEWER_STDOUT_SELF)); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 1629566063dSJacob Faibussowitsch PetscCall(ISEqual(is1[i], is2[i], &flg)); 163c4762a1bSJed Brown if (!flg) { 164dd400576SPatrick Sanan if (rank == 0) { 1659566063dSJacob Faibussowitsch PetscCall(ISSort(is1[i])); 1669566063dSJacob Faibussowitsch PetscCall(ISSort(is2[i])); 167c4762a1bSJed Brown } 16898921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "i=%" PetscInt_FMT ", is1 != is2", i); 169c4762a1bSJed Brown } 170c4762a1bSJed Brown } 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* Test MatCreateSubmatrices */ 174c4762a1bSJed Brown if (TestSubMat) { 175c4762a1bSJed Brown if (test_sorted) { 17648a46eb9SPierre Jolivet for (i = 0; i < nd; ++i) PetscCall(ISSort(is1[i])); 177c4762a1bSJed Brown } 1789566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_INITIAL_MATRIX, &submatA)); 1799566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(sA, nd, is1, is1, MAT_INITIAL_MATRIX, &submatsA)); 180c4762a1bSJed Brown 1819566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A, sA, 10, &flg)); 18228b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "A != sA"); 183c4762a1bSJed Brown 184c4762a1bSJed Brown /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */ 1859566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_REUSE_MATRIX, &submatA)); 1869566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(sA, nd, is1, is1, MAT_REUSE_MATRIX, &submatsA)); 1879566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A, sA, 10, &flg)); 18828b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "MatCreateSubmatrices(): A != sA"); 189c4762a1bSJed Brown 1909566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd, &submatA)); 1919566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd, &submatsA)); 192c4762a1bSJed Brown } 193c4762a1bSJed Brown 194c4762a1bSJed Brown /* Free allocated memory */ 195c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 1969566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1[i])); 1979566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2[i])); 198c4762a1bSJed Brown } 1999566063dSJacob Faibussowitsch PetscCall(PetscFree(is1)); 2009566063dSJacob Faibussowitsch PetscCall(PetscFree(is2)); 2019566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 2029566063dSJacob Faibussowitsch PetscCall(PetscFree(rows)); 2039566063dSJacob Faibussowitsch PetscCall(PetscFree(cols)); 2049566063dSJacob Faibussowitsch PetscCall(PetscFree(vals)); 2059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sA)); 2079566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 2089566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 209b122ec5aSJacob Faibussowitsch return 0; 210c4762a1bSJed Brown } 211c4762a1bSJed Brown 212c4762a1bSJed Brown /*TEST 213c4762a1bSJed Brown 214c4762a1bSJed Brown test: 215c4762a1bSJed Brown args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_submat 2163886731fSPierre Jolivet output_file: output/empty.out 217c4762a1bSJed Brown 218c4762a1bSJed Brown test: 219c4762a1bSJed Brown suffix: 2 220c4762a1bSJed Brown nsize: {{3 4}} 221c4762a1bSJed Brown args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_submat 2223886731fSPierre Jolivet output_file: output/empty.out 223c4762a1bSJed Brown 224c4762a1bSJed Brown test: 225c4762a1bSJed Brown suffix: 3 226c4762a1bSJed Brown nsize: {{3 4}} 227c4762a1bSJed Brown args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_allcols 2283886731fSPierre Jolivet output_file: output/empty.out 229c4762a1bSJed Brown 230c4762a1bSJed Brown test: 231c4762a1bSJed Brown suffix: 3_sorted 232c4762a1bSJed Brown nsize: {{3 4}} 233c4762a1bSJed Brown args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_allcols -test_sorted 2343886731fSPierre Jolivet output_file: output/empty.out 235c4762a1bSJed Brown 236c4762a1bSJed Brown test: 237c4762a1bSJed Brown suffix: 4 238c4762a1bSJed Brown nsize: {{3 4}} 239c4762a1bSJed Brown args: -ov {{1 3}} -mat_block_size {{2 8}} -test_submat -test_allcols 2403886731fSPierre Jolivet output_file: output/empty.out 241c4762a1bSJed Brown 242c4762a1bSJed Brown TEST*/ 243