1c4762a1bSJed Brown static char help[] = "Tests for parallel mesh loading and parallel topological interpolation\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petsc/private/dmpleximpl.h> 4c4762a1bSJed Brown /* List of test meshes 5c4762a1bSJed Brown 6c4762a1bSJed Brown Network 7c4762a1bSJed Brown ------- 8c4762a1bSJed Brown Test 0 (2 ranks): 9c4762a1bSJed Brown 10c4762a1bSJed Brown network=0: 11c4762a1bSJed Brown --------- 12c4762a1bSJed Brown cell 0 cell 1 cell 2 nCells-1 (edge) 13c4762a1bSJed Brown 0 ------ 1 ------ 2 ------ 3 -- -- v -- -- nCells (vertex) 14c4762a1bSJed Brown 15c4762a1bSJed Brown vertex distribution: 16c4762a1bSJed Brown rank 0: 0 1 17c4762a1bSJed Brown rank 1: 2 3 ... nCells 18c4762a1bSJed Brown cell(edge) distribution: 19c4762a1bSJed Brown rank 0: 0 1 20c4762a1bSJed Brown rank 1: 2 ... nCells-1 21c4762a1bSJed Brown 22c4762a1bSJed Brown network=1: 23c4762a1bSJed Brown --------- 24c4762a1bSJed Brown v2 25c4762a1bSJed Brown ^ 26c4762a1bSJed Brown | 27c4762a1bSJed Brown cell 2 28c4762a1bSJed Brown | 29c4762a1bSJed Brown v0 --cell 0--> v3--cell 1--> v1 30c4762a1bSJed Brown 31c4762a1bSJed Brown vertex distribution: 32c4762a1bSJed Brown rank 0: 0 1 3 33c4762a1bSJed Brown rank 1: 2 34c4762a1bSJed Brown cell(edge) distribution: 35c4762a1bSJed Brown rank 0: 0 1 36c4762a1bSJed Brown rank 1: 2 37c4762a1bSJed Brown 38c4762a1bSJed Brown example: 39c4762a1bSJed Brown mpiexec -n 2 ./ex18 -distribute 1 -dim 1 -orig_dm_view -dist_dm_view -dist_dm_view -petscpartitioner_type parmetis -ncells 50 40c4762a1bSJed Brown 41c4762a1bSJed Brown Triangle 42c4762a1bSJed Brown -------- 43c4762a1bSJed Brown Test 0 (2 ranks): 44c4762a1bSJed Brown Two triangles sharing a face 45c4762a1bSJed Brown 46c4762a1bSJed Brown 2 47c4762a1bSJed Brown / | \ 48c4762a1bSJed Brown / | \ 49c4762a1bSJed Brown / | \ 50c4762a1bSJed Brown 0 0 | 1 3 51c4762a1bSJed Brown \ | / 52c4762a1bSJed Brown \ | / 53c4762a1bSJed Brown \ | / 54c4762a1bSJed Brown 1 55c4762a1bSJed Brown 56c4762a1bSJed Brown vertex distribution: 57c4762a1bSJed Brown rank 0: 0 1 58c4762a1bSJed Brown rank 1: 2 3 59c4762a1bSJed Brown cell distribution: 60c4762a1bSJed Brown rank 0: 0 61c4762a1bSJed Brown rank 1: 1 62c4762a1bSJed Brown 63c4762a1bSJed Brown Test 1 (3 ranks): 64c4762a1bSJed Brown Four triangles partitioned across 3 ranks 65c4762a1bSJed Brown 66c4762a1bSJed Brown 0 _______ 3 67c4762a1bSJed Brown | \ / | 68c4762a1bSJed Brown | \ 1 / | 69c4762a1bSJed Brown | \ / | 70c4762a1bSJed Brown | 0 2 2 | 71c4762a1bSJed Brown | / \ | 72c4762a1bSJed Brown | / 3 \ | 73c4762a1bSJed Brown | / \ | 74c4762a1bSJed Brown 1 ------- 4 75c4762a1bSJed Brown 76c4762a1bSJed Brown vertex distribution: 77c4762a1bSJed Brown rank 0: 0 1 78c4762a1bSJed Brown rank 1: 2 3 79c4762a1bSJed Brown rank 2: 4 80c4762a1bSJed Brown cell distribution: 81c4762a1bSJed Brown rank 0: 0 82c4762a1bSJed Brown rank 1: 1 83c4762a1bSJed Brown rank 2: 2 3 84c4762a1bSJed Brown 85c4762a1bSJed Brown Test 2 (3 ranks): 86c4762a1bSJed Brown Four triangles partitioned across 3 ranks 87c4762a1bSJed Brown 88c4762a1bSJed Brown 1 _______ 3 89c4762a1bSJed Brown | \ / | 90c4762a1bSJed Brown | \ 1 / | 91c4762a1bSJed Brown | \ / | 92c4762a1bSJed Brown | 0 0 2 | 93c4762a1bSJed Brown | / \ | 94c4762a1bSJed Brown | / 3 \ | 95c4762a1bSJed Brown | / \ | 96c4762a1bSJed Brown 2 ------- 4 97c4762a1bSJed Brown 98c4762a1bSJed Brown vertex distribution: 99c4762a1bSJed Brown rank 0: 0 1 100c4762a1bSJed Brown rank 1: 2 3 101c4762a1bSJed Brown rank 2: 4 102c4762a1bSJed Brown cell distribution: 103c4762a1bSJed Brown rank 0: 0 104c4762a1bSJed Brown rank 1: 1 105c4762a1bSJed Brown rank 2: 2 3 106c4762a1bSJed Brown 107c4762a1bSJed Brown Tetrahedron 108c4762a1bSJed Brown ----------- 109c4762a1bSJed Brown Test 0: 110c4762a1bSJed Brown Two tets sharing a face 111c4762a1bSJed Brown 112c4762a1bSJed Brown cell 3 _______ cell 113c4762a1bSJed Brown 0 / | \ \ 1 114c4762a1bSJed Brown / | \ \ 115c4762a1bSJed Brown / | \ \ 116c4762a1bSJed Brown 0----|----4-----2 117c4762a1bSJed Brown \ | / / 118c4762a1bSJed Brown \ | / / 119c4762a1bSJed Brown \ | / / 120c4762a1bSJed Brown 1------- 121c4762a1bSJed Brown y 122c4762a1bSJed Brown | x 123c4762a1bSJed Brown |/ 124c4762a1bSJed Brown *----z 125c4762a1bSJed Brown 126c4762a1bSJed Brown vertex distribution: 127c4762a1bSJed Brown rank 0: 0 1 128c4762a1bSJed Brown rank 1: 2 3 4 129c4762a1bSJed Brown cell distribution: 130c4762a1bSJed Brown rank 0: 0 131c4762a1bSJed Brown rank 1: 1 132c4762a1bSJed Brown 133c4762a1bSJed Brown Quadrilateral 134c4762a1bSJed Brown ------------- 135c4762a1bSJed Brown Test 0 (2 ranks): 136c4762a1bSJed Brown Two quads sharing a face 137c4762a1bSJed Brown 138c4762a1bSJed Brown 3-------2-------5 139c4762a1bSJed Brown | | | 140c4762a1bSJed Brown | 0 | 1 | 141c4762a1bSJed Brown | | | 142c4762a1bSJed Brown 0-------1-------4 143c4762a1bSJed Brown 144c4762a1bSJed Brown vertex distribution: 145c4762a1bSJed Brown rank 0: 0 1 2 146c4762a1bSJed Brown rank 1: 3 4 5 147c4762a1bSJed Brown cell distribution: 148c4762a1bSJed Brown rank 0: 0 149c4762a1bSJed Brown rank 1: 1 150c4762a1bSJed Brown 151c4762a1bSJed Brown TODO Test 1: 152c4762a1bSJed Brown A quad and a triangle sharing a face 153c4762a1bSJed Brown 154c4762a1bSJed Brown 5-------4 155c4762a1bSJed Brown | | \ 156c4762a1bSJed Brown | 0 | \ 157c4762a1bSJed Brown | | 1 \ 158c4762a1bSJed Brown 2-------3----6 159c4762a1bSJed Brown 160c4762a1bSJed Brown Hexahedron 161c4762a1bSJed Brown ---------- 162c4762a1bSJed Brown Test 0 (2 ranks): 163c4762a1bSJed Brown Two hexes sharing a face 164c4762a1bSJed Brown 165c4762a1bSJed Brown cell 7-------------6-------------11 cell 166c4762a1bSJed Brown 0 /| /| /| 1 167c4762a1bSJed Brown / | F1 / | F7 / | 168c4762a1bSJed Brown / | / | / | 169c4762a1bSJed Brown 4-------------5-------------10 | 170c4762a1bSJed Brown | | F4 | | F10 | | 171c4762a1bSJed Brown | | | | | | 172c4762a1bSJed Brown |F5 | |F3 | |F9 | 173c4762a1bSJed Brown | | F2 | | F8 | | 174c4762a1bSJed Brown | 3---------|---2---------|---9 175c4762a1bSJed Brown | / | / | / 176c4762a1bSJed Brown | / F0 | / F6 | / 177c4762a1bSJed Brown |/ |/ |/ 178c4762a1bSJed Brown 0-------------1-------------8 179c4762a1bSJed Brown 180c4762a1bSJed Brown vertex distribution: 181c4762a1bSJed Brown rank 0: 0 1 2 3 4 5 182c4762a1bSJed Brown rank 1: 6 7 8 9 10 11 183c4762a1bSJed Brown cell distribution: 184c4762a1bSJed Brown rank 0: 0 185c4762a1bSJed Brown rank 1: 1 186c4762a1bSJed Brown 187c4762a1bSJed Brown */ 188c4762a1bSJed Brown 1899371c9d4SSatish Balay typedef enum { 1909371c9d4SSatish Balay NONE, 1919371c9d4SSatish Balay CREATE, 1929371c9d4SSatish Balay AFTER_CREATE, 1939371c9d4SSatish Balay AFTER_DISTRIBUTE 1949371c9d4SSatish Balay } InterpType; 195c4762a1bSJed Brown 196c4762a1bSJed Brown typedef struct { 197c4762a1bSJed Brown PetscInt debug; /* The debugging level */ 198c4762a1bSJed Brown PetscInt testNum; /* Indicates the mesh to create */ 199c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 200c4762a1bSJed Brown PetscBool cellSimplex; /* Use simplices or hexes */ 201c4762a1bSJed Brown PetscBool distribute; /* Distribute the mesh */ 202c4762a1bSJed Brown InterpType interpolate; /* Interpolate the mesh before or after DMPlexDistribute() */ 203c4762a1bSJed Brown PetscBool useGenerator; /* Construct mesh with a mesh generator */ 204c4762a1bSJed Brown PetscBool testOrientIF; /* Test for different original interface orientations */ 205c4762a1bSJed Brown PetscBool testHeavy; /* Run the heavy PointSF test */ 206c4762a1bSJed Brown PetscBool customView; /* Show results of DMPlexIsInterpolated() etc. */ 207c4762a1bSJed Brown PetscInt ornt[2]; /* Orientation of interface on rank 0 and rank 1 */ 208c4762a1bSJed Brown PetscInt faces[3]; /* Number of faces per dimension for generator */ 209d3e1f4ccSVaclav Hapla PetscScalar coords[128]; 210c4762a1bSJed Brown PetscReal coordsTol; 211c4762a1bSJed Brown PetscInt ncoords; 212c4762a1bSJed Brown PetscInt pointsToExpand[128]; 213c4762a1bSJed Brown PetscInt nPointsToExpand; 214c4762a1bSJed Brown PetscBool testExpandPointsEmpty; 215c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */ 216c4762a1bSJed Brown } AppCtx; 217c4762a1bSJed Brown 218c4762a1bSJed Brown struct _n_PortableBoundary { 219c4762a1bSJed Brown Vec coordinates; 220c4762a1bSJed Brown PetscInt depth; 221c4762a1bSJed Brown PetscSection *sections; 222c4762a1bSJed Brown }; 223c4762a1bSJed Brown typedef struct _n_PortableBoundary *PortableBoundary; 224c4762a1bSJed Brown 225c4762a1bSJed Brown static PetscLogStage stage[3]; 226c4762a1bSJed Brown 227c4762a1bSJed Brown static PetscErrorCode DMPlexCheckPointSFHeavy(DM, PortableBoundary); 228c4762a1bSJed Brown static PetscErrorCode DMPlexSetOrientInterface_Private(DM, PetscBool); 229c4762a1bSJed Brown static PetscErrorCode DMPlexGetExpandedBoundary_Private(DM, PortableBoundary *); 230c4762a1bSJed Brown static PetscErrorCode DMPlexExpandedConesToFaces_Private(DM, IS, PetscSection, IS *); 231c4762a1bSJed Brown 232d71ae5a4SJacob Faibussowitsch static PetscErrorCode PortableBoundaryDestroy(PortableBoundary *bnd) 233d71ae5a4SJacob Faibussowitsch { 234c4762a1bSJed Brown PetscInt d; 235c4762a1bSJed Brown 236c4762a1bSJed Brown PetscFunctionBegin; 2373ba16761SJacob Faibussowitsch if (!*bnd) PetscFunctionReturn(PETSC_SUCCESS); 2389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*bnd)->coordinates)); 23948a46eb9SPierre Jolivet for (d = 0; d < (*bnd)->depth; d++) PetscCall(PetscSectionDestroy(&(*bnd)->sections[d])); 2409566063dSJacob Faibussowitsch PetscCall(PetscFree((*bnd)->sections)); 2419566063dSJacob Faibussowitsch PetscCall(PetscFree(*bnd)); 2423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 243c4762a1bSJed Brown } 244c4762a1bSJed Brown 245d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 246d71ae5a4SJacob Faibussowitsch { 247c4762a1bSJed Brown const char *interpTypes[4] = {"none", "create", "after_create", "after_distribute"}; 248c4762a1bSJed Brown PetscInt interp = NONE, dim; 249c4762a1bSJed Brown PetscBool flg1, flg2; 250c4762a1bSJed Brown 251c4762a1bSJed Brown PetscFunctionBegin; 252c4762a1bSJed Brown options->debug = 0; 253c4762a1bSJed Brown options->testNum = 0; 254c4762a1bSJed Brown options->dim = 2; 255c4762a1bSJed Brown options->cellSimplex = PETSC_TRUE; 256c4762a1bSJed Brown options->distribute = PETSC_FALSE; 257c4762a1bSJed Brown options->interpolate = NONE; 258c4762a1bSJed Brown options->useGenerator = PETSC_FALSE; 259c4762a1bSJed Brown options->testOrientIF = PETSC_FALSE; 260c4762a1bSJed Brown options->testHeavy = PETSC_TRUE; 261c4762a1bSJed Brown options->customView = PETSC_FALSE; 262c4762a1bSJed Brown options->testExpandPointsEmpty = PETSC_FALSE; 263c4762a1bSJed Brown options->ornt[0] = 0; 264c4762a1bSJed Brown options->ornt[1] = 0; 265c4762a1bSJed Brown options->faces[0] = 2; 266c4762a1bSJed Brown options->faces[1] = 2; 267c4762a1bSJed Brown options->faces[2] = 2; 268c4762a1bSJed Brown options->filename[0] = '\0'; 269c4762a1bSJed Brown options->coordsTol = PETSC_DEFAULT; 270c4762a1bSJed Brown 271d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX"); 2729566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-debug", "The debugging level", "ex18.c", options->debug, &options->debug, NULL, 0)); 2739566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-testnum", "The mesh to create", "ex18.c", options->testNum, &options->testNum, NULL, 0)); 2749566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-cell_simplex", "Generate simplices if true, otherwise hexes", "ex18.c", options->cellSimplex, &options->cellSimplex, NULL)); 2759566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-distribute", "Distribute the mesh", "ex18.c", options->distribute, &options->distribute, NULL)); 2769566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-interpolate", "Type of mesh interpolation (none, create, after_create, after_distribute)", "ex18.c", interpTypes, 4, interpTypes[options->interpolate], &interp, NULL)); 277c4762a1bSJed Brown options->interpolate = (InterpType)interp; 2781dca8a05SBarry Smith PetscCheck(options->distribute || options->interpolate != AFTER_DISTRIBUTE, comm, PETSC_ERR_SUP, "-interpolate after_distribute needs -distribute 1"); 2799566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_generator", "Use a mesh generator to build the mesh", "ex18.c", options->useGenerator, &options->useGenerator, NULL)); 280c4762a1bSJed Brown options->ncoords = 128; 2819566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalarArray("-view_vertices_from_coords", "Print DAG points corresponding to vertices with given coordinates", "ex18.c", options->coords, &options->ncoords, NULL)); 2829566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-view_vertices_from_coords_tol", "Tolerance for -view_vertices_from_coords", "ex18.c", options->coordsTol, &options->coordsTol, NULL)); 283c4762a1bSJed Brown options->nPointsToExpand = 128; 2849566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-test_expand_points", "Expand given array of DAG point using DMPlexGetConeRecursive() and print results", "ex18.c", options->pointsToExpand, &options->nPointsToExpand, NULL)); 28548a46eb9SPierre Jolivet if (options->nPointsToExpand) PetscCall(PetscOptionsBool("-test_expand_points_empty", "For -test_expand_points, rank 0 will have empty input array", "ex18.c", options->testExpandPointsEmpty, &options->testExpandPointsEmpty, NULL)); 2869566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_heavy", "Run the heavy PointSF test", "ex18.c", options->testHeavy, &options->testHeavy, NULL)); 2879566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-custom_view", "Custom DMPlex view", "ex18.c", options->customView, &options->customView, NULL)); 2889566063dSJacob Faibussowitsch PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex18.c", options->dim, &options->dim, &flg1, 1, 3)); 289c4762a1bSJed Brown dim = 3; 2909566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-faces", "Number of faces per dimension", "ex18.c", options->faces, &dim, &flg2)); 291c4762a1bSJed Brown if (flg2) { 2921dca8a05SBarry Smith PetscCheck(!flg1 || dim == options->dim, comm, PETSC_ERR_ARG_OUTOFRANGE, "specified -dim %" PetscInt_FMT " is not equal to length %" PetscInt_FMT " of -faces (note that -dim can be omitted)", options->dim, dim); 293c4762a1bSJed Brown options->dim = dim; 294c4762a1bSJed Brown } 2959566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-filename", "The mesh file", "ex18.c", options->filename, options->filename, sizeof(options->filename), NULL)); 2969566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-rotate_interface_0", "Rotation (relative orientation) of interface on rank 0; implies -interpolate create -distribute 0", "ex18.c", options->ornt[0], &options->ornt[0], &options->testOrientIF, 0)); 2979566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-rotate_interface_1", "Rotation (relative orientation) of interface on rank 1; implies -interpolate create -distribute 0", "ex18.c", options->ornt[1], &options->ornt[1], &flg2, 0)); 29808401ef6SPierre Jolivet PetscCheck(flg2 == options->testOrientIF, comm, PETSC_ERR_ARG_OUTOFRANGE, "neither or both -rotate_interface_0 -rotate_interface_1 must be set"); 299c4762a1bSJed Brown if (options->testOrientIF) { 300c4762a1bSJed Brown PetscInt i; 301c4762a1bSJed Brown for (i = 0; i < 2; i++) { 302c4762a1bSJed Brown if (options->ornt[i] >= 10) options->ornt[i] = -(options->ornt[i] - 10); /* 11 12 13 become -1 -2 -3 */ 303c4762a1bSJed Brown } 304c4762a1bSJed Brown options->filename[0] = 0; 305c4762a1bSJed Brown options->useGenerator = PETSC_FALSE; 306c4762a1bSJed Brown options->dim = 3; 307c4762a1bSJed Brown options->cellSimplex = PETSC_TRUE; 308c4762a1bSJed Brown options->interpolate = CREATE; 309c4762a1bSJed Brown options->distribute = PETSC_FALSE; 310c4762a1bSJed Brown } 311d0609cedSBarry Smith PetscOptionsEnd(); 3123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 313c4762a1bSJed Brown } 314c4762a1bSJed Brown 315d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh_1D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm) 316d71ae5a4SJacob Faibussowitsch { 317c4762a1bSJed Brown PetscInt testNum = user->testNum; 318c4762a1bSJed Brown PetscMPIInt rank, size; 319c4762a1bSJed Brown PetscInt numCorners = 2, i; 320c4762a1bSJed Brown PetscInt numCells, numVertices, network; 321a4a685f2SJacob Faibussowitsch PetscInt *cells; 322c4762a1bSJed Brown PetscReal *coords; 323c4762a1bSJed Brown 324c4762a1bSJed Brown PetscFunctionBegin; 3259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 3269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 32763a3b9bcSJacob Faibussowitsch PetscCheck(size <= 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for <=2 processes", testNum); 328c4762a1bSJed Brown 329c4762a1bSJed Brown numCells = 3; 3309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ncells", &numCells, NULL)); 33163a3b9bcSJacob Faibussowitsch PetscCheck(numCells >= 3, comm, PETSC_ERR_ARG_OUTOFRANGE, "Test ncells %" PetscInt_FMT " must >=3", numCells); 332c4762a1bSJed Brown 333c4762a1bSJed Brown if (size == 1) { 334c4762a1bSJed Brown numVertices = numCells + 1; 3359566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2 * numCells, &cells, 2 * numVertices, &coords)); 336c4762a1bSJed Brown for (i = 0; i < numCells; i++) { 3379371c9d4SSatish Balay cells[2 * i] = i; 3389371c9d4SSatish Balay cells[2 * i + 1] = i + 1; 3399371c9d4SSatish Balay coords[2 * i] = i; 3409371c9d4SSatish Balay coords[2 * i + 1] = i + 1; 341c4762a1bSJed Brown } 342c4762a1bSJed Brown 3439566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListPetsc(comm, user->dim, numCells, numVertices, numCorners, PETSC_FALSE, cells, user->dim, coords, dm)); 3449566063dSJacob Faibussowitsch PetscCall(PetscFree2(cells, coords)); 3453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 346c4762a1bSJed Brown } 347c4762a1bSJed Brown 348c4762a1bSJed Brown network = 0; 3499566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-network_case", &network, NULL)); 350c4762a1bSJed Brown if (network == 0) { 351c4762a1bSJed Brown switch (rank) { 3529371c9d4SSatish Balay case 0: { 353c4762a1bSJed Brown numCells = 2; 354c4762a1bSJed Brown numVertices = numCells; 3559566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2 * numCells, &cells, 2 * numCells, &coords)); 3569371c9d4SSatish Balay cells[0] = 0; 3579371c9d4SSatish Balay cells[1] = 1; 3589371c9d4SSatish Balay cells[2] = 1; 3599371c9d4SSatish Balay cells[3] = 2; 3609371c9d4SSatish Balay coords[0] = 0.; 3619371c9d4SSatish Balay coords[1] = 1.; 3629371c9d4SSatish Balay coords[2] = 1.; 3639371c9d4SSatish Balay coords[3] = 2.; 3649371c9d4SSatish Balay } break; 3659371c9d4SSatish Balay case 1: { 366c4762a1bSJed Brown numCells -= 2; 367c4762a1bSJed Brown numVertices = numCells + 1; 3689566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2 * numCells, &cells, 2 * numCells, &coords)); 369c4762a1bSJed Brown for (i = 0; i < numCells; i++) { 3709371c9d4SSatish Balay cells[2 * i] = 2 + i; 3719371c9d4SSatish Balay cells[2 * i + 1] = 2 + i + 1; 3729371c9d4SSatish Balay coords[2 * i] = 2 + i; 3739371c9d4SSatish Balay coords[2 * i + 1] = 2 + i + 1; 374c4762a1bSJed Brown } 3759371c9d4SSatish Balay } break; 376d71ae5a4SJacob Faibussowitsch default: 377d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank); 378c4762a1bSJed Brown } 379c4762a1bSJed Brown } else { /* network_case = 1 */ 380c4762a1bSJed Brown /* ----------------------- */ 381c4762a1bSJed Brown switch (rank) { 3829371c9d4SSatish Balay case 0: { 383c4762a1bSJed Brown numCells = 2; 384c4762a1bSJed Brown numVertices = 3; 3859566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2 * numCells, &cells, 2 * numCells, &coords)); 3869371c9d4SSatish Balay cells[0] = 0; 3879371c9d4SSatish Balay cells[1] = 3; 3889371c9d4SSatish Balay cells[2] = 3; 3899371c9d4SSatish Balay cells[3] = 1; 3909371c9d4SSatish Balay } break; 3919371c9d4SSatish Balay case 1: { 392c4762a1bSJed Brown numCells = 1; 393c4762a1bSJed Brown numVertices = 1; 3949566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2 * numCells, &cells, 2 * numCells, &coords)); 3959371c9d4SSatish Balay cells[0] = 3; 3969371c9d4SSatish Balay cells[1] = 2; 3979371c9d4SSatish Balay } break; 398d71ae5a4SJacob Faibussowitsch default: 399d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank); 400c4762a1bSJed Brown } 401c4762a1bSJed Brown } 4029566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, PETSC_FALSE, cells, user->dim, coords, NULL, NULL, dm)); 4039566063dSJacob Faibussowitsch PetscCall(PetscFree2(cells, coords)); 4043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 405c4762a1bSJed Brown } 406c4762a1bSJed Brown 407d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSimplex_2D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm) 408d71ae5a4SJacob Faibussowitsch { 409c4762a1bSJed Brown PetscInt testNum = user->testNum, p; 410c4762a1bSJed Brown PetscMPIInt rank, size; 411c4762a1bSJed Brown 412c4762a1bSJed Brown PetscFunctionBegin; 4139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 4149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 415c4762a1bSJed Brown switch (testNum) { 416c4762a1bSJed Brown case 0: 41763a3b9bcSJacob Faibussowitsch PetscCheck(size == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 2 processes", testNum); 418c4762a1bSJed Brown switch (rank) { 4199371c9d4SSatish Balay case 0: { 420c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 2, numCorners = 3; 421a4a685f2SJacob Faibussowitsch const PetscInt cells[3] = {0, 1, 2}; 422c4762a1bSJed Brown PetscReal coords[4] = {-0.5, 0.5, 0.0, 0.0}; 423c4762a1bSJed Brown PetscInt markerPoints[6] = {1, 1, 2, 1, 3, 1}; 424c4762a1bSJed Brown 4259566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 4269566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 4279371c9d4SSatish Balay } break; 4289371c9d4SSatish Balay case 1: { 429c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 2, numCorners = 3; 430a4a685f2SJacob Faibussowitsch const PetscInt cells[3] = {1, 3, 2}; 431c4762a1bSJed Brown PetscReal coords[4] = {0.0, 1.0, 0.5, 0.5}; 432c4762a1bSJed Brown PetscInt markerPoints[6] = {1, 1, 2, 1, 3, 1}; 433c4762a1bSJed Brown 4349566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 4359566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 4369371c9d4SSatish Balay } break; 437d71ae5a4SJacob Faibussowitsch default: 438d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank); 439c4762a1bSJed Brown } 440c4762a1bSJed Brown break; 441c4762a1bSJed Brown case 1: 44263a3b9bcSJacob Faibussowitsch PetscCheck(size == 3, comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 3 processes", testNum); 443c4762a1bSJed Brown switch (rank) { 4449371c9d4SSatish Balay case 0: { 445c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 2, numCorners = 3; 446a4a685f2SJacob Faibussowitsch const PetscInt cells[3] = {0, 1, 2}; 447c4762a1bSJed Brown PetscReal coords[4] = {0.0, 1.0, 0.0, 0.0}; 448c4762a1bSJed Brown PetscInt markerPoints[6] = {1, 1, 2, 1, 3, 1}; 449c4762a1bSJed Brown 4509566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 4519566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 4529371c9d4SSatish Balay } break; 4539371c9d4SSatish Balay case 1: { 454c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 2, numCorners = 3; 455a4a685f2SJacob Faibussowitsch const PetscInt cells[3] = {0, 2, 3}; 456c4762a1bSJed Brown PetscReal coords[4] = {0.5, 0.5, 1.0, 1.0}; 457c4762a1bSJed Brown PetscInt markerPoints[6] = {1, 1, 2, 1, 3, 1}; 458c4762a1bSJed Brown 4599566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 4609566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 4619371c9d4SSatish Balay } break; 4629371c9d4SSatish Balay case 2: { 463c4762a1bSJed Brown const PetscInt numCells = 2, numVertices = 1, numCorners = 3; 464a4a685f2SJacob Faibussowitsch const PetscInt cells[6] = {2, 4, 3, 2, 1, 4}; 465c4762a1bSJed Brown PetscReal coords[2] = {1.0, 0.0}; 466c4762a1bSJed Brown PetscInt markerPoints[10] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1}; 467c4762a1bSJed Brown 4689566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 4699566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 4709371c9d4SSatish Balay } break; 471d71ae5a4SJacob Faibussowitsch default: 472d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank); 473c4762a1bSJed Brown } 474c4762a1bSJed Brown break; 475c4762a1bSJed Brown case 2: 47663a3b9bcSJacob Faibussowitsch PetscCheck(size == 3, comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 3 processes", testNum); 477c4762a1bSJed Brown switch (rank) { 4789371c9d4SSatish Balay case 0: { 479c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 2, numCorners = 3; 480a4a685f2SJacob Faibussowitsch const PetscInt cells[3] = {1, 2, 0}; 481c4762a1bSJed Brown PetscReal coords[4] = {0.5, 0.5, 0.0, 1.0}; 482c4762a1bSJed Brown PetscInt markerPoints[6] = {1, 1, 2, 1, 3, 1}; 483c4762a1bSJed Brown 4849566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 4859566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 4869371c9d4SSatish Balay } break; 4879371c9d4SSatish Balay case 1: { 488c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 2, numCorners = 3; 489a4a685f2SJacob Faibussowitsch const PetscInt cells[3] = {1, 0, 3}; 490c4762a1bSJed Brown PetscReal coords[4] = {0.0, 0.0, 1.0, 1.0}; 491c4762a1bSJed Brown PetscInt markerPoints[6] = {1, 1, 2, 1, 3, 1}; 492c4762a1bSJed Brown 4939566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 4949566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 4959371c9d4SSatish Balay } break; 4969371c9d4SSatish Balay case 2: { 497c4762a1bSJed Brown const PetscInt numCells = 2, numVertices = 1, numCorners = 3; 498a4a685f2SJacob Faibussowitsch const PetscInt cells[6] = {0, 4, 3, 0, 2, 4}; 499c4762a1bSJed Brown PetscReal coords[2] = {1.0, 0.0}; 500c4762a1bSJed Brown PetscInt markerPoints[10] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1}; 501c4762a1bSJed Brown 5029566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 5039566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 5049371c9d4SSatish Balay } break; 505d71ae5a4SJacob Faibussowitsch default: 506d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank); 507c4762a1bSJed Brown } 508c4762a1bSJed Brown break; 509d71ae5a4SJacob Faibussowitsch default: 510d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 511c4762a1bSJed Brown } 5123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 513c4762a1bSJed Brown } 514c4762a1bSJed Brown 515d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSimplex_3D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm) 516d71ae5a4SJacob Faibussowitsch { 517c4762a1bSJed Brown PetscInt testNum = user->testNum, p; 518c4762a1bSJed Brown PetscMPIInt rank, size; 519c4762a1bSJed Brown 520c4762a1bSJed Brown PetscFunctionBegin; 5219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 5229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 523c4762a1bSJed Brown switch (testNum) { 524c4762a1bSJed Brown case 0: 52563a3b9bcSJacob Faibussowitsch PetscCheck(size == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 2 processes", testNum); 526c4762a1bSJed Brown switch (rank) { 5279371c9d4SSatish Balay case 0: { 528c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 2, numCorners = 4; 529a4a685f2SJacob Faibussowitsch const PetscInt cells[4] = {0, 2, 1, 3}; 530c4762a1bSJed Brown PetscReal coords[6] = {0.0, 0.0, -0.5, 0.0, -0.5, 0.0}; 531c4762a1bSJed Brown PetscInt markerPoints[8] = {1, 1, 2, 1, 3, 1, 4, 1}; 532c4762a1bSJed Brown 5339566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 5349566063dSJacob Faibussowitsch for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 5359371c9d4SSatish Balay } break; 5369371c9d4SSatish Balay case 1: { 537c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 3, numCorners = 4; 538a4a685f2SJacob Faibussowitsch const PetscInt cells[4] = {1, 2, 4, 3}; 539c4762a1bSJed Brown PetscReal coords[9] = {1.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.5}; 540c4762a1bSJed Brown PetscInt markerPoints[8] = {1, 1, 2, 1, 3, 1, 4, 1}; 541c4762a1bSJed Brown 5429566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 5439566063dSJacob Faibussowitsch for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 5449371c9d4SSatish Balay } break; 545d71ae5a4SJacob Faibussowitsch default: 546d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank); 547c4762a1bSJed Brown } 548c4762a1bSJed Brown break; 549d71ae5a4SJacob Faibussowitsch default: 550d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 551c4762a1bSJed Brown } 552c4762a1bSJed Brown if (user->testOrientIF) { 553c4762a1bSJed Brown PetscInt ifp[] = {8, 6}; 554c4762a1bSJed Brown 5559566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Mesh before orientation")); 5569566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-before_orientation_dm_view")); 557c4762a1bSJed Brown /* rotate interface face ifp[rank] by given orientation ornt[rank] */ 5589566063dSJacob Faibussowitsch PetscCall(DMPlexOrientPoint(*dm, ifp[rank], user->ornt[rank])); 5599566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-before_orientation_dm_view")); 5609566063dSJacob Faibussowitsch PetscCall(DMPlexCheckFaces(*dm, 0)); 5619566063dSJacob Faibussowitsch PetscCall(DMPlexOrientInterface_Internal(*dm)); 5629566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Orientation test PASSED\n")); 563c4762a1bSJed Brown } 5643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 565c4762a1bSJed Brown } 566c4762a1bSJed Brown 567d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateQuad_2D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm) 568d71ae5a4SJacob Faibussowitsch { 569c4762a1bSJed Brown PetscInt testNum = user->testNum, p; 570c4762a1bSJed Brown PetscMPIInt rank, size; 571c4762a1bSJed Brown 572c4762a1bSJed Brown PetscFunctionBegin; 5739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 5749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 575c4762a1bSJed Brown switch (testNum) { 576c4762a1bSJed Brown case 0: 57763a3b9bcSJacob Faibussowitsch PetscCheck(size == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 2 processes", testNum); 578c4762a1bSJed Brown switch (rank) { 5799371c9d4SSatish Balay case 0: { 580c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 3, numCorners = 4; 581a4a685f2SJacob Faibussowitsch const PetscInt cells[4] = {0, 1, 2, 3}; 582c4762a1bSJed Brown PetscReal coords[6] = {-0.5, 0.0, 0.0, 0.0, 0.0, 1.0}; 583c4762a1bSJed Brown PetscInt markerPoints[4 * 2] = {1, 1, 2, 1, 3, 1, 4, 1}; 584c4762a1bSJed Brown 5859566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 5869566063dSJacob Faibussowitsch for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 5879371c9d4SSatish Balay } break; 5889371c9d4SSatish Balay case 1: { 589c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 3, numCorners = 4; 590a4a685f2SJacob Faibussowitsch const PetscInt cells[4] = {1, 4, 5, 2}; 591c4762a1bSJed Brown PetscReal coords[6] = {-0.5, 1.0, 0.5, 0.0, 0.5, 1.0}; 592c4762a1bSJed Brown PetscInt markerPoints[4 * 2] = {1, 1, 2, 1, 3, 1, 4, 1}; 593c4762a1bSJed Brown 5949566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 5959566063dSJacob Faibussowitsch for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 5969371c9d4SSatish Balay } break; 597d71ae5a4SJacob Faibussowitsch default: 598d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank); 599c4762a1bSJed Brown } 600c4762a1bSJed Brown break; 601d71ae5a4SJacob Faibussowitsch default: 602d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 603c4762a1bSJed Brown } 6043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 605c4762a1bSJed Brown } 606c4762a1bSJed Brown 607d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateHex_3D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm) 608d71ae5a4SJacob Faibussowitsch { 609c4762a1bSJed Brown PetscInt testNum = user->testNum, p; 610c4762a1bSJed Brown PetscMPIInt rank, size; 611c4762a1bSJed Brown 612c4762a1bSJed Brown PetscFunctionBegin; 6139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 6149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 615c4762a1bSJed Brown switch (testNum) { 616c4762a1bSJed Brown case 0: 61763a3b9bcSJacob Faibussowitsch PetscCheck(size == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 2 processes", testNum); 618c4762a1bSJed Brown switch (rank) { 6199371c9d4SSatish Balay case 0: { 620c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 6, numCorners = 8; 621a4a685f2SJacob Faibussowitsch const PetscInt cells[8] = {0, 3, 2, 1, 4, 5, 6, 7}; 622c4762a1bSJed Brown PetscReal coords[6 * 3] = {-0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, -0.5, 1.0, 0.0, -0.5, 0.0, 1.0, 0.0, 0.0, 1.0}; 623c4762a1bSJed Brown PetscInt markerPoints[8 * 2] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1, 8, 1, 9, 1}; 624c4762a1bSJed Brown 6259566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 6269566063dSJacob Faibussowitsch for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 6279371c9d4SSatish Balay } break; 6289371c9d4SSatish Balay case 1: { 629c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 6, numCorners = 8; 630a4a685f2SJacob Faibussowitsch const PetscInt cells[8] = {1, 2, 9, 8, 5, 10, 11, 6}; 631c4762a1bSJed Brown PetscReal coords[6 * 3] = {0.0, 1.0, 1.0, -0.5, 1.0, 1.0, 0.5, 0.0, 0.0, 0.5, 1.0, 0.0, 0.5, 0.0, 1.0, 0.5, 1.0, 1.0}; 632c4762a1bSJed Brown PetscInt markerPoints[8 * 2] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1, 8, 1, 9, 1}; 633c4762a1bSJed Brown 6349566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 6359566063dSJacob Faibussowitsch for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 6369371c9d4SSatish Balay } break; 637d71ae5a4SJacob Faibussowitsch default: 638d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank); 639c4762a1bSJed Brown } 640c4762a1bSJed Brown break; 641d71ae5a4SJacob Faibussowitsch default: 642d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 643c4762a1bSJed Brown } 6443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 645c4762a1bSJed Brown } 646c4762a1bSJed Brown 647d71ae5a4SJacob Faibussowitsch static PetscErrorCode CustomView(DM dm, PetscViewer v) 648d71ae5a4SJacob Faibussowitsch { 649c4762a1bSJed Brown DMPlexInterpolatedFlag interpolated; 650c4762a1bSJed Brown PetscBool distributed; 651c4762a1bSJed Brown 652c4762a1bSJed Brown PetscFunctionBegin; 6539566063dSJacob Faibussowitsch PetscCall(DMPlexIsDistributed(dm, &distributed)); 6549566063dSJacob Faibussowitsch PetscCall(DMPlexIsInterpolatedCollective(dm, &interpolated)); 6559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, "DMPlexIsDistributed: %s\n", PetscBools[distributed])); 6569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, "DMPlexIsInterpolatedCollective: %s\n", DMPlexInterpolatedFlags[interpolated])); 6573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 658c4762a1bSJed Brown } 659c4762a1bSJed Brown 660d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMeshFromFile(MPI_Comm comm, AppCtx *user, DM *dm, DM *serialDM) 661d71ae5a4SJacob Faibussowitsch { 662c4762a1bSJed Brown const char *filename = user->filename; 663c4762a1bSJed Brown PetscBool testHeavy = user->testHeavy; 664c4762a1bSJed Brown PetscBool interpCreate = user->interpolate == CREATE ? PETSC_TRUE : PETSC_FALSE; 665c4762a1bSJed Brown PetscBool distributed = PETSC_FALSE; 666c4762a1bSJed Brown 667c4762a1bSJed Brown PetscFunctionBegin; 668c4762a1bSJed Brown *serialDM = NULL; 6699566063dSJacob Faibussowitsch if (testHeavy && interpCreate) PetscCall(DMPlexSetOrientInterface_Private(NULL, PETSC_FALSE)); 6709566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[0])); 6719566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromFile(comm, filename, "ex18_plex", interpCreate, dm)); /* with DMPlexOrientInterface_Internal() call skipped so that PointSF issues are left to DMPlexCheckPointSFHeavy() */ 6729566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 6739566063dSJacob Faibussowitsch if (testHeavy && interpCreate) PetscCall(DMPlexSetOrientInterface_Private(NULL, PETSC_TRUE)); 6749566063dSJacob Faibussowitsch PetscCall(DMPlexIsDistributed(*dm, &distributed)); 6759566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "DMPlexCreateFromFile produced %s mesh.\n", distributed ? "distributed" : "serial")); 676c4762a1bSJed Brown if (testHeavy && distributed) { 6779566063dSJacob Faibussowitsch PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_hdf5_force_sequential", NULL)); 6789566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromFile(comm, filename, "ex18_plex", interpCreate, serialDM)); 6799566063dSJacob Faibussowitsch PetscCall(DMPlexIsDistributed(*serialDM, &distributed)); 68028b400f6SJacob Faibussowitsch PetscCheck(!distributed, comm, PETSC_ERR_PLIB, "unable to create a serial DM from file"); 681c4762a1bSJed Brown } 6829566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &user->dim)); 6833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 684c4762a1bSJed Brown } 685c4762a1bSJed Brown 686d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 687d71ae5a4SJacob Faibussowitsch { 688c4762a1bSJed Brown PetscPartitioner part; 689c4762a1bSJed Brown PortableBoundary boundary = NULL; 690c4762a1bSJed Brown DM serialDM = NULL; 691c4762a1bSJed Brown PetscBool cellSimplex = user->cellSimplex; 692c4762a1bSJed Brown PetscBool useGenerator = user->useGenerator; 693c4762a1bSJed Brown PetscBool interpCreate = user->interpolate == CREATE ? PETSC_TRUE : PETSC_FALSE; 694c4762a1bSJed Brown PetscBool interpSerial = user->interpolate == AFTER_CREATE ? PETSC_TRUE : PETSC_FALSE; 695c4762a1bSJed Brown PetscBool interpParallel = user->interpolate == AFTER_DISTRIBUTE ? PETSC_TRUE : PETSC_FALSE; 696c4762a1bSJed Brown PetscBool testHeavy = user->testHeavy; 697c4762a1bSJed Brown PetscMPIInt rank; 698c4762a1bSJed Brown 699c4762a1bSJed Brown PetscFunctionBegin; 7009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 701c4762a1bSJed Brown if (user->filename[0]) { 7029566063dSJacob Faibussowitsch PetscCall(CreateMeshFromFile(comm, user, dm, &serialDM)); 703c4762a1bSJed Brown } else if (useGenerator) { 7049566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[0])); 7059566063dSJacob Faibussowitsch PetscCall(DMPlexCreateBoxMesh(comm, user->dim, cellSimplex, user->faces, NULL, NULL, NULL, interpCreate, dm)); 7069566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 707c4762a1bSJed Brown } else { 7089566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[0])); 709c4762a1bSJed Brown switch (user->dim) { 710d71ae5a4SJacob Faibussowitsch case 1: 711d71ae5a4SJacob Faibussowitsch PetscCall(CreateMesh_1D(comm, interpCreate, user, dm)); 712d71ae5a4SJacob Faibussowitsch break; 713c4762a1bSJed Brown case 2: 714c4762a1bSJed Brown if (cellSimplex) { 7159566063dSJacob Faibussowitsch PetscCall(CreateSimplex_2D(comm, interpCreate, user, dm)); 716c4762a1bSJed Brown } else { 7179566063dSJacob Faibussowitsch PetscCall(CreateQuad_2D(comm, interpCreate, user, dm)); 718c4762a1bSJed Brown } 719c4762a1bSJed Brown break; 720c4762a1bSJed Brown case 3: 721c4762a1bSJed Brown if (cellSimplex) { 7229566063dSJacob Faibussowitsch PetscCall(CreateSimplex_3D(comm, interpCreate, user, dm)); 723c4762a1bSJed Brown } else { 7249566063dSJacob Faibussowitsch PetscCall(CreateHex_3D(comm, interpCreate, user, dm)); 725c4762a1bSJed Brown } 726c4762a1bSJed Brown break; 727d71ae5a4SJacob Faibussowitsch default: 728d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %" PetscInt_FMT, user->dim); 729c4762a1bSJed Brown } 7309566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 731c4762a1bSJed Brown } 732da81f932SPierre Jolivet PetscCheck(user->ncoords % user->dim == 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "length of coordinates array %" PetscInt_FMT " must be divisible by spatial dimension %" PetscInt_FMT, user->ncoords, user->dim); 7339566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Original Mesh")); 7349566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view")); 735c4762a1bSJed Brown 736c4762a1bSJed Brown if (interpSerial) { 737c4762a1bSJed Brown DM idm; 738c4762a1bSJed Brown 7399566063dSJacob Faibussowitsch if (testHeavy) PetscCall(DMPlexSetOrientInterface_Private(*dm, PETSC_FALSE)); 7409566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[2])); 7419566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); /* with DMPlexOrientInterface_Internal() call skipped so that PointSF issues are left to DMPlexCheckPointSFHeavy() */ 7429566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 7439566063dSJacob Faibussowitsch if (testHeavy) PetscCall(DMPlexSetOrientInterface_Private(*dm, PETSC_TRUE)); 7449566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 745c4762a1bSJed Brown *dm = idm; 7469566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Interpolated Mesh")); 7479566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-intp_dm_view")); 748c4762a1bSJed Brown } 749c4762a1bSJed Brown 750c4762a1bSJed Brown /* Set partitioner options */ 7519566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(*dm, &part)); 752c4762a1bSJed Brown if (part) { 7539566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSIMPLE)); 7549566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part)); 755c4762a1bSJed Brown } 756c4762a1bSJed Brown 7579566063dSJacob Faibussowitsch if (user->customView) PetscCall(CustomView(*dm, PETSC_VIEWER_STDOUT_(comm))); 758c4762a1bSJed Brown if (testHeavy) { 759c4762a1bSJed Brown PetscBool distributed; 760c4762a1bSJed Brown 7619566063dSJacob Faibussowitsch PetscCall(DMPlexIsDistributed(*dm, &distributed)); 762c4762a1bSJed Brown if (!serialDM && !distributed) { 763c4762a1bSJed Brown serialDM = *dm; 7649566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)*dm)); 765c4762a1bSJed Brown } 76648a46eb9SPierre Jolivet if (serialDM) PetscCall(DMPlexGetExpandedBoundary_Private(serialDM, &boundary)); 767c4762a1bSJed Brown if (boundary) { 768c4762a1bSJed Brown /* check DM which has been created in parallel and already interpolated */ 7699566063dSJacob Faibussowitsch PetscCall(DMPlexCheckPointSFHeavy(*dm, boundary)); 770c4762a1bSJed Brown } 771c4762a1bSJed Brown /* Orient interface because it could be deliberately skipped above. It is idempotent. */ 7729566063dSJacob Faibussowitsch PetscCall(DMPlexOrientInterface_Internal(*dm)); 773c4762a1bSJed Brown } 774c4762a1bSJed Brown if (user->distribute) { 775c4762a1bSJed Brown DM pdm = NULL; 776c4762a1bSJed Brown 777c4762a1bSJed Brown /* Redistribute mesh over processes using that partitioner */ 7789566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[1])); 7799566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(*dm, 0, NULL, &pdm)); 7809566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 781c4762a1bSJed Brown if (pdm) { 7829566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 783c4762a1bSJed Brown *dm = pdm; 7849566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Redistributed Mesh")); 7859566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dist_dm_view")); 786c4762a1bSJed Brown } 787c4762a1bSJed Brown 788c4762a1bSJed Brown if (interpParallel) { 789c4762a1bSJed Brown DM idm; 790c4762a1bSJed Brown 7919566063dSJacob Faibussowitsch if (testHeavy) PetscCall(DMPlexSetOrientInterface_Private(*dm, PETSC_FALSE)); 7929566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[2])); 7939566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); /* with DMPlexOrientInterface_Internal() call skipped so that PointSF issues are left to DMPlexCheckPointSFHeavy() */ 7949566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 7959566063dSJacob Faibussowitsch if (testHeavy) PetscCall(DMPlexSetOrientInterface_Private(*dm, PETSC_TRUE)); 7969566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 797c4762a1bSJed Brown *dm = idm; 7989566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Interpolated Redistributed Mesh")); 7999566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-intp_dm_view")); 800c4762a1bSJed Brown } 801c4762a1bSJed Brown } 802c4762a1bSJed Brown if (testHeavy) { 8031baa6e33SBarry Smith if (boundary) PetscCall(DMPlexCheckPointSFHeavy(*dm, boundary)); 804c4762a1bSJed Brown /* Orient interface because it could be deliberately skipped above. It is idempotent. */ 8059566063dSJacob Faibussowitsch PetscCall(DMPlexOrientInterface_Internal(*dm)); 806c4762a1bSJed Brown } 807c4762a1bSJed Brown 8089566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Parallel Mesh")); 8099566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 8109566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 8119566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 812c4762a1bSJed Brown 8139566063dSJacob Faibussowitsch if (user->customView) PetscCall(CustomView(*dm, PETSC_VIEWER_STDOUT_(comm))); 8149566063dSJacob Faibussowitsch PetscCall(DMDestroy(&serialDM)); 8159566063dSJacob Faibussowitsch PetscCall(PortableBoundaryDestroy(&boundary)); 8163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 817c4762a1bSJed Brown } 818c4762a1bSJed Brown 819d3e1f4ccSVaclav Hapla #define ps2d(number) ((double)PetscRealPart(number)) 820d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode coord2str(char buf[], size_t len, PetscInt dim, const PetscScalar coords[], PetscReal tol) 821d71ae5a4SJacob Faibussowitsch { 822c4762a1bSJed Brown PetscFunctionBegin; 82308401ef6SPierre Jolivet PetscCheck(dim <= 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "dim must be less than or equal 3"); 824c4762a1bSJed Brown if (tol >= 1e-3) { 825c4762a1bSJed Brown switch (dim) { 826d71ae5a4SJacob Faibussowitsch case 1: 827d71ae5a4SJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, len, "(%12.3f)", ps2d(coords[0]))); 828d71ae5a4SJacob Faibussowitsch case 2: 829d71ae5a4SJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, len, "(%12.3f, %12.3f)", ps2d(coords[0]), ps2d(coords[1]))); 830d71ae5a4SJacob Faibussowitsch case 3: 831d71ae5a4SJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, len, "(%12.3f, %12.3f, %12.3f)", ps2d(coords[0]), ps2d(coords[1]), ps2d(coords[2]))); 832c4762a1bSJed Brown } 833c4762a1bSJed Brown } else { 834c4762a1bSJed Brown switch (dim) { 835d71ae5a4SJacob Faibussowitsch case 1: 836d71ae5a4SJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, len, "(%12.6f)", ps2d(coords[0]))); 837d71ae5a4SJacob Faibussowitsch case 2: 838d71ae5a4SJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, len, "(%12.6f, %12.6f)", ps2d(coords[0]), ps2d(coords[1]))); 839d71ae5a4SJacob Faibussowitsch case 3: 840d71ae5a4SJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, len, "(%12.6f, %12.6f, %12.6f)", ps2d(coords[0]), ps2d(coords[1]), ps2d(coords[2]))); 841c4762a1bSJed Brown } 842c4762a1bSJed Brown } 8433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 844c4762a1bSJed Brown } 845c4762a1bSJed Brown 846d71ae5a4SJacob Faibussowitsch static PetscErrorCode ViewVerticesFromCoords(DM dm, Vec coordsVec, PetscReal tol, PetscViewer viewer) 847d71ae5a4SJacob Faibussowitsch { 848d3e1f4ccSVaclav Hapla PetscInt dim, i, npoints; 849d3e1f4ccSVaclav Hapla IS pointsIS; 850d3e1f4ccSVaclav Hapla const PetscInt *points; 851d3e1f4ccSVaclav Hapla const PetscScalar *coords; 852c4762a1bSJed Brown char coordstr[128]; 853c4762a1bSJed Brown MPI_Comm comm; 854c4762a1bSJed Brown PetscMPIInt rank; 855c4762a1bSJed Brown 856c4762a1bSJed Brown PetscFunctionBegin; 8579566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 8589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 8599566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 8609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 8619566063dSJacob Faibussowitsch PetscCall(DMPlexFindVertices(dm, coordsVec, tol, &pointsIS)); 8629566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pointsIS, &points)); 8639566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pointsIS, &npoints)); 8649566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coordsVec, &coords)); 865c4762a1bSJed Brown for (i = 0; i < npoints; i++) { 8669566063dSJacob Faibussowitsch PetscCall(coord2str(coordstr, sizeof(coordstr), dim, &coords[i * dim], tol)); 8679566063dSJacob Faibussowitsch if (rank == 0 && i) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "-----\n")); 86863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %s --> points[%" PetscInt_FMT "] = %" PetscInt_FMT "\n", rank, coordstr, i, points[i])); 8699566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 870c4762a1bSJed Brown } 8719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 8729566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coordsVec, &coords)); 8739566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pointsIS, &points)); 8749566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointsIS)); 8753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 876c4762a1bSJed Brown } 877c4762a1bSJed Brown 878d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestExpandPoints(DM dm, AppCtx *user) 879d71ae5a4SJacob Faibussowitsch { 880c4762a1bSJed Brown IS is; 881c4762a1bSJed Brown PetscSection *sects; 882c4762a1bSJed Brown IS *iss; 883c4762a1bSJed Brown PetscInt d, depth; 884c4762a1bSJed Brown PetscMPIInt rank; 885c4762a1bSJed Brown PetscViewer viewer = PETSC_VIEWER_STDOUT_WORLD, sviewer; 886c4762a1bSJed Brown 887c4762a1bSJed Brown PetscFunctionBegin; 8889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 889dd400576SPatrick Sanan if (user->testExpandPointsEmpty && rank == 0) { 8909566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_USE_POINTER, &is)); 891c4762a1bSJed Brown } else { 8929566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, user->nPointsToExpand, user->pointsToExpand, PETSC_USE_POINTER, &is)); 893c4762a1bSJed Brown } 8949566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeRecursive(dm, is, &depth, &iss, §s)); 8959566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 8969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(sviewer, "[%d] ==========================\n", rank)); 897c4762a1bSJed Brown for (d = depth - 1; d >= 0; d--) { 898c4762a1bSJed Brown IS checkIS; 899c4762a1bSJed Brown PetscBool flg; 900c4762a1bSJed Brown 90163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(sviewer, "depth %" PetscInt_FMT " ---------------\n", d)); 9029566063dSJacob Faibussowitsch PetscCall(PetscSectionView(sects[d], sviewer)); 9039566063dSJacob Faibussowitsch PetscCall(ISView(iss[d], sviewer)); 904c4762a1bSJed Brown /* check reverse operation */ 905c4762a1bSJed Brown if (d < depth - 1) { 9069566063dSJacob Faibussowitsch PetscCall(DMPlexExpandedConesToFaces_Private(dm, iss[d], sects[d], &checkIS)); 9079566063dSJacob Faibussowitsch PetscCall(ISEqualUnsorted(checkIS, iss[d + 1], &flg)); 90828b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)checkIS), PETSC_ERR_PLIB, "DMPlexExpandedConesToFaces_Private produced wrong IS"); 9099566063dSJacob Faibussowitsch PetscCall(ISDestroy(&checkIS)); 910c4762a1bSJed Brown } 911c4762a1bSJed Brown } 9129566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 9139566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 9149566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreConeRecursive(dm, is, &depth, &iss, §s)); 9159566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 9163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 917c4762a1bSJed Brown } 918c4762a1bSJed Brown 919d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexExpandedConesToFaces_Private(DM dm, IS is, PetscSection section, IS *newis) 920d71ae5a4SJacob Faibussowitsch { 921c4762a1bSJed Brown PetscInt n, n1, ncone, numCoveredPoints, o, p, q, start, end; 922c4762a1bSJed Brown const PetscInt *coveredPoints; 923c4762a1bSJed Brown const PetscInt *arr, *cone; 924c4762a1bSJed Brown PetscInt *newarr; 925c4762a1bSJed Brown 926c4762a1bSJed Brown PetscFunctionBegin; 9279566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &n)); 9289566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(section, &n1)); 9299566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &start, &end)); 93063a3b9bcSJacob Faibussowitsch PetscCheck(n == n1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "IS size = %" PetscInt_FMT " != %" PetscInt_FMT " = section storage size", n, n1); 9319566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &arr)); 9329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(end - start, &newarr)); 933c4762a1bSJed Brown for (q = start; q < end; q++) { 9349566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, q, &ncone)); 9359566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, q, &o)); 936c4762a1bSJed Brown cone = &arr[o]; 937c4762a1bSJed Brown if (ncone == 1) { 938c4762a1bSJed Brown numCoveredPoints = 1; 939c4762a1bSJed Brown p = cone[0]; 940c4762a1bSJed Brown } else { 941c4762a1bSJed Brown PetscInt i; 942c4762a1bSJed Brown p = PETSC_MAX_INT; 9439371c9d4SSatish Balay for (i = 0; i < ncone; i++) 9449371c9d4SSatish Balay if (cone[i] < 0) { 9459371c9d4SSatish Balay p = -1; 9469371c9d4SSatish Balay break; 9479371c9d4SSatish Balay } 948c4762a1bSJed Brown if (p >= 0) { 9499566063dSJacob Faibussowitsch PetscCall(DMPlexGetJoin(dm, ncone, cone, &numCoveredPoints, &coveredPoints)); 95063a3b9bcSJacob Faibussowitsch PetscCheck(numCoveredPoints <= 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "more than one covered points for section point %" PetscInt_FMT, q); 951c4762a1bSJed Brown if (numCoveredPoints) p = coveredPoints[0]; 952c4762a1bSJed Brown else p = -2; 9539566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(dm, ncone, cone, &numCoveredPoints, &coveredPoints)); 954c4762a1bSJed Brown } 955c4762a1bSJed Brown } 956c4762a1bSJed Brown newarr[q - start] = p; 957c4762a1bSJed Brown } 9589566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &arr)); 9599566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, end - start, newarr, PETSC_OWN_POINTER, newis)); 9603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 961c4762a1bSJed Brown } 962c4762a1bSJed Brown 963d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexExpandedVerticesToFaces_Private(DM dm, IS boundary_expanded_is, PetscInt depth, PetscSection sections[], IS *boundary_is) 964d71ae5a4SJacob Faibussowitsch { 965c4762a1bSJed Brown PetscInt d; 966c4762a1bSJed Brown IS is, newis; 967c4762a1bSJed Brown 968c4762a1bSJed Brown PetscFunctionBegin; 969c4762a1bSJed Brown is = boundary_expanded_is; 9709566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)is)); 971c4762a1bSJed Brown for (d = 0; d < depth - 1; ++d) { 9729566063dSJacob Faibussowitsch PetscCall(DMPlexExpandedConesToFaces_Private(dm, is, sections[d], &newis)); 9739566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 974c4762a1bSJed Brown is = newis; 975c4762a1bSJed Brown } 976c4762a1bSJed Brown *boundary_is = is; 9773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 978c4762a1bSJed Brown } 979c4762a1bSJed Brown 9809371c9d4SSatish Balay #define CHKERRQI(incall, ierr) \ 981ad540459SPierre Jolivet if (ierr) incall = PETSC_FALSE; 982c4762a1bSJed Brown 983d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMLabelViewFromOptionsOnComm_Private(DMLabel label, const char optionname[], MPI_Comm comm) 984d71ae5a4SJacob Faibussowitsch { 985c4762a1bSJed Brown PetscViewer viewer; 986c4762a1bSJed Brown PetscBool flg; 987c4762a1bSJed Brown static PetscBool incall = PETSC_FALSE; 988c4762a1bSJed Brown PetscViewerFormat format; 989c4762a1bSJed Brown 990c4762a1bSJed Brown PetscFunctionBegin; 9913ba16761SJacob Faibussowitsch if (incall) PetscFunctionReturn(PETSC_SUCCESS); 992c4762a1bSJed Brown incall = PETSC_TRUE; 9935f80ce2aSJacob Faibussowitsch CHKERRQI(incall, PetscOptionsGetViewer(comm, ((PetscObject)label)->options, ((PetscObject)label)->prefix, optionname, &viewer, &format, &flg)); 994c4762a1bSJed Brown if (flg) { 9955f80ce2aSJacob Faibussowitsch CHKERRQI(incall, PetscViewerPushFormat(viewer, format)); 9965f80ce2aSJacob Faibussowitsch CHKERRQI(incall, DMLabelView(label, viewer)); 9975f80ce2aSJacob Faibussowitsch CHKERRQI(incall, PetscViewerPopFormat(viewer)); 9985f80ce2aSJacob Faibussowitsch CHKERRQI(incall, PetscViewerDestroy(&viewer)); 999c4762a1bSJed Brown } 1000c4762a1bSJed Brown incall = PETSC_FALSE; 10013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1002c4762a1bSJed Brown } 1003c4762a1bSJed Brown 1004c4762a1bSJed Brown /* TODO: this is hotfixing DMLabelGetStratumIS() - it should be fixed systematically instead */ 1005d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode DMLabelGetStratumISOnComm_Private(DMLabel label, PetscInt value, MPI_Comm comm, IS *is) 1006d71ae5a4SJacob Faibussowitsch { 1007c4762a1bSJed Brown IS tmpis; 1008c4762a1bSJed Brown 1009c4762a1bSJed Brown PetscFunctionBegin; 10109566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, value, &tmpis)); 10119566063dSJacob Faibussowitsch if (!tmpis) PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_USE_POINTER, &tmpis)); 10129566063dSJacob Faibussowitsch PetscCall(ISOnComm(tmpis, comm, PETSC_COPY_VALUES, is)); 10139566063dSJacob Faibussowitsch PetscCall(ISDestroy(&tmpis)); 10143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1015c4762a1bSJed Brown } 1016c4762a1bSJed Brown 1017c4762a1bSJed Brown /* currently only for simple PetscSection without fields or constraints */ 1018d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSectionReplicate_Private(MPI_Comm comm, PetscMPIInt rootrank, PetscSection sec0, PetscSection *secout) 1019d71ae5a4SJacob Faibussowitsch { 1020c4762a1bSJed Brown PetscSection sec; 1021c4762a1bSJed Brown PetscInt chart[2], p; 1022c4762a1bSJed Brown PetscInt *dofarr; 1023c4762a1bSJed Brown PetscMPIInt rank; 1024c4762a1bSJed Brown 1025c4762a1bSJed Brown PetscFunctionBegin; 10269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 102748a46eb9SPierre Jolivet if (rank == rootrank) PetscCall(PetscSectionGetChart(sec0, &chart[0], &chart[1])); 10289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(chart, 2, MPIU_INT, rootrank, comm)); 10299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(chart[1] - chart[0], &dofarr)); 1030c4762a1bSJed Brown if (rank == rootrank) { 103148a46eb9SPierre Jolivet for (p = chart[0]; p < chart[1]; p++) PetscCall(PetscSectionGetDof(sec0, p, &dofarr[p - chart[0]])); 1032c4762a1bSJed Brown } 10339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(dofarr, chart[1] - chart[0], MPIU_INT, rootrank, comm)); 10349566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, &sec)); 10359566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(sec, chart[0], chart[1])); 103648a46eb9SPierre Jolivet for (p = chart[0]; p < chart[1]; p++) PetscCall(PetscSectionSetDof(sec, p, dofarr[p - chart[0]])); 10379566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(sec)); 10389566063dSJacob Faibussowitsch PetscCall(PetscFree(dofarr)); 1039c4762a1bSJed Brown *secout = sec; 10403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1041c4762a1bSJed Brown } 1042c4762a1bSJed Brown 1043d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexExpandedVerticesCoordinatesToFaces_Private(DM ipdm, PortableBoundary bnd, IS *face_is) 1044d71ae5a4SJacob Faibussowitsch { 1045c4762a1bSJed Brown IS faces_expanded_is; 1046c4762a1bSJed Brown 1047c4762a1bSJed Brown PetscFunctionBegin; 10489566063dSJacob Faibussowitsch PetscCall(DMPlexFindVertices(ipdm, bnd->coordinates, 0.0, &faces_expanded_is)); 10499566063dSJacob Faibussowitsch PetscCall(DMPlexExpandedVerticesToFaces_Private(ipdm, faces_expanded_is, bnd->depth, bnd->sections, face_is)); 10509566063dSJacob Faibussowitsch PetscCall(ISDestroy(&faces_expanded_is)); 10513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1052c4762a1bSJed Brown } 1053c4762a1bSJed Brown 1054c4762a1bSJed Brown /* hack disabling DMPlexOrientInterface() call in DMPlexInterpolate() via -dm_plex_interpolate_orient_interfaces option */ 1055d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexSetOrientInterface_Private(DM dm, PetscBool enable) 1056d71ae5a4SJacob Faibussowitsch { 1057c4762a1bSJed Brown PetscOptions options = NULL; 1058c4762a1bSJed Brown const char *prefix = NULL; 1059c4762a1bSJed Brown const char opt[] = "-dm_plex_interpolate_orient_interfaces"; 1060c4762a1bSJed Brown char prefix_opt[512]; 1061c4762a1bSJed Brown PetscBool flg, set; 1062c4762a1bSJed Brown static PetscBool wasSetTrue = PETSC_FALSE; 1063c4762a1bSJed Brown 1064c4762a1bSJed Brown PetscFunctionBegin; 1065c4762a1bSJed Brown if (dm) { 10669566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 1067c4762a1bSJed Brown options = ((PetscObject)dm)->options; 1068c4762a1bSJed Brown } 1069c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(prefix_opt, "-", sizeof(prefix_opt))); 10709566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(prefix_opt, prefix, sizeof(prefix_opt))); 10719566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(prefix_opt, &opt[1], sizeof(prefix_opt))); 10729566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(options, prefix, opt, &flg, &set)); 1073c4762a1bSJed Brown if (!enable) { 1074c4762a1bSJed Brown if (set && flg) wasSetTrue = PETSC_TRUE; 10759566063dSJacob Faibussowitsch PetscCall(PetscOptionsSetValue(options, prefix_opt, "0")); 1076c4762a1bSJed Brown } else if (set && !flg) { 1077c4762a1bSJed Brown if (wasSetTrue) { 10789566063dSJacob Faibussowitsch PetscCall(PetscOptionsSetValue(options, prefix_opt, "1")); 1079c4762a1bSJed Brown } else { 1080c4762a1bSJed Brown /* default is PETSC_TRUE */ 10819566063dSJacob Faibussowitsch PetscCall(PetscOptionsClearValue(options, prefix_opt)); 1082c4762a1bSJed Brown } 1083c4762a1bSJed Brown wasSetTrue = PETSC_FALSE; 1084c4762a1bSJed Brown } 108576bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 10869566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(options, prefix, opt, &flg, &set)); 10871dca8a05SBarry Smith PetscCheck(!set || flg == enable, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "PetscOptionsSetValue did not have the desired effect"); 1088c4762a1bSJed Brown } 10893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1090c4762a1bSJed Brown } 1091c4762a1bSJed Brown 1092c4762a1bSJed Brown /* get coordinate description of the whole-domain boundary */ 1093d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexGetExpandedBoundary_Private(DM dm, PortableBoundary *boundary) 1094d71ae5a4SJacob Faibussowitsch { 1095c4762a1bSJed Brown PortableBoundary bnd0, bnd; 1096c4762a1bSJed Brown MPI_Comm comm; 1097c4762a1bSJed Brown DM idm; 1098c4762a1bSJed Brown DMLabel label; 1099c4762a1bSJed Brown PetscInt d; 1100c4762a1bSJed Brown const char boundaryName[] = "DMPlexDistributeInterpolateMarkInterface_boundary"; 1101c4762a1bSJed Brown IS boundary_is; 1102c4762a1bSJed Brown IS *boundary_expanded_iss; 1103c4762a1bSJed Brown PetscMPIInt rootrank = 0; 1104c4762a1bSJed Brown PetscMPIInt rank, size; 1105c4762a1bSJed Brown PetscInt value = 1; 1106c4762a1bSJed Brown DMPlexInterpolatedFlag intp; 1107c4762a1bSJed Brown PetscBool flg; 1108c4762a1bSJed Brown 1109c4762a1bSJed Brown PetscFunctionBegin; 11109566063dSJacob Faibussowitsch PetscCall(PetscNew(&bnd)); 11119566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 11129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 11139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 11149566063dSJacob Faibussowitsch PetscCall(DMPlexIsDistributed(dm, &flg)); 111528b400f6SJacob Faibussowitsch PetscCheck(!flg, comm, PETSC_ERR_ARG_WRONG, "serial DM (all points on one rank) needed"); 1116c4762a1bSJed Brown 1117c4762a1bSJed Brown /* interpolate serial DM if not yet interpolated */ 11189566063dSJacob Faibussowitsch PetscCall(DMPlexIsInterpolatedCollective(dm, &intp)); 1119c4762a1bSJed Brown if (intp == DMPLEX_INTERPOLATED_FULL) { 1120c4762a1bSJed Brown idm = dm; 11219566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 1122c4762a1bSJed Brown } else { 11239566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(dm, &idm)); 11249566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(idm, NULL, "-idm_view")); 1125c4762a1bSJed Brown } 1126c4762a1bSJed Brown 1127c4762a1bSJed Brown /* mark whole-domain boundary of the serial DM */ 11289566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, boundaryName, &label)); 11299566063dSJacob Faibussowitsch PetscCall(DMAddLabel(idm, label)); 11309566063dSJacob Faibussowitsch PetscCall(DMPlexMarkBoundaryFaces(idm, value, label)); 11319566063dSJacob Faibussowitsch PetscCall(DMLabelViewFromOptionsOnComm_Private(label, "-idm_boundary_view", comm)); 11329566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, value, &boundary_is)); 1133c4762a1bSJed Brown 1134c4762a1bSJed Brown /* translate to coordinates */ 11359566063dSJacob Faibussowitsch PetscCall(PetscNew(&bnd0)); 11369566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalSetUp(idm)); 1137c4762a1bSJed Brown if (rank == rootrank) { 11389566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeRecursive(idm, boundary_is, &bnd0->depth, &boundary_expanded_iss, &bnd0->sections)); 11399566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalTuple(dm, boundary_expanded_iss[0], NULL, &bnd0->coordinates)); 1140c4762a1bSJed Brown /* self-check */ 1141c4762a1bSJed Brown { 1142c4762a1bSJed Brown IS is0; 11439566063dSJacob Faibussowitsch PetscCall(DMPlexExpandedVerticesCoordinatesToFaces_Private(idm, bnd0, &is0)); 11449566063dSJacob Faibussowitsch PetscCall(ISEqual(is0, boundary_is, &flg)); 114528b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "DMPlexExpandedVerticesCoordinatesToFaces_Private produced a wrong IS"); 11469566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is0)); 1147c4762a1bSJed Brown } 1148c4762a1bSJed Brown } else { 1149*77433607SBarry Smith PetscCall(VecCreateFromOptions(PETSC_COMM_SELF, NULL, 1, 0, 0, &bnd0->coordinates)); 1150c4762a1bSJed Brown } 1151c4762a1bSJed Brown 1152c4762a1bSJed Brown { 1153c4762a1bSJed Brown Vec tmp; 1154c4762a1bSJed Brown VecScatter sc; 1155c4762a1bSJed Brown IS xis; 1156c4762a1bSJed Brown PetscInt n; 1157c4762a1bSJed Brown 1158c4762a1bSJed Brown /* just convert seq vectors to mpi vector */ 11599566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(bnd0->coordinates, &n)); 11609566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&n, 1, MPIU_INT, rootrank, comm)); 1161c4762a1bSJed Brown if (rank == rootrank) { 1162*77433607SBarry Smith PetscCall(VecCreateFromOptions(comm, NULL, 1, n, n, &tmp)); 1163c4762a1bSJed Brown } else { 1164*77433607SBarry Smith PetscCall(VecCreateFromOptions(comm, NULL, 1, 0, n, &tmp)); 1165c4762a1bSJed Brown } 11669566063dSJacob Faibussowitsch PetscCall(VecCopy(bnd0->coordinates, tmp)); 11679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnd0->coordinates)); 1168c4762a1bSJed Brown bnd0->coordinates = tmp; 1169c4762a1bSJed Brown 1170c4762a1bSJed Brown /* replicate coordinates from root rank to all ranks */ 1171*77433607SBarry Smith PetscCall(VecCreateFromOptions(comm, NULL, 1, n, n * size, &bnd->coordinates)); 11729566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, n, 0, 1, &xis)); 11739566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(bnd0->coordinates, xis, bnd->coordinates, NULL, &sc)); 11749566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sc, bnd0->coordinates, bnd->coordinates, INSERT_VALUES, SCATTER_FORWARD)); 11759566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sc, bnd0->coordinates, bnd->coordinates, INSERT_VALUES, SCATTER_FORWARD)); 11769566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&sc)); 11779566063dSJacob Faibussowitsch PetscCall(ISDestroy(&xis)); 1178c4762a1bSJed Brown } 1179c4762a1bSJed Brown bnd->depth = bnd0->depth; 11809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&bnd->depth, 1, MPIU_INT, rootrank, comm)); 11819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bnd->depth, &bnd->sections)); 118248a46eb9SPierre Jolivet for (d = 0; d < bnd->depth; d++) PetscCall(PetscSectionReplicate_Private(comm, rootrank, (rank == rootrank) ? bnd0->sections[d] : NULL, &bnd->sections[d])); 1183c4762a1bSJed Brown 118448a46eb9SPierre Jolivet if (rank == rootrank) PetscCall(DMPlexRestoreConeRecursive(idm, boundary_is, &bnd0->depth, &boundary_expanded_iss, &bnd0->sections)); 11859566063dSJacob Faibussowitsch PetscCall(PortableBoundaryDestroy(&bnd0)); 11869566063dSJacob Faibussowitsch PetscCall(DMRemoveLabelBySelf(idm, &label, PETSC_TRUE)); 11879566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&label)); 11889566063dSJacob Faibussowitsch PetscCall(ISDestroy(&boundary_is)); 11899566063dSJacob Faibussowitsch PetscCall(DMDestroy(&idm)); 1190c4762a1bSJed Brown *boundary = bnd; 11913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1192c4762a1bSJed Brown } 1193c4762a1bSJed Brown 1194c4762a1bSJed Brown /* get faces of inter-partition interface */ 1195d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexGetInterfaceFaces_Private(DM ipdm, IS boundary_faces_is, IS *interface_faces_is) 1196d71ae5a4SJacob Faibussowitsch { 1197c4762a1bSJed Brown MPI_Comm comm; 1198c4762a1bSJed Brown DMLabel label; 1199c4762a1bSJed Brown IS part_boundary_faces_is; 1200c4762a1bSJed Brown const char partBoundaryName[] = "DMPlexDistributeInterpolateMarkInterface_partBoundary"; 1201c4762a1bSJed Brown PetscInt value = 1; 1202c4762a1bSJed Brown DMPlexInterpolatedFlag intp; 1203c4762a1bSJed Brown 1204c4762a1bSJed Brown PetscFunctionBegin; 12059566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ipdm, &comm)); 12069566063dSJacob Faibussowitsch PetscCall(DMPlexIsInterpolatedCollective(ipdm, &intp)); 120708401ef6SPierre Jolivet PetscCheck(intp == DMPLEX_INTERPOLATED_FULL, comm, PETSC_ERR_ARG_WRONG, "only for fully interpolated DMPlex"); 1208c4762a1bSJed Brown 1209c4762a1bSJed Brown /* get ipdm partition boundary (partBoundary) */ 1210429fa399SMatthew G. Knepley { 1211429fa399SMatthew G. Knepley PetscSF sf; 1212429fa399SMatthew G. Knepley 12139566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, partBoundaryName, &label)); 12149566063dSJacob Faibussowitsch PetscCall(DMAddLabel(ipdm, label)); 1215429fa399SMatthew G. Knepley PetscCall(DMGetPointSF(ipdm, &sf)); 1216429fa399SMatthew G. Knepley PetscCall(DMSetPointSF(ipdm, NULL)); 12179566063dSJacob Faibussowitsch PetscCall(DMPlexMarkBoundaryFaces(ipdm, value, label)); 1218429fa399SMatthew G. Knepley PetscCall(DMSetPointSF(ipdm, sf)); 12199566063dSJacob Faibussowitsch PetscCall(DMLabelViewFromOptionsOnComm_Private(label, "-ipdm_part_boundary_view", comm)); 12209566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumISOnComm_Private(label, value, comm, &part_boundary_faces_is)); 12219566063dSJacob Faibussowitsch PetscCall(DMRemoveLabelBySelf(ipdm, &label, PETSC_TRUE)); 12229566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&label)); 1223429fa399SMatthew G. Knepley } 1224c4762a1bSJed Brown 1225c4762a1bSJed Brown /* remove ipdm whole-domain boundary (boundary_faces_is) from ipdm partition boundary (part_boundary_faces_is), resulting just in inter-partition interface */ 12269566063dSJacob Faibussowitsch PetscCall(ISDifference(part_boundary_faces_is, boundary_faces_is, interface_faces_is)); 12279566063dSJacob Faibussowitsch PetscCall(ISDestroy(&part_boundary_faces_is)); 12283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1229c4762a1bSJed Brown } 1230c4762a1bSJed Brown 1231c4762a1bSJed Brown /* compute inter-partition interface including edges and vertices */ 1232d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexComputeCompleteInterface_Private(DM ipdm, IS interface_faces_is, IS *interface_is) 1233d71ae5a4SJacob Faibussowitsch { 1234c4762a1bSJed Brown DMLabel label; 1235c4762a1bSJed Brown PetscInt value = 1; 1236c4762a1bSJed Brown const char interfaceName[] = "DMPlexDistributeInterpolateMarkInterface_interface"; 1237c4762a1bSJed Brown DMPlexInterpolatedFlag intp; 1238c4762a1bSJed Brown MPI_Comm comm; 1239c4762a1bSJed Brown 1240c4762a1bSJed Brown PetscFunctionBegin; 12419566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ipdm, &comm)); 12429566063dSJacob Faibussowitsch PetscCall(DMPlexIsInterpolatedCollective(ipdm, &intp)); 124308401ef6SPierre Jolivet PetscCheck(intp == DMPLEX_INTERPOLATED_FULL, comm, PETSC_ERR_ARG_WRONG, "only for fully interpolated DMPlex"); 1244c4762a1bSJed Brown 12459566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, interfaceName, &label)); 12469566063dSJacob Faibussowitsch PetscCall(DMAddLabel(ipdm, label)); 12479566063dSJacob Faibussowitsch PetscCall(DMLabelSetStratumIS(label, value, interface_faces_is)); 12489566063dSJacob Faibussowitsch PetscCall(DMLabelViewFromOptionsOnComm_Private(label, "-interface_faces_view", comm)); 12499566063dSJacob Faibussowitsch PetscCall(DMPlexLabelComplete(ipdm, label)); 12509566063dSJacob Faibussowitsch PetscCall(DMLabelViewFromOptionsOnComm_Private(label, "-interface_view", comm)); 12519566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumISOnComm_Private(label, value, comm, interface_is)); 12529566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*interface_is, "interface_is")); 12539566063dSJacob Faibussowitsch PetscCall(ISViewFromOptions(*interface_is, NULL, "-interface_is_view")); 12549566063dSJacob Faibussowitsch PetscCall(DMRemoveLabelBySelf(ipdm, &label, PETSC_TRUE)); 12559566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&label)); 12563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1257c4762a1bSJed Brown } 1258c4762a1bSJed Brown 1259d71ae5a4SJacob Faibussowitsch static PetscErrorCode PointSFGetOutwardInterfacePoints(PetscSF sf, IS *is) 1260d71ae5a4SJacob Faibussowitsch { 1261c4762a1bSJed Brown PetscInt n; 1262c4762a1bSJed Brown const PetscInt *arr; 1263c4762a1bSJed Brown 1264c4762a1bSJed Brown PetscFunctionBegin; 12659566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, NULL, &n, &arr, NULL)); 12669566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)sf), n, arr, PETSC_USE_POINTER, is)); 12673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1268c4762a1bSJed Brown } 1269c4762a1bSJed Brown 1270d71ae5a4SJacob Faibussowitsch static PetscErrorCode PointSFGetInwardInterfacePoints(PetscSF sf, IS *is) 1271d71ae5a4SJacob Faibussowitsch { 1272c4762a1bSJed Brown PetscInt n; 1273c4762a1bSJed Brown const PetscInt *rootdegree; 1274c4762a1bSJed Brown PetscInt *arr; 1275c4762a1bSJed Brown 1276c4762a1bSJed Brown PetscFunctionBegin; 12779566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sf)); 12789566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sf, &rootdegree)); 12799566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sf, &rootdegree)); 12809566063dSJacob Faibussowitsch PetscCall(PetscSFComputeMultiRootOriginalNumbering(sf, rootdegree, &n, &arr)); 12819566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)sf), n, arr, PETSC_OWN_POINTER, is)); 12823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1283c4762a1bSJed Brown } 1284c4762a1bSJed Brown 1285d71ae5a4SJacob Faibussowitsch static PetscErrorCode PointSFGetInterfacePoints_Private(PetscSF pointSF, IS *is) 1286d71ae5a4SJacob Faibussowitsch { 1287c4762a1bSJed Brown IS pointSF_out_is, pointSF_in_is; 1288c4762a1bSJed Brown 1289c4762a1bSJed Brown PetscFunctionBegin; 12909566063dSJacob Faibussowitsch PetscCall(PointSFGetOutwardInterfacePoints(pointSF, &pointSF_out_is)); 12919566063dSJacob Faibussowitsch PetscCall(PointSFGetInwardInterfacePoints(pointSF, &pointSF_in_is)); 12929566063dSJacob Faibussowitsch PetscCall(ISExpand(pointSF_out_is, pointSF_in_is, is)); 12939566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointSF_out_is)); 12949566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointSF_in_is)); 12953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1296c4762a1bSJed Brown } 1297c4762a1bSJed Brown 12985f80ce2aSJacob Faibussowitsch #define CHKERRMY(ierr) PetscCheck(!ierr, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PointSF is wrong. Unable to show details!") 1299c4762a1bSJed Brown 1300d71ae5a4SJacob Faibussowitsch static PetscErrorCode ViewPointsWithType_Internal(DM dm, IS pointsIS, PetscViewer v) 1301d71ae5a4SJacob Faibussowitsch { 1302c4762a1bSJed Brown DMLabel label; 1303c4762a1bSJed Brown PetscSection coordsSection; 1304c4762a1bSJed Brown Vec coordsVec; 1305c4762a1bSJed Brown PetscScalar *coordsScalar; 1306c4762a1bSJed Brown PetscInt coneSize, depth, dim, i, p, npoints; 1307c4762a1bSJed Brown const PetscInt *points; 1308c4762a1bSJed Brown 1309c4762a1bSJed Brown PetscFunctionBegin; 13109566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 13119566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordsSection)); 13129566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordsVec)); 13139566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordsVec, &coordsScalar)); 13149566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pointsIS, &npoints)); 13159566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pointsIS, &points)); 13169566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(dm, &label)); 13179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(v)); 1318c4762a1bSJed Brown for (i = 0; i < npoints; i++) { 1319c4762a1bSJed Brown p = points[i]; 13209566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, p, &depth)); 1321c4762a1bSJed Brown if (!depth) { 1322d3e1f4ccSVaclav Hapla PetscInt n, o; 1323c4762a1bSJed Brown char coordstr[128]; 1324c4762a1bSJed Brown 13259566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(coordsSection, p, &n)); 13269566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordsSection, p, &o)); 13279566063dSJacob Faibussowitsch PetscCall(coord2str(coordstr, sizeof(coordstr), n, &coordsScalar[o], 1.0)); 132863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(v, "vertex %" PetscInt_FMT " w/ coordinates %s\n", p, coordstr)); 1329c4762a1bSJed Brown } else { 1330c4762a1bSJed Brown char entityType[16]; 1331c4762a1bSJed Brown 1332c4762a1bSJed Brown switch (depth) { 1333d71ae5a4SJacob Faibussowitsch case 1: 1334c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(entityType, "edge", sizeof(entityType))); 1335d71ae5a4SJacob Faibussowitsch break; 1336d71ae5a4SJacob Faibussowitsch case 2: 1337c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(entityType, "face", sizeof(entityType))); 1338d71ae5a4SJacob Faibussowitsch break; 1339d71ae5a4SJacob Faibussowitsch case 3: 1340c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(entityType, "cell", sizeof(entityType))); 1341d71ae5a4SJacob Faibussowitsch break; 1342d71ae5a4SJacob Faibussowitsch default: 1343d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Only for depth <= 3"); 1344c4762a1bSJed Brown } 134548a46eb9SPierre Jolivet if (depth == dim && dim < 3) PetscCall(PetscStrlcat(entityType, " (cell)", sizeof(entityType))); 134663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(v, "%s %" PetscInt_FMT "\n", entityType, p)); 1347c4762a1bSJed Brown } 13489566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, p, &coneSize)); 1349c4762a1bSJed Brown if (coneSize) { 1350c4762a1bSJed Brown const PetscInt *cone; 1351c4762a1bSJed Brown IS coneIS; 1352c4762a1bSJed Brown 13539566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, p, &cone)); 13549566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, coneSize, cone, PETSC_USE_POINTER, &coneIS)); 13559566063dSJacob Faibussowitsch PetscCall(ViewPointsWithType_Internal(dm, coneIS, v)); 13569566063dSJacob Faibussowitsch PetscCall(ISDestroy(&coneIS)); 1357c4762a1bSJed Brown } 1358c4762a1bSJed Brown } 13599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(v)); 13609566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordsVec, &coordsScalar)); 13619566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pointsIS, &points)); 13623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1363c4762a1bSJed Brown } 1364c4762a1bSJed Brown 1365d71ae5a4SJacob Faibussowitsch static PetscErrorCode ViewPointsWithType(DM dm, IS points, PetscViewer v) 1366d71ae5a4SJacob Faibussowitsch { 1367c4762a1bSJed Brown PetscBool flg; 1368c4762a1bSJed Brown PetscInt npoints; 1369c4762a1bSJed Brown PetscMPIInt rank; 1370c4762a1bSJed Brown 1371c4762a1bSJed Brown PetscFunctionBegin; 13729566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &flg)); 137328b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Only for ASCII viewer"); 13749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)v), &rank)); 13759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(v)); 13769566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(points, &npoints)); 1377c4762a1bSJed Brown if (npoints) { 13789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d] --------\n", rank)); 13799566063dSJacob Faibussowitsch PetscCall(ViewPointsWithType_Internal(dm, points, v)); 1380c4762a1bSJed Brown } 13819566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(v)); 13829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(v)); 13833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1384c4762a1bSJed Brown } 1385c4762a1bSJed Brown 1386d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexComparePointSFWithInterface_Private(DM ipdm, IS interface_is) 1387d71ae5a4SJacob Faibussowitsch { 1388c4762a1bSJed Brown PetscSF pointsf; 1389c4762a1bSJed Brown IS pointsf_is; 1390c4762a1bSJed Brown PetscBool flg; 1391c4762a1bSJed Brown MPI_Comm comm; 1392c4762a1bSJed Brown PetscMPIInt size; 1393c4762a1bSJed Brown 1394c4762a1bSJed Brown PetscFunctionBegin; 13959566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ipdm, &comm)); 13969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 13979566063dSJacob Faibussowitsch PetscCall(DMGetPointSF(ipdm, &pointsf)); 1398c4762a1bSJed Brown if (pointsf) { 1399c4762a1bSJed Brown PetscInt nroots; 14009566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(pointsf, &nroots, NULL, NULL, NULL)); 1401c4762a1bSJed Brown if (nroots < 0) pointsf = NULL; /* uninitialized SF */ 1402c4762a1bSJed Brown } 1403c4762a1bSJed Brown if (!pointsf) { 1404c4762a1bSJed Brown PetscInt N = 0; 14059566063dSJacob Faibussowitsch if (interface_is) PetscCall(ISGetSize(interface_is, &N)); 140628b400f6SJacob Faibussowitsch PetscCheck(!N, comm, PETSC_ERR_PLIB, "interface_is should be NULL or empty for PointSF being NULL"); 14073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1408c4762a1bSJed Brown } 1409c4762a1bSJed Brown 1410c4762a1bSJed Brown /* get PointSF points as IS pointsf_is */ 14119566063dSJacob Faibussowitsch PetscCall(PointSFGetInterfacePoints_Private(pointsf, &pointsf_is)); 1412c4762a1bSJed Brown 1413c4762a1bSJed Brown /* compare pointsf_is with interface_is */ 14149566063dSJacob Faibussowitsch PetscCall(ISEqual(interface_is, pointsf_is, &flg)); 1415712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &flg, 1, MPIU_BOOL, MPI_LAND, comm)); 1416c4762a1bSJed Brown if (!flg) { 1417c4762a1bSJed Brown IS pointsf_extra_is, pointsf_missing_is; 1418c4762a1bSJed Brown PetscViewer errv = PETSC_VIEWER_STDERR_(comm); 14195f80ce2aSJacob Faibussowitsch CHKERRMY(ISDifference(interface_is, pointsf_is, &pointsf_missing_is)); 14205f80ce2aSJacob Faibussowitsch CHKERRMY(ISDifference(pointsf_is, interface_is, &pointsf_extra_is)); 14215f80ce2aSJacob Faibussowitsch CHKERRMY(PetscViewerASCIIPrintf(errv, "Points missing in PointSF:\n")); 14225f80ce2aSJacob Faibussowitsch CHKERRMY(ViewPointsWithType(ipdm, pointsf_missing_is, errv)); 14235f80ce2aSJacob Faibussowitsch CHKERRMY(PetscViewerASCIIPrintf(errv, "Extra points in PointSF:\n")); 14245f80ce2aSJacob Faibussowitsch CHKERRMY(ViewPointsWithType(ipdm, pointsf_extra_is, errv)); 14255f80ce2aSJacob Faibussowitsch CHKERRMY(ISDestroy(&pointsf_extra_is)); 14265f80ce2aSJacob Faibussowitsch CHKERRMY(ISDestroy(&pointsf_missing_is)); 1427c4762a1bSJed Brown SETERRQ(comm, PETSC_ERR_PLIB, "PointSF is wrong! See details above."); 1428c4762a1bSJed Brown } 14299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointsf_is)); 14303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1431c4762a1bSJed Brown } 1432c4762a1bSJed Brown 1433c4762a1bSJed Brown /* remove faces & edges from label, leave just vertices */ 1434d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexISFilterVertices_Private(DM dm, IS points) 1435d71ae5a4SJacob Faibussowitsch { 1436c4762a1bSJed Brown PetscInt vStart, vEnd; 1437c4762a1bSJed Brown MPI_Comm comm; 1438c4762a1bSJed Brown 1439c4762a1bSJed Brown PetscFunctionBegin; 14409566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 14419566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 14429566063dSJacob Faibussowitsch PetscCall(ISGeneralFilter(points, vStart, vEnd)); 14433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1444c4762a1bSJed Brown } 1445c4762a1bSJed Brown 1446c4762a1bSJed Brown /* 14472e58f0efSBarry Smith DMPlexCheckPointSFHeavy - Thoroughly test that the PointSF after parallel DMPlexInterpolate() includes exactly all interface points. 1448c4762a1bSJed Brown 1449c4762a1bSJed Brown Collective 1450c4762a1bSJed Brown 14512fe279fdSBarry Smith Input Parameter: 1452c4762a1bSJed Brown . dm - The DMPlex object 1453c4762a1bSJed Brown 1454c4762a1bSJed Brown Notes: 1455c4762a1bSJed Brown The input DMPlex must be serial (one partition has all points, the other partitions have no points). 1456c4762a1bSJed Brown This is a heavy test which involves DMPlexInterpolate() if the input DM is not interpolated yet, and depends on having a representation of the whole-domain boundary (PortableBoundary), which can be obtained only by DMPlexGetExpandedBoundary_Private() (which involves DMPlexInterpolate() of a sequential DM). 1457c4762a1bSJed Brown This is mainly intended for debugging/testing purposes. 1458c4762a1bSJed Brown 1459c4762a1bSJed Brown Algorithm: 1460c4762a1bSJed Brown 1. boundary faces of the serial version of the whole mesh are found using DMPlexMarkBoundaryFaces() 1461c4762a1bSJed Brown 2. boundary faces are translated into vertices using DMPlexGetConeRecursive() and these are translated into coordinates - this description (aka PortableBoundary) is completely independent of partitioning and point numbering 1462c4762a1bSJed Brown 3. the mesh is distributed or loaded in parallel 1463c4762a1bSJed Brown 4. boundary faces of the distributed mesh are reconstructed from PortableBoundary using DMPlexFindVertices() 1464c4762a1bSJed Brown 5. partition boundary faces of the parallel mesh are found using DMPlexMarkBoundaryFaces() 1465c4762a1bSJed Brown 6. partition interfaces are computed as set difference of partition boundary faces minus the reconstructed boundary 1466c4762a1bSJed Brown 7. check that interface covered by PointSF (union of inward and outward points) is equal to the partition interface for each rank, otherwise print the difference and throw an error 1467c4762a1bSJed Brown 1468c4762a1bSJed Brown Level: developer 1469c4762a1bSJed Brown 1470c4762a1bSJed Brown .seealso: DMGetPointSF(), DMPlexCheckSymmetry(), DMPlexCheckSkeleton(), DMPlexCheckFaces() 1471c4762a1bSJed Brown */ 1472d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCheckPointSFHeavy(DM dm, PortableBoundary bnd) 1473d71ae5a4SJacob Faibussowitsch { 1474c4762a1bSJed Brown DM ipdm = NULL; 1475c4762a1bSJed Brown IS boundary_faces_is, interface_faces_is, interface_is; 1476c4762a1bSJed Brown DMPlexInterpolatedFlag intp; 1477c4762a1bSJed Brown MPI_Comm comm; 1478c4762a1bSJed Brown 1479c4762a1bSJed Brown PetscFunctionBegin; 14809566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 1481c4762a1bSJed Brown 14829566063dSJacob Faibussowitsch PetscCall(DMPlexIsInterpolatedCollective(dm, &intp)); 1483c4762a1bSJed Brown if (intp == DMPLEX_INTERPOLATED_FULL) { 1484c4762a1bSJed Brown ipdm = dm; 1485c4762a1bSJed Brown } else { 1486c4762a1bSJed Brown /* create temporary interpolated DM if input DM is not interpolated */ 14879566063dSJacob Faibussowitsch PetscCall(DMPlexSetOrientInterface_Private(dm, PETSC_FALSE)); 14889566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(dm, &ipdm)); /* with DMPlexOrientInterface_Internal() call skipped so that PointSF issues are left to DMPlexComparePointSFWithInterface_Private() below */ 14899566063dSJacob Faibussowitsch PetscCall(DMPlexSetOrientInterface_Private(dm, PETSC_TRUE)); 1490c4762a1bSJed Brown } 14919566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(ipdm, NULL, "-ipdm_view")); 1492c4762a1bSJed Brown 1493c4762a1bSJed Brown /* recover ipdm whole-domain boundary faces from the expanded vertices coordinates */ 14949566063dSJacob Faibussowitsch PetscCall(DMPlexExpandedVerticesCoordinatesToFaces_Private(ipdm, bnd, &boundary_faces_is)); 1495c4762a1bSJed Brown /* get inter-partition interface faces (interface_faces_is)*/ 14969566063dSJacob Faibussowitsch PetscCall(DMPlexGetInterfaceFaces_Private(ipdm, boundary_faces_is, &interface_faces_is)); 1497c4762a1bSJed Brown /* compute inter-partition interface including edges and vertices (interface_is) */ 14989566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCompleteInterface_Private(ipdm, interface_faces_is, &interface_is)); 1499c4762a1bSJed Brown /* destroy immediate ISs */ 15009566063dSJacob Faibussowitsch PetscCall(ISDestroy(&boundary_faces_is)); 15019566063dSJacob Faibussowitsch PetscCall(ISDestroy(&interface_faces_is)); 1502c4762a1bSJed Brown 1503c4762a1bSJed Brown /* for uninterpolated case, keep just vertices in interface */ 1504c4762a1bSJed Brown if (!intp) { 15059566063dSJacob Faibussowitsch PetscCall(DMPlexISFilterVertices_Private(ipdm, interface_is)); 15069566063dSJacob Faibussowitsch PetscCall(DMDestroy(&ipdm)); 1507c4762a1bSJed Brown } 1508c4762a1bSJed Brown 1509c4762a1bSJed Brown /* compare PointSF with the boundary reconstructed from coordinates */ 15109566063dSJacob Faibussowitsch PetscCall(DMPlexComparePointSFWithInterface_Private(dm, interface_is)); 15119566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "DMPlexCheckPointSFHeavy PASSED\n")); 15129566063dSJacob Faibussowitsch PetscCall(ISDestroy(&interface_is)); 15133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1514c4762a1bSJed Brown } 1515c4762a1bSJed Brown 1516d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 1517d71ae5a4SJacob Faibussowitsch { 1518c4762a1bSJed Brown DM dm; 1519c4762a1bSJed Brown AppCtx user; 1520c4762a1bSJed Brown 1521327415f7SBarry Smith PetscFunctionBeginUser; 15229566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 15239566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("create", &stage[0])); 15249566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("distribute", &stage[1])); 15259566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("interpolate", &stage[2])); 15269566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 15279566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 152848a46eb9SPierre Jolivet if (user.nPointsToExpand) PetscCall(TestExpandPoints(dm, &user)); 1529c4762a1bSJed Brown if (user.ncoords) { 1530d3e1f4ccSVaclav Hapla Vec coords; 1531d3e1f4ccSVaclav Hapla 15329566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, user.ncoords, user.ncoords, user.coords, &coords)); 15339566063dSJacob Faibussowitsch PetscCall(ViewVerticesFromCoords(dm, coords, user.coordsTol, PETSC_VIEWER_STDOUT_WORLD)); 15349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coords)); 1535c4762a1bSJed Brown } 15369566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 15379566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 1538b122ec5aSJacob Faibussowitsch return 0; 1539c4762a1bSJed Brown } 1540c4762a1bSJed Brown 1541c4762a1bSJed Brown /*TEST 1542c4762a1bSJed Brown 1543c4762a1bSJed Brown testset: 1544c4762a1bSJed Brown nsize: 2 1545c4762a1bSJed Brown args: -dm_view ascii::ascii_info_detail 1546c4762a1bSJed Brown args: -dm_plex_check_all 1547c4762a1bSJed Brown test: 1548c4762a1bSJed Brown suffix: 1_tri_dist0 1549c4762a1bSJed Brown args: -distribute 0 -interpolate {{none create}separate output} 1550c4762a1bSJed Brown test: 1551c4762a1bSJed Brown suffix: 1_tri_dist1 1552c4762a1bSJed Brown args: -distribute 1 -interpolate {{none create after_distribute}separate output} 1553c4762a1bSJed Brown test: 1554c4762a1bSJed Brown suffix: 1_quad_dist0 1555c4762a1bSJed Brown args: -cell_simplex 0 -distribute 0 -interpolate {{none create}separate output} 1556c4762a1bSJed Brown test: 1557c4762a1bSJed Brown suffix: 1_quad_dist1 1558c4762a1bSJed Brown args: -cell_simplex 0 -distribute 1 -interpolate {{none create after_distribute}separate output} 1559c4762a1bSJed Brown test: 1560c4762a1bSJed Brown suffix: 1_1d_dist1 1561c4762a1bSJed Brown args: -dim 1 -distribute 1 1562c4762a1bSJed Brown 1563c4762a1bSJed Brown testset: 1564c4762a1bSJed Brown nsize: 3 1565c4762a1bSJed Brown args: -testnum 1 -interpolate create 1566c4762a1bSJed Brown args: -dm_plex_check_all 1567c4762a1bSJed Brown test: 1568c4762a1bSJed Brown suffix: 2 1569c4762a1bSJed Brown args: -dm_view ascii::ascii_info_detail 1570c4762a1bSJed Brown test: 1571c4762a1bSJed Brown suffix: 2a 1572c4762a1bSJed Brown args: -dm_plex_check_cones_conform_on_interfaces_verbose 1573c4762a1bSJed Brown test: 1574c4762a1bSJed Brown suffix: 2b 1575c4762a1bSJed Brown args: -test_expand_points 0,1,2,5,6 1576c4762a1bSJed Brown test: 1577c4762a1bSJed Brown suffix: 2c 1578c4762a1bSJed Brown args: -test_expand_points 0,1,2,5,6 -test_expand_points_empty 1579c4762a1bSJed Brown 1580c4762a1bSJed Brown testset: 1581c4762a1bSJed Brown # the same as 1% for 3D 1582c4762a1bSJed Brown nsize: 2 1583c4762a1bSJed Brown args: -dim 3 -dm_view ascii::ascii_info_detail 1584c4762a1bSJed Brown args: -dm_plex_check_all 1585c4762a1bSJed Brown test: 1586c4762a1bSJed Brown suffix: 4_tet_dist0 1587c4762a1bSJed Brown args: -distribute 0 -interpolate {{none create}separate output} 1588c4762a1bSJed Brown test: 1589c4762a1bSJed Brown suffix: 4_tet_dist1 1590c4762a1bSJed Brown args: -distribute 1 -interpolate {{none create after_distribute}separate output} 1591c4762a1bSJed Brown test: 1592c4762a1bSJed Brown suffix: 4_hex_dist0 1593c4762a1bSJed Brown args: -cell_simplex 0 -distribute 0 -interpolate {{none create}separate output} 1594c4762a1bSJed Brown test: 1595c4762a1bSJed Brown suffix: 4_hex_dist1 1596c4762a1bSJed Brown args: -cell_simplex 0 -distribute 1 -interpolate {{none create after_distribute}separate output} 1597c4762a1bSJed Brown 1598c4762a1bSJed Brown test: 1599c4762a1bSJed Brown # the same as 4_tet_dist0 but test different initial orientations 1600c4762a1bSJed Brown suffix: 4_tet_test_orient 1601c4762a1bSJed Brown nsize: 2 1602c4762a1bSJed Brown args: -dim 3 -distribute 0 1603c4762a1bSJed Brown args: -dm_plex_check_all 1604c4762a1bSJed Brown args: -rotate_interface_0 {{0 1 2 11 12 13}} 1605c4762a1bSJed Brown args: -rotate_interface_1 {{0 1 2 11 12 13}} 1606c4762a1bSJed Brown 1607c4762a1bSJed Brown testset: 1608c4762a1bSJed Brown requires: exodusii 1609c4762a1bSJed Brown args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/TwoQuads.exo 1610c4762a1bSJed Brown args: -dm_view ascii::ascii_info_detail 1611c4762a1bSJed Brown args: -dm_plex_check_all 1612c4762a1bSJed Brown args: -custom_view 1613c4762a1bSJed Brown test: 1614c4762a1bSJed Brown suffix: 5_seq 1615c4762a1bSJed Brown nsize: 1 1616c4762a1bSJed Brown args: -distribute 0 -interpolate {{none create}separate output} 1617c4762a1bSJed Brown test: 161830602db0SMatthew G. Knepley # Detail viewing in a non-distributed mesh is broken because the DMLabelView() is collective, but the label is not shared 1619c4762a1bSJed Brown suffix: 5_dist0 1620c4762a1bSJed Brown nsize: 2 162130602db0SMatthew G. Knepley args: -distribute 0 -interpolate {{none create}separate output} -dm_view 1622c4762a1bSJed Brown test: 1623c4762a1bSJed Brown suffix: 5_dist1 1624c4762a1bSJed Brown nsize: 2 1625c4762a1bSJed Brown args: -distribute 1 -interpolate {{none create after_distribute}separate output} 1626c4762a1bSJed Brown 1627c4762a1bSJed Brown testset: 1628c4762a1bSJed Brown nsize: {{1 2 4}} 1629c4762a1bSJed Brown args: -use_generator 1630c4762a1bSJed Brown args: -dm_plex_check_all 1631c4762a1bSJed Brown args: -distribute -interpolate none 1632c4762a1bSJed Brown test: 1633c4762a1bSJed Brown suffix: 6_tri 1634c4762a1bSJed Brown requires: triangle 1635c0517cd5SMatthew G. Knepley args: -faces {{2,2 1,3 7,4}} -cell_simplex 1 -dm_generator triangle 1636c4762a1bSJed Brown test: 1637c4762a1bSJed Brown suffix: 6_quad 1638c4762a1bSJed Brown args: -faces {{2,2 1,3 7,4}} -cell_simplex 0 1639c4762a1bSJed Brown test: 1640c4762a1bSJed Brown suffix: 6_tet 1641c4762a1bSJed Brown requires: ctetgen 1642c0517cd5SMatthew G. Knepley args: -faces {{2,2,2 1,3,5 3,4,7}} -cell_simplex 1 -dm_generator ctetgen 1643c4762a1bSJed Brown test: 1644c4762a1bSJed Brown suffix: 6_hex 1645c4762a1bSJed Brown args: -faces {{2,2,2 1,3,5 3,4,7}} -cell_simplex 0 1646c4762a1bSJed Brown testset: 1647c4762a1bSJed Brown nsize: {{1 2 4}} 1648c4762a1bSJed Brown args: -use_generator 1649c4762a1bSJed Brown args: -dm_plex_check_all 1650c4762a1bSJed Brown args: -distribute -interpolate create 1651c4762a1bSJed Brown test: 1652c4762a1bSJed Brown suffix: 6_int_tri 1653c4762a1bSJed Brown requires: triangle 1654c0517cd5SMatthew G. Knepley args: -faces {{2,2 1,3 7,4}} -cell_simplex 1 -dm_generator triangle 1655c4762a1bSJed Brown test: 1656c4762a1bSJed Brown suffix: 6_int_quad 1657c4762a1bSJed Brown args: -faces {{2,2 1,3 7,4}} -cell_simplex 0 1658c4762a1bSJed Brown test: 1659c4762a1bSJed Brown suffix: 6_int_tet 1660c4762a1bSJed Brown requires: ctetgen 1661c0517cd5SMatthew G. Knepley args: -faces {{2,2,2 1,3,5 3,4,7}} -cell_simplex 1 -dm_generator ctetgen 1662c4762a1bSJed Brown test: 1663c4762a1bSJed Brown suffix: 6_int_hex 1664c4762a1bSJed Brown args: -faces {{2,2,2 1,3,5 3,4,7}} -cell_simplex 0 1665c4762a1bSJed Brown testset: 1666c4762a1bSJed Brown nsize: {{2 4}} 1667c4762a1bSJed Brown args: -use_generator 1668c4762a1bSJed Brown args: -dm_plex_check_all 1669c4762a1bSJed Brown args: -distribute -interpolate after_distribute 1670c4762a1bSJed Brown test: 1671c4762a1bSJed Brown suffix: 6_parint_tri 1672c4762a1bSJed Brown requires: triangle 1673c0517cd5SMatthew G. Knepley args: -faces {{2,2 1,3 7,4}} -cell_simplex 1 -dm_generator triangle 1674c4762a1bSJed Brown test: 1675c4762a1bSJed Brown suffix: 6_parint_quad 1676c4762a1bSJed Brown args: -faces {{2,2 1,3 7,4}} -cell_simplex 0 1677c4762a1bSJed Brown test: 1678c4762a1bSJed Brown suffix: 6_parint_tet 1679c4762a1bSJed Brown requires: ctetgen 1680c0517cd5SMatthew G. Knepley args: -faces {{2,2,2 1,3,5 3,4,7}} -cell_simplex 1 -dm_generator ctetgen 1681c4762a1bSJed Brown test: 1682c4762a1bSJed Brown suffix: 6_parint_hex 1683c4762a1bSJed Brown args: -faces {{2,2,2 1,3,5 3,4,7}} -cell_simplex 0 1684c4762a1bSJed Brown 1685c4762a1bSJed Brown testset: # 7 EXODUS 1686c4762a1bSJed Brown requires: exodusii 1687c4762a1bSJed Brown args: -dm_plex_check_all 1688c4762a1bSJed Brown args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo 1689c4762a1bSJed Brown args: -distribute 1690c4762a1bSJed Brown test: # seq load, simple partitioner 1691c4762a1bSJed Brown suffix: 7_exo 1692c4762a1bSJed Brown nsize: {{1 2 4 5}} 1693c4762a1bSJed Brown args: -interpolate none 1694c4762a1bSJed Brown test: # seq load, seq interpolation, simple partitioner 1695c4762a1bSJed Brown suffix: 7_exo_int_simple 1696c4762a1bSJed Brown nsize: {{1 2 4 5}} 1697c4762a1bSJed Brown args: -interpolate create 1698c4762a1bSJed Brown test: # seq load, seq interpolation, metis partitioner 1699c4762a1bSJed Brown suffix: 7_exo_int_metis 1700c4762a1bSJed Brown requires: parmetis 1701c4762a1bSJed Brown nsize: {{2 4 5}} 1702c4762a1bSJed Brown args: -interpolate create 1703c4762a1bSJed Brown args: -petscpartitioner_type parmetis 1704c4762a1bSJed Brown test: # seq load, simple partitioner, par interpolation 1705c4762a1bSJed Brown suffix: 7_exo_simple_int 1706c4762a1bSJed Brown nsize: {{2 4 5}} 1707c4762a1bSJed Brown args: -interpolate after_distribute 1708c4762a1bSJed Brown test: # seq load, metis partitioner, par interpolation 1709c4762a1bSJed Brown suffix: 7_exo_metis_int 1710c4762a1bSJed Brown requires: parmetis 1711c4762a1bSJed Brown nsize: {{2 4 5}} 1712c4762a1bSJed Brown args: -interpolate after_distribute 1713c4762a1bSJed Brown args: -petscpartitioner_type parmetis 1714c4762a1bSJed Brown 1715c4762a1bSJed Brown testset: # 7 HDF5 SEQUANTIAL LOAD 1716c4762a1bSJed Brown requires: hdf5 !complex 1717c4762a1bSJed Brown args: -dm_plex_check_all 1718c4762a1bSJed Brown args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf 1719c4762a1bSJed Brown args: -dm_plex_hdf5_force_sequential 1720c4762a1bSJed Brown args: -distribute 1721c4762a1bSJed Brown test: # seq load, simple partitioner 1722c4762a1bSJed Brown suffix: 7_seq_hdf5_simple 1723c4762a1bSJed Brown nsize: {{1 2 4 5}} 1724c4762a1bSJed Brown args: -interpolate none 1725c4762a1bSJed Brown test: # seq load, seq interpolation, simple partitioner 1726c4762a1bSJed Brown suffix: 7_seq_hdf5_int_simple 1727c4762a1bSJed Brown nsize: {{1 2 4 5}} 1728c4762a1bSJed Brown args: -interpolate after_create 1729c4762a1bSJed Brown test: # seq load, seq interpolation, metis partitioner 1730c4762a1bSJed Brown nsize: {{2 4 5}} 1731c4762a1bSJed Brown suffix: 7_seq_hdf5_int_metis 1732c4762a1bSJed Brown requires: parmetis 1733c4762a1bSJed Brown args: -interpolate after_create 1734c4762a1bSJed Brown args: -petscpartitioner_type parmetis 1735c4762a1bSJed Brown test: # seq load, simple partitioner, par interpolation 1736c4762a1bSJed Brown suffix: 7_seq_hdf5_simple_int 1737c4762a1bSJed Brown nsize: {{2 4 5}} 1738c4762a1bSJed Brown args: -interpolate after_distribute 1739c4762a1bSJed Brown test: # seq load, metis partitioner, par interpolation 1740c4762a1bSJed Brown nsize: {{2 4 5}} 1741c4762a1bSJed Brown suffix: 7_seq_hdf5_metis_int 1742c4762a1bSJed Brown requires: parmetis 1743c4762a1bSJed Brown args: -interpolate after_distribute 1744c4762a1bSJed Brown args: -petscpartitioner_type parmetis 1745c4762a1bSJed Brown 1746c4762a1bSJed Brown testset: # 7 HDF5 PARALLEL LOAD 1747c4762a1bSJed Brown requires: hdf5 !complex 1748c4762a1bSJed Brown nsize: {{2 4 5}} 1749c4762a1bSJed Brown args: -dm_plex_check_all 1750c4762a1bSJed Brown args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf 1751c4762a1bSJed Brown test: # par load 1752c4762a1bSJed Brown suffix: 7_par_hdf5 1753c4762a1bSJed Brown args: -interpolate none 1754c4762a1bSJed Brown test: # par load, par interpolation 1755c4762a1bSJed Brown suffix: 7_par_hdf5_int 1756c4762a1bSJed Brown args: -interpolate after_create 1757c4762a1bSJed Brown test: # par load, parmetis repartitioner 1758c4762a1bSJed Brown TODO: Parallel partitioning of uninterpolated meshes not supported 1759c4762a1bSJed Brown suffix: 7_par_hdf5_parmetis 1760c4762a1bSJed Brown requires: parmetis 1761c4762a1bSJed Brown args: -distribute -petscpartitioner_type parmetis 1762c4762a1bSJed Brown args: -interpolate none 1763c4762a1bSJed Brown test: # par load, par interpolation, parmetis repartitioner 1764c4762a1bSJed Brown suffix: 7_par_hdf5_int_parmetis 1765c4762a1bSJed Brown requires: parmetis 1766c4762a1bSJed Brown args: -distribute -petscpartitioner_type parmetis 1767c4762a1bSJed Brown args: -interpolate after_create 1768c4762a1bSJed Brown test: # par load, parmetis partitioner, par interpolation 1769c4762a1bSJed Brown TODO: Parallel partitioning of uninterpolated meshes not supported 1770c4762a1bSJed Brown suffix: 7_par_hdf5_parmetis_int 1771c4762a1bSJed Brown requires: parmetis 1772c4762a1bSJed Brown args: -distribute -petscpartitioner_type parmetis 1773c4762a1bSJed Brown args: -interpolate after_distribute 1774c4762a1bSJed Brown 1775c4762a1bSJed Brown test: 1776c4762a1bSJed Brown suffix: 7_hdf5_hierarch 1777c4762a1bSJed Brown requires: hdf5 ptscotch !complex 1778c4762a1bSJed Brown nsize: {{2 3 4}separate output} 1779c4762a1bSJed Brown args: -distribute 1780c4762a1bSJed Brown args: -interpolate after_create 1781c4762a1bSJed Brown args: -petscpartitioner_type matpartitioning -petscpartitioner_view ::ascii_info 1782c4762a1bSJed Brown args: -mat_partitioning_type hierarch -mat_partitioning_hierarchical_nfineparts 2 1783c4762a1bSJed Brown args: -mat_partitioning_hierarchical_coarseparttype ptscotch -mat_partitioning_hierarchical_fineparttype ptscotch 1784c4762a1bSJed Brown 1785c4762a1bSJed Brown test: 1786c4762a1bSJed Brown suffix: 8 1787c4762a1bSJed Brown requires: hdf5 !complex 1788c4762a1bSJed Brown nsize: 4 1789c4762a1bSJed Brown args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf 1790c4762a1bSJed Brown args: -distribute 0 -interpolate after_create 1791c4762a1bSJed Brown args: -view_vertices_from_coords 0.,1.,0.,-0.5,1.,0.,0.583,-0.644,0.,-2.,-2.,-2. -view_vertices_from_coords_tol 1e-3 1792c4762a1bSJed Brown args: -dm_plex_check_all 1793c4762a1bSJed Brown args: -custom_view 1794c4762a1bSJed Brown 1795c4762a1bSJed Brown testset: # 9 HDF5 SEQUANTIAL LOAD 1796c4762a1bSJed Brown requires: hdf5 !complex datafilespath 1797c4762a1bSJed Brown args: -dm_plex_check_all 1798c4762a1bSJed Brown args: -filename ${DATAFILESPATH}/meshes/cube-hexahedra-refined.h5 -dm_plex_create_from_hdf5_xdmf -dm_plex_hdf5_topology_path /cells -dm_plex_hdf5_geometry_path /coordinates 1799c4762a1bSJed Brown args: -dm_plex_hdf5_force_sequential 1800c4762a1bSJed Brown args: -distribute 1801c4762a1bSJed Brown test: # seq load, simple partitioner 1802c4762a1bSJed Brown suffix: 9_seq_hdf5_simple 1803c4762a1bSJed Brown nsize: {{1 2 4 5}} 1804c4762a1bSJed Brown args: -interpolate none 1805c4762a1bSJed Brown test: # seq load, seq interpolation, simple partitioner 1806c4762a1bSJed Brown suffix: 9_seq_hdf5_int_simple 1807c4762a1bSJed Brown nsize: {{1 2 4 5}} 1808c4762a1bSJed Brown args: -interpolate after_create 1809c4762a1bSJed Brown test: # seq load, seq interpolation, metis partitioner 1810c4762a1bSJed Brown nsize: {{2 4 5}} 1811c4762a1bSJed Brown suffix: 9_seq_hdf5_int_metis 1812c4762a1bSJed Brown requires: parmetis 1813c4762a1bSJed Brown args: -interpolate after_create 1814c4762a1bSJed Brown args: -petscpartitioner_type parmetis 1815c4762a1bSJed Brown test: # seq load, simple partitioner, par interpolation 1816c4762a1bSJed Brown suffix: 9_seq_hdf5_simple_int 1817c4762a1bSJed Brown nsize: {{2 4 5}} 1818c4762a1bSJed Brown args: -interpolate after_distribute 1819c4762a1bSJed Brown test: # seq load, simple partitioner, par interpolation 1820c4762a1bSJed Brown # This is like 9_seq_hdf5_simple_int but testing error output of DMPlexCheckPointSFHeavy(). 1821c4762a1bSJed Brown # Once 9_seq_hdf5_simple_int gets fixed, this one gets broken. 1822c4762a1bSJed Brown # We can then provide an intentionally broken mesh instead. 1823c4762a1bSJed Brown TODO: This test is broken because PointSF is fixed. 1824c4762a1bSJed Brown suffix: 9_seq_hdf5_simple_int_err 1825c4762a1bSJed Brown nsize: 4 1826c4762a1bSJed Brown args: -interpolate after_distribute 1827c4762a1bSJed Brown filter: sed -e "/PETSC ERROR/,$$d" 1828c4762a1bSJed Brown test: # seq load, metis partitioner, par interpolation 1829c4762a1bSJed Brown nsize: {{2 4 5}} 1830c4762a1bSJed Brown suffix: 9_seq_hdf5_metis_int 1831c4762a1bSJed Brown requires: parmetis 1832c4762a1bSJed Brown args: -interpolate after_distribute 1833c4762a1bSJed Brown args: -petscpartitioner_type parmetis 1834c4762a1bSJed Brown 1835c4762a1bSJed Brown testset: # 9 HDF5 PARALLEL LOAD 1836c4762a1bSJed Brown requires: hdf5 !complex datafilespath 1837c4762a1bSJed Brown nsize: {{2 4 5}} 1838c4762a1bSJed Brown args: -dm_plex_check_all 1839c4762a1bSJed Brown args: -filename ${DATAFILESPATH}/meshes/cube-hexahedra-refined.h5 -dm_plex_create_from_hdf5_xdmf -dm_plex_hdf5_topology_path /cells -dm_plex_hdf5_geometry_path /coordinates 1840c4762a1bSJed Brown test: # par load 1841c4762a1bSJed Brown suffix: 9_par_hdf5 1842c4762a1bSJed Brown args: -interpolate none 1843c4762a1bSJed Brown test: # par load, par interpolation 1844c4762a1bSJed Brown suffix: 9_par_hdf5_int 1845c4762a1bSJed Brown args: -interpolate after_create 1846c4762a1bSJed Brown test: # par load, parmetis repartitioner 1847c4762a1bSJed Brown TODO: Parallel partitioning of uninterpolated meshes not supported 1848c4762a1bSJed Brown suffix: 9_par_hdf5_parmetis 1849c4762a1bSJed Brown requires: parmetis 1850c4762a1bSJed Brown args: -distribute -petscpartitioner_type parmetis 1851c4762a1bSJed Brown args: -interpolate none 1852c4762a1bSJed Brown test: # par load, par interpolation, parmetis repartitioner 1853c4762a1bSJed Brown suffix: 9_par_hdf5_int_parmetis 1854c4762a1bSJed Brown requires: parmetis 1855c4762a1bSJed Brown args: -distribute -petscpartitioner_type parmetis 1856c4762a1bSJed Brown args: -interpolate after_create 1857c4762a1bSJed Brown test: # par load, parmetis partitioner, par interpolation 1858c4762a1bSJed Brown TODO: Parallel partitioning of uninterpolated meshes not supported 1859c4762a1bSJed Brown suffix: 9_par_hdf5_parmetis_int 1860c4762a1bSJed Brown requires: parmetis 1861c4762a1bSJed Brown args: -distribute -petscpartitioner_type parmetis 1862c4762a1bSJed Brown args: -interpolate after_distribute 1863c4762a1bSJed Brown 1864c4762a1bSJed Brown testset: # 10 HDF5 PARALLEL LOAD 1865c4762a1bSJed Brown requires: hdf5 !complex datafilespath 1866c4762a1bSJed Brown nsize: {{2 4 7}} 1867c4762a1bSJed Brown args: -dm_plex_check_all 1868c4762a1bSJed Brown args: -filename ${DATAFILESPATH}/meshes/cube-hexahedra-refined2.h5 -dm_plex_create_from_hdf5_xdmf -dm_plex_hdf5_topology_path /topo -dm_plex_hdf5_geometry_path /geom 1869c4762a1bSJed Brown test: # par load, par interpolation 1870c4762a1bSJed Brown suffix: 10_par_hdf5_int 1871c4762a1bSJed Brown args: -interpolate after_create 1872c4762a1bSJed Brown TEST*/ 1873