Actual source code: dmmbio.cxx
petsc-3.8.4 2018-03-24
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: {
7: char *wopts;
8: char wopts_par[PETSC_MAX_PATH_LEN];
9: char wopts_parid[PETSC_MAX_PATH_LEN];
10: char wopts_dbg[PETSC_MAX_PATH_LEN];
13: PetscMalloc1(PETSC_MAX_PATH_LEN, &wopts);
14: PetscMemzero(&wopts_par, PETSC_MAX_PATH_LEN);
15: PetscMemzero(&wopts_parid, PETSC_MAX_PATH_LEN);
16: PetscMemzero(&wopts_dbg, PETSC_MAX_PATH_LEN);
18: // do parallel read unless only one processor
19: if (numproc > 1) {
20: PetscSNPrintf(wopts_par, PETSC_MAX_PATH_LEN, "PARALLEL=%s;", MoabWriteModes[mode]);
21: if (fsetid >= 0) {
22: PetscSNPrintf(wopts_parid, PETSC_MAX_PATH_LEN, "PARALLEL_COMM=%d;", fsetid);
23: }
24: }
26: if (dbglevel) {
27: PetscSNPrintf(wopts_dbg, PETSC_MAX_PATH_LEN, "CPUTIME;DEBUG_IO=%d;", dbglevel);
28: }
30: 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 : ""));
31: *write_opts = wopts;
32: return(0);
33: }
36: /*@C
37: DMMoabOutput - Output the solution vectors that are stored in the DMMoab object as tags
38: along with the complete mesh data structure in the native H5M or VTK format. The H5M output file
39: can be visualized directly with Paraview (if compiled with appropriate plugin) or converted
40: with MOAB/tools/mbconvert to a VTK or Exodus file.
42: This routine can also be used for check-pointing purposes to store a complete history of
43: the solution along with any other necessary data to restart computations.
45: Collective
47: Input Parameters:
48: + dm - the discretization manager object containing solution in MOAB tags.
49: . filename - the name of the output file: e.g., poisson.h5m
50: - usrwriteopts - the parallel write options needed for serializing a MOAB mesh database. Can be NULL.
51: Reference (Parallel Mesh Initialization: http://ftp.mcs.anl.gov/pub/fathom/moab-docs/contents.html#fivetwo)
53: Level: intermediate
55: .keywords: discretization manager, set, component solution
57: .seealso: DMMoabLoadFromFile(), DMMoabSetGlobalFieldVector()
58: @*/
59: PetscErrorCode DMMoabOutput(DM dm, const char* filename, const char* usrwriteopts)
60: {
61: DM_Moab *dmmoab;
62: const char *writeopts;
63: PetscBool isftype;
64: PetscErrorCode ierr;
65: moab::ErrorCode merr;
69: dmmoab = (DM_Moab*)(dm)->data;
71: PetscStrendswith(filename, "h5m", &isftype);
73: /* add mesh loading options specific to the DM */
74: if (isftype) {
75: #ifdef MOAB_HAVE_MPI
76: DMMoab_GetWriteOptions_Private(dmmoab->pcomm->get_id(), dmmoab->pcomm->size(), dmmoab->dim, dmmoab->write_mode,
77: dmmoab->rw_dbglevel, dmmoab->extra_write_options, usrwriteopts, &writeopts);
78: #else
79: DMMoab_GetWriteOptions_Private(0, 1, dmmoab->dim, dmmoab->write_mode,
80: dmmoab->rw_dbglevel, dmmoab->extra_write_options, usrwriteopts, &writeopts);
81: #endif
82: PetscInfo2(dm, "Writing file %s with options: %s\n", filename, writeopts);
83: }
84: else {
85: writeopts = NULL;
86: }
88: /* output file, using parallel write */
89: merr = dmmoab->mbiface->write_file(filename, NULL, writeopts, &dmmoab->fileset, 1); MBERRVM(dmmoab->mbiface, "Writing output of DMMoab failed.", merr);
90: PetscFree(writeopts);
91: return(0);
92: }