Actual source code: dmmbio.cxx

  1: #include <petsc/private/dmmbimpl.h>
  2: #include <petscdmmoab.h>

  4: static PetscErrorCode DMMoab_GetWriteOptions_Private(PetscInt fsetid, PetscInt numproc, PetscInt dim, MoabWriteMode mode, PetscInt dbglevel, const char *dm_opts, const char *extra_opts, const char **write_opts)
  5: {
  6:   char *wopts;
  7:   char  wopts_par[PETSC_MAX_PATH_LEN];
  8:   char  wopts_parid[PETSC_MAX_PATH_LEN];
  9:   char  wopts_dbg[PETSC_MAX_PATH_LEN];
 10:   PetscFunctionBegin;

 12:   PetscCall(PetscMalloc1(PETSC_MAX_PATH_LEN, &wopts));
 13:   PetscCall(PetscMemzero(&wopts_par, PETSC_MAX_PATH_LEN));
 14:   PetscCall(PetscMemzero(&wopts_parid, PETSC_MAX_PATH_LEN));
 15:   PetscCall(PetscMemzero(&wopts_dbg, PETSC_MAX_PATH_LEN));

 17:   // do parallel read unless only one processor
 18:   if (numproc > 1) {
 19:     PetscCall(PetscSNPrintf(wopts_par, PETSC_MAX_PATH_LEN, "PARALLEL=%s;", MoabWriteModes[mode]));
 20:     if (fsetid >= 0) PetscCall(PetscSNPrintf(wopts_parid, PETSC_MAX_PATH_LEN, "PARALLEL_COMM=%d;", fsetid));
 21:   }

 23:   if (dbglevel) PetscCall(PetscSNPrintf(wopts_dbg, PETSC_MAX_PATH_LEN, "CPUTIME;DEBUG_IO=%d;", dbglevel));

 25:   PetscCall(PetscSNPrintf(wopts, PETSC_MAX_PATH_LEN, "%s%s%s%s%s", wopts_par, wopts_parid, wopts_dbg, (extra_opts ? extra_opts : ""), (dm_opts ? dm_opts : "")));
 26:   *write_opts = wopts;
 27:   PetscFunctionReturn(PETSC_SUCCESS);
 28: }

 30: /*@C
 31:   DMMoabOutput - Output the solution vectors that are stored in the DMMoab object as tags
 32:   along with the complete mesh data structure in the native H5M or VTK format.

 34:   Collective

 36:   Input Parameters:
 37: + dm           - the discretization manager object containing solution in MOAB tags.
 38: . filename     - the name of the output file: e.g., poisson.h5m
 39: - usrwriteopts - the parallel write options needed for serializing a MOAB mesh database. Can be NULL.

 41:   Level: intermediate

 43:   Notes:
 44:   The H5M output file can be visualized directly with Paraview (if compiled with appropriate
 45:   plugin) or converted with MOAB/tools/mbconvert to a VTK or Exodus file.

 47:   This routine can also be used for check-pointing purposes to store a complete history of the
 48:   solution along with any other necessary data to restart computations.

 50:   References:
 51: . * - Parallel Mesh Initialization: http://ftp.mcs.anl.gov/pub/fathom/moab-docs/contents.html#fivetwo

 53: .seealso: `DMMoabLoadFromFile()`, `DMMoabSetGlobalFieldVector()`
 54: @*/
 55: PetscErrorCode DMMoabOutput(DM dm, const char *filename, const char *usrwriteopts)
 56: {
 57:   DM_Moab        *dmmoab;
 58:   const char     *writeopts;
 59:   PetscBool       isftype;
 60:   moab::ErrorCode merr;

 62:   PetscFunctionBegin;
 64:   dmmoab = (DM_Moab *)(dm)->data;

 66:   PetscCall(PetscStrendswith(filename, "h5m", &isftype));

 68:   /* add mesh loading options specific to the DM */
 69:   if (isftype) {
 70: #ifdef MOAB_HAVE_MPI
 71:     PetscCall(DMMoab_GetWriteOptions_Private(dmmoab->pcomm->get_id(), dmmoab->pcomm->size(), dmmoab->dim, dmmoab->write_mode, dmmoab->rw_dbglevel, dmmoab->extra_write_options, usrwriteopts, &writeopts));
 72: #else
 73:     PetscCall(DMMoab_GetWriteOptions_Private(0, 1, dmmoab->dim, dmmoab->write_mode, dmmoab->rw_dbglevel, dmmoab->extra_write_options, usrwriteopts, &writeopts));
 74: #endif
 75:     PetscCall(PetscInfo(dm, "Writing file %s with options: %s\n", filename, writeopts));
 76:   } else {
 77:     writeopts = NULL;
 78:   }

 80:   /* output file, using parallel write */
 81:   merr = dmmoab->mbiface->write_file(filename, NULL, writeopts, &dmmoab->fileset, 1);
 82:   MBERRVM(dmmoab->mbiface, "Writing output of DMMoab failed.", merr);
 83:   PetscCall(PetscFree(writeopts));
 84:   PetscFunctionReturn(PETSC_SUCCESS);
 85: }