1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests the parallel case for MatIncreaseOverlap(). Input arguments are:\n\ 3c4762a1bSJed Brown -f <input_file> : file to load. For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\ 4c4762a1bSJed Brown -nd <size> : > 0 number of domains per processor \n\ 5c4762a1bSJed Brown -ov <overlap> : >=0 amount of overlap between domains\n\n"; 6c4762a1bSJed Brown 7c4762a1bSJed Brown #include <petscmat.h> 8c4762a1bSJed Brown 9*9371c9d4SSatish Balay PetscErrorCode ISAllGatherDisjoint(IS iis, IS **ois) { 10c4762a1bSJed Brown IS *is2, is; 11c4762a1bSJed Brown const PetscInt *idxs; 12c4762a1bSJed Brown PetscInt i, ls, *sizes; 13c4762a1bSJed Brown PetscMPIInt size; 14c4762a1bSJed Brown 15c4762a1bSJed Brown PetscFunctionBeginUser; 169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)iis), &size)); 179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &is2)); 189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &sizes)); 199566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iis, &ls)); 20c4762a1bSJed Brown /* we don't have a public ISGetLayout */ 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&ls, 1, MPIU_INT, sizes, 1, MPIU_INT, PetscObjectComm((PetscObject)iis))); 229566063dSJacob Faibussowitsch PetscCall(ISAllGather(iis, &is)); 239566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &idxs)); 24c4762a1bSJed Brown for (i = 0, ls = 0; i < size; i++) { 259566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, sizes[i], idxs + ls, PETSC_COPY_VALUES, &is2[i])); 26c4762a1bSJed Brown ls += sizes[i]; 27c4762a1bSJed Brown } 289566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &idxs)); 299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 309566063dSJacob Faibussowitsch PetscCall(PetscFree(sizes)); 31c4762a1bSJed Brown *ois = is2; 32c4762a1bSJed Brown PetscFunctionReturn(0); 33c4762a1bSJed Brown } 34c4762a1bSJed Brown 35*9371c9d4SSatish Balay int main(int argc, char **args) { 36c4762a1bSJed Brown PetscInt nd = 2, ov = 1, ndpar, i, start, m, n, end, lsize; 37c4762a1bSJed Brown PetscMPIInt rank; 38c4762a1bSJed Brown PetscBool flg, useND = PETSC_FALSE; 39c4762a1bSJed Brown Mat A, B; 40c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; 41c4762a1bSJed Brown PetscViewer fd; 42c4762a1bSJed Brown IS *is1, *is2; 43c4762a1bSJed Brown PetscRandom r; 44c4762a1bSJed Brown PetscScalar rand; 45c4762a1bSJed Brown 46327415f7SBarry Smith PetscFunctionBeginUser; 479566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 499566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &flg)); 5028b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must use -f filename to indicate a file containing a PETSc binary matrix"); 519566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL)); 529566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL)); 539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-nested_dissection", &useND, NULL)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Read matrix */ 569566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd)); 579566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 589566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATMPIAIJ)); 599566063dSJacob Faibussowitsch PetscCall(MatLoad(A, fd)); 609566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 619566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* Read the matrix again as a sequential matrix */ 649566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, file, FILE_MODE_READ, &fd)); 659566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &B)); 669566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATSEQAIJ)); 679566063dSJacob Faibussowitsch PetscCall(MatLoad(B, fd)); 689566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 699566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* Create the IS corresponding to subdomains */ 72c4762a1bSJed Brown if (useND) { 73c4762a1bSJed Brown MatPartitioning part; 74c4762a1bSJed Brown IS ndmap; 75c4762a1bSJed Brown PetscMPIInt size; 76c4762a1bSJed Brown 77c4762a1bSJed Brown ndpar = 1; 789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 79c4762a1bSJed Brown nd = (PetscInt)size; 809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ndpar, &is1)); 819566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(PETSC_COMM_WORLD, &part)); 829566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(part, A)); 839566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(part)); 849566063dSJacob Faibussowitsch PetscCall(MatPartitioningApplyND(part, &ndmap)); 859566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&part)); 869566063dSJacob Faibussowitsch PetscCall(ISBuildTwoSided(ndmap, NULL, &is1[0])); 879566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ndmap)); 889566063dSJacob Faibussowitsch PetscCall(ISAllGatherDisjoint(is1[0], &is2)); 89c4762a1bSJed Brown } else { 90c4762a1bSJed Brown /* Create the random Index Sets */ 919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is1)); 929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is2)); 93c4762a1bSJed Brown 949566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &m, &n)); 959566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r)); 969566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(r)); 97c4762a1bSJed Brown for (i = 0; i < nd; i++) { 989566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r, &rand)); 99c4762a1bSJed Brown start = (PetscInt)(rand * m); 1009566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r, &rand)); 101c4762a1bSJed Brown end = (PetscInt)(rand * m); 102c4762a1bSJed Brown lsize = end - start; 103*9371c9d4SSatish Balay if (start > end) { 104*9371c9d4SSatish Balay start = end; 105*9371c9d4SSatish Balay lsize = -lsize; 106*9371c9d4SSatish Balay } 1079566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, lsize, start, 1, is1 + i)); 1089566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, lsize, start, 1, is2 + i)); 109c4762a1bSJed Brown } 110c4762a1bSJed Brown ndpar = nd; 1119566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&r)); 112c4762a1bSJed Brown } 1139566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A, ndpar, is1, ov)); 1149566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(B, nd, is2, ov)); 115c4762a1bSJed Brown if (useND) { 116c4762a1bSJed Brown IS *is; 117c4762a1bSJed Brown 1189566063dSJacob Faibussowitsch PetscCall(ISAllGatherDisjoint(is1[0], &is)); 1199566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1[0])); 1209566063dSJacob Faibussowitsch PetscCall(PetscFree(is1)); 121c4762a1bSJed Brown is1 = is; 122c4762a1bSJed Brown } 123c4762a1bSJed Brown /* Now see if the serial and parallel case have the same answers */ 124c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 1259566063dSJacob Faibussowitsch PetscCall(ISEqual(is1[i], is2[i], &flg)); 126c4762a1bSJed Brown if (!flg) { 1279566063dSJacob Faibussowitsch PetscCall(ISViewFromOptions(is1[i], NULL, "-err_view")); 1289566063dSJacob Faibussowitsch PetscCall(ISViewFromOptions(is2[i], NULL, "-err_view")); 12998921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "proc:[%d], i=%" PetscInt_FMT ", flg =%d", rank, i, (int)flg); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown } 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* Free allocated memory */ 134c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 1359566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1[i])); 1369566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2[i])); 137c4762a1bSJed Brown } 1389566063dSJacob Faibussowitsch PetscCall(PetscFree(is1)); 1399566063dSJacob Faibussowitsch PetscCall(PetscFree(is2)); 1409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1429566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 143b122ec5aSJacob Faibussowitsch return 0; 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146c4762a1bSJed Brown /*TEST 147c4762a1bSJed Brown 148c4762a1bSJed Brown build: 149c4762a1bSJed Brown requires: !complex 150c4762a1bSJed Brown 151c4762a1bSJed Brown testset: 152c4762a1bSJed Brown nsize: 5 153dfd57a17SPierre Jolivet requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex 154c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/arco1 -viewer_binary_skip_info -ov 2 155c4762a1bSJed Brown output_file: output/ex40_1.out 156c4762a1bSJed Brown test: 157c4762a1bSJed Brown suffix: 1 158c4762a1bSJed Brown args: -nd 7 159c4762a1bSJed Brown test: 160c4762a1bSJed Brown requires: parmetis 161c4762a1bSJed Brown suffix: 1_nd 162c4762a1bSJed Brown args: -nested_dissection -mat_partitioning_type parmetis 163c4762a1bSJed Brown 164c4762a1bSJed Brown testset: 165c4762a1bSJed Brown nsize: 3 166dfd57a17SPierre Jolivet requires: double !defined(PETSC_USE_64BIT_INDICES) !complex 167c4762a1bSJed Brown args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64 -mat_increase_overlap_scalable 1 -ov 2 168c4762a1bSJed Brown output_file: output/ex40_1.out 169c4762a1bSJed Brown test: 170c4762a1bSJed Brown suffix: 2 171c4762a1bSJed Brown args: -nd 7 172c4762a1bSJed Brown test: 173c4762a1bSJed Brown requires: parmetis 174c4762a1bSJed Brown suffix: 2_nd 175c4762a1bSJed Brown args: -nested_dissection -mat_partitioning_type parmetis 176c4762a1bSJed Brown 177c4762a1bSJed Brown TEST*/ 178