15cca0b87SVaclav Hapla static char help[] = "Load and save the mesh and fields to HDF5 and ExodusII\n\n"; 25cca0b87SVaclav Hapla 3bd412c90SVaclav Hapla #include <petsc/private/dmpleximpl.h> 45cca0b87SVaclav Hapla #include <petscviewerhdf5.h> 55cca0b87SVaclav Hapla #include <petscsf.h> 65cca0b87SVaclav Hapla 75cca0b87SVaclav Hapla typedef struct { 85cca0b87SVaclav Hapla PetscBool compare; /* Compare the meshes using DMPlexEqual() */ 938883f1bSVaclav Hapla PetscBool compare_labels; /* Compare labels in the meshes using DMCompareLabels() */ 105cca0b87SVaclav Hapla PetscBool distribute; /* Distribute the mesh */ 11*bcde9acaSMatthew G. Knepley PetscBool field; /* Layout a field over the mesh */ 12*bcde9acaSMatthew G. Knepley PetscBool reorder; /* Reorder the points in the section */ 13bd412c90SVaclav Hapla char ofname[PETSC_MAX_PATH_LEN]; /* Output mesh filename */ 145cca0b87SVaclav Hapla PetscViewerFormat format; /* Format to write and read */ 155cca0b87SVaclav Hapla PetscBool second_write_read; /* Write and read for the 2nd time */ 16806c51deSksagiyam PetscBool use_low_level_functions; /* Use low level functions for viewing and loading */ 175cca0b87SVaclav Hapla } AppCtx; 185cca0b87SVaclav Hapla 19d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 20d71ae5a4SJacob Faibussowitsch { 215cca0b87SVaclav Hapla PetscFunctionBeginUser; 225cca0b87SVaclav Hapla options->compare = PETSC_FALSE; 2338883f1bSVaclav Hapla options->compare_labels = PETSC_FALSE; 245cca0b87SVaclav Hapla options->distribute = PETSC_TRUE; 25*bcde9acaSMatthew G. Knepley options->field = PETSC_FALSE; 26*bcde9acaSMatthew G. Knepley options->reorder = PETSC_FALSE; 275cca0b87SVaclav Hapla options->format = PETSC_VIEWER_DEFAULT; 285cca0b87SVaclav Hapla options->second_write_read = PETSC_FALSE; 29806c51deSksagiyam options->use_low_level_functions = PETSC_FALSE; 30c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(options->ofname, "ex55.h5", sizeof(options->ofname))); 315cca0b87SVaclav Hapla 32d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX"); 339566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-compare", "Compare the meshes using DMPlexEqual()", "ex55.c", options->compare, &options->compare, NULL)); 349566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-compare_labels", "Compare labels in the meshes using DMCompareLabels()", "ex55.c", options->compare_labels, &options->compare_labels, NULL)); 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-distribute", "Distribute the mesh", "ex55.c", options->distribute, &options->distribute, NULL)); 36*bcde9acaSMatthew G. Knepley PetscCall(PetscOptionsBool("-field", "Layout a field over the mesh", "ex55.c", options->field, &options->field, NULL)); 37*bcde9acaSMatthew G. Knepley PetscCall(PetscOptionsBool("-reorder", "Reorder the points in the section", "ex55.c", options->reorder, &options->reorder, NULL)); 38bd412c90SVaclav Hapla PetscCall(PetscOptionsString("-ofilename", "The output mesh file", "ex55.c", options->ofname, options->ofname, sizeof(options->ofname), NULL)); 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-format", "Format to write and read", "ex55.c", PetscViewerFormats, (PetscEnum)options->format, (PetscEnum *)&options->format, NULL)); 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-second_write_read", "Write and read for the 2nd time", "ex55.c", options->second_write_read, &options->second_write_read, NULL)); 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_low_level_functions", "Use low level functions for viewing and loading", "ex55.c", options->use_low_level_functions, &options->use_low_level_functions, NULL)); 42d0609cedSBarry Smith PetscOptionsEnd(); 433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44bd412c90SVaclav Hapla } 455cca0b87SVaclav Hapla 46*bcde9acaSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) 47*bcde9acaSMatthew G. Knepley { 48*bcde9acaSMatthew G. Knepley PetscFunctionBeginUser; 49*bcde9acaSMatthew G. Knepley PetscCall(DMCreate(comm, dm)); 50*bcde9acaSMatthew G. Knepley PetscCall(DMSetType(*dm, DMPLEX)); 51*bcde9acaSMatthew G. Knepley PetscCall(DMSetOptionsPrefix(*dm, "orig_")); 52*bcde9acaSMatthew G. Knepley PetscCall(DMSetFromOptions(*dm)); 53*bcde9acaSMatthew G. Knepley PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 54*bcde9acaSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 55*bcde9acaSMatthew G. Knepley } 56*bcde9acaSMatthew G. Knepley 57*bcde9acaSMatthew G. Knepley static PetscErrorCode CreateDiscretization(DM dm) 58*bcde9acaSMatthew G. Knepley { 59*bcde9acaSMatthew G. Knepley PetscFE fe; 60*bcde9acaSMatthew G. Knepley DMPolytopeType ct; 61*bcde9acaSMatthew G. Knepley PetscInt dim, cStart; 62*bcde9acaSMatthew G. Knepley 63*bcde9acaSMatthew G. Knepley PetscFunctionBeginUser; 64*bcde9acaSMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim)); 65*bcde9acaSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 66*bcde9acaSMatthew G. Knepley PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 67*bcde9acaSMatthew G. Knepley PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, 1, ct, 1, PETSC_DETERMINE, &fe)); 68*bcde9acaSMatthew G. Knepley PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 69*bcde9acaSMatthew G. Knepley PetscCall(PetscFEDestroy(&fe)); 70*bcde9acaSMatthew G. Knepley PetscCall(DMCreateDS(dm)); 71*bcde9acaSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 72*bcde9acaSMatthew G. Knepley } 73*bcde9acaSMatthew G. Knepley 74bd412c90SVaclav Hapla static PetscErrorCode CheckDistributed(DM dm, PetscBool expectedDistributed) 75bd412c90SVaclav Hapla { 76bd412c90SVaclav Hapla PetscMPIInt size; 77bd412c90SVaclav Hapla PetscBool distributed; 78bd412c90SVaclav Hapla const char YES[] = "DISTRIBUTED"; 79bd412c90SVaclav Hapla const char NO[] = "NOT DISTRIBUTED"; 80bd412c90SVaclav Hapla 81bd412c90SVaclav Hapla PetscFunctionBeginUser; 82bd412c90SVaclav Hapla PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 833ba16761SJacob Faibussowitsch if (size < 2) PetscFunctionReturn(PETSC_SUCCESS); 84bd412c90SVaclav Hapla PetscCall(DMPlexIsDistributed(dm, &distributed)); 85bd412c90SVaclav Hapla PetscCheck(distributed == expectedDistributed, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Expected DM being %s but actually is %s", expectedDistributed ? YES : NO, distributed ? YES : NO); 863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 87bd412c90SVaclav Hapla } 88bd412c90SVaclav Hapla 89bd412c90SVaclav Hapla static PetscErrorCode CheckInterpolated(DM dm, PetscBool expectedInterpolated) 90bd412c90SVaclav Hapla { 91bd412c90SVaclav Hapla DMPlexInterpolatedFlag iflg; 92bd412c90SVaclav Hapla PetscBool interpolated; 93bd412c90SVaclav Hapla const char YES[] = "INTERPOLATED"; 94bd412c90SVaclav Hapla const char NO[] = "NOT INTERPOLATED"; 95bd412c90SVaclav Hapla 96bd412c90SVaclav Hapla PetscFunctionBeginUser; 97bd412c90SVaclav Hapla PetscCall(DMPlexIsInterpolatedCollective(dm, &iflg)); 98bd412c90SVaclav Hapla interpolated = (PetscBool)(iflg == DMPLEX_INTERPOLATED_FULL); 99bd412c90SVaclav Hapla PetscCheck(interpolated == expectedInterpolated, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Expected DM being %s but actually is %s", expectedInterpolated ? YES : NO, interpolated ? YES : NO); 1003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 101bd412c90SVaclav Hapla } 102bd412c90SVaclav Hapla 103*bcde9acaSMatthew G. Knepley static PetscErrorCode CheckDistributedInterpolated(DM dm, PetscBool expectedInterpolated, PetscViewer v, AppCtx *user) 104bd412c90SVaclav Hapla { 105bd412c90SVaclav Hapla PetscViewerFormat format; 106*bcde9acaSMatthew G. Knepley PetscBool distributed, interpolated = expectedInterpolated; 107bd412c90SVaclav Hapla 108bd412c90SVaclav Hapla PetscFunctionBeginUser; 109bd412c90SVaclav Hapla PetscCall(PetscViewerGetFormat(v, &format)); 110bd412c90SVaclav Hapla switch (format) { 111bd412c90SVaclav Hapla case PETSC_VIEWER_HDF5_XDMF: 112bd412c90SVaclav Hapla case PETSC_VIEWER_HDF5_VIZ: { 113bd412c90SVaclav Hapla distributed = PETSC_TRUE; 114bd412c90SVaclav Hapla interpolated = PETSC_FALSE; 115bd412c90SVaclav Hapla }; break; 116bd412c90SVaclav Hapla case PETSC_VIEWER_HDF5_PETSC: 117bd412c90SVaclav Hapla case PETSC_VIEWER_DEFAULT: 118bd412c90SVaclav Hapla case PETSC_VIEWER_NATIVE: { 119bd412c90SVaclav Hapla DMPlexStorageVersion version; 120bd412c90SVaclav Hapla 121bd412c90SVaclav Hapla PetscCall(PetscViewerHDF5GetDMPlexStorageVersionReading(v, &version)); 122bd412c90SVaclav Hapla distributed = (PetscBool)(version->major >= 3); 123bd412c90SVaclav Hapla }; break; 124*bcde9acaSMatthew G. Knepley default: 125bd412c90SVaclav Hapla distributed = PETSC_FALSE; 126bd412c90SVaclav Hapla } 127bd412c90SVaclav Hapla PetscCall(CheckDistributed(dm, distributed)); 128bd412c90SVaclav Hapla PetscCall(CheckInterpolated(dm, interpolated)); 1293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 130bd412c90SVaclav Hapla } 131bd412c90SVaclav Hapla 132*bcde9acaSMatthew G. Knepley static PetscErrorCode DMPlexWriteAndReadHDF5(DM dm, Vec vec, const char filename[], const char prefix[], PetscBool expectedInterpolated, AppCtx *user, DM *dm_new, Vec *v_new) 133d71ae5a4SJacob Faibussowitsch { 1345cca0b87SVaclav Hapla DM dmnew; 135*bcde9acaSMatthew G. Knepley Vec vnew = NULL; 136ed5e4e85SVaclav Hapla const char savedName[] = "Mesh"; 137ed5e4e85SVaclav Hapla const char loadedName[] = "Mesh_new"; 1385cca0b87SVaclav Hapla PetscViewer v; 1395cca0b87SVaclav Hapla 1405cca0b87SVaclav Hapla PetscFunctionBeginUser; 1419566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Open(PetscObjectComm((PetscObject)dm), filename, FILE_MODE_WRITE, &v)); 142bd412c90SVaclav Hapla PetscCall(PetscViewerPushFormat(v, user->format)); 1439566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dm, savedName)); 144bd412c90SVaclav Hapla if (user->use_low_level_functions) { 1459566063dSJacob Faibussowitsch PetscCall(DMPlexTopologyView(dm, v)); 1469566063dSJacob Faibussowitsch PetscCall(DMPlexCoordinatesView(dm, v)); 1479566063dSJacob Faibussowitsch PetscCall(DMPlexLabelsView(dm, v)); 148806c51deSksagiyam } else { 1499566063dSJacob Faibussowitsch PetscCall(DMView(dm, v)); 150*bcde9acaSMatthew G. Knepley if (vec) PetscCall(VecView(vec, v)); 151806c51deSksagiyam } 1529566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(v, FILE_MODE_READ)); 1539566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dmnew)); 1549566063dSJacob Faibussowitsch PetscCall(DMSetType(dmnew, DMPLEX)); 155bd412c90SVaclav Hapla PetscCall(DMPlexDistributeSetDefault(dmnew, PETSC_FALSE)); 1569566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dmnew, savedName)); 1579566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(dmnew, prefix)); 158bd412c90SVaclav Hapla if (user->use_low_level_functions) { 159c9ad657eSksagiyam PetscSF sfXC; 160c9ad657eSksagiyam 1619566063dSJacob Faibussowitsch PetscCall(DMPlexTopologyLoad(dmnew, v, &sfXC)); 1629566063dSJacob Faibussowitsch PetscCall(DMPlexCoordinatesLoad(dmnew, v, sfXC)); 1639566063dSJacob Faibussowitsch PetscCall(DMPlexLabelsLoad(dmnew, v, sfXC)); 1649566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfXC)); 165806c51deSksagiyam } else { 1669566063dSJacob Faibussowitsch PetscCall(DMLoad(dmnew, v)); 167*bcde9acaSMatthew G. Knepley if (vec) { 168*bcde9acaSMatthew G. Knepley PetscCall(CreateDiscretization(dmnew)); 169*bcde9acaSMatthew G. Knepley PetscCall(DMCreateGlobalVector(dmnew, &vnew)); 170*bcde9acaSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)vnew, "solution")); 171*bcde9acaSMatthew G. Knepley PetscCall(VecLoad(vnew, v)); 172806c51deSksagiyam } 173*bcde9acaSMatthew G. Knepley } 174*bcde9acaSMatthew G. Knepley PetscCall(CheckDistributedInterpolated(dmnew, expectedInterpolated, v, user)); 1759566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dmnew, loadedName)); 1769566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(v)); 1779566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&v)); 1785cca0b87SVaclav Hapla *dm_new = dmnew; 179*bcde9acaSMatthew G. Knepley *v_new = vnew; 1803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1815cca0b87SVaclav Hapla } 1825cca0b87SVaclav Hapla 183d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 184d71ae5a4SJacob Faibussowitsch { 1855cca0b87SVaclav Hapla DM dm, dmnew; 186*bcde9acaSMatthew G. Knepley Vec v = NULL, vnew = NULL; 1875cca0b87SVaclav Hapla PetscPartitioner part; 1885cca0b87SVaclav Hapla AppCtx user; 189*bcde9acaSMatthew G. Knepley PetscBool interpolated = PETSC_TRUE, flg; 1905cca0b87SVaclav Hapla 191327415f7SBarry Smith PetscFunctionBeginUser; 1929566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1939566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 194*bcde9acaSMatthew G. Knepley PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm)); 195*bcde9acaSMatthew G. Knepley PetscCall(PetscOptionsGetBool(NULL, "orig_", "-dm_plex_interpolate", &interpolated, NULL)); 196*bcde9acaSMatthew G. Knepley PetscCall(CheckInterpolated(dm, interpolated)); 1975cca0b87SVaclav Hapla 1985cca0b87SVaclav Hapla if (user.distribute) { 1995cca0b87SVaclav Hapla DM dmdist; 2005cca0b87SVaclav Hapla 2019566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(dm, &part)); 202bd412c90SVaclav Hapla PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSIMPLE)); 2039566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part)); 2049566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(dm, 0, NULL, &dmdist)); 2055cca0b87SVaclav Hapla if (dmdist) { 2069566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2075cca0b87SVaclav Hapla dm = dmdist; 208bd412c90SVaclav Hapla PetscCall(CheckDistributed(dm, PETSC_TRUE)); 209*bcde9acaSMatthew G. Knepley PetscCall(CheckInterpolated(dm, interpolated)); 2105cca0b87SVaclav Hapla } 2115cca0b87SVaclav Hapla } 212*bcde9acaSMatthew G. Knepley if (user.field) { 213*bcde9acaSMatthew G. Knepley PetscSection gs; 214*bcde9acaSMatthew G. Knepley PetscScalar *a; 215*bcde9acaSMatthew G. Knepley PetscInt pStart, pEnd, rStart; 216*bcde9acaSMatthew G. Knepley 217*bcde9acaSMatthew G. Knepley PetscCall(CreateDiscretization(dm)); 218*bcde9acaSMatthew G. Knepley 219*bcde9acaSMatthew G. Knepley PetscCall(DMCreateGlobalVector(dm, &v)); 220*bcde9acaSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)v, "solution")); 221*bcde9acaSMatthew G. Knepley PetscCall(DMGetGlobalSection(dm, &gs)); 222*bcde9acaSMatthew G. Knepley PetscCall(PetscSectionGetChart(gs, &pStart, &pEnd)); 223*bcde9acaSMatthew G. Knepley PetscCall(VecGetOwnershipRange(v, &rStart, NULL)); 224*bcde9acaSMatthew G. Knepley PetscCall(VecGetArrayWrite(v, &a)); 225*bcde9acaSMatthew G. Knepley for (PetscInt p = pStart; p < pEnd; ++p) { 226*bcde9acaSMatthew G. Knepley PetscInt dof, off; 227*bcde9acaSMatthew G. Knepley 228*bcde9acaSMatthew G. Knepley PetscCall(PetscSectionGetDof(gs, p, &dof)); 229*bcde9acaSMatthew G. Knepley PetscCall(PetscSectionGetOffset(gs, p, &off)); 230*bcde9acaSMatthew G. Knepley if (off < 0) continue; 231*bcde9acaSMatthew G. Knepley for (PetscInt d = 0; d < dof; ++d) a[off + d] = p; 232*bcde9acaSMatthew G. Knepley } 233*bcde9acaSMatthew G. Knepley PetscCall(VecRestoreArrayWrite(v, &a)); 234*bcde9acaSMatthew G. Knepley } 2355cca0b87SVaclav Hapla 2369566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(dm, NULL)); 2379566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE)); 2389566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 2399566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 2405cca0b87SVaclav Hapla 241*bcde9acaSMatthew G. Knepley PetscCall(DMPlexWriteAndReadHDF5(dm, v, user.ofname, "new_", interpolated, &user, &dmnew, &vnew)); 2425cca0b87SVaclav Hapla 2435cca0b87SVaclav Hapla if (user.second_write_read) { 2449566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2455cca0b87SVaclav Hapla dm = dmnew; 246*bcde9acaSMatthew G. Knepley PetscCall(DMPlexWriteAndReadHDF5(dm, v, user.ofname, "new_", interpolated, &user, &dmnew, &vnew)); 2475cca0b87SVaclav Hapla } 2485cca0b87SVaclav Hapla 2499566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dmnew, NULL, "-dm_view")); 2505cca0b87SVaclav Hapla 2515cca0b87SVaclav Hapla /* This currently makes sense only for sequential meshes. */ 2525cca0b87SVaclav Hapla if (user.compare) { 2539566063dSJacob Faibussowitsch PetscCall(DMPlexEqual(dmnew, dm, &flg)); 25428b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "DMs are not equal"); 2559566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "DMs equal\n")); 256*bcde9acaSMatthew G. Knepley if (v) { 257*bcde9acaSMatthew G. Knepley PetscCall(VecEqual(vnew, v, &flg)); 258*bcde9acaSMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Vecs %s\n", flg ? "equal" : "are not equal")); 259*bcde9acaSMatthew G. Knepley PetscCall(VecViewFromOptions(v, NULL, "-old_vec_view")); 260*bcde9acaSMatthew G. Knepley PetscCall(VecViewFromOptions(vnew, NULL, "-new_vec_view")); 261*bcde9acaSMatthew G. Knepley } 26238883f1bSVaclav Hapla } 26338883f1bSVaclav Hapla if (user.compare_labels) { 2649566063dSJacob Faibussowitsch PetscCall(DMCompareLabels(dmnew, dm, NULL, NULL)); 2659566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "DMLabels equal\n")); 2665cca0b87SVaclav Hapla } 2675cca0b87SVaclav Hapla 268*bcde9acaSMatthew G. Knepley PetscCall(VecDestroy(&v)); 2699566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 270*bcde9acaSMatthew G. Knepley PetscCall(VecDestroy(&vnew)); 2719566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmnew)); 2729566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 273b122ec5aSJacob Faibussowitsch return 0; 2745cca0b87SVaclav Hapla } 2755cca0b87SVaclav Hapla 2765cca0b87SVaclav Hapla /*TEST 2775cca0b87SVaclav Hapla build: 2785cca0b87SVaclav Hapla requires: hdf5 2795cca0b87SVaclav Hapla # Idempotence of saving/loading 2805cca0b87SVaclav Hapla # Have to replace Exodus file, which is creating uninterpolated edges 2815cca0b87SVaclav Hapla test: 2825cca0b87SVaclav Hapla suffix: 0 28338883f1bSVaclav Hapla TODO: broken 28438883f1bSVaclav Hapla requires: exodusii 2855cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail 2865cca0b87SVaclav Hapla args: -format hdf5_petsc -compare 2875cca0b87SVaclav Hapla test: 2885cca0b87SVaclav Hapla suffix: 1 28938883f1bSVaclav Hapla TODO: broken 29038883f1bSVaclav Hapla requires: exodusii parmetis !defined(PETSC_USE_64BIT_INDICES) 2915cca0b87SVaclav Hapla nsize: 2 2925cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail 2935cca0b87SVaclav Hapla args: -petscpartitioner_type parmetis 2945cca0b87SVaclav Hapla args: -format hdf5_petsc -new_dm_view ascii::ascii_info_detail 295bd412c90SVaclav Hapla 2965cca0b87SVaclav Hapla testset: 2975cca0b87SVaclav Hapla requires: exodusii 298*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo 2995cca0b87SVaclav Hapla test: 3005cca0b87SVaclav Hapla suffix: 2 3015cca0b87SVaclav Hapla nsize: {{1 2 4 8}separate output} 3025cca0b87SVaclav Hapla args: -format {{default hdf5_petsc}separate output} 3038886822aSVaclav Hapla args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0 3.0.0}} 304*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_interpolate {{0 1}separate output} 3055cca0b87SVaclav Hapla test: 3065cca0b87SVaclav Hapla suffix: 2a 3075cca0b87SVaclav Hapla nsize: {{1 2 4 8}separate output} 3085cca0b87SVaclav Hapla args: -format {{hdf5_xdmf hdf5_viz}separate output} 309bd412c90SVaclav Hapla 3105cca0b87SVaclav Hapla test: 3115cca0b87SVaclav Hapla suffix: 3 3125cca0b87SVaclav Hapla requires: exodusii 313*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -compare -compare_labels 3148886822aSVaclav Hapla args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0 3.0.0}} 3155cca0b87SVaclav Hapla 3165cca0b87SVaclav Hapla # Load HDF5 file in XDMF format in parallel, write, read dm1, write, read dm2, and compare dm1 and dm2 31738883f1bSVaclav Hapla testset: 3185cca0b87SVaclav Hapla suffix: 4 3195cca0b87SVaclav Hapla requires: !complex 320*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf 32138883f1bSVaclav Hapla args: -distribute 0 -second_write_read -compare 32238883f1bSVaclav Hapla test: 32338883f1bSVaclav Hapla suffix: hdf5_petsc 32438883f1bSVaclav Hapla nsize: {{1 2}} 32538883f1bSVaclav Hapla args: -format hdf5_petsc -compare_labels 3268886822aSVaclav Hapla args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0 3.0.0}} 32738883f1bSVaclav Hapla test: 32838883f1bSVaclav Hapla suffix: hdf5_xdmf 32938883f1bSVaclav Hapla nsize: {{1 3 8}} 33038883f1bSVaclav Hapla args: -format hdf5_xdmf 3315cca0b87SVaclav Hapla 332806c51deSksagiyam # Use low level functions, DMPlexTopologyView()/Load(), DMPlexCoordinatesView()/Load(), and DMPlexLabelsView()/Load() 333bd412c90SVaclav Hapla # TODO: The output is very long so keeping just 1.0.0 version. This test should be redesigned or removed. 334806c51deSksagiyam test: 335806c51deSksagiyam suffix: 5 336806c51deSksagiyam requires: exodusii 337806c51deSksagiyam nsize: 2 338*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo 339*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_interpolate 0 -orig_dm_distribute 0 340806c51deSksagiyam args: -dm_view ascii::ascii_info_detail 341806c51deSksagiyam args: -new_dm_view ascii::ascii_info_detail 342bd412c90SVaclav Hapla args: -format hdf5_petsc -use_low_level_functions {{0 1}} 343bd412c90SVaclav Hapla args: -dm_plex_view_hdf5_storage_version 1.0.0 344806c51deSksagiyam 3455cca0b87SVaclav Hapla testset: 34638883f1bSVaclav Hapla suffix: 6 34738883f1bSVaclav Hapla requires: hdf5 !complex datafilespath 34838883f1bSVaclav Hapla nsize: {{1 3}} 34938883f1bSVaclav Hapla args: -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_geometry 350*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_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 351*bcde9acaSMatthew G. Knepley args: -orig_dm_distribute 0 35238883f1bSVaclav Hapla args: -format hdf5_petsc -second_write_read -compare -compare_labels 353*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_interpolate {{0 1}} -distribute {{0 1}} 3548886822aSVaclav Hapla args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0 3.0.0}} 35538883f1bSVaclav Hapla 35638883f1bSVaclav Hapla testset: 3575cca0b87SVaclav Hapla # the same data and settings as dm_impls_plex_tests-ex18_9% 35838883f1bSVaclav Hapla suffix: 9 3595cca0b87SVaclav Hapla requires: hdf5 !complex datafilespath 3604d056bb6SVaclav Hapla nsize: {{1 2 4}} 3615cca0b87SVaclav Hapla args: -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_geometry 362*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_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 363*bcde9acaSMatthew G. Knepley args: -orig_dm_distribute 0 364bd412c90SVaclav Hapla args: -format {{hdf5_petsc hdf5_xdmf}} -second_write_read -compare 3658886822aSVaclav Hapla args: -dm_plex_view_hdf5_storage_version 3.0.0 3665cca0b87SVaclav Hapla test: 36738883f1bSVaclav Hapla suffix: hdf5_seqload 368bd412c90SVaclav Hapla args: -distribute 369*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_interpolate {{0 1}} 3705cca0b87SVaclav Hapla args: -dm_plex_hdf5_force_sequential 3715cca0b87SVaclav Hapla test: 37238883f1bSVaclav Hapla suffix: hdf5_seqload_metis 3735cca0b87SVaclav Hapla requires: parmetis 3745cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis 375*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_interpolate 1 3765cca0b87SVaclav Hapla args: -dm_plex_hdf5_force_sequential 3775cca0b87SVaclav Hapla test: 37838883f1bSVaclav Hapla suffix: hdf5 379*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_interpolate 1 3805cca0b87SVaclav Hapla test: 38138883f1bSVaclav Hapla suffix: hdf5_repart 3825cca0b87SVaclav Hapla requires: parmetis 3835cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis 384*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_interpolate 1 3855cca0b87SVaclav Hapla test: 3865cca0b87SVaclav Hapla TODO: Parallel partitioning of uninterpolated meshes not supported 38738883f1bSVaclav Hapla suffix: hdf5_repart_ppu 3885cca0b87SVaclav Hapla requires: parmetis 3895cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis 390*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_interpolate 0 3915cca0b87SVaclav Hapla 3925cca0b87SVaclav Hapla # reproduce PetscSFView() crash - fixed, left as regression test 3935cca0b87SVaclav Hapla test: 3945cca0b87SVaclav Hapla suffix: new_dm_view 3955cca0b87SVaclav Hapla requires: exodusii 3965cca0b87SVaclav Hapla nsize: 2 397*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/TwoQuads.exo -new_dm_view ascii:ex5_new_dm_view.log:ascii_info_detail 3988886822aSVaclav Hapla args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0 3.0.0}} 39996b7d01aSVaclav Hapla 400bd412c90SVaclav Hapla # test backward compatibility with petsc_hdf5 format version 1.0.0, serial idempotence 40196b7d01aSVaclav Hapla testset: 40296b7d01aSVaclav Hapla suffix: 10-v3.16.0-v1.0.0 40396b7d01aSVaclav Hapla requires: hdf5 !complex datafilespath 404e6ef0f20SVaclav Hapla args: -dm_plex_check_all -compare -compare_labels 4059ed1248dSVaclav Hapla args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0 3.0.0}} -use_low_level_functions {{0 1}} 40696b7d01aSVaclav Hapla test: 40796b7d01aSVaclav Hapla suffix: a 408*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5 40996b7d01aSVaclav Hapla test: 41096b7d01aSVaclav Hapla suffix: b 41196b7d01aSVaclav Hapla TODO: broken 412*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5 41396b7d01aSVaclav Hapla test: 41496b7d01aSVaclav Hapla suffix: c 415*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5 41696b7d01aSVaclav Hapla test: 41796b7d01aSVaclav Hapla suffix: d 418*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 41996b7d01aSVaclav Hapla test: 42096b7d01aSVaclav Hapla suffix: e 421*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5 42296b7d01aSVaclav Hapla test: 42396b7d01aSVaclav Hapla suffix: f 424*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5 425*bcde9acaSMatthew G. Knepley 426*bcde9acaSMatthew G. Knepley # test permuted sections with petsc_hdf5 format version 1.0.0 427*bcde9acaSMatthew G. Knepley testset: 428*bcde9acaSMatthew G. Knepley suffix: 11 429*bcde9acaSMatthew G. Knepley requires: hdf5 triangle 430*bcde9acaSMatthew G. Knepley args: -field 431*bcde9acaSMatthew G. Knepley args: -dm_plex_check_all -compare -compare_labels 432*bcde9acaSMatthew G. Knepley args: -orig_dm_plex_box_faces 3,3 -dm_plex_view_hdf5_storage_version 1.0.0 433*bcde9acaSMatthew G. Knepley args: -orig_dm_reorder_section -orig_dm_reorder_section_type reverse 434*bcde9acaSMatthew G. Knepley 435*bcde9acaSMatthew G. Knepley test: 436*bcde9acaSMatthew G. Knepley suffix: serial 437*bcde9acaSMatthew G. Knepley test: 438*bcde9acaSMatthew G. Knepley suffix: serial_no_perm 439*bcde9acaSMatthew G. Knepley args: -orig_dm_ignore_perm_output 44096b7d01aSVaclav Hapla 4415cca0b87SVaclav Hapla TEST*/ 442