Actual source code: ex5.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1: static char help[] = "Load and save the mesh and fields to HDF5 and ExodusII\n\n";

  3:  #include <petscdmplex.h>
  4: #include <petscviewerhdf5.h>
  5:  #include <petscsf.h>

  7: typedef struct {
  8:   PetscBool interpolate;                  /* Generate intermediate mesh elements */
  9:   char      filename[PETSC_MAX_PATH_LEN]; /* Mesh filename */
 10: } AppCtx;

 12: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
 13: {

 17:   options->interpolate = PETSC_FALSE;
 18:   options->filename[0] = '\0';

 20:   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
 21:   PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex5.c", options->interpolate, &options->interpolate, NULL);
 22:   PetscOptionsString("-filename", "The mesh file", "ex5.c", options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);
 23:   PetscOptionsEnd();
 24:   return(0);
 25: };

 27: int main(int argc, char **argv)
 28: {
 29:   DM             dm, dmdist, dmnew;
 30:   PetscPartitioner part;
 31:   AppCtx         user;
 32:   PetscViewer    v;
 33:   PetscSF        pointSF;
 34:   PetscBool      flg;

 37:   PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
 38:   ProcessOptions(PETSC_COMM_WORLD, &user);
 39:   DMPlexCreateFromFile(PETSC_COMM_WORLD, user.filename, user.interpolate, &dm);
 40:   DMViewFromOptions(dm, NULL, "-orig_dm_view");
 41:   DMPlexGetPartitioner(dm, &part);
 42:   PetscObjectSetOptionsPrefix((PetscObject)part, "orig_");
 43:   PetscPartitionerSetFromOptions(part);
 44:   DMPlexDistribute(dm, 0, &pointSF, &dmdist);
 45:   if (dmdist) {
 46:     DMDestroy(&dm);
 47:     dm   = dmdist;
 48:   }
 49:   PetscSFDestroy(&pointSF);
 50:   DMSetFromOptions(dm);
 51:   DMViewFromOptions(dm, NULL, "-dm_view");

 53:   PetscViewerHDF5Open(PetscObjectComm((PetscObject) dm), "dmdist.h5", FILE_MODE_WRITE, &v);
 54:   DMView(dm, v);
 55:   PetscViewerDestroy(&v);

 57:   DMCreate(PetscObjectComm((PetscObject) dm), &dmnew);
 58:   DMSetType(dmnew, DMPLEX);
 59:   PetscViewerHDF5Open(PETSC_COMM_WORLD, "dmdist.h5", FILE_MODE_READ, &v);
 60:   DMLoad(dmnew, v);
 61:   PetscViewerDestroy(&v);
 62:   DMViewFromOptions(dmnew, NULL, "-new_dm_view");

 64:   /* The NATIVE format for coordiante viewing is killing parallel output, since we have a local vector. Map it to global, and it will work. */
 65:   DMPlexEqual(dmnew, dm, &flg);
 66:   if (flg) {PetscPrintf(PETSC_COMM_WORLD,"DMs equal\n");}
 67:   else     {PetscPrintf(PETSC_COMM_WORLD,"DMs are not equal\n");}

 69:   DMDestroy(&dm);
 70:   DMDestroy(&dmnew);
 71:   PetscFinalize();
 72:   return ierr;
 73: }

 75: /*TEST
 76:   build:
 77:     requires: exodusii
 78:   # Idempotence of saving/loading
 79:   test:
 80:     suffix: 0
 81:     requires: exodusii broken
 82:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail
 83:   test:
 84:     suffix: 1
 85:     requires: exodusii broken
 86:     nsize: 2
 87:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail
 88:   # reproduce PetscSFView() crash - fixed, left as regression test
 89:   test:
 90:     suffix: new_dm_view
 91:     requires: exodusii
 92:     nsize: 2
 93:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/TwoQuads.exo -new_dm_view ascii::ascii_info_detail

 95: TEST*/