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 9d71ae5a4SJacob Faibussowitsch PetscErrorCode ISAllGatherDisjoint(IS iis, IS **ois) 10d71ae5a4SJacob Faibussowitsch { 11c4762a1bSJed Brown IS *is2, is; 12c4762a1bSJed Brown const PetscInt *idxs; 13c4762a1bSJed Brown PetscInt i, ls, *sizes; 14c4762a1bSJed Brown PetscMPIInt size; 15c4762a1bSJed Brown 16c4762a1bSJed Brown PetscFunctionBeginUser; 179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)iis), &size)); 189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &is2)); 199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &sizes)); 209566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iis, &ls)); 21c4762a1bSJed Brown /* we don't have a public ISGetLayout */ 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&ls, 1, MPIU_INT, sizes, 1, MPIU_INT, PetscObjectComm((PetscObject)iis))); 239566063dSJacob Faibussowitsch PetscCall(ISAllGather(iis, &is)); 249566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &idxs)); 25c4762a1bSJed Brown for (i = 0, ls = 0; i < size; i++) { 269566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, sizes[i], idxs + ls, PETSC_COPY_VALUES, &is2[i])); 27c4762a1bSJed Brown ls += sizes[i]; 28c4762a1bSJed Brown } 299566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &idxs)); 309566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 319566063dSJacob Faibussowitsch PetscCall(PetscFree(sizes)); 32c4762a1bSJed Brown *ois = is2; 33*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34c4762a1bSJed Brown } 35c4762a1bSJed Brown 36d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 37d71ae5a4SJacob Faibussowitsch { 38c4762a1bSJed Brown PetscInt nd = 2, ov = 1, ndpar, i, start, m, n, end, lsize; 39c4762a1bSJed Brown PetscMPIInt rank; 40c4762a1bSJed Brown PetscBool flg, useND = PETSC_FALSE; 41c4762a1bSJed Brown Mat A, B; 42c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; 43c4762a1bSJed Brown PetscViewer fd; 44c4762a1bSJed Brown IS *is1, *is2; 45c4762a1bSJed Brown PetscRandom r; 46c4762a1bSJed Brown PetscScalar rand; 47c4762a1bSJed Brown 48327415f7SBarry Smith PetscFunctionBeginUser; 499566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 509566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 519566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &flg)); 5228b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must use -f filename to indicate a file containing a PETSc binary matrix"); 539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL)); 549566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL)); 559566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-nested_dissection", &useND, NULL)); 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* Read matrix */ 589566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd)); 599566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 609566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATMPIAIJ)); 619566063dSJacob Faibussowitsch PetscCall(MatLoad(A, fd)); 629566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 639566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* Read the matrix again as a sequential matrix */ 669566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, file, FILE_MODE_READ, &fd)); 679566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &B)); 689566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATSEQAIJ)); 699566063dSJacob Faibussowitsch PetscCall(MatLoad(B, fd)); 709566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 719566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* Create the IS corresponding to subdomains */ 74c4762a1bSJed Brown if (useND) { 75c4762a1bSJed Brown MatPartitioning part; 76c4762a1bSJed Brown IS ndmap; 77c4762a1bSJed Brown PetscMPIInt size; 78c4762a1bSJed Brown 79c4762a1bSJed Brown ndpar = 1; 809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 81c4762a1bSJed Brown nd = (PetscInt)size; 829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ndpar, &is1)); 839566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(PETSC_COMM_WORLD, &part)); 849566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(part, A)); 859566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(part)); 869566063dSJacob Faibussowitsch PetscCall(MatPartitioningApplyND(part, &ndmap)); 879566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&part)); 889566063dSJacob Faibussowitsch PetscCall(ISBuildTwoSided(ndmap, NULL, &is1[0])); 899566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ndmap)); 909566063dSJacob Faibussowitsch PetscCall(ISAllGatherDisjoint(is1[0], &is2)); 91c4762a1bSJed Brown } else { 92c4762a1bSJed Brown /* Create the random Index Sets */ 939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is1)); 949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is2)); 95c4762a1bSJed Brown 969566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &m, &n)); 979566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r)); 989566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(r)); 99c4762a1bSJed Brown for (i = 0; i < nd; i++) { 1009566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r, &rand)); 101c4762a1bSJed Brown start = (PetscInt)(rand * m); 1029566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r, &rand)); 103c4762a1bSJed Brown end = (PetscInt)(rand * m); 104c4762a1bSJed Brown lsize = end - start; 1059371c9d4SSatish Balay if (start > end) { 1069371c9d4SSatish Balay start = end; 1079371c9d4SSatish Balay lsize = -lsize; 1089371c9d4SSatish Balay } 1099566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, lsize, start, 1, is1 + i)); 1109566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, lsize, start, 1, is2 + i)); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown ndpar = nd; 1139566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&r)); 114c4762a1bSJed Brown } 1159566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A, ndpar, is1, ov)); 1169566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(B, nd, is2, ov)); 117c4762a1bSJed Brown if (useND) { 118c4762a1bSJed Brown IS *is; 119c4762a1bSJed Brown 1209566063dSJacob Faibussowitsch PetscCall(ISAllGatherDisjoint(is1[0], &is)); 1219566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1[0])); 1229566063dSJacob Faibussowitsch PetscCall(PetscFree(is1)); 123c4762a1bSJed Brown is1 = is; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown /* Now see if the serial and parallel case have the same answers */ 126c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 1279566063dSJacob Faibussowitsch PetscCall(ISEqual(is1[i], is2[i], &flg)); 128c4762a1bSJed Brown if (!flg) { 1299566063dSJacob Faibussowitsch PetscCall(ISViewFromOptions(is1[i], NULL, "-err_view")); 1309566063dSJacob Faibussowitsch PetscCall(ISViewFromOptions(is2[i], NULL, "-err_view")); 13198921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "proc:[%d], i=%" PetscInt_FMT ", flg =%d", rank, i, (int)flg); 132c4762a1bSJed Brown } 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* Free allocated memory */ 136c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 1379566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1[i])); 1389566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2[i])); 139c4762a1bSJed Brown } 1409566063dSJacob Faibussowitsch PetscCall(PetscFree(is1)); 1419566063dSJacob Faibussowitsch PetscCall(PetscFree(is2)); 1429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1449566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 145b122ec5aSJacob Faibussowitsch return 0; 146c4762a1bSJed Brown } 147c4762a1bSJed Brown 148c4762a1bSJed Brown /*TEST 149c4762a1bSJed Brown 150c4762a1bSJed Brown build: 151c4762a1bSJed Brown requires: !complex 152c4762a1bSJed Brown 153c4762a1bSJed Brown testset: 154c4762a1bSJed Brown nsize: 5 155dfd57a17SPierre Jolivet requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex 156c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/arco1 -viewer_binary_skip_info -ov 2 157c4762a1bSJed Brown output_file: output/ex40_1.out 158c4762a1bSJed Brown test: 159c4762a1bSJed Brown suffix: 1 160c4762a1bSJed Brown args: -nd 7 161c4762a1bSJed Brown test: 162c4762a1bSJed Brown requires: parmetis 163c4762a1bSJed Brown suffix: 1_nd 164c4762a1bSJed Brown args: -nested_dissection -mat_partitioning_type parmetis 165c4762a1bSJed Brown 166c4762a1bSJed Brown testset: 167c4762a1bSJed Brown nsize: 3 168dfd57a17SPierre Jolivet requires: double !defined(PETSC_USE_64BIT_INDICES) !complex 169c4762a1bSJed Brown args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64 -mat_increase_overlap_scalable 1 -ov 2 170c4762a1bSJed Brown output_file: output/ex40_1.out 171c4762a1bSJed Brown test: 172c4762a1bSJed Brown suffix: 2 173c4762a1bSJed Brown args: -nd 7 174c4762a1bSJed Brown test: 175c4762a1bSJed Brown requires: parmetis 176c4762a1bSJed Brown suffix: 2_nd 177c4762a1bSJed Brown args: -nested_dissection -mat_partitioning_type parmetis 178c4762a1bSJed Brown 179c4762a1bSJed Brown TEST*/ 180