1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for parallel MatSBAIJ format.\n"; 3c4762a1bSJed Brown /* Example of usage: 4c4762a1bSJed Brown mpiexec -n 2 ./ex92 -nd 2 -ov 3 -mat_block_size 2 -view_id 0 -test_overlap -test_submat 5c4762a1bSJed Brown */ 6c4762a1bSJed Brown #include <petscmat.h> 7c4762a1bSJed Brown 89371c9d4SSatish Balay int main(int argc, char **args) { 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 #if defined(PETSC_USE_LOG) 18c4762a1bSJed Brown PetscLogStage stages[2]; 19c4762a1bSJed Brown #endif 20c4762a1bSJed Brown 21327415f7SBarry Smith PetscFunctionBeginUser; 229566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 25c4762a1bSJed Brown 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL)); 279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_mbs", &mbs, NULL)); 289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL)); 299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL)); 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-view_id", &vid, NULL)); 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-test_overlap", &TestOverlap)); 329566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-test_submat", &TestSubMat)); 339566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-test_allcols", &TestAllcols)); 349566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_sorted", &test_sorted, NULL)); 35c4762a1bSJed Brown 369566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 379566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, mbs * bs, mbs * bs, PETSC_DECIDE, PETSC_DECIDE)); 389566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATBAIJ)); 399566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(A, bs, PETSC_DEFAULT, NULL)); 409566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(A, bs, PETSC_DEFAULT, NULL, PETSC_DEFAULT, NULL)); 41c4762a1bSJed Brown 429566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand)); 439566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 44c4762a1bSJed Brown 459566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 469566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 47c4762a1bSJed Brown Mbs = M / bs; 48c4762a1bSJed Brown 499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs, &rows)); 509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs, &cols)); 519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * bs, &vals)); 529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M, &idx)); 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* Now set blocks of values */ 55c4762a1bSJed Brown for (j = 0; j < bs * bs; j++) vals[j] = 0.0; 56c4762a1bSJed Brown for (i = 0; i < Mbs; i++) { 579371c9d4SSatish Balay cols[0] = i * bs; 589371c9d4SSatish Balay rows[0] = i * bs; 59c4762a1bSJed Brown for (j = 1; j < bs; j++) { 60c4762a1bSJed Brown rows[j] = rows[j - 1] + 1; 61c4762a1bSJed Brown cols[j] = cols[j - 1] + 1; 62c4762a1bSJed Brown } 639566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, bs, rows, bs, cols, vals, ADD_VALUES)); 64c4762a1bSJed Brown } 65c4762a1bSJed Brown /* second, add random blocks */ 66c4762a1bSJed Brown for (i = 0; i < 20 * bs; i++) { 679566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &rval)); 68c4762a1bSJed Brown cols[0] = bs * (PetscInt)(PetscRealPart(rval) * Mbs); 699566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &rval)); 70c4762a1bSJed Brown rows[0] = rstart + bs * (PetscInt)(PetscRealPart(rval) * mbs); 71c4762a1bSJed Brown for (j = 1; j < bs; j++) { 72c4762a1bSJed Brown rows[j] = rows[j - 1] + 1; 73c4762a1bSJed Brown cols[j] = cols[j - 1] + 1; 74c4762a1bSJed Brown } 75c4762a1bSJed Brown 76c4762a1bSJed Brown for (j = 0; j < bs * bs; j++) { 779566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &rval)); 78c4762a1bSJed Brown vals[j] = rval; 79c4762a1bSJed Brown } 809566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, bs, rows, bs, cols, vals, ADD_VALUES)); 81c4762a1bSJed Brown } 82c4762a1bSJed Brown 839566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 849566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* make A a symmetric matrix: A <- A^T + A */ 879566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &Atrans)); 889566063dSJacob Faibussowitsch PetscCall(MatAXPY(A, one, Atrans, DIFFERENT_NONZERO_PATTERN)); 899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Atrans)); 909566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &Atrans)); 919566063dSJacob Faibussowitsch PetscCall(MatEqual(A, Atrans, &flg)); 92c4762a1bSJed Brown if (flg) { 939566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE)); 94c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "A+A^T is non-symmetric"); 959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Atrans)); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* create a SeqSBAIJ matrix sA (= A) */ 989566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSBAIJ, MAT_INITIAL_MATRIX, &sA)); 99c4762a1bSJed Brown if (vid >= 0 && vid < size) { 1009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "A:\n")); 1019566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 1029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "sA:\n")); 1039566063dSJacob Faibussowitsch PetscCall(MatView(sA, PETSC_VIEWER_STDOUT_WORLD)); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* Test sA==A through MatMult() */ 1079566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A, sA, 10, &flg)); 10828b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Error in MatConvert(): A != sA"); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* Test MatIncreaseOverlap() */ 1119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is1)); 1129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is2)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown for (i = 0; i < nd; i++) { 115c4762a1bSJed Brown if (!TestAllcols) { 1169566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &rval)); 117c4762a1bSJed Brown sz = (PetscInt)((0.5 + 0.2 * PetscRealPart(rval)) * mbs); /* 0.5*mbs < sz < 0.7*mbs */ 118c4762a1bSJed Brown 119c4762a1bSJed Brown for (j = 0; j < sz; j++) { 1209566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &rval)); 121c4762a1bSJed Brown idx[j * bs] = bs * (PetscInt)(PetscRealPart(rval) * Mbs); 122c4762a1bSJed Brown for (k = 1; k < bs; k++) idx[j * bs + k] = idx[j * bs] + k; 123c4762a1bSJed Brown } 1249566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, sz * bs, idx, PETSC_COPY_VALUES, is1 + i)); 1259566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, sz * bs, idx, PETSC_COPY_VALUES, is2 + i)); 126c4762a1bSJed Brown if (rank == vid) { 1279566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " [%d] IS sz[%" PetscInt_FMT "]: %" PetscInt_FMT "\n", rank, i, sz)); 1289566063dSJacob Faibussowitsch PetscCall(ISView(is2[i], PETSC_VIEWER_STDOUT_SELF)); 129c4762a1bSJed Brown } 130a5b23f4aSJose E. Roman } else { /* Test all rows and columns */ 131c4762a1bSJed Brown sz = M; 1329566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, sz, 0, 1, is1 + i)); 1339566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, sz, 0, 1, is2 + i)); 134c4762a1bSJed Brown 135c4762a1bSJed Brown if (rank == vid) { 136c4762a1bSJed Brown PetscBool colflag; 1379566063dSJacob Faibussowitsch PetscCall(ISIdentity(is2[i], &colflag)); 1389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] is2[%" PetscInt_FMT "], colflag %d\n", rank, i, colflag)); 1399566063dSJacob Faibussowitsch PetscCall(ISView(is2[i], PETSC_VIEWER_STDOUT_SELF)); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown } 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 1449566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("MatOv_SBAIJ", &stages[0])); 1459566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("MatOv_BAIJ", &stages[1])); 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* Test MatIncreaseOverlap */ 148c4762a1bSJed Brown if (TestOverlap) { 1499566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stages[0])); 1509566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(sA, nd, is2, ov)); 1519566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 152c4762a1bSJed Brown 1539566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stages[1])); 1549566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A, nd, is1, ov)); 1559566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 156c4762a1bSJed Brown 157c4762a1bSJed Brown if (rank == vid) { 1589566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n[%d] IS from BAIJ:\n", rank)); 1599566063dSJacob Faibussowitsch PetscCall(ISView(is1[0], PETSC_VIEWER_STDOUT_SELF)); 1609566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n[%d] IS from SBAIJ:\n", rank)); 1619566063dSJacob Faibussowitsch PetscCall(ISView(is2[0], PETSC_VIEWER_STDOUT_SELF)); 162c4762a1bSJed Brown } 163c4762a1bSJed Brown 164c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 1659566063dSJacob Faibussowitsch PetscCall(ISEqual(is1[i], is2[i], &flg)); 166c4762a1bSJed Brown if (!flg) { 167dd400576SPatrick Sanan if (rank == 0) { 1689566063dSJacob Faibussowitsch PetscCall(ISSort(is1[i])); 1699566063dSJacob Faibussowitsch PetscCall(ISSort(is2[i])); 170c4762a1bSJed Brown } 17198921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "i=%" PetscInt_FMT ", is1 != is2", i); 172c4762a1bSJed Brown } 173c4762a1bSJed Brown } 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176c4762a1bSJed Brown /* Test MatCreateSubmatrices */ 177c4762a1bSJed Brown if (TestSubMat) { 178c4762a1bSJed Brown if (test_sorted) { 179*48a46eb9SPierre Jolivet for (i = 0; i < nd; ++i) PetscCall(ISSort(is1[i])); 180c4762a1bSJed Brown } 1819566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_INITIAL_MATRIX, &submatA)); 1829566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(sA, nd, is1, is1, MAT_INITIAL_MATRIX, &submatsA)); 183c4762a1bSJed Brown 1849566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A, sA, 10, &flg)); 18528b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "A != sA"); 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */ 1889566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_REUSE_MATRIX, &submatA)); 1899566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(sA, nd, is1, is1, MAT_REUSE_MATRIX, &submatsA)); 1909566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A, sA, 10, &flg)); 19128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "MatCreateSubmatrices(): A != sA"); 192c4762a1bSJed Brown 1939566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd, &submatA)); 1949566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd, &submatsA)); 195c4762a1bSJed Brown } 196c4762a1bSJed Brown 197c4762a1bSJed Brown /* Free allocated memory */ 198c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 1999566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1[i])); 2009566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2[i])); 201c4762a1bSJed Brown } 2029566063dSJacob Faibussowitsch PetscCall(PetscFree(is1)); 2039566063dSJacob Faibussowitsch PetscCall(PetscFree(is2)); 2049566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 2059566063dSJacob Faibussowitsch PetscCall(PetscFree(rows)); 2069566063dSJacob Faibussowitsch PetscCall(PetscFree(cols)); 2079566063dSJacob Faibussowitsch PetscCall(PetscFree(vals)); 2089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sA)); 2109566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 2119566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 212b122ec5aSJacob Faibussowitsch return 0; 213c4762a1bSJed Brown } 214c4762a1bSJed Brown 215c4762a1bSJed Brown /*TEST 216c4762a1bSJed Brown 217c4762a1bSJed Brown test: 218c4762a1bSJed Brown args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_submat 219c4762a1bSJed Brown output_file: output/ex92_1.out 220c4762a1bSJed Brown 221c4762a1bSJed Brown test: 222c4762a1bSJed Brown suffix: 2 223c4762a1bSJed Brown nsize: {{3 4}} 224c4762a1bSJed Brown args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_submat 225c4762a1bSJed Brown output_file: output/ex92_1.out 226c4762a1bSJed Brown 227c4762a1bSJed Brown test: 228c4762a1bSJed Brown suffix: 3 229c4762a1bSJed Brown nsize: {{3 4}} 230c4762a1bSJed Brown args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_allcols 231c4762a1bSJed Brown output_file: output/ex92_1.out 232c4762a1bSJed Brown 233c4762a1bSJed Brown test: 234c4762a1bSJed Brown suffix: 3_sorted 235c4762a1bSJed Brown nsize: {{3 4}} 236c4762a1bSJed Brown args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_allcols -test_sorted 237c4762a1bSJed Brown output_file: output/ex92_1.out 238c4762a1bSJed Brown 239c4762a1bSJed Brown test: 240c4762a1bSJed Brown suffix: 4 241c4762a1bSJed Brown nsize: {{3 4}} 242c4762a1bSJed Brown args: -ov {{1 3}} -mat_block_size {{2 8}} -test_submat -test_allcols 243c4762a1bSJed Brown output_file: output/ex92_1.out 244c4762a1bSJed Brown 245c4762a1bSJed Brown TEST*/ 246