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 225956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 226c4762a1bSJed Brown static PetscLogStage stage[3]; 227956f8c0dSBarry Smith #endif 228c4762a1bSJed Brown 229c4762a1bSJed Brown static PetscErrorCode DMPlexCheckPointSFHeavy(DM, PortableBoundary); 230c4762a1bSJed Brown static PetscErrorCode DMPlexSetOrientInterface_Private(DM, PetscBool); 231c4762a1bSJed Brown static PetscErrorCode DMPlexGetExpandedBoundary_Private(DM, PortableBoundary *); 232c4762a1bSJed Brown static PetscErrorCode DMPlexExpandedConesToFaces_Private(DM, IS, PetscSection, IS *); 233c4762a1bSJed Brown 234d71ae5a4SJacob Faibussowitsch static PetscErrorCode PortableBoundaryDestroy(PortableBoundary *bnd) 235d71ae5a4SJacob Faibussowitsch { 236c4762a1bSJed Brown PetscInt d; 237c4762a1bSJed Brown 238c4762a1bSJed Brown PetscFunctionBegin; 2393ba16761SJacob Faibussowitsch if (!*bnd) PetscFunctionReturn(PETSC_SUCCESS); 2409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*bnd)->coordinates)); 24148a46eb9SPierre Jolivet for (d = 0; d < (*bnd)->depth; d++) PetscCall(PetscSectionDestroy(&(*bnd)->sections[d])); 2429566063dSJacob Faibussowitsch PetscCall(PetscFree((*bnd)->sections)); 2439566063dSJacob Faibussowitsch PetscCall(PetscFree(*bnd)); 2443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 245c4762a1bSJed Brown } 246c4762a1bSJed Brown 247d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 248d71ae5a4SJacob Faibussowitsch { 249c4762a1bSJed Brown const char *interpTypes[4] = {"none", "create", "after_create", "after_distribute"}; 250c4762a1bSJed Brown PetscInt interp = NONE, dim; 251c4762a1bSJed Brown PetscBool flg1, flg2; 252c4762a1bSJed Brown 253c4762a1bSJed Brown PetscFunctionBegin; 254c4762a1bSJed Brown options->debug = 0; 255c4762a1bSJed Brown options->testNum = 0; 256c4762a1bSJed Brown options->dim = 2; 257c4762a1bSJed Brown options->cellSimplex = PETSC_TRUE; 258c4762a1bSJed Brown options->distribute = PETSC_FALSE; 259c4762a1bSJed Brown options->interpolate = NONE; 260c4762a1bSJed Brown options->useGenerator = PETSC_FALSE; 261c4762a1bSJed Brown options->testOrientIF = PETSC_FALSE; 262c4762a1bSJed Brown options->testHeavy = PETSC_TRUE; 263c4762a1bSJed Brown options->customView = PETSC_FALSE; 264c4762a1bSJed Brown options->testExpandPointsEmpty = PETSC_FALSE; 265c4762a1bSJed Brown options->ornt[0] = 0; 266c4762a1bSJed Brown options->ornt[1] = 0; 267c4762a1bSJed Brown options->faces[0] = 2; 268c4762a1bSJed Brown options->faces[1] = 2; 269c4762a1bSJed Brown options->faces[2] = 2; 270c4762a1bSJed Brown options->filename[0] = '\0'; 271c4762a1bSJed Brown options->coordsTol = PETSC_DEFAULT; 272c4762a1bSJed Brown 273d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX"); 2749566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-debug", "The debugging level", "ex18.c", options->debug, &options->debug, NULL, 0)); 2759566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-testnum", "The mesh to create", "ex18.c", options->testNum, &options->testNum, NULL, 0)); 2769566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-cell_simplex", "Generate simplices if true, otherwise hexes", "ex18.c", options->cellSimplex, &options->cellSimplex, NULL)); 2779566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-distribute", "Distribute the mesh", "ex18.c", options->distribute, &options->distribute, NULL)); 2789566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-interpolate", "Type of mesh interpolation (none, create, after_create, after_distribute)", "ex18.c", interpTypes, 4, interpTypes[options->interpolate], &interp, NULL)); 279c4762a1bSJed Brown options->interpolate = (InterpType)interp; 2801dca8a05SBarry Smith PetscCheck(options->distribute || options->interpolate != AFTER_DISTRIBUTE, comm, PETSC_ERR_SUP, "-interpolate after_distribute needs -distribute 1"); 2819566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_generator", "Use a mesh generator to build the mesh", "ex18.c", options->useGenerator, &options->useGenerator, NULL)); 282c4762a1bSJed Brown options->ncoords = 128; 2839566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalarArray("-view_vertices_from_coords", "Print DAG points corresponding to vertices with given coordinates", "ex18.c", options->coords, &options->ncoords, NULL)); 2849566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-view_vertices_from_coords_tol", "Tolerance for -view_vertices_from_coords", "ex18.c", options->coordsTol, &options->coordsTol, NULL)); 285c4762a1bSJed Brown options->nPointsToExpand = 128; 2869566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-test_expand_points", "Expand given array of DAG point using DMPlexGetConeRecursive() and print results", "ex18.c", options->pointsToExpand, &options->nPointsToExpand, NULL)); 28748a46eb9SPierre 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)); 2889566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_heavy", "Run the heavy PointSF test", "ex18.c", options->testHeavy, &options->testHeavy, NULL)); 2899566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-custom_view", "Custom DMPlex view", "ex18.c", options->customView, &options->customView, NULL)); 2909566063dSJacob Faibussowitsch PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex18.c", options->dim, &options->dim, &flg1, 1, 3)); 291c4762a1bSJed Brown dim = 3; 2929566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-faces", "Number of faces per dimension", "ex18.c", options->faces, &dim, &flg2)); 293c4762a1bSJed Brown if (flg2) { 2941dca8a05SBarry 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); 295c4762a1bSJed Brown options->dim = dim; 296c4762a1bSJed Brown } 2979566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-filename", "The mesh file", "ex18.c", options->filename, options->filename, sizeof(options->filename), NULL)); 2989566063dSJacob 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)); 2999566063dSJacob 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)); 30008401ef6SPierre Jolivet PetscCheck(flg2 == options->testOrientIF, comm, PETSC_ERR_ARG_OUTOFRANGE, "neither or both -rotate_interface_0 -rotate_interface_1 must be set"); 301c4762a1bSJed Brown if (options->testOrientIF) { 302c4762a1bSJed Brown PetscInt i; 303c4762a1bSJed Brown for (i = 0; i < 2; i++) { 304c4762a1bSJed Brown if (options->ornt[i] >= 10) options->ornt[i] = -(options->ornt[i] - 10); /* 11 12 13 become -1 -2 -3 */ 305c4762a1bSJed Brown } 306c4762a1bSJed Brown options->filename[0] = 0; 307c4762a1bSJed Brown options->useGenerator = PETSC_FALSE; 308c4762a1bSJed Brown options->dim = 3; 309c4762a1bSJed Brown options->cellSimplex = PETSC_TRUE; 310c4762a1bSJed Brown options->interpolate = CREATE; 311c4762a1bSJed Brown options->distribute = PETSC_FALSE; 312c4762a1bSJed Brown } 313d0609cedSBarry Smith PetscOptionsEnd(); 3143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 315c4762a1bSJed Brown } 316c4762a1bSJed Brown 317d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh_1D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm) 318d71ae5a4SJacob Faibussowitsch { 319c4762a1bSJed Brown PetscInt testNum = user->testNum; 320c4762a1bSJed Brown PetscMPIInt rank, size; 321c4762a1bSJed Brown PetscInt numCorners = 2, i; 322c4762a1bSJed Brown PetscInt numCells, numVertices, network; 323a4a685f2SJacob Faibussowitsch PetscInt *cells; 324c4762a1bSJed Brown PetscReal *coords; 325c4762a1bSJed Brown 326c4762a1bSJed Brown PetscFunctionBegin; 3279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 3289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 32963a3b9bcSJacob Faibussowitsch PetscCheck(size <= 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for <=2 processes", testNum); 330c4762a1bSJed Brown 331c4762a1bSJed Brown numCells = 3; 3329566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ncells", &numCells, NULL)); 33363a3b9bcSJacob Faibussowitsch PetscCheck(numCells >= 3, comm, PETSC_ERR_ARG_OUTOFRANGE, "Test ncells %" PetscInt_FMT " must >=3", numCells); 334c4762a1bSJed Brown 335c4762a1bSJed Brown if (size == 1) { 336c4762a1bSJed Brown numVertices = numCells + 1; 3379566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2 * numCells, &cells, 2 * numVertices, &coords)); 338c4762a1bSJed Brown for (i = 0; i < numCells; i++) { 3399371c9d4SSatish Balay cells[2 * i] = i; 3409371c9d4SSatish Balay cells[2 * i + 1] = i + 1; 3419371c9d4SSatish Balay coords[2 * i] = i; 3429371c9d4SSatish Balay coords[2 * i + 1] = i + 1; 343c4762a1bSJed Brown } 344c4762a1bSJed Brown 3459566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListPetsc(comm, user->dim, numCells, numVertices, numCorners, PETSC_FALSE, cells, user->dim, coords, dm)); 3469566063dSJacob Faibussowitsch PetscCall(PetscFree2(cells, coords)); 3473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 348c4762a1bSJed Brown } 349c4762a1bSJed Brown 350c4762a1bSJed Brown network = 0; 3519566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-network_case", &network, NULL)); 352c4762a1bSJed Brown if (network == 0) { 353c4762a1bSJed Brown switch (rank) { 3549371c9d4SSatish Balay case 0: { 355c4762a1bSJed Brown numCells = 2; 356c4762a1bSJed Brown numVertices = numCells; 3579566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2 * numCells, &cells, 2 * numCells, &coords)); 3589371c9d4SSatish Balay cells[0] = 0; 3599371c9d4SSatish Balay cells[1] = 1; 3609371c9d4SSatish Balay cells[2] = 1; 3619371c9d4SSatish Balay cells[3] = 2; 3629371c9d4SSatish Balay coords[0] = 0.; 3639371c9d4SSatish Balay coords[1] = 1.; 3649371c9d4SSatish Balay coords[2] = 1.; 3659371c9d4SSatish Balay coords[3] = 2.; 3669371c9d4SSatish Balay } break; 3679371c9d4SSatish Balay case 1: { 368c4762a1bSJed Brown numCells -= 2; 369c4762a1bSJed Brown numVertices = numCells + 1; 3709566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2 * numCells, &cells, 2 * numCells, &coords)); 371c4762a1bSJed Brown for (i = 0; i < numCells; i++) { 3729371c9d4SSatish Balay cells[2 * i] = 2 + i; 3739371c9d4SSatish Balay cells[2 * i + 1] = 2 + i + 1; 3749371c9d4SSatish Balay coords[2 * i] = 2 + i; 3759371c9d4SSatish Balay coords[2 * i + 1] = 2 + i + 1; 376c4762a1bSJed Brown } 3779371c9d4SSatish Balay } break; 378d71ae5a4SJacob Faibussowitsch default: 379d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank); 380c4762a1bSJed Brown } 381c4762a1bSJed Brown } else { /* network_case = 1 */ 382c4762a1bSJed Brown /* ----------------------- */ 383c4762a1bSJed Brown switch (rank) { 3849371c9d4SSatish Balay case 0: { 385c4762a1bSJed Brown numCells = 2; 386c4762a1bSJed Brown numVertices = 3; 3879566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2 * numCells, &cells, 2 * numCells, &coords)); 3889371c9d4SSatish Balay cells[0] = 0; 3899371c9d4SSatish Balay cells[1] = 3; 3909371c9d4SSatish Balay cells[2] = 3; 3919371c9d4SSatish Balay cells[3] = 1; 3929371c9d4SSatish Balay } break; 3939371c9d4SSatish Balay case 1: { 394c4762a1bSJed Brown numCells = 1; 395c4762a1bSJed Brown numVertices = 1; 3969566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(2 * numCells, &cells, 2 * numCells, &coords)); 3979371c9d4SSatish Balay cells[0] = 3; 3989371c9d4SSatish Balay cells[1] = 2; 3999371c9d4SSatish Balay } break; 400d71ae5a4SJacob Faibussowitsch default: 401d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank); 402c4762a1bSJed Brown } 403c4762a1bSJed Brown } 4049566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, PETSC_FALSE, cells, user->dim, coords, NULL, NULL, dm)); 4059566063dSJacob Faibussowitsch PetscCall(PetscFree2(cells, coords)); 4063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 407c4762a1bSJed Brown } 408c4762a1bSJed Brown 409d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSimplex_2D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm) 410d71ae5a4SJacob Faibussowitsch { 411c4762a1bSJed Brown PetscInt testNum = user->testNum, p; 412c4762a1bSJed Brown PetscMPIInt rank, size; 413c4762a1bSJed Brown 414c4762a1bSJed Brown PetscFunctionBegin; 4159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 4169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 417c4762a1bSJed Brown switch (testNum) { 418c4762a1bSJed Brown case 0: 41963a3b9bcSJacob Faibussowitsch PetscCheck(size == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 2 processes", testNum); 420c4762a1bSJed Brown switch (rank) { 4219371c9d4SSatish Balay case 0: { 422c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 2, numCorners = 3; 423a4a685f2SJacob Faibussowitsch const PetscInt cells[3] = {0, 1, 2}; 424c4762a1bSJed Brown PetscReal coords[4] = {-0.5, 0.5, 0.0, 0.0}; 425c4762a1bSJed Brown PetscInt markerPoints[6] = {1, 1, 2, 1, 3, 1}; 426c4762a1bSJed Brown 4279566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 4289566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 4299371c9d4SSatish Balay } break; 4309371c9d4SSatish Balay case 1: { 431c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 2, numCorners = 3; 432a4a685f2SJacob Faibussowitsch const PetscInt cells[3] = {1, 3, 2}; 433c4762a1bSJed Brown PetscReal coords[4] = {0.0, 1.0, 0.5, 0.5}; 434c4762a1bSJed Brown PetscInt markerPoints[6] = {1, 1, 2, 1, 3, 1}; 435c4762a1bSJed Brown 4369566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 4379566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 4389371c9d4SSatish Balay } break; 439d71ae5a4SJacob Faibussowitsch default: 440d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank); 441c4762a1bSJed Brown } 442c4762a1bSJed Brown break; 443c4762a1bSJed Brown case 1: 44463a3b9bcSJacob Faibussowitsch PetscCheck(size == 3, comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 3 processes", testNum); 445c4762a1bSJed Brown switch (rank) { 4469371c9d4SSatish Balay case 0: { 447c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 2, numCorners = 3; 448a4a685f2SJacob Faibussowitsch const PetscInt cells[3] = {0, 1, 2}; 449c4762a1bSJed Brown PetscReal coords[4] = {0.0, 1.0, 0.0, 0.0}; 450c4762a1bSJed Brown PetscInt markerPoints[6] = {1, 1, 2, 1, 3, 1}; 451c4762a1bSJed Brown 4529566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 4539566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 4549371c9d4SSatish Balay } break; 4559371c9d4SSatish Balay case 1: { 456c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 2, numCorners = 3; 457a4a685f2SJacob Faibussowitsch const PetscInt cells[3] = {0, 2, 3}; 458c4762a1bSJed Brown PetscReal coords[4] = {0.5, 0.5, 1.0, 1.0}; 459c4762a1bSJed Brown PetscInt markerPoints[6] = {1, 1, 2, 1, 3, 1}; 460c4762a1bSJed Brown 4619566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 4629566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 4639371c9d4SSatish Balay } break; 4649371c9d4SSatish Balay case 2: { 465c4762a1bSJed Brown const PetscInt numCells = 2, numVertices = 1, numCorners = 3; 466a4a685f2SJacob Faibussowitsch const PetscInt cells[6] = {2, 4, 3, 2, 1, 4}; 467c4762a1bSJed Brown PetscReal coords[2] = {1.0, 0.0}; 468c4762a1bSJed Brown PetscInt markerPoints[10] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1}; 469c4762a1bSJed Brown 4709566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 4719566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 4729371c9d4SSatish Balay } break; 473d71ae5a4SJacob Faibussowitsch default: 474d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank); 475c4762a1bSJed Brown } 476c4762a1bSJed Brown break; 477c4762a1bSJed Brown case 2: 47863a3b9bcSJacob Faibussowitsch PetscCheck(size == 3, comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 3 processes", testNum); 479c4762a1bSJed Brown switch (rank) { 4809371c9d4SSatish Balay case 0: { 481c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 2, numCorners = 3; 482a4a685f2SJacob Faibussowitsch const PetscInt cells[3] = {1, 2, 0}; 483c4762a1bSJed Brown PetscReal coords[4] = {0.5, 0.5, 0.0, 1.0}; 484c4762a1bSJed Brown PetscInt markerPoints[6] = {1, 1, 2, 1, 3, 1}; 485c4762a1bSJed Brown 4869566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 4879566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 4889371c9d4SSatish Balay } break; 4899371c9d4SSatish Balay case 1: { 490c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 2, numCorners = 3; 491a4a685f2SJacob Faibussowitsch const PetscInt cells[3] = {1, 0, 3}; 492c4762a1bSJed Brown PetscReal coords[4] = {0.0, 0.0, 1.0, 1.0}; 493c4762a1bSJed Brown PetscInt markerPoints[6] = {1, 1, 2, 1, 3, 1}; 494c4762a1bSJed Brown 4959566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 4969566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 4979371c9d4SSatish Balay } break; 4989371c9d4SSatish Balay case 2: { 499c4762a1bSJed Brown const PetscInt numCells = 2, numVertices = 1, numCorners = 3; 500a4a685f2SJacob Faibussowitsch const PetscInt cells[6] = {0, 4, 3, 0, 2, 4}; 501c4762a1bSJed Brown PetscReal coords[2] = {1.0, 0.0}; 502c4762a1bSJed Brown PetscInt markerPoints[10] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1}; 503c4762a1bSJed Brown 5049566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 5059566063dSJacob Faibussowitsch for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 5069371c9d4SSatish Balay } break; 507d71ae5a4SJacob Faibussowitsch default: 508d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank); 509c4762a1bSJed Brown } 510c4762a1bSJed Brown break; 511d71ae5a4SJacob Faibussowitsch default: 512d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 513c4762a1bSJed Brown } 5143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 515c4762a1bSJed Brown } 516c4762a1bSJed Brown 517d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSimplex_3D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm) 518d71ae5a4SJacob Faibussowitsch { 519c4762a1bSJed Brown PetscInt testNum = user->testNum, p; 520c4762a1bSJed Brown PetscMPIInt rank, size; 521c4762a1bSJed Brown 522c4762a1bSJed Brown PetscFunctionBegin; 5239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 5249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 525c4762a1bSJed Brown switch (testNum) { 526c4762a1bSJed Brown case 0: 52763a3b9bcSJacob Faibussowitsch PetscCheck(size == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 2 processes", testNum); 528c4762a1bSJed Brown switch (rank) { 5299371c9d4SSatish Balay case 0: { 530c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 2, numCorners = 4; 531a4a685f2SJacob Faibussowitsch const PetscInt cells[4] = {0, 2, 1, 3}; 532c4762a1bSJed Brown PetscReal coords[6] = {0.0, 0.0, -0.5, 0.0, -0.5, 0.0}; 533c4762a1bSJed Brown PetscInt markerPoints[8] = {1, 1, 2, 1, 3, 1, 4, 1}; 534c4762a1bSJed Brown 5359566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 5369566063dSJacob Faibussowitsch for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 5379371c9d4SSatish Balay } break; 5389371c9d4SSatish Balay case 1: { 539c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 3, numCorners = 4; 540a4a685f2SJacob Faibussowitsch const PetscInt cells[4] = {1, 2, 4, 3}; 541c4762a1bSJed Brown PetscReal coords[9] = {1.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.5}; 542c4762a1bSJed Brown PetscInt markerPoints[8] = {1, 1, 2, 1, 3, 1, 4, 1}; 543c4762a1bSJed Brown 5449566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 5459566063dSJacob Faibussowitsch for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 5469371c9d4SSatish Balay } break; 547d71ae5a4SJacob Faibussowitsch default: 548d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank); 549c4762a1bSJed Brown } 550c4762a1bSJed Brown break; 551d71ae5a4SJacob Faibussowitsch default: 552d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 553c4762a1bSJed Brown } 554c4762a1bSJed Brown if (user->testOrientIF) { 555c4762a1bSJed Brown PetscInt ifp[] = {8, 6}; 556c4762a1bSJed Brown 5579566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Mesh before orientation")); 5589566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-before_orientation_dm_view")); 559c4762a1bSJed Brown /* rotate interface face ifp[rank] by given orientation ornt[rank] */ 5609566063dSJacob Faibussowitsch PetscCall(DMPlexOrientPoint(*dm, ifp[rank], user->ornt[rank])); 5619566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-before_orientation_dm_view")); 5629566063dSJacob Faibussowitsch PetscCall(DMPlexCheckFaces(*dm, 0)); 5639566063dSJacob Faibussowitsch PetscCall(DMPlexOrientInterface_Internal(*dm)); 5649566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Orientation test PASSED\n")); 565c4762a1bSJed Brown } 5663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 567c4762a1bSJed Brown } 568c4762a1bSJed Brown 569d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateQuad_2D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm) 570d71ae5a4SJacob Faibussowitsch { 571c4762a1bSJed Brown PetscInt testNum = user->testNum, p; 572c4762a1bSJed Brown PetscMPIInt rank, size; 573c4762a1bSJed Brown 574c4762a1bSJed Brown PetscFunctionBegin; 5759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 5769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 577c4762a1bSJed Brown switch (testNum) { 578c4762a1bSJed Brown case 0: 57963a3b9bcSJacob Faibussowitsch PetscCheck(size == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 2 processes", testNum); 580c4762a1bSJed Brown switch (rank) { 5819371c9d4SSatish Balay case 0: { 582c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 3, numCorners = 4; 583a4a685f2SJacob Faibussowitsch const PetscInt cells[4] = {0, 1, 2, 3}; 584c4762a1bSJed Brown PetscReal coords[6] = {-0.5, 0.0, 0.0, 0.0, 0.0, 1.0}; 585c4762a1bSJed Brown PetscInt markerPoints[4 * 2] = {1, 1, 2, 1, 3, 1, 4, 1}; 586c4762a1bSJed Brown 5879566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 5889566063dSJacob Faibussowitsch for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 5899371c9d4SSatish Balay } break; 5909371c9d4SSatish Balay case 1: { 591c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 3, numCorners = 4; 592a4a685f2SJacob Faibussowitsch const PetscInt cells[4] = {1, 4, 5, 2}; 593c4762a1bSJed Brown PetscReal coords[6] = {-0.5, 1.0, 0.5, 0.0, 0.5, 1.0}; 594c4762a1bSJed Brown PetscInt markerPoints[4 * 2] = {1, 1, 2, 1, 3, 1, 4, 1}; 595c4762a1bSJed Brown 5969566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 5979566063dSJacob Faibussowitsch for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 5989371c9d4SSatish Balay } break; 599d71ae5a4SJacob Faibussowitsch default: 600d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank); 601c4762a1bSJed Brown } 602c4762a1bSJed Brown break; 603d71ae5a4SJacob Faibussowitsch default: 604d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 605c4762a1bSJed Brown } 6063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 607c4762a1bSJed Brown } 608c4762a1bSJed Brown 609d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateHex_3D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm) 610d71ae5a4SJacob Faibussowitsch { 611c4762a1bSJed Brown PetscInt testNum = user->testNum, p; 612c4762a1bSJed Brown PetscMPIInt rank, size; 613c4762a1bSJed Brown 614c4762a1bSJed Brown PetscFunctionBegin; 6159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 6169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 617c4762a1bSJed Brown switch (testNum) { 618c4762a1bSJed Brown case 0: 61963a3b9bcSJacob Faibussowitsch PetscCheck(size == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 2 processes", testNum); 620c4762a1bSJed Brown switch (rank) { 6219371c9d4SSatish Balay case 0: { 622c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 6, numCorners = 8; 623a4a685f2SJacob Faibussowitsch const PetscInt cells[8] = {0, 3, 2, 1, 4, 5, 6, 7}; 624c4762a1bSJed 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}; 625c4762a1bSJed Brown PetscInt markerPoints[8 * 2] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1, 8, 1, 9, 1}; 626c4762a1bSJed Brown 6279566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 6289566063dSJacob Faibussowitsch for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 6299371c9d4SSatish Balay } break; 6309371c9d4SSatish Balay case 1: { 631c4762a1bSJed Brown const PetscInt numCells = 1, numVertices = 6, numCorners = 8; 632a4a685f2SJacob Faibussowitsch const PetscInt cells[8] = {1, 2, 9, 8, 5, 10, 11, 6}; 633c4762a1bSJed 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}; 634c4762a1bSJed Brown PetscInt markerPoints[8 * 2] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1, 8, 1, 9, 1}; 635c4762a1bSJed Brown 6369566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm)); 6379566063dSJacob Faibussowitsch for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); 6389371c9d4SSatish Balay } break; 639d71ae5a4SJacob Faibussowitsch default: 640d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank); 641c4762a1bSJed Brown } 642c4762a1bSJed Brown break; 643d71ae5a4SJacob Faibussowitsch default: 644d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum); 645c4762a1bSJed Brown } 6463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 647c4762a1bSJed Brown } 648c4762a1bSJed Brown 649d71ae5a4SJacob Faibussowitsch static PetscErrorCode CustomView(DM dm, PetscViewer v) 650d71ae5a4SJacob Faibussowitsch { 651c4762a1bSJed Brown DMPlexInterpolatedFlag interpolated; 652c4762a1bSJed Brown PetscBool distributed; 653c4762a1bSJed Brown 654c4762a1bSJed Brown PetscFunctionBegin; 6559566063dSJacob Faibussowitsch PetscCall(DMPlexIsDistributed(dm, &distributed)); 6569566063dSJacob Faibussowitsch PetscCall(DMPlexIsInterpolatedCollective(dm, &interpolated)); 6579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, "DMPlexIsDistributed: %s\n", PetscBools[distributed])); 6589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, "DMPlexIsInterpolatedCollective: %s\n", DMPlexInterpolatedFlags[interpolated])); 6593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 660c4762a1bSJed Brown } 661c4762a1bSJed Brown 662d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMeshFromFile(MPI_Comm comm, AppCtx *user, DM *dm, DM *serialDM) 663d71ae5a4SJacob Faibussowitsch { 664c4762a1bSJed Brown const char *filename = user->filename; 665c4762a1bSJed Brown PetscBool testHeavy = user->testHeavy; 666c4762a1bSJed Brown PetscBool interpCreate = user->interpolate == CREATE ? PETSC_TRUE : PETSC_FALSE; 667c4762a1bSJed Brown PetscBool distributed = PETSC_FALSE; 668c4762a1bSJed Brown 669c4762a1bSJed Brown PetscFunctionBegin; 670c4762a1bSJed Brown *serialDM = NULL; 6719566063dSJacob Faibussowitsch if (testHeavy && interpCreate) PetscCall(DMPlexSetOrientInterface_Private(NULL, PETSC_FALSE)); 6729566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[0])); 6739566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromFile(comm, filename, "ex18_plex", interpCreate, dm)); /* with DMPlexOrientInterface_Internal() call skipped so that PointSF issues are left to DMPlexCheckPointSFHeavy() */ 6749566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 6759566063dSJacob Faibussowitsch if (testHeavy && interpCreate) PetscCall(DMPlexSetOrientInterface_Private(NULL, PETSC_TRUE)); 6769566063dSJacob Faibussowitsch PetscCall(DMPlexIsDistributed(*dm, &distributed)); 6779566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "DMPlexCreateFromFile produced %s mesh.\n", distributed ? "distributed" : "serial")); 678c4762a1bSJed Brown if (testHeavy && distributed) { 6799566063dSJacob Faibussowitsch PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_hdf5_force_sequential", NULL)); 6809566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromFile(comm, filename, "ex18_plex", interpCreate, serialDM)); 6819566063dSJacob Faibussowitsch PetscCall(DMPlexIsDistributed(*serialDM, &distributed)); 68228b400f6SJacob Faibussowitsch PetscCheck(!distributed, comm, PETSC_ERR_PLIB, "unable to create a serial DM from file"); 683c4762a1bSJed Brown } 6849566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &user->dim)); 6853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 686c4762a1bSJed Brown } 687c4762a1bSJed Brown 688d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 689d71ae5a4SJacob Faibussowitsch { 690c4762a1bSJed Brown PetscPartitioner part; 691c4762a1bSJed Brown PortableBoundary boundary = NULL; 692c4762a1bSJed Brown DM serialDM = NULL; 693c4762a1bSJed Brown PetscBool cellSimplex = user->cellSimplex; 694c4762a1bSJed Brown PetscBool useGenerator = user->useGenerator; 695c4762a1bSJed Brown PetscBool interpCreate = user->interpolate == CREATE ? PETSC_TRUE : PETSC_FALSE; 696c4762a1bSJed Brown PetscBool interpSerial = user->interpolate == AFTER_CREATE ? PETSC_TRUE : PETSC_FALSE; 697c4762a1bSJed Brown PetscBool interpParallel = user->interpolate == AFTER_DISTRIBUTE ? PETSC_TRUE : PETSC_FALSE; 698c4762a1bSJed Brown PetscBool testHeavy = user->testHeavy; 699c4762a1bSJed Brown PetscMPIInt rank; 700c4762a1bSJed Brown 701c4762a1bSJed Brown PetscFunctionBegin; 7029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 703c4762a1bSJed Brown if (user->filename[0]) { 7049566063dSJacob Faibussowitsch PetscCall(CreateMeshFromFile(comm, user, dm, &serialDM)); 705c4762a1bSJed Brown } else if (useGenerator) { 7069566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[0])); 7079566063dSJacob Faibussowitsch PetscCall(DMPlexCreateBoxMesh(comm, user->dim, cellSimplex, user->faces, NULL, NULL, NULL, interpCreate, dm)); 7089566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 709c4762a1bSJed Brown } else { 7109566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[0])); 711c4762a1bSJed Brown switch (user->dim) { 712d71ae5a4SJacob Faibussowitsch case 1: 713d71ae5a4SJacob Faibussowitsch PetscCall(CreateMesh_1D(comm, interpCreate, user, dm)); 714d71ae5a4SJacob Faibussowitsch break; 715c4762a1bSJed Brown case 2: 716c4762a1bSJed Brown if (cellSimplex) { 7179566063dSJacob Faibussowitsch PetscCall(CreateSimplex_2D(comm, interpCreate, user, dm)); 718c4762a1bSJed Brown } else { 7199566063dSJacob Faibussowitsch PetscCall(CreateQuad_2D(comm, interpCreate, user, dm)); 720c4762a1bSJed Brown } 721c4762a1bSJed Brown break; 722c4762a1bSJed Brown case 3: 723c4762a1bSJed Brown if (cellSimplex) { 7249566063dSJacob Faibussowitsch PetscCall(CreateSimplex_3D(comm, interpCreate, user, dm)); 725c4762a1bSJed Brown } else { 7269566063dSJacob Faibussowitsch PetscCall(CreateHex_3D(comm, interpCreate, user, dm)); 727c4762a1bSJed Brown } 728c4762a1bSJed Brown break; 729d71ae5a4SJacob Faibussowitsch default: 730d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %" PetscInt_FMT, user->dim); 731c4762a1bSJed Brown } 7329566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 733c4762a1bSJed Brown } 734da81f932SPierre 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); 7359566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Original Mesh")); 7369566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view")); 737c4762a1bSJed Brown 738c4762a1bSJed Brown if (interpSerial) { 739c4762a1bSJed Brown DM idm; 740c4762a1bSJed Brown 7419566063dSJacob Faibussowitsch if (testHeavy) PetscCall(DMPlexSetOrientInterface_Private(*dm, PETSC_FALSE)); 7429566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[2])); 7439566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); /* with DMPlexOrientInterface_Internal() call skipped so that PointSF issues are left to DMPlexCheckPointSFHeavy() */ 7449566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 7459566063dSJacob Faibussowitsch if (testHeavy) PetscCall(DMPlexSetOrientInterface_Private(*dm, PETSC_TRUE)); 7469566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 747c4762a1bSJed Brown *dm = idm; 7489566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Interpolated Mesh")); 7499566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-intp_dm_view")); 750c4762a1bSJed Brown } 751c4762a1bSJed Brown 752c4762a1bSJed Brown /* Set partitioner options */ 7539566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(*dm, &part)); 754c4762a1bSJed Brown if (part) { 7559566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSIMPLE)); 7569566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part)); 757c4762a1bSJed Brown } 758c4762a1bSJed Brown 7599566063dSJacob Faibussowitsch if (user->customView) PetscCall(CustomView(*dm, PETSC_VIEWER_STDOUT_(comm))); 760c4762a1bSJed Brown if (testHeavy) { 761c4762a1bSJed Brown PetscBool distributed; 762c4762a1bSJed Brown 7639566063dSJacob Faibussowitsch PetscCall(DMPlexIsDistributed(*dm, &distributed)); 764c4762a1bSJed Brown if (!serialDM && !distributed) { 765c4762a1bSJed Brown serialDM = *dm; 7669566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)*dm)); 767c4762a1bSJed Brown } 76848a46eb9SPierre Jolivet if (serialDM) PetscCall(DMPlexGetExpandedBoundary_Private(serialDM, &boundary)); 769c4762a1bSJed Brown if (boundary) { 770c4762a1bSJed Brown /* check DM which has been created in parallel and already interpolated */ 7719566063dSJacob Faibussowitsch PetscCall(DMPlexCheckPointSFHeavy(*dm, boundary)); 772c4762a1bSJed Brown } 773c4762a1bSJed Brown /* Orient interface because it could be deliberately skipped above. It is idempotent. */ 7749566063dSJacob Faibussowitsch PetscCall(DMPlexOrientInterface_Internal(*dm)); 775c4762a1bSJed Brown } 776c4762a1bSJed Brown if (user->distribute) { 777c4762a1bSJed Brown DM pdm = NULL; 778c4762a1bSJed Brown 779c4762a1bSJed Brown /* Redistribute mesh over processes using that partitioner */ 7809566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[1])); 7819566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(*dm, 0, NULL, &pdm)); 7829566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 783c4762a1bSJed Brown if (pdm) { 7849566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 785c4762a1bSJed Brown *dm = pdm; 7869566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Redistributed Mesh")); 7879566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dist_dm_view")); 788c4762a1bSJed Brown } 789c4762a1bSJed Brown 790c4762a1bSJed Brown if (interpParallel) { 791c4762a1bSJed Brown DM idm; 792c4762a1bSJed Brown 7939566063dSJacob Faibussowitsch if (testHeavy) PetscCall(DMPlexSetOrientInterface_Private(*dm, PETSC_FALSE)); 7949566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[2])); 7959566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); /* with DMPlexOrientInterface_Internal() call skipped so that PointSF issues are left to DMPlexCheckPointSFHeavy() */ 7969566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 7979566063dSJacob Faibussowitsch if (testHeavy) PetscCall(DMPlexSetOrientInterface_Private(*dm, PETSC_TRUE)); 7989566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 799c4762a1bSJed Brown *dm = idm; 8009566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Interpolated Redistributed Mesh")); 8019566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-intp_dm_view")); 802c4762a1bSJed Brown } 803c4762a1bSJed Brown } 804c4762a1bSJed Brown if (testHeavy) { 8051baa6e33SBarry Smith if (boundary) PetscCall(DMPlexCheckPointSFHeavy(*dm, boundary)); 806c4762a1bSJed Brown /* Orient interface because it could be deliberately skipped above. It is idempotent. */ 8079566063dSJacob Faibussowitsch PetscCall(DMPlexOrientInterface_Internal(*dm)); 808c4762a1bSJed Brown } 809c4762a1bSJed Brown 8109566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Parallel Mesh")); 8119566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 8129566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 8139566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 814c4762a1bSJed Brown 8159566063dSJacob Faibussowitsch if (user->customView) PetscCall(CustomView(*dm, PETSC_VIEWER_STDOUT_(comm))); 8169566063dSJacob Faibussowitsch PetscCall(DMDestroy(&serialDM)); 8179566063dSJacob Faibussowitsch PetscCall(PortableBoundaryDestroy(&boundary)); 8183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 819c4762a1bSJed Brown } 820c4762a1bSJed Brown 821d3e1f4ccSVaclav Hapla #define ps2d(number) ((double)PetscRealPart(number)) 822d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode coord2str(char buf[], size_t len, PetscInt dim, const PetscScalar coords[], PetscReal tol) 823d71ae5a4SJacob Faibussowitsch { 824c4762a1bSJed Brown PetscFunctionBegin; 82508401ef6SPierre Jolivet PetscCheck(dim <= 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "dim must be less than or equal 3"); 826c4762a1bSJed Brown if (tol >= 1e-3) { 827c4762a1bSJed Brown switch (dim) { 828d71ae5a4SJacob Faibussowitsch case 1: 829d71ae5a4SJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, len, "(%12.3f)", ps2d(coords[0]))); 830d71ae5a4SJacob Faibussowitsch case 2: 831d71ae5a4SJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, len, "(%12.3f, %12.3f)", ps2d(coords[0]), ps2d(coords[1]))); 832d71ae5a4SJacob Faibussowitsch case 3: 833d71ae5a4SJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, len, "(%12.3f, %12.3f, %12.3f)", ps2d(coords[0]), ps2d(coords[1]), ps2d(coords[2]))); 834c4762a1bSJed Brown } 835c4762a1bSJed Brown } else { 836c4762a1bSJed Brown switch (dim) { 837d71ae5a4SJacob Faibussowitsch case 1: 838d71ae5a4SJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, len, "(%12.6f)", ps2d(coords[0]))); 839d71ae5a4SJacob Faibussowitsch case 2: 840d71ae5a4SJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, len, "(%12.6f, %12.6f)", ps2d(coords[0]), ps2d(coords[1]))); 841d71ae5a4SJacob Faibussowitsch case 3: 842d71ae5a4SJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, len, "(%12.6f, %12.6f, %12.6f)", ps2d(coords[0]), ps2d(coords[1]), ps2d(coords[2]))); 843c4762a1bSJed Brown } 844c4762a1bSJed Brown } 8453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 846c4762a1bSJed Brown } 847c4762a1bSJed Brown 848d71ae5a4SJacob Faibussowitsch static PetscErrorCode ViewVerticesFromCoords(DM dm, Vec coordsVec, PetscReal tol, PetscViewer viewer) 849d71ae5a4SJacob Faibussowitsch { 850d3e1f4ccSVaclav Hapla PetscInt dim, i, npoints; 851d3e1f4ccSVaclav Hapla IS pointsIS; 852d3e1f4ccSVaclav Hapla const PetscInt *points; 853d3e1f4ccSVaclav Hapla const PetscScalar *coords; 854c4762a1bSJed Brown char coordstr[128]; 855c4762a1bSJed Brown MPI_Comm comm; 856c4762a1bSJed Brown PetscMPIInt rank; 857c4762a1bSJed Brown 858c4762a1bSJed Brown PetscFunctionBegin; 8599566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 8609566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 8619566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 8629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 8639566063dSJacob Faibussowitsch PetscCall(DMPlexFindVertices(dm, coordsVec, tol, &pointsIS)); 8649566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pointsIS, &points)); 8659566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pointsIS, &npoints)); 8669566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coordsVec, &coords)); 867c4762a1bSJed Brown for (i = 0; i < npoints; i++) { 8689566063dSJacob Faibussowitsch PetscCall(coord2str(coordstr, sizeof(coordstr), dim, &coords[i * dim], tol)); 8699566063dSJacob Faibussowitsch if (rank == 0 && i) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "-----\n")); 87063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %s --> points[%" PetscInt_FMT "] = %" PetscInt_FMT "\n", rank, coordstr, i, points[i])); 8719566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 872c4762a1bSJed Brown } 8739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 8749566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coordsVec, &coords)); 8759566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pointsIS, &points)); 8769566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointsIS)); 8773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 878c4762a1bSJed Brown } 879c4762a1bSJed Brown 880d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestExpandPoints(DM dm, AppCtx *user) 881d71ae5a4SJacob Faibussowitsch { 882c4762a1bSJed Brown IS is; 883c4762a1bSJed Brown PetscSection *sects; 884c4762a1bSJed Brown IS *iss; 885c4762a1bSJed Brown PetscInt d, depth; 886c4762a1bSJed Brown PetscMPIInt rank; 887c4762a1bSJed Brown PetscViewer viewer = PETSC_VIEWER_STDOUT_WORLD, sviewer; 888c4762a1bSJed Brown 889c4762a1bSJed Brown PetscFunctionBegin; 8909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 891dd400576SPatrick Sanan if (user->testExpandPointsEmpty && rank == 0) { 8929566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_USE_POINTER, &is)); 893c4762a1bSJed Brown } else { 8949566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, user->nPointsToExpand, user->pointsToExpand, PETSC_USE_POINTER, &is)); 895c4762a1bSJed Brown } 8969566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeRecursive(dm, is, &depth, &iss, §s)); 8979566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 8989566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(sviewer, "[%d] ==========================\n", rank)); 899c4762a1bSJed Brown for (d = depth - 1; d >= 0; d--) { 900c4762a1bSJed Brown IS checkIS; 901c4762a1bSJed Brown PetscBool flg; 902c4762a1bSJed Brown 90363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(sviewer, "depth %" PetscInt_FMT " ---------------\n", d)); 9049566063dSJacob Faibussowitsch PetscCall(PetscSectionView(sects[d], sviewer)); 9059566063dSJacob Faibussowitsch PetscCall(ISView(iss[d], sviewer)); 906c4762a1bSJed Brown /* check reverse operation */ 907c4762a1bSJed Brown if (d < depth - 1) { 9089566063dSJacob Faibussowitsch PetscCall(DMPlexExpandedConesToFaces_Private(dm, iss[d], sects[d], &checkIS)); 9099566063dSJacob Faibussowitsch PetscCall(ISEqualUnsorted(checkIS, iss[d + 1], &flg)); 91028b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)checkIS), PETSC_ERR_PLIB, "DMPlexExpandedConesToFaces_Private produced wrong IS"); 9119566063dSJacob Faibussowitsch PetscCall(ISDestroy(&checkIS)); 912c4762a1bSJed Brown } 913c4762a1bSJed Brown } 9149566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 9159566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 9169566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreConeRecursive(dm, is, &depth, &iss, §s)); 9179566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 9183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 919c4762a1bSJed Brown } 920c4762a1bSJed Brown 921d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexExpandedConesToFaces_Private(DM dm, IS is, PetscSection section, IS *newis) 922d71ae5a4SJacob Faibussowitsch { 923c4762a1bSJed Brown PetscInt n, n1, ncone, numCoveredPoints, o, p, q, start, end; 924c4762a1bSJed Brown const PetscInt *coveredPoints; 925c4762a1bSJed Brown const PetscInt *arr, *cone; 926c4762a1bSJed Brown PetscInt *newarr; 927c4762a1bSJed Brown 928c4762a1bSJed Brown PetscFunctionBegin; 9299566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &n)); 9309566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(section, &n1)); 9319566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &start, &end)); 93263a3b9bcSJacob Faibussowitsch PetscCheck(n == n1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "IS size = %" PetscInt_FMT " != %" PetscInt_FMT " = section storage size", n, n1); 9339566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &arr)); 9349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(end - start, &newarr)); 935c4762a1bSJed Brown for (q = start; q < end; q++) { 9369566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, q, &ncone)); 9379566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, q, &o)); 938c4762a1bSJed Brown cone = &arr[o]; 939c4762a1bSJed Brown if (ncone == 1) { 940c4762a1bSJed Brown numCoveredPoints = 1; 941c4762a1bSJed Brown p = cone[0]; 942c4762a1bSJed Brown } else { 943c4762a1bSJed Brown PetscInt i; 944c4762a1bSJed Brown p = PETSC_MAX_INT; 9459371c9d4SSatish Balay for (i = 0; i < ncone; i++) 9469371c9d4SSatish Balay if (cone[i] < 0) { 9479371c9d4SSatish Balay p = -1; 9489371c9d4SSatish Balay break; 9499371c9d4SSatish Balay } 950c4762a1bSJed Brown if (p >= 0) { 9519566063dSJacob Faibussowitsch PetscCall(DMPlexGetJoin(dm, ncone, cone, &numCoveredPoints, &coveredPoints)); 95263a3b9bcSJacob Faibussowitsch PetscCheck(numCoveredPoints <= 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "more than one covered points for section point %" PetscInt_FMT, q); 953c4762a1bSJed Brown if (numCoveredPoints) p = coveredPoints[0]; 954c4762a1bSJed Brown else p = -2; 9559566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(dm, ncone, cone, &numCoveredPoints, &coveredPoints)); 956c4762a1bSJed Brown } 957c4762a1bSJed Brown } 958c4762a1bSJed Brown newarr[q - start] = p; 959c4762a1bSJed Brown } 9609566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &arr)); 9619566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, end - start, newarr, PETSC_OWN_POINTER, newis)); 9623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 963c4762a1bSJed Brown } 964c4762a1bSJed Brown 965d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexExpandedVerticesToFaces_Private(DM dm, IS boundary_expanded_is, PetscInt depth, PetscSection sections[], IS *boundary_is) 966d71ae5a4SJacob Faibussowitsch { 967c4762a1bSJed Brown PetscInt d; 968c4762a1bSJed Brown IS is, newis; 969c4762a1bSJed Brown 970c4762a1bSJed Brown PetscFunctionBegin; 971c4762a1bSJed Brown is = boundary_expanded_is; 9729566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)is)); 973c4762a1bSJed Brown for (d = 0; d < depth - 1; ++d) { 9749566063dSJacob Faibussowitsch PetscCall(DMPlexExpandedConesToFaces_Private(dm, is, sections[d], &newis)); 9759566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 976c4762a1bSJed Brown is = newis; 977c4762a1bSJed Brown } 978c4762a1bSJed Brown *boundary_is = is; 9793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 980c4762a1bSJed Brown } 981c4762a1bSJed Brown 9829371c9d4SSatish Balay #define CHKERRQI(incall, ierr) \ 983ad540459SPierre Jolivet if (ierr) incall = PETSC_FALSE; 984c4762a1bSJed Brown 985d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMLabelViewFromOptionsOnComm_Private(DMLabel label, const char optionname[], MPI_Comm comm) 986d71ae5a4SJacob Faibussowitsch { 987c4762a1bSJed Brown PetscViewer viewer; 988c4762a1bSJed Brown PetscBool flg; 989c4762a1bSJed Brown static PetscBool incall = PETSC_FALSE; 990c4762a1bSJed Brown PetscViewerFormat format; 991c4762a1bSJed Brown 992c4762a1bSJed Brown PetscFunctionBegin; 9933ba16761SJacob Faibussowitsch if (incall) PetscFunctionReturn(PETSC_SUCCESS); 994c4762a1bSJed Brown incall = PETSC_TRUE; 9955f80ce2aSJacob Faibussowitsch CHKERRQI(incall, PetscOptionsGetViewer(comm, ((PetscObject)label)->options, ((PetscObject)label)->prefix, optionname, &viewer, &format, &flg)); 996c4762a1bSJed Brown if (flg) { 9975f80ce2aSJacob Faibussowitsch CHKERRQI(incall, PetscViewerPushFormat(viewer, format)); 9985f80ce2aSJacob Faibussowitsch CHKERRQI(incall, DMLabelView(label, viewer)); 9995f80ce2aSJacob Faibussowitsch CHKERRQI(incall, PetscViewerPopFormat(viewer)); 10005f80ce2aSJacob Faibussowitsch CHKERRQI(incall, PetscViewerDestroy(&viewer)); 1001c4762a1bSJed Brown } 1002c4762a1bSJed Brown incall = PETSC_FALSE; 10033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1004c4762a1bSJed Brown } 1005c4762a1bSJed Brown 1006c4762a1bSJed Brown /* TODO: this is hotfixing DMLabelGetStratumIS() - it should be fixed systematically instead */ 1007d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode DMLabelGetStratumISOnComm_Private(DMLabel label, PetscInt value, MPI_Comm comm, IS *is) 1008d71ae5a4SJacob Faibussowitsch { 1009c4762a1bSJed Brown IS tmpis; 1010c4762a1bSJed Brown 1011c4762a1bSJed Brown PetscFunctionBegin; 10129566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, value, &tmpis)); 10139566063dSJacob Faibussowitsch if (!tmpis) PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_USE_POINTER, &tmpis)); 10149566063dSJacob Faibussowitsch PetscCall(ISOnComm(tmpis, comm, PETSC_COPY_VALUES, is)); 10159566063dSJacob Faibussowitsch PetscCall(ISDestroy(&tmpis)); 10163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1017c4762a1bSJed Brown } 1018c4762a1bSJed Brown 1019c4762a1bSJed Brown /* currently only for simple PetscSection without fields or constraints */ 1020d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSectionReplicate_Private(MPI_Comm comm, PetscMPIInt rootrank, PetscSection sec0, PetscSection *secout) 1021d71ae5a4SJacob Faibussowitsch { 1022c4762a1bSJed Brown PetscSection sec; 1023c4762a1bSJed Brown PetscInt chart[2], p; 1024c4762a1bSJed Brown PetscInt *dofarr; 1025c4762a1bSJed Brown PetscMPIInt rank; 1026c4762a1bSJed Brown 1027c4762a1bSJed Brown PetscFunctionBegin; 10289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 102948a46eb9SPierre Jolivet if (rank == rootrank) PetscCall(PetscSectionGetChart(sec0, &chart[0], &chart[1])); 10309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(chart, 2, MPIU_INT, rootrank, comm)); 10319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(chart[1] - chart[0], &dofarr)); 1032c4762a1bSJed Brown if (rank == rootrank) { 103348a46eb9SPierre Jolivet for (p = chart[0]; p < chart[1]; p++) PetscCall(PetscSectionGetDof(sec0, p, &dofarr[p - chart[0]])); 1034c4762a1bSJed Brown } 10359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(dofarr, chart[1] - chart[0], MPIU_INT, rootrank, comm)); 10369566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, &sec)); 10379566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(sec, chart[0], chart[1])); 103848a46eb9SPierre Jolivet for (p = chart[0]; p < chart[1]; p++) PetscCall(PetscSectionSetDof(sec, p, dofarr[p - chart[0]])); 10399566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(sec)); 10409566063dSJacob Faibussowitsch PetscCall(PetscFree(dofarr)); 1041c4762a1bSJed Brown *secout = sec; 10423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1043c4762a1bSJed Brown } 1044c4762a1bSJed Brown 1045d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexExpandedVerticesCoordinatesToFaces_Private(DM ipdm, PortableBoundary bnd, IS *face_is) 1046d71ae5a4SJacob Faibussowitsch { 1047c4762a1bSJed Brown IS faces_expanded_is; 1048c4762a1bSJed Brown 1049c4762a1bSJed Brown PetscFunctionBegin; 10509566063dSJacob Faibussowitsch PetscCall(DMPlexFindVertices(ipdm, bnd->coordinates, 0.0, &faces_expanded_is)); 10519566063dSJacob Faibussowitsch PetscCall(DMPlexExpandedVerticesToFaces_Private(ipdm, faces_expanded_is, bnd->depth, bnd->sections, face_is)); 10529566063dSJacob Faibussowitsch PetscCall(ISDestroy(&faces_expanded_is)); 10533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1054c4762a1bSJed Brown } 1055c4762a1bSJed Brown 1056c4762a1bSJed Brown /* hack disabling DMPlexOrientInterface() call in DMPlexInterpolate() via -dm_plex_interpolate_orient_interfaces option */ 1057d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexSetOrientInterface_Private(DM dm, PetscBool enable) 1058d71ae5a4SJacob Faibussowitsch { 1059c4762a1bSJed Brown PetscOptions options = NULL; 1060c4762a1bSJed Brown const char *prefix = NULL; 1061c4762a1bSJed Brown const char opt[] = "-dm_plex_interpolate_orient_interfaces"; 1062c4762a1bSJed Brown char prefix_opt[512]; 1063c4762a1bSJed Brown PetscBool flg, set; 1064c4762a1bSJed Brown static PetscBool wasSetTrue = PETSC_FALSE; 1065c4762a1bSJed Brown 1066c4762a1bSJed Brown PetscFunctionBegin; 1067c4762a1bSJed Brown if (dm) { 10689566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 1069c4762a1bSJed Brown options = ((PetscObject)dm)->options; 1070c4762a1bSJed Brown } 1071c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(prefix_opt, "-", sizeof(prefix_opt))); 10729566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(prefix_opt, prefix, sizeof(prefix_opt))); 10739566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(prefix_opt, &opt[1], sizeof(prefix_opt))); 10749566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(options, prefix, opt, &flg, &set)); 1075c4762a1bSJed Brown if (!enable) { 1076c4762a1bSJed Brown if (set && flg) wasSetTrue = PETSC_TRUE; 10779566063dSJacob Faibussowitsch PetscCall(PetscOptionsSetValue(options, prefix_opt, "0")); 1078c4762a1bSJed Brown } else if (set && !flg) { 1079c4762a1bSJed Brown if (wasSetTrue) { 10809566063dSJacob Faibussowitsch PetscCall(PetscOptionsSetValue(options, prefix_opt, "1")); 1081c4762a1bSJed Brown } else { 1082c4762a1bSJed Brown /* default is PETSC_TRUE */ 10839566063dSJacob Faibussowitsch PetscCall(PetscOptionsClearValue(options, prefix_opt)); 1084c4762a1bSJed Brown } 1085c4762a1bSJed Brown wasSetTrue = PETSC_FALSE; 1086c4762a1bSJed Brown } 108776bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 10889566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(options, prefix, opt, &flg, &set)); 10891dca8a05SBarry Smith PetscCheck(!set || flg == enable, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "PetscOptionsSetValue did not have the desired effect"); 1090c4762a1bSJed Brown } 10913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1092c4762a1bSJed Brown } 1093c4762a1bSJed Brown 1094c4762a1bSJed Brown /* get coordinate description of the whole-domain boundary */ 1095d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexGetExpandedBoundary_Private(DM dm, PortableBoundary *boundary) 1096d71ae5a4SJacob Faibussowitsch { 1097c4762a1bSJed Brown PortableBoundary bnd0, bnd; 1098c4762a1bSJed Brown MPI_Comm comm; 1099c4762a1bSJed Brown DM idm; 1100c4762a1bSJed Brown DMLabel label; 1101c4762a1bSJed Brown PetscInt d; 1102c4762a1bSJed Brown const char boundaryName[] = "DMPlexDistributeInterpolateMarkInterface_boundary"; 1103c4762a1bSJed Brown IS boundary_is; 1104c4762a1bSJed Brown IS *boundary_expanded_iss; 1105c4762a1bSJed Brown PetscMPIInt rootrank = 0; 1106c4762a1bSJed Brown PetscMPIInt rank, size; 1107c4762a1bSJed Brown PetscInt value = 1; 1108c4762a1bSJed Brown DMPlexInterpolatedFlag intp; 1109c4762a1bSJed Brown PetscBool flg; 1110c4762a1bSJed Brown 1111c4762a1bSJed Brown PetscFunctionBegin; 11129566063dSJacob Faibussowitsch PetscCall(PetscNew(&bnd)); 11139566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 11149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 11159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 11169566063dSJacob Faibussowitsch PetscCall(DMPlexIsDistributed(dm, &flg)); 111728b400f6SJacob Faibussowitsch PetscCheck(!flg, comm, PETSC_ERR_ARG_WRONG, "serial DM (all points on one rank) needed"); 1118c4762a1bSJed Brown 1119c4762a1bSJed Brown /* interpolate serial DM if not yet interpolated */ 11209566063dSJacob Faibussowitsch PetscCall(DMPlexIsInterpolatedCollective(dm, &intp)); 1121c4762a1bSJed Brown if (intp == DMPLEX_INTERPOLATED_FULL) { 1122c4762a1bSJed Brown idm = dm; 11239566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 1124c4762a1bSJed Brown } else { 11259566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(dm, &idm)); 11269566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(idm, NULL, "-idm_view")); 1127c4762a1bSJed Brown } 1128c4762a1bSJed Brown 1129c4762a1bSJed Brown /* mark whole-domain boundary of the serial DM */ 11309566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, boundaryName, &label)); 11319566063dSJacob Faibussowitsch PetscCall(DMAddLabel(idm, label)); 11329566063dSJacob Faibussowitsch PetscCall(DMPlexMarkBoundaryFaces(idm, value, label)); 11339566063dSJacob Faibussowitsch PetscCall(DMLabelViewFromOptionsOnComm_Private(label, "-idm_boundary_view", comm)); 11349566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, value, &boundary_is)); 1135c4762a1bSJed Brown 1136c4762a1bSJed Brown /* translate to coordinates */ 11379566063dSJacob Faibussowitsch PetscCall(PetscNew(&bnd0)); 11389566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalSetUp(idm)); 1139c4762a1bSJed Brown if (rank == rootrank) { 11409566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeRecursive(idm, boundary_is, &bnd0->depth, &boundary_expanded_iss, &bnd0->sections)); 11419566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalTuple(dm, boundary_expanded_iss[0], NULL, &bnd0->coordinates)); 1142c4762a1bSJed Brown /* self-check */ 1143c4762a1bSJed Brown { 1144c4762a1bSJed Brown IS is0; 11459566063dSJacob Faibussowitsch PetscCall(DMPlexExpandedVerticesCoordinatesToFaces_Private(idm, bnd0, &is0)); 11469566063dSJacob Faibussowitsch PetscCall(ISEqual(is0, boundary_is, &flg)); 114728b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "DMPlexExpandedVerticesCoordinatesToFaces_Private produced a wrong IS"); 11489566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is0)); 1149c4762a1bSJed Brown } 1150c4762a1bSJed Brown } else { 11519566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, 0, &bnd0->coordinates)); 1152c4762a1bSJed Brown } 1153c4762a1bSJed Brown 1154c4762a1bSJed Brown { 1155c4762a1bSJed Brown Vec tmp; 1156c4762a1bSJed Brown VecScatter sc; 1157c4762a1bSJed Brown IS xis; 1158c4762a1bSJed Brown PetscInt n; 1159c4762a1bSJed Brown 1160c4762a1bSJed Brown /* just convert seq vectors to mpi vector */ 11619566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(bnd0->coordinates, &n)); 11629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&n, 1, MPIU_INT, rootrank, comm)); 1163c4762a1bSJed Brown if (rank == rootrank) { 11649566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(comm, n, n, &tmp)); 1165c4762a1bSJed Brown } else { 11669566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(comm, 0, n, &tmp)); 1167c4762a1bSJed Brown } 11689566063dSJacob Faibussowitsch PetscCall(VecCopy(bnd0->coordinates, tmp)); 11699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnd0->coordinates)); 1170c4762a1bSJed Brown bnd0->coordinates = tmp; 1171c4762a1bSJed Brown 1172c4762a1bSJed Brown /* replicate coordinates from root rank to all ranks */ 11739566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(comm, n, n * size, &bnd->coordinates)); 11749566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, n, 0, 1, &xis)); 11759566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(bnd0->coordinates, xis, bnd->coordinates, NULL, &sc)); 11769566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sc, bnd0->coordinates, bnd->coordinates, INSERT_VALUES, SCATTER_FORWARD)); 11779566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sc, bnd0->coordinates, bnd->coordinates, INSERT_VALUES, SCATTER_FORWARD)); 11789566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&sc)); 11799566063dSJacob Faibussowitsch PetscCall(ISDestroy(&xis)); 1180c4762a1bSJed Brown } 1181c4762a1bSJed Brown bnd->depth = bnd0->depth; 11829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&bnd->depth, 1, MPIU_INT, rootrank, comm)); 11839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bnd->depth, &bnd->sections)); 118448a46eb9SPierre Jolivet for (d = 0; d < bnd->depth; d++) PetscCall(PetscSectionReplicate_Private(comm, rootrank, (rank == rootrank) ? bnd0->sections[d] : NULL, &bnd->sections[d])); 1185c4762a1bSJed Brown 118648a46eb9SPierre Jolivet if (rank == rootrank) PetscCall(DMPlexRestoreConeRecursive(idm, boundary_is, &bnd0->depth, &boundary_expanded_iss, &bnd0->sections)); 11879566063dSJacob Faibussowitsch PetscCall(PortableBoundaryDestroy(&bnd0)); 11889566063dSJacob Faibussowitsch PetscCall(DMRemoveLabelBySelf(idm, &label, PETSC_TRUE)); 11899566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&label)); 11909566063dSJacob Faibussowitsch PetscCall(ISDestroy(&boundary_is)); 11919566063dSJacob Faibussowitsch PetscCall(DMDestroy(&idm)); 1192c4762a1bSJed Brown *boundary = bnd; 11933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1194c4762a1bSJed Brown } 1195c4762a1bSJed Brown 1196c4762a1bSJed Brown /* get faces of inter-partition interface */ 1197d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexGetInterfaceFaces_Private(DM ipdm, IS boundary_faces_is, IS *interface_faces_is) 1198d71ae5a4SJacob Faibussowitsch { 1199c4762a1bSJed Brown MPI_Comm comm; 1200c4762a1bSJed Brown DMLabel label; 1201c4762a1bSJed Brown IS part_boundary_faces_is; 1202c4762a1bSJed Brown const char partBoundaryName[] = "DMPlexDistributeInterpolateMarkInterface_partBoundary"; 1203c4762a1bSJed Brown PetscInt value = 1; 1204c4762a1bSJed Brown DMPlexInterpolatedFlag intp; 1205c4762a1bSJed Brown 1206c4762a1bSJed Brown PetscFunctionBegin; 12079566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ipdm, &comm)); 12089566063dSJacob Faibussowitsch PetscCall(DMPlexIsInterpolatedCollective(ipdm, &intp)); 120908401ef6SPierre Jolivet PetscCheck(intp == DMPLEX_INTERPOLATED_FULL, comm, PETSC_ERR_ARG_WRONG, "only for fully interpolated DMPlex"); 1210c4762a1bSJed Brown 1211c4762a1bSJed Brown /* get ipdm partition boundary (partBoundary) */ 1212429fa399SMatthew G. Knepley { 1213429fa399SMatthew G. Knepley PetscSF sf; 1214429fa399SMatthew G. Knepley 12159566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, partBoundaryName, &label)); 12169566063dSJacob Faibussowitsch PetscCall(DMAddLabel(ipdm, label)); 1217429fa399SMatthew G. Knepley PetscCall(DMGetPointSF(ipdm, &sf)); 1218429fa399SMatthew G. Knepley PetscCall(DMSetPointSF(ipdm, NULL)); 12199566063dSJacob Faibussowitsch PetscCall(DMPlexMarkBoundaryFaces(ipdm, value, label)); 1220429fa399SMatthew G. Knepley PetscCall(DMSetPointSF(ipdm, sf)); 12219566063dSJacob Faibussowitsch PetscCall(DMLabelViewFromOptionsOnComm_Private(label, "-ipdm_part_boundary_view", comm)); 12229566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumISOnComm_Private(label, value, comm, &part_boundary_faces_is)); 12239566063dSJacob Faibussowitsch PetscCall(DMRemoveLabelBySelf(ipdm, &label, PETSC_TRUE)); 12249566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&label)); 1225429fa399SMatthew G. Knepley } 1226c4762a1bSJed Brown 1227c4762a1bSJed Brown /* remove ipdm whole-domain boundary (boundary_faces_is) from ipdm partition boundary (part_boundary_faces_is), resulting just in inter-partition interface */ 12289566063dSJacob Faibussowitsch PetscCall(ISDifference(part_boundary_faces_is, boundary_faces_is, interface_faces_is)); 12299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&part_boundary_faces_is)); 12303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1231c4762a1bSJed Brown } 1232c4762a1bSJed Brown 1233c4762a1bSJed Brown /* compute inter-partition interface including edges and vertices */ 1234d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexComputeCompleteInterface_Private(DM ipdm, IS interface_faces_is, IS *interface_is) 1235d71ae5a4SJacob Faibussowitsch { 1236c4762a1bSJed Brown DMLabel label; 1237c4762a1bSJed Brown PetscInt value = 1; 1238c4762a1bSJed Brown const char interfaceName[] = "DMPlexDistributeInterpolateMarkInterface_interface"; 1239c4762a1bSJed Brown DMPlexInterpolatedFlag intp; 1240c4762a1bSJed Brown MPI_Comm comm; 1241c4762a1bSJed Brown 1242c4762a1bSJed Brown PetscFunctionBegin; 12439566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ipdm, &comm)); 12449566063dSJacob Faibussowitsch PetscCall(DMPlexIsInterpolatedCollective(ipdm, &intp)); 124508401ef6SPierre Jolivet PetscCheck(intp == DMPLEX_INTERPOLATED_FULL, comm, PETSC_ERR_ARG_WRONG, "only for fully interpolated DMPlex"); 1246c4762a1bSJed Brown 12479566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, interfaceName, &label)); 12489566063dSJacob Faibussowitsch PetscCall(DMAddLabel(ipdm, label)); 12499566063dSJacob Faibussowitsch PetscCall(DMLabelSetStratumIS(label, value, interface_faces_is)); 12509566063dSJacob Faibussowitsch PetscCall(DMLabelViewFromOptionsOnComm_Private(label, "-interface_faces_view", comm)); 12519566063dSJacob Faibussowitsch PetscCall(DMPlexLabelComplete(ipdm, label)); 12529566063dSJacob Faibussowitsch PetscCall(DMLabelViewFromOptionsOnComm_Private(label, "-interface_view", comm)); 12539566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumISOnComm_Private(label, value, comm, interface_is)); 12549566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*interface_is, "interface_is")); 12559566063dSJacob Faibussowitsch PetscCall(ISViewFromOptions(*interface_is, NULL, "-interface_is_view")); 12569566063dSJacob Faibussowitsch PetscCall(DMRemoveLabelBySelf(ipdm, &label, PETSC_TRUE)); 12579566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&label)); 12583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1259c4762a1bSJed Brown } 1260c4762a1bSJed Brown 1261d71ae5a4SJacob Faibussowitsch static PetscErrorCode PointSFGetOutwardInterfacePoints(PetscSF sf, IS *is) 1262d71ae5a4SJacob Faibussowitsch { 1263c4762a1bSJed Brown PetscInt n; 1264c4762a1bSJed Brown const PetscInt *arr; 1265c4762a1bSJed Brown 1266c4762a1bSJed Brown PetscFunctionBegin; 12679566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, NULL, &n, &arr, NULL)); 12689566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)sf), n, arr, PETSC_USE_POINTER, is)); 12693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1270c4762a1bSJed Brown } 1271c4762a1bSJed Brown 1272d71ae5a4SJacob Faibussowitsch static PetscErrorCode PointSFGetInwardInterfacePoints(PetscSF sf, IS *is) 1273d71ae5a4SJacob Faibussowitsch { 1274c4762a1bSJed Brown PetscInt n; 1275c4762a1bSJed Brown const PetscInt *rootdegree; 1276c4762a1bSJed Brown PetscInt *arr; 1277c4762a1bSJed Brown 1278c4762a1bSJed Brown PetscFunctionBegin; 12799566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sf)); 12809566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sf, &rootdegree)); 12819566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sf, &rootdegree)); 12829566063dSJacob Faibussowitsch PetscCall(PetscSFComputeMultiRootOriginalNumbering(sf, rootdegree, &n, &arr)); 12839566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)sf), n, arr, PETSC_OWN_POINTER, is)); 12843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1285c4762a1bSJed Brown } 1286c4762a1bSJed Brown 1287d71ae5a4SJacob Faibussowitsch static PetscErrorCode PointSFGetInterfacePoints_Private(PetscSF pointSF, IS *is) 1288d71ae5a4SJacob Faibussowitsch { 1289c4762a1bSJed Brown IS pointSF_out_is, pointSF_in_is; 1290c4762a1bSJed Brown 1291c4762a1bSJed Brown PetscFunctionBegin; 12929566063dSJacob Faibussowitsch PetscCall(PointSFGetOutwardInterfacePoints(pointSF, &pointSF_out_is)); 12939566063dSJacob Faibussowitsch PetscCall(PointSFGetInwardInterfacePoints(pointSF, &pointSF_in_is)); 12949566063dSJacob Faibussowitsch PetscCall(ISExpand(pointSF_out_is, pointSF_in_is, is)); 12959566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointSF_out_is)); 12969566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointSF_in_is)); 12973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1298c4762a1bSJed Brown } 1299c4762a1bSJed Brown 13005f80ce2aSJacob Faibussowitsch #define CHKERRMY(ierr) PetscCheck(!ierr, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PointSF is wrong. Unable to show details!") 1301c4762a1bSJed Brown 1302d71ae5a4SJacob Faibussowitsch static PetscErrorCode ViewPointsWithType_Internal(DM dm, IS pointsIS, PetscViewer v) 1303d71ae5a4SJacob Faibussowitsch { 1304c4762a1bSJed Brown DMLabel label; 1305c4762a1bSJed Brown PetscSection coordsSection; 1306c4762a1bSJed Brown Vec coordsVec; 1307c4762a1bSJed Brown PetscScalar *coordsScalar; 1308c4762a1bSJed Brown PetscInt coneSize, depth, dim, i, p, npoints; 1309c4762a1bSJed Brown const PetscInt *points; 1310c4762a1bSJed Brown 1311c4762a1bSJed Brown PetscFunctionBegin; 13129566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 13139566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordsSection)); 13149566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordsVec)); 13159566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordsVec, &coordsScalar)); 13169566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pointsIS, &npoints)); 13179566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pointsIS, &points)); 13189566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(dm, &label)); 13199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(v)); 1320c4762a1bSJed Brown for (i = 0; i < npoints; i++) { 1321c4762a1bSJed Brown p = points[i]; 13229566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, p, &depth)); 1323c4762a1bSJed Brown if (!depth) { 1324d3e1f4ccSVaclav Hapla PetscInt n, o; 1325c4762a1bSJed Brown char coordstr[128]; 1326c4762a1bSJed Brown 13279566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(coordsSection, p, &n)); 13289566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordsSection, p, &o)); 13299566063dSJacob Faibussowitsch PetscCall(coord2str(coordstr, sizeof(coordstr), n, &coordsScalar[o], 1.0)); 133063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(v, "vertex %" PetscInt_FMT " w/ coordinates %s\n", p, coordstr)); 1331c4762a1bSJed Brown } else { 1332c4762a1bSJed Brown char entityType[16]; 1333c4762a1bSJed Brown 1334c4762a1bSJed Brown switch (depth) { 1335d71ae5a4SJacob Faibussowitsch case 1: 1336c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(entityType, "edge", sizeof(entityType))); 1337d71ae5a4SJacob Faibussowitsch break; 1338d71ae5a4SJacob Faibussowitsch case 2: 1339c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(entityType, "face", sizeof(entityType))); 1340d71ae5a4SJacob Faibussowitsch break; 1341d71ae5a4SJacob Faibussowitsch case 3: 1342c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(entityType, "cell", sizeof(entityType))); 1343d71ae5a4SJacob Faibussowitsch break; 1344d71ae5a4SJacob Faibussowitsch default: 1345d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Only for depth <= 3"); 1346c4762a1bSJed Brown } 134748a46eb9SPierre Jolivet if (depth == dim && dim < 3) PetscCall(PetscStrlcat(entityType, " (cell)", sizeof(entityType))); 134863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(v, "%s %" PetscInt_FMT "\n", entityType, p)); 1349c4762a1bSJed Brown } 13509566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, p, &coneSize)); 1351c4762a1bSJed Brown if (coneSize) { 1352c4762a1bSJed Brown const PetscInt *cone; 1353c4762a1bSJed Brown IS coneIS; 1354c4762a1bSJed Brown 13559566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, p, &cone)); 13569566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, coneSize, cone, PETSC_USE_POINTER, &coneIS)); 13579566063dSJacob Faibussowitsch PetscCall(ViewPointsWithType_Internal(dm, coneIS, v)); 13589566063dSJacob Faibussowitsch PetscCall(ISDestroy(&coneIS)); 1359c4762a1bSJed Brown } 1360c4762a1bSJed Brown } 13619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(v)); 13629566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordsVec, &coordsScalar)); 13639566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pointsIS, &points)); 13643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1365c4762a1bSJed Brown } 1366c4762a1bSJed Brown 1367d71ae5a4SJacob Faibussowitsch static PetscErrorCode ViewPointsWithType(DM dm, IS points, PetscViewer v) 1368d71ae5a4SJacob Faibussowitsch { 1369c4762a1bSJed Brown PetscBool flg; 1370c4762a1bSJed Brown PetscInt npoints; 1371c4762a1bSJed Brown PetscMPIInt rank; 1372c4762a1bSJed Brown 1373c4762a1bSJed Brown PetscFunctionBegin; 13749566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &flg)); 137528b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Only for ASCII viewer"); 13769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)v), &rank)); 13779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(v)); 13789566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(points, &npoints)); 1379c4762a1bSJed Brown if (npoints) { 13809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d] --------\n", rank)); 13819566063dSJacob Faibussowitsch PetscCall(ViewPointsWithType_Internal(dm, points, v)); 1382c4762a1bSJed Brown } 13839566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(v)); 13849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(v)); 13853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1386c4762a1bSJed Brown } 1387c4762a1bSJed Brown 1388d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexComparePointSFWithInterface_Private(DM ipdm, IS interface_is) 1389d71ae5a4SJacob Faibussowitsch { 1390c4762a1bSJed Brown PetscSF pointsf; 1391c4762a1bSJed Brown IS pointsf_is; 1392c4762a1bSJed Brown PetscBool flg; 1393c4762a1bSJed Brown MPI_Comm comm; 1394c4762a1bSJed Brown PetscMPIInt size; 1395c4762a1bSJed Brown 1396c4762a1bSJed Brown PetscFunctionBegin; 13979566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ipdm, &comm)); 13989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 13999566063dSJacob Faibussowitsch PetscCall(DMGetPointSF(ipdm, &pointsf)); 1400c4762a1bSJed Brown if (pointsf) { 1401c4762a1bSJed Brown PetscInt nroots; 14029566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(pointsf, &nroots, NULL, NULL, NULL)); 1403c4762a1bSJed Brown if (nroots < 0) pointsf = NULL; /* uninitialized SF */ 1404c4762a1bSJed Brown } 1405c4762a1bSJed Brown if (!pointsf) { 1406c4762a1bSJed Brown PetscInt N = 0; 14079566063dSJacob Faibussowitsch if (interface_is) PetscCall(ISGetSize(interface_is, &N)); 140828b400f6SJacob Faibussowitsch PetscCheck(!N, comm, PETSC_ERR_PLIB, "interface_is should be NULL or empty for PointSF being NULL"); 14093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1410c4762a1bSJed Brown } 1411c4762a1bSJed Brown 1412c4762a1bSJed Brown /* get PointSF points as IS pointsf_is */ 14139566063dSJacob Faibussowitsch PetscCall(PointSFGetInterfacePoints_Private(pointsf, &pointsf_is)); 1414c4762a1bSJed Brown 1415c4762a1bSJed Brown /* compare pointsf_is with interface_is */ 14169566063dSJacob Faibussowitsch PetscCall(ISEqual(interface_is, pointsf_is, &flg)); 14179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &flg, 1, MPIU_BOOL, MPI_LAND, comm)); 1418c4762a1bSJed Brown if (!flg) { 1419c4762a1bSJed Brown IS pointsf_extra_is, pointsf_missing_is; 1420c4762a1bSJed Brown PetscViewer errv = PETSC_VIEWER_STDERR_(comm); 14215f80ce2aSJacob Faibussowitsch CHKERRMY(ISDifference(interface_is, pointsf_is, &pointsf_missing_is)); 14225f80ce2aSJacob Faibussowitsch CHKERRMY(ISDifference(pointsf_is, interface_is, &pointsf_extra_is)); 14235f80ce2aSJacob Faibussowitsch CHKERRMY(PetscViewerASCIIPrintf(errv, "Points missing in PointSF:\n")); 14245f80ce2aSJacob Faibussowitsch CHKERRMY(ViewPointsWithType(ipdm, pointsf_missing_is, errv)); 14255f80ce2aSJacob Faibussowitsch CHKERRMY(PetscViewerASCIIPrintf(errv, "Extra points in PointSF:\n")); 14265f80ce2aSJacob Faibussowitsch CHKERRMY(ViewPointsWithType(ipdm, pointsf_extra_is, errv)); 14275f80ce2aSJacob Faibussowitsch CHKERRMY(ISDestroy(&pointsf_extra_is)); 14285f80ce2aSJacob Faibussowitsch CHKERRMY(ISDestroy(&pointsf_missing_is)); 1429c4762a1bSJed Brown SETERRQ(comm, PETSC_ERR_PLIB, "PointSF is wrong! See details above."); 1430c4762a1bSJed Brown } 14319566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointsf_is)); 14323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1433c4762a1bSJed Brown } 1434c4762a1bSJed Brown 1435c4762a1bSJed Brown /* remove faces & edges from label, leave just vertices */ 1436d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexISFilterVertices_Private(DM dm, IS points) 1437d71ae5a4SJacob Faibussowitsch { 1438c4762a1bSJed Brown PetscInt vStart, vEnd; 1439c4762a1bSJed Brown MPI_Comm comm; 1440c4762a1bSJed Brown 1441c4762a1bSJed Brown PetscFunctionBegin; 14429566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 14439566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 14449566063dSJacob Faibussowitsch PetscCall(ISGeneralFilter(points, vStart, vEnd)); 14453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1446c4762a1bSJed Brown } 1447c4762a1bSJed Brown 1448c4762a1bSJed Brown /* 14492e58f0efSBarry Smith DMPlexCheckPointSFHeavy - Thoroughly test that the PointSF after parallel DMPlexInterpolate() includes exactly all interface points. 1450c4762a1bSJed Brown 1451c4762a1bSJed Brown Collective 1452c4762a1bSJed Brown 1453*2fe279fdSBarry Smith Input Parameter: 1454c4762a1bSJed Brown . dm - The DMPlex object 1455c4762a1bSJed Brown 1456c4762a1bSJed Brown Notes: 1457c4762a1bSJed Brown The input DMPlex must be serial (one partition has all points, the other partitions have no points). 1458c4762a1bSJed 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). 1459c4762a1bSJed Brown This is mainly intended for debugging/testing purposes. 1460c4762a1bSJed Brown 1461c4762a1bSJed Brown Algorithm: 1462c4762a1bSJed Brown 1. boundary faces of the serial version of the whole mesh are found using DMPlexMarkBoundaryFaces() 1463c4762a1bSJed 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 1464c4762a1bSJed Brown 3. the mesh is distributed or loaded in parallel 1465c4762a1bSJed Brown 4. boundary faces of the distributed mesh are reconstructed from PortableBoundary using DMPlexFindVertices() 1466c4762a1bSJed Brown 5. partition boundary faces of the parallel mesh are found using DMPlexMarkBoundaryFaces() 1467c4762a1bSJed Brown 6. partition interfaces are computed as set difference of partition boundary faces minus the reconstructed boundary 1468c4762a1bSJed 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 1469c4762a1bSJed Brown 1470c4762a1bSJed Brown Level: developer 1471c4762a1bSJed Brown 1472c4762a1bSJed Brown .seealso: DMGetPointSF(), DMPlexCheckSymmetry(), DMPlexCheckSkeleton(), DMPlexCheckFaces() 1473c4762a1bSJed Brown */ 1474d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCheckPointSFHeavy(DM dm, PortableBoundary bnd) 1475d71ae5a4SJacob Faibussowitsch { 1476c4762a1bSJed Brown DM ipdm = NULL; 1477c4762a1bSJed Brown IS boundary_faces_is, interface_faces_is, interface_is; 1478c4762a1bSJed Brown DMPlexInterpolatedFlag intp; 1479c4762a1bSJed Brown MPI_Comm comm; 1480c4762a1bSJed Brown 1481c4762a1bSJed Brown PetscFunctionBegin; 14829566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 1483c4762a1bSJed Brown 14849566063dSJacob Faibussowitsch PetscCall(DMPlexIsInterpolatedCollective(dm, &intp)); 1485c4762a1bSJed Brown if (intp == DMPLEX_INTERPOLATED_FULL) { 1486c4762a1bSJed Brown ipdm = dm; 1487c4762a1bSJed Brown } else { 1488c4762a1bSJed Brown /* create temporary interpolated DM if input DM is not interpolated */ 14899566063dSJacob Faibussowitsch PetscCall(DMPlexSetOrientInterface_Private(dm, PETSC_FALSE)); 14909566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(dm, &ipdm)); /* with DMPlexOrientInterface_Internal() call skipped so that PointSF issues are left to DMPlexComparePointSFWithInterface_Private() below */ 14919566063dSJacob Faibussowitsch PetscCall(DMPlexSetOrientInterface_Private(dm, PETSC_TRUE)); 1492c4762a1bSJed Brown } 14939566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(ipdm, NULL, "-ipdm_view")); 1494c4762a1bSJed Brown 1495c4762a1bSJed Brown /* recover ipdm whole-domain boundary faces from the expanded vertices coordinates */ 14969566063dSJacob Faibussowitsch PetscCall(DMPlexExpandedVerticesCoordinatesToFaces_Private(ipdm, bnd, &boundary_faces_is)); 1497c4762a1bSJed Brown /* get inter-partition interface faces (interface_faces_is)*/ 14989566063dSJacob Faibussowitsch PetscCall(DMPlexGetInterfaceFaces_Private(ipdm, boundary_faces_is, &interface_faces_is)); 1499c4762a1bSJed Brown /* compute inter-partition interface including edges and vertices (interface_is) */ 15009566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCompleteInterface_Private(ipdm, interface_faces_is, &interface_is)); 1501c4762a1bSJed Brown /* destroy immediate ISs */ 15029566063dSJacob Faibussowitsch PetscCall(ISDestroy(&boundary_faces_is)); 15039566063dSJacob Faibussowitsch PetscCall(ISDestroy(&interface_faces_is)); 1504c4762a1bSJed Brown 1505c4762a1bSJed Brown /* for uninterpolated case, keep just vertices in interface */ 1506c4762a1bSJed Brown if (!intp) { 15079566063dSJacob Faibussowitsch PetscCall(DMPlexISFilterVertices_Private(ipdm, interface_is)); 15089566063dSJacob Faibussowitsch PetscCall(DMDestroy(&ipdm)); 1509c4762a1bSJed Brown } 1510c4762a1bSJed Brown 1511c4762a1bSJed Brown /* compare PointSF with the boundary reconstructed from coordinates */ 15129566063dSJacob Faibussowitsch PetscCall(DMPlexComparePointSFWithInterface_Private(dm, interface_is)); 15139566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "DMPlexCheckPointSFHeavy PASSED\n")); 15149566063dSJacob Faibussowitsch PetscCall(ISDestroy(&interface_is)); 15153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1516c4762a1bSJed Brown } 1517c4762a1bSJed Brown 1518d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 1519d71ae5a4SJacob Faibussowitsch { 1520c4762a1bSJed Brown DM dm; 1521c4762a1bSJed Brown AppCtx user; 1522c4762a1bSJed Brown 1523327415f7SBarry Smith PetscFunctionBeginUser; 15249566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 15259566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("create", &stage[0])); 15269566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("distribute", &stage[1])); 15279566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("interpolate", &stage[2])); 15289566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 15299566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 153048a46eb9SPierre Jolivet if (user.nPointsToExpand) PetscCall(TestExpandPoints(dm, &user)); 1531c4762a1bSJed Brown if (user.ncoords) { 1532d3e1f4ccSVaclav Hapla Vec coords; 1533d3e1f4ccSVaclav Hapla 15349566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, user.ncoords, user.ncoords, user.coords, &coords)); 15359566063dSJacob Faibussowitsch PetscCall(ViewVerticesFromCoords(dm, coords, user.coordsTol, PETSC_VIEWER_STDOUT_WORLD)); 15369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coords)); 1537c4762a1bSJed Brown } 15389566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 15399566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 1540b122ec5aSJacob Faibussowitsch return 0; 1541c4762a1bSJed Brown } 1542c4762a1bSJed Brown 1543c4762a1bSJed Brown /*TEST 1544c4762a1bSJed Brown 1545c4762a1bSJed Brown testset: 1546c4762a1bSJed Brown nsize: 2 1547c4762a1bSJed Brown args: -dm_view ascii::ascii_info_detail 1548c4762a1bSJed Brown args: -dm_plex_check_all 1549c4762a1bSJed Brown test: 1550c4762a1bSJed Brown suffix: 1_tri_dist0 1551c4762a1bSJed Brown args: -distribute 0 -interpolate {{none create}separate output} 1552c4762a1bSJed Brown test: 1553c4762a1bSJed Brown suffix: 1_tri_dist1 1554c4762a1bSJed Brown args: -distribute 1 -interpolate {{none create after_distribute}separate output} 1555c4762a1bSJed Brown test: 1556c4762a1bSJed Brown suffix: 1_quad_dist0 1557c4762a1bSJed Brown args: -cell_simplex 0 -distribute 0 -interpolate {{none create}separate output} 1558c4762a1bSJed Brown test: 1559c4762a1bSJed Brown suffix: 1_quad_dist1 1560c4762a1bSJed Brown args: -cell_simplex 0 -distribute 1 -interpolate {{none create after_distribute}separate output} 1561c4762a1bSJed Brown test: 1562c4762a1bSJed Brown suffix: 1_1d_dist1 1563c4762a1bSJed Brown args: -dim 1 -distribute 1 1564c4762a1bSJed Brown 1565c4762a1bSJed Brown testset: 1566c4762a1bSJed Brown nsize: 3 1567c4762a1bSJed Brown args: -testnum 1 -interpolate create 1568c4762a1bSJed Brown args: -dm_plex_check_all 1569c4762a1bSJed Brown test: 1570c4762a1bSJed Brown suffix: 2 1571c4762a1bSJed Brown args: -dm_view ascii::ascii_info_detail 1572c4762a1bSJed Brown test: 1573c4762a1bSJed Brown suffix: 2a 1574c4762a1bSJed Brown args: -dm_plex_check_cones_conform_on_interfaces_verbose 1575c4762a1bSJed Brown test: 1576c4762a1bSJed Brown suffix: 2b 1577c4762a1bSJed Brown args: -test_expand_points 0,1,2,5,6 1578c4762a1bSJed Brown test: 1579c4762a1bSJed Brown suffix: 2c 1580c4762a1bSJed Brown args: -test_expand_points 0,1,2,5,6 -test_expand_points_empty 1581c4762a1bSJed Brown 1582c4762a1bSJed Brown testset: 1583c4762a1bSJed Brown # the same as 1% for 3D 1584c4762a1bSJed Brown nsize: 2 1585c4762a1bSJed Brown args: -dim 3 -dm_view ascii::ascii_info_detail 1586c4762a1bSJed Brown args: -dm_plex_check_all 1587c4762a1bSJed Brown test: 1588c4762a1bSJed Brown suffix: 4_tet_dist0 1589c4762a1bSJed Brown args: -distribute 0 -interpolate {{none create}separate output} 1590c4762a1bSJed Brown test: 1591c4762a1bSJed Brown suffix: 4_tet_dist1 1592c4762a1bSJed Brown args: -distribute 1 -interpolate {{none create after_distribute}separate output} 1593c4762a1bSJed Brown test: 1594c4762a1bSJed Brown suffix: 4_hex_dist0 1595c4762a1bSJed Brown args: -cell_simplex 0 -distribute 0 -interpolate {{none create}separate output} 1596c4762a1bSJed Brown test: 1597c4762a1bSJed Brown suffix: 4_hex_dist1 1598c4762a1bSJed Brown args: -cell_simplex 0 -distribute 1 -interpolate {{none create after_distribute}separate output} 1599c4762a1bSJed Brown 1600c4762a1bSJed Brown test: 1601c4762a1bSJed Brown # the same as 4_tet_dist0 but test different initial orientations 1602c4762a1bSJed Brown suffix: 4_tet_test_orient 1603c4762a1bSJed Brown nsize: 2 1604c4762a1bSJed Brown args: -dim 3 -distribute 0 1605c4762a1bSJed Brown args: -dm_plex_check_all 1606c4762a1bSJed Brown args: -rotate_interface_0 {{0 1 2 11 12 13}} 1607c4762a1bSJed Brown args: -rotate_interface_1 {{0 1 2 11 12 13}} 1608c4762a1bSJed Brown 1609c4762a1bSJed Brown testset: 1610c4762a1bSJed Brown requires: exodusii 1611c4762a1bSJed Brown args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/TwoQuads.exo 1612c4762a1bSJed Brown args: -dm_view ascii::ascii_info_detail 1613c4762a1bSJed Brown args: -dm_plex_check_all 1614c4762a1bSJed Brown args: -custom_view 1615c4762a1bSJed Brown test: 1616c4762a1bSJed Brown suffix: 5_seq 1617c4762a1bSJed Brown nsize: 1 1618c4762a1bSJed Brown args: -distribute 0 -interpolate {{none create}separate output} 1619c4762a1bSJed Brown test: 162030602db0SMatthew G. Knepley # Detail viewing in a non-distributed mesh is broken because the DMLabelView() is collective, but the label is not shared 1621c4762a1bSJed Brown suffix: 5_dist0 1622c4762a1bSJed Brown nsize: 2 162330602db0SMatthew G. Knepley args: -distribute 0 -interpolate {{none create}separate output} -dm_view 1624c4762a1bSJed Brown test: 1625c4762a1bSJed Brown suffix: 5_dist1 1626c4762a1bSJed Brown nsize: 2 1627c4762a1bSJed Brown args: -distribute 1 -interpolate {{none create after_distribute}separate output} 1628c4762a1bSJed Brown 1629c4762a1bSJed Brown testset: 1630c4762a1bSJed Brown nsize: {{1 2 4}} 1631c4762a1bSJed Brown args: -use_generator 1632c4762a1bSJed Brown args: -dm_plex_check_all 1633c4762a1bSJed Brown args: -distribute -interpolate none 1634c4762a1bSJed Brown test: 1635c4762a1bSJed Brown suffix: 6_tri 1636c4762a1bSJed Brown requires: triangle 1637c0517cd5SMatthew G. Knepley args: -faces {{2,2 1,3 7,4}} -cell_simplex 1 -dm_generator triangle 1638c4762a1bSJed Brown test: 1639c4762a1bSJed Brown suffix: 6_quad 1640c4762a1bSJed Brown args: -faces {{2,2 1,3 7,4}} -cell_simplex 0 1641c4762a1bSJed Brown test: 1642c4762a1bSJed Brown suffix: 6_tet 1643c4762a1bSJed Brown requires: ctetgen 1644c0517cd5SMatthew G. Knepley args: -faces {{2,2,2 1,3,5 3,4,7}} -cell_simplex 1 -dm_generator ctetgen 1645c4762a1bSJed Brown test: 1646c4762a1bSJed Brown suffix: 6_hex 1647c4762a1bSJed Brown args: -faces {{2,2,2 1,3,5 3,4,7}} -cell_simplex 0 1648c4762a1bSJed Brown testset: 1649c4762a1bSJed Brown nsize: {{1 2 4}} 1650c4762a1bSJed Brown args: -use_generator 1651c4762a1bSJed Brown args: -dm_plex_check_all 1652c4762a1bSJed Brown args: -distribute -interpolate create 1653c4762a1bSJed Brown test: 1654c4762a1bSJed Brown suffix: 6_int_tri 1655c4762a1bSJed Brown requires: triangle 1656c0517cd5SMatthew G. Knepley args: -faces {{2,2 1,3 7,4}} -cell_simplex 1 -dm_generator triangle 1657c4762a1bSJed Brown test: 1658c4762a1bSJed Brown suffix: 6_int_quad 1659c4762a1bSJed Brown args: -faces {{2,2 1,3 7,4}} -cell_simplex 0 1660c4762a1bSJed Brown test: 1661c4762a1bSJed Brown suffix: 6_int_tet 1662c4762a1bSJed Brown requires: ctetgen 1663c0517cd5SMatthew G. Knepley args: -faces {{2,2,2 1,3,5 3,4,7}} -cell_simplex 1 -dm_generator ctetgen 1664c4762a1bSJed Brown test: 1665c4762a1bSJed Brown suffix: 6_int_hex 1666c4762a1bSJed Brown args: -faces {{2,2,2 1,3,5 3,4,7}} -cell_simplex 0 1667c4762a1bSJed Brown testset: 1668c4762a1bSJed Brown nsize: {{2 4}} 1669c4762a1bSJed Brown args: -use_generator 1670c4762a1bSJed Brown args: -dm_plex_check_all 1671c4762a1bSJed Brown args: -distribute -interpolate after_distribute 1672c4762a1bSJed Brown test: 1673c4762a1bSJed Brown suffix: 6_parint_tri 1674c4762a1bSJed Brown requires: triangle 1675c0517cd5SMatthew G. Knepley args: -faces {{2,2 1,3 7,4}} -cell_simplex 1 -dm_generator triangle 1676c4762a1bSJed Brown test: 1677c4762a1bSJed Brown suffix: 6_parint_quad 1678c4762a1bSJed Brown args: -faces {{2,2 1,3 7,4}} -cell_simplex 0 1679c4762a1bSJed Brown test: 1680c4762a1bSJed Brown suffix: 6_parint_tet 1681c4762a1bSJed Brown requires: ctetgen 1682c0517cd5SMatthew G. Knepley args: -faces {{2,2,2 1,3,5 3,4,7}} -cell_simplex 1 -dm_generator ctetgen 1683c4762a1bSJed Brown test: 1684c4762a1bSJed Brown suffix: 6_parint_hex 1685c4762a1bSJed Brown args: -faces {{2,2,2 1,3,5 3,4,7}} -cell_simplex 0 1686c4762a1bSJed Brown 1687c4762a1bSJed Brown testset: # 7 EXODUS 1688c4762a1bSJed Brown requires: exodusii 1689c4762a1bSJed Brown args: -dm_plex_check_all 1690c4762a1bSJed Brown args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo 1691c4762a1bSJed Brown args: -distribute 1692c4762a1bSJed Brown test: # seq load, simple partitioner 1693c4762a1bSJed Brown suffix: 7_exo 1694c4762a1bSJed Brown nsize: {{1 2 4 5}} 1695c4762a1bSJed Brown args: -interpolate none 1696c4762a1bSJed Brown test: # seq load, seq interpolation, simple partitioner 1697c4762a1bSJed Brown suffix: 7_exo_int_simple 1698c4762a1bSJed Brown nsize: {{1 2 4 5}} 1699c4762a1bSJed Brown args: -interpolate create 1700c4762a1bSJed Brown test: # seq load, seq interpolation, metis partitioner 1701c4762a1bSJed Brown suffix: 7_exo_int_metis 1702c4762a1bSJed Brown requires: parmetis 1703c4762a1bSJed Brown nsize: {{2 4 5}} 1704c4762a1bSJed Brown args: -interpolate create 1705c4762a1bSJed Brown args: -petscpartitioner_type parmetis 1706c4762a1bSJed Brown test: # seq load, simple partitioner, par interpolation 1707c4762a1bSJed Brown suffix: 7_exo_simple_int 1708c4762a1bSJed Brown nsize: {{2 4 5}} 1709c4762a1bSJed Brown args: -interpolate after_distribute 1710c4762a1bSJed Brown test: # seq load, metis partitioner, par interpolation 1711c4762a1bSJed Brown suffix: 7_exo_metis_int 1712c4762a1bSJed Brown requires: parmetis 1713c4762a1bSJed Brown nsize: {{2 4 5}} 1714c4762a1bSJed Brown args: -interpolate after_distribute 1715c4762a1bSJed Brown args: -petscpartitioner_type parmetis 1716c4762a1bSJed Brown 1717c4762a1bSJed Brown testset: # 7 HDF5 SEQUANTIAL LOAD 1718c4762a1bSJed Brown requires: hdf5 !complex 1719c4762a1bSJed Brown args: -dm_plex_check_all 1720c4762a1bSJed Brown args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf 1721c4762a1bSJed Brown args: -dm_plex_hdf5_force_sequential 1722c4762a1bSJed Brown args: -distribute 1723c4762a1bSJed Brown test: # seq load, simple partitioner 1724c4762a1bSJed Brown suffix: 7_seq_hdf5_simple 1725c4762a1bSJed Brown nsize: {{1 2 4 5}} 1726c4762a1bSJed Brown args: -interpolate none 1727c4762a1bSJed Brown test: # seq load, seq interpolation, simple partitioner 1728c4762a1bSJed Brown suffix: 7_seq_hdf5_int_simple 1729c4762a1bSJed Brown nsize: {{1 2 4 5}} 1730c4762a1bSJed Brown args: -interpolate after_create 1731c4762a1bSJed Brown test: # seq load, seq interpolation, metis partitioner 1732c4762a1bSJed Brown nsize: {{2 4 5}} 1733c4762a1bSJed Brown suffix: 7_seq_hdf5_int_metis 1734c4762a1bSJed Brown requires: parmetis 1735c4762a1bSJed Brown args: -interpolate after_create 1736c4762a1bSJed Brown args: -petscpartitioner_type parmetis 1737c4762a1bSJed Brown test: # seq load, simple partitioner, par interpolation 1738c4762a1bSJed Brown suffix: 7_seq_hdf5_simple_int 1739c4762a1bSJed Brown nsize: {{2 4 5}} 1740c4762a1bSJed Brown args: -interpolate after_distribute 1741c4762a1bSJed Brown test: # seq load, metis partitioner, par interpolation 1742c4762a1bSJed Brown nsize: {{2 4 5}} 1743c4762a1bSJed Brown suffix: 7_seq_hdf5_metis_int 1744c4762a1bSJed Brown requires: parmetis 1745c4762a1bSJed Brown args: -interpolate after_distribute 1746c4762a1bSJed Brown args: -petscpartitioner_type parmetis 1747c4762a1bSJed Brown 1748c4762a1bSJed Brown testset: # 7 HDF5 PARALLEL LOAD 1749c4762a1bSJed Brown requires: hdf5 !complex 1750c4762a1bSJed Brown nsize: {{2 4 5}} 1751c4762a1bSJed Brown args: -dm_plex_check_all 1752c4762a1bSJed Brown args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf 1753c4762a1bSJed Brown test: # par load 1754c4762a1bSJed Brown suffix: 7_par_hdf5 1755c4762a1bSJed Brown args: -interpolate none 1756c4762a1bSJed Brown test: # par load, par interpolation 1757c4762a1bSJed Brown suffix: 7_par_hdf5_int 1758c4762a1bSJed Brown args: -interpolate after_create 1759c4762a1bSJed Brown test: # par load, parmetis repartitioner 1760c4762a1bSJed Brown TODO: Parallel partitioning of uninterpolated meshes not supported 1761c4762a1bSJed Brown suffix: 7_par_hdf5_parmetis 1762c4762a1bSJed Brown requires: parmetis 1763c4762a1bSJed Brown args: -distribute -petscpartitioner_type parmetis 1764c4762a1bSJed Brown args: -interpolate none 1765c4762a1bSJed Brown test: # par load, par interpolation, parmetis repartitioner 1766c4762a1bSJed Brown suffix: 7_par_hdf5_int_parmetis 1767c4762a1bSJed Brown requires: parmetis 1768c4762a1bSJed Brown args: -distribute -petscpartitioner_type parmetis 1769c4762a1bSJed Brown args: -interpolate after_create 1770c4762a1bSJed Brown test: # par load, parmetis partitioner, par interpolation 1771c4762a1bSJed Brown TODO: Parallel partitioning of uninterpolated meshes not supported 1772c4762a1bSJed Brown suffix: 7_par_hdf5_parmetis_int 1773c4762a1bSJed Brown requires: parmetis 1774c4762a1bSJed Brown args: -distribute -petscpartitioner_type parmetis 1775c4762a1bSJed Brown args: -interpolate after_distribute 1776c4762a1bSJed Brown 1777c4762a1bSJed Brown test: 1778c4762a1bSJed Brown suffix: 7_hdf5_hierarch 1779c4762a1bSJed Brown requires: hdf5 ptscotch !complex 1780c4762a1bSJed Brown nsize: {{2 3 4}separate output} 1781c4762a1bSJed Brown args: -distribute 1782c4762a1bSJed Brown args: -interpolate after_create 1783c4762a1bSJed Brown args: -petscpartitioner_type matpartitioning -petscpartitioner_view ::ascii_info 1784c4762a1bSJed Brown args: -mat_partitioning_type hierarch -mat_partitioning_hierarchical_nfineparts 2 1785c4762a1bSJed Brown args: -mat_partitioning_hierarchical_coarseparttype ptscotch -mat_partitioning_hierarchical_fineparttype ptscotch 1786c4762a1bSJed Brown 1787c4762a1bSJed Brown test: 1788c4762a1bSJed Brown suffix: 8 1789c4762a1bSJed Brown requires: hdf5 !complex 1790c4762a1bSJed Brown nsize: 4 1791c4762a1bSJed Brown args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf 1792c4762a1bSJed Brown args: -distribute 0 -interpolate after_create 1793c4762a1bSJed 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 1794c4762a1bSJed Brown args: -dm_plex_check_all 1795c4762a1bSJed Brown args: -custom_view 1796c4762a1bSJed Brown 1797c4762a1bSJed Brown testset: # 9 HDF5 SEQUANTIAL LOAD 1798c4762a1bSJed Brown requires: hdf5 !complex datafilespath 1799c4762a1bSJed Brown args: -dm_plex_check_all 1800c4762a1bSJed 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 1801c4762a1bSJed Brown args: -dm_plex_hdf5_force_sequential 1802c4762a1bSJed Brown args: -distribute 1803c4762a1bSJed Brown test: # seq load, simple partitioner 1804c4762a1bSJed Brown suffix: 9_seq_hdf5_simple 1805c4762a1bSJed Brown nsize: {{1 2 4 5}} 1806c4762a1bSJed Brown args: -interpolate none 1807c4762a1bSJed Brown test: # seq load, seq interpolation, simple partitioner 1808c4762a1bSJed Brown suffix: 9_seq_hdf5_int_simple 1809c4762a1bSJed Brown nsize: {{1 2 4 5}} 1810c4762a1bSJed Brown args: -interpolate after_create 1811c4762a1bSJed Brown test: # seq load, seq interpolation, metis partitioner 1812c4762a1bSJed Brown nsize: {{2 4 5}} 1813c4762a1bSJed Brown suffix: 9_seq_hdf5_int_metis 1814c4762a1bSJed Brown requires: parmetis 1815c4762a1bSJed Brown args: -interpolate after_create 1816c4762a1bSJed Brown args: -petscpartitioner_type parmetis 1817c4762a1bSJed Brown test: # seq load, simple partitioner, par interpolation 1818c4762a1bSJed Brown suffix: 9_seq_hdf5_simple_int 1819c4762a1bSJed Brown nsize: {{2 4 5}} 1820c4762a1bSJed Brown args: -interpolate after_distribute 1821c4762a1bSJed Brown test: # seq load, simple partitioner, par interpolation 1822c4762a1bSJed Brown # This is like 9_seq_hdf5_simple_int but testing error output of DMPlexCheckPointSFHeavy(). 1823c4762a1bSJed Brown # Once 9_seq_hdf5_simple_int gets fixed, this one gets broken. 1824c4762a1bSJed Brown # We can then provide an intentionally broken mesh instead. 1825c4762a1bSJed Brown TODO: This test is broken because PointSF is fixed. 1826c4762a1bSJed Brown suffix: 9_seq_hdf5_simple_int_err 1827c4762a1bSJed Brown nsize: 4 1828c4762a1bSJed Brown args: -interpolate after_distribute 1829c4762a1bSJed Brown filter: sed -e "/PETSC ERROR/,$$d" 1830c4762a1bSJed Brown test: # seq load, metis partitioner, par interpolation 1831c4762a1bSJed Brown nsize: {{2 4 5}} 1832c4762a1bSJed Brown suffix: 9_seq_hdf5_metis_int 1833c4762a1bSJed Brown requires: parmetis 1834c4762a1bSJed Brown args: -interpolate after_distribute 1835c4762a1bSJed Brown args: -petscpartitioner_type parmetis 1836c4762a1bSJed Brown 1837c4762a1bSJed Brown testset: # 9 HDF5 PARALLEL LOAD 1838c4762a1bSJed Brown requires: hdf5 !complex datafilespath 1839c4762a1bSJed Brown nsize: {{2 4 5}} 1840c4762a1bSJed Brown args: -dm_plex_check_all 1841c4762a1bSJed 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 1842c4762a1bSJed Brown test: # par load 1843c4762a1bSJed Brown suffix: 9_par_hdf5 1844c4762a1bSJed Brown args: -interpolate none 1845c4762a1bSJed Brown test: # par load, par interpolation 1846c4762a1bSJed Brown suffix: 9_par_hdf5_int 1847c4762a1bSJed Brown args: -interpolate after_create 1848c4762a1bSJed Brown test: # par load, parmetis repartitioner 1849c4762a1bSJed Brown TODO: Parallel partitioning of uninterpolated meshes not supported 1850c4762a1bSJed Brown suffix: 9_par_hdf5_parmetis 1851c4762a1bSJed Brown requires: parmetis 1852c4762a1bSJed Brown args: -distribute -petscpartitioner_type parmetis 1853c4762a1bSJed Brown args: -interpolate none 1854c4762a1bSJed Brown test: # par load, par interpolation, parmetis repartitioner 1855c4762a1bSJed Brown suffix: 9_par_hdf5_int_parmetis 1856c4762a1bSJed Brown requires: parmetis 1857c4762a1bSJed Brown args: -distribute -petscpartitioner_type parmetis 1858c4762a1bSJed Brown args: -interpolate after_create 1859c4762a1bSJed Brown test: # par load, parmetis partitioner, par interpolation 1860c4762a1bSJed Brown TODO: Parallel partitioning of uninterpolated meshes not supported 1861c4762a1bSJed Brown suffix: 9_par_hdf5_parmetis_int 1862c4762a1bSJed Brown requires: parmetis 1863c4762a1bSJed Brown args: -distribute -petscpartitioner_type parmetis 1864c4762a1bSJed Brown args: -interpolate after_distribute 1865c4762a1bSJed Brown 1866c4762a1bSJed Brown testset: # 10 HDF5 PARALLEL LOAD 1867c4762a1bSJed Brown requires: hdf5 !complex datafilespath 1868c4762a1bSJed Brown nsize: {{2 4 7}} 1869c4762a1bSJed Brown args: -dm_plex_check_all 1870c4762a1bSJed 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 1871c4762a1bSJed Brown test: # par load, par interpolation 1872c4762a1bSJed Brown suffix: 10_par_hdf5_int 1873c4762a1bSJed Brown args: -interpolate after_create 1874c4762a1bSJed Brown TEST*/ 1875