1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap() and MatCreateSubmatrices() for the parallel case.\n\ 3c4762a1bSJed Brown This example is similar to ex40.c; here the index sets used are random.\n\ 4c4762a1bSJed Brown Input arguments are:\n\ 5c4762a1bSJed Brown -f <input_file> : file to load. For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\ 6c4762a1bSJed Brown -nd <size> : > 0 no of domains per processor \n\ 7c4762a1bSJed Brown -ov <overlap> : >=0 amount of overlap between domains\n\n"; 8c4762a1bSJed Brown 9c4762a1bSJed Brown #include <petscmat.h> 10c4762a1bSJed Brown 119371c9d4SSatish Balay int main(int argc, char **args) { 12c4762a1bSJed Brown PetscInt nd = 2, ov = 1, i, j, lsize, m, n, *idx, bs; 13c4762a1bSJed Brown PetscMPIInt rank, size; 14c4762a1bSJed Brown PetscBool flg; 15c4762a1bSJed Brown Mat A, B, *submatA, *submatB; 16c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; 17c4762a1bSJed Brown PetscViewer fd; 18c4762a1bSJed Brown IS *is1, *is2; 19c4762a1bSJed Brown PetscRandom r; 20c4762a1bSJed Brown PetscBool test_unsorted = PETSC_FALSE; 21c4762a1bSJed Brown PetscScalar rand; 22c4762a1bSJed Brown 23327415f7SBarry Smith PetscFunctionBeginUser; 249566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), NULL)); 289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL)); 299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL)); 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_unsorted", &test_unsorted, NULL)); 31c4762a1bSJed Brown 32c4762a1bSJed Brown /* Read matrix A and RHS */ 339566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd)); 349566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 359566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATAIJ)); 369566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 379566063dSJacob Faibussowitsch PetscCall(MatLoad(A, fd)); 389566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* Read the same matrix as a seq matrix B */ 419566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, file, FILE_MODE_READ, &fd)); 429566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &B)); 439566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATSEQAIJ)); 449566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 459566063dSJacob Faibussowitsch PetscCall(MatLoad(B, fd)); 469566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 47c4762a1bSJed Brown 489566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &bs)); 49c4762a1bSJed Brown 50c4762a1bSJed Brown /* Create the Random no generator */ 519566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &m, &n)); 529566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r)); 539566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(r)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Create the IS corresponding to subdomains */ 569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is1)); 579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is2)); 589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &idx)); 592f613bf5SBarry Smith for (i = 0; i < m; i++) { idx[i] = i; } 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* Create the random Index Sets */ 62c4762a1bSJed Brown for (i = 0; i < nd; i++) { 63c4762a1bSJed Brown /* Skip a few,so that the IS on different procs are diffeent*/ 64*48a46eb9SPierre Jolivet for (j = 0; j < rank; j++) PetscCall(PetscRandomGetValue(r, &rand)); 659566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r, &rand)); 66c4762a1bSJed Brown lsize = (PetscInt)(rand * (m / bs)); 67c4762a1bSJed Brown /* shuffle */ 68c4762a1bSJed Brown for (j = 0; j < lsize; j++) { 69c4762a1bSJed Brown PetscInt k, swap, l; 70c4762a1bSJed Brown 719566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r, &rand)); 72c4762a1bSJed Brown k = j + (PetscInt)(rand * ((m / bs) - j)); 73c4762a1bSJed Brown for (l = 0; l < bs; l++) { 74c4762a1bSJed Brown swap = idx[bs * j + l]; 75c4762a1bSJed Brown idx[bs * j + l] = idx[bs * k + l]; 76c4762a1bSJed Brown idx[bs * k + l] = swap; 77c4762a1bSJed Brown } 78c4762a1bSJed Brown } 799566063dSJacob Faibussowitsch if (!test_unsorted) PetscCall(PetscSortInt(lsize * bs, idx)); 809566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsize * bs, idx, PETSC_COPY_VALUES, is1 + i)); 819566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsize * bs, idx, PETSC_COPY_VALUES, is2 + i)); 829566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(is1[i], bs)); 839566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(is2[i], bs)); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown 86c4762a1bSJed Brown if (!test_unsorted) { 879566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A, nd, is1, ov)); 889566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(B, nd, is2, ov)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 919566063dSJacob Faibussowitsch PetscCall(ISSort(is1[i])); 929566063dSJacob Faibussowitsch PetscCall(ISSort(is2[i])); 93c4762a1bSJed Brown } 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 969566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_INITIAL_MATRIX, &submatA)); 979566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_INITIAL_MATRIX, &submatB)); 98c4762a1bSJed Brown 99c4762a1bSJed Brown /* Now see if the serial and parallel case have the same answers */ 100c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 1019566063dSJacob Faibussowitsch PetscCall(MatEqual(submatA[i], submatB[i], &flg)); 1026aad120cSJose E. Roman PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "%" PetscInt_FMT "-th parallel submatA != seq submatB", i); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* Free Allocated Memory */ 106c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 1079566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1[i])); 1089566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2[i])); 109c4762a1bSJed Brown } 1109566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd, &submatA)); 1119566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd, &submatB)); 112c4762a1bSJed Brown 1139566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&r)); 1149566063dSJacob Faibussowitsch PetscCall(PetscFree(is1)); 1159566063dSJacob Faibussowitsch PetscCall(PetscFree(is2)); 1169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1189566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 1199566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 120b122ec5aSJacob Faibussowitsch return 0; 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 123c4762a1bSJed Brown /*TEST 124c4762a1bSJed Brown 125c4762a1bSJed Brown build: 126c4762a1bSJed Brown requires: !complex 127c4762a1bSJed Brown 128c4762a1bSJed Brown test: 129c4762a1bSJed Brown nsize: 3 130dfd57a17SPierre Jolivet requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex 131c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/arco1 -nd 5 -ov 2 132c4762a1bSJed Brown 133c4762a1bSJed Brown test: 134c4762a1bSJed Brown suffix: 2 135c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/arco1 -nd 8 -ov 2 136dfd57a17SPierre Jolivet requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex 137c4762a1bSJed Brown 138c4762a1bSJed Brown test: 139c4762a1bSJed Brown suffix: unsorted_baij_mpi 140c4762a1bSJed Brown nsize: 3 141dfd57a17SPierre Jolivet requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex 142c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/cfd.1.10 -nd 8 -mat_type baij -test_unsorted 143c4762a1bSJed Brown 144c4762a1bSJed Brown test: 145c4762a1bSJed Brown suffix: unsorted_baij_seq 146dfd57a17SPierre Jolivet requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex 147c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/cfd.1.10 -nd 8 -mat_type baij -test_unsorted 148c4762a1bSJed Brown 149c4762a1bSJed Brown test: 150c4762a1bSJed Brown suffix: unsorted_mpi 151c4762a1bSJed Brown nsize: 3 152dfd57a17SPierre Jolivet requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex 153c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/arco1 -nd 8 -test_unsorted 154c4762a1bSJed Brown 155c4762a1bSJed Brown test: 156c4762a1bSJed Brown suffix: unsorted_seq 157dfd57a17SPierre Jolivet requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex 158c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/arco1 -nd 8 -test_unsorted 159c4762a1bSJed Brown 160c4762a1bSJed Brown TEST*/ 161