xref: /petsc/src/dm/impls/plex/tutorials/ex16.c (revision ff018026d6c90175dcf01f9821e7a548ec4226a8)
1*ff018026SKenneth E. Jansen #include "petscsf.h"
2*ff018026SKenneth E. Jansen static char help[] = "Simple demonstration of CGNS parallel load-save including data\n\n";
3*ff018026SKenneth E. Jansen // As this is a tutorial that is intended to be an easy starting point feel free to make new
4*ff018026SKenneth E. Jansen // example files that extend this but please keep this one simple.
5*ff018026SKenneth E. Jansen // In subsequent examples we will also provide tools to generate an arbitrary size initital
6*ff018026SKenneth E. Jansen // CGNS file to support performance benchmarking.
7*ff018026SKenneth E. Jansen 
8*ff018026SKenneth E. Jansen #include <petscdmplex.h>
9*ff018026SKenneth E. Jansen #include <petscviewerhdf5.h>
10*ff018026SKenneth E. Jansen #define EX "ex16.c"
11*ff018026SKenneth E. Jansen 
12*ff018026SKenneth E. Jansen typedef struct {
13*ff018026SKenneth E. Jansen   char infile[PETSC_MAX_PATH_LEN]; /* Input mesh filename */
14*ff018026SKenneth E. Jansen } AppCtx;
15*ff018026SKenneth E. Jansen 
16*ff018026SKenneth E. Jansen static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
17*ff018026SKenneth E. Jansen {
18*ff018026SKenneth E. Jansen   PetscFunctionBeginUser;
19*ff018026SKenneth E. Jansen   options->infile[0] = '\0';
20*ff018026SKenneth E. Jansen   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
21*ff018026SKenneth E. Jansen   PetscCall(PetscOptionsString("-infile", "The input CGNS file", EX, options->infile, options->infile, sizeof(options->infile), NULL));
22*ff018026SKenneth E. Jansen   PetscOptionsEnd();
23*ff018026SKenneth E. Jansen   PetscCheck(options->infile[0], comm, PETSC_ERR_USER_INPUT, "-infile needs to be specified");
24*ff018026SKenneth E. Jansen   PetscFunctionReturn(PETSC_SUCCESS);
25*ff018026SKenneth E. Jansen }
26*ff018026SKenneth E. Jansen 
27*ff018026SKenneth E. Jansen // @brief Create DM from CGNS file and setup PetscFE to VecLoad solution from that file
28*ff018026SKenneth E. Jansen PetscErrorCode ReadCGNSDM(MPI_Comm comm, const char filename[], DM *dm)
29*ff018026SKenneth E. Jansen {
30*ff018026SKenneth E. Jansen   PetscInt degree;
31*ff018026SKenneth E. Jansen 
32*ff018026SKenneth E. Jansen   PetscFunctionBeginUser;
33*ff018026SKenneth E. Jansen   PetscCall(DMPlexCreateFromFile(comm, filename, "ex16_plex", PETSC_TRUE, dm));
34*ff018026SKenneth E. Jansen   PetscCall(DMSetFromOptions(*dm));
35*ff018026SKenneth E. Jansen   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
36*ff018026SKenneth E. Jansen 
37*ff018026SKenneth E. Jansen   { // Get degree of the natural section
38*ff018026SKenneth E. Jansen     PetscFE        fe_natural;
39*ff018026SKenneth E. Jansen     PetscDualSpace dual_space_natural;
40*ff018026SKenneth E. Jansen 
41*ff018026SKenneth E. Jansen     PetscCall(DMGetField(*dm, 0, NULL, (PetscObject *)&fe_natural));
42*ff018026SKenneth E. Jansen     PetscCall(PetscFEGetDualSpace(fe_natural, &dual_space_natural));
43*ff018026SKenneth E. Jansen     PetscCall(PetscDualSpaceGetOrder(dual_space_natural, &degree));
44*ff018026SKenneth E. Jansen     PetscCall(DMClearFields(*dm));
45*ff018026SKenneth E. Jansen     PetscCall(DMSetLocalSection(*dm, NULL));
46*ff018026SKenneth E. Jansen   }
47*ff018026SKenneth E. Jansen 
48*ff018026SKenneth E. Jansen   { // Setup fe to load in the initial condition data
49*ff018026SKenneth E. Jansen     PetscFE        fe;
50*ff018026SKenneth E. Jansen     PetscInt       dim, cStart, cEnd;
51*ff018026SKenneth E. Jansen     PetscInt       ctInt, mincti, maxcti;
52*ff018026SKenneth E. Jansen     DMPolytopeType dm_polytope, cti;
53*ff018026SKenneth E. Jansen 
54*ff018026SKenneth E. Jansen     PetscCall(DMGetDimension(*dm, &dim));
55*ff018026SKenneth E. Jansen     // Limiting to single topology in this simple example
56*ff018026SKenneth E. Jansen     PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd));
57*ff018026SKenneth E. Jansen     PetscCall(DMPlexGetCellType(*dm, cStart, &dm_polytope));
58*ff018026SKenneth E. Jansen     for (PetscInt i = cStart + 1; i < cEnd; i++) {
59*ff018026SKenneth E. Jansen       PetscCall(DMPlexGetCellType(*dm, i, &cti));
60*ff018026SKenneth E. Jansen       PetscCheck(cti == dm_polytope, comm, PETSC_ERR_RETURN, "Multi-topology not yet supported in this example!");
61*ff018026SKenneth E. Jansen     }
62*ff018026SKenneth E. Jansen     ctInt = cti;
63*ff018026SKenneth E. Jansen     PetscCallMPI(MPIU_Allreduce(&ctInt, &maxcti, 1, MPIU_INT, MPI_MAX, comm));
64*ff018026SKenneth E. Jansen     PetscCallMPI(MPIU_Allreduce(&ctInt, &mincti, 1, MPIU_INT, MPI_MIN, comm));
65*ff018026SKenneth E. Jansen     PetscCheck(mincti == maxcti, comm, PETSC_ERR_RETURN, "Multi-topology not yet supported in this example!");
66*ff018026SKenneth E. Jansen     PetscCall(PetscPrintf(comm, "Mesh confirmed to be single topology degree %" PetscInt_FMT " %s\n", degree, DMPolytopeTypes[cti]));
67*ff018026SKenneth E. Jansen     PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, 5, dm_polytope, degree, PETSC_DETERMINE, &fe));
68*ff018026SKenneth E. Jansen     PetscCall(PetscObjectSetName((PetscObject)fe, "FE for VecLoad"));
69*ff018026SKenneth E. Jansen     PetscCall(DMAddField(*dm, NULL, (PetscObject)fe));
70*ff018026SKenneth E. Jansen     PetscCall(DMCreateDS(*dm));
71*ff018026SKenneth E. Jansen     PetscCall(PetscFEDestroy(&fe));
72*ff018026SKenneth E. Jansen   }
73*ff018026SKenneth E. Jansen 
74*ff018026SKenneth E. Jansen   // Set section component names, used when writing out CGNS files
75*ff018026SKenneth E. Jansen   PetscSection section;
76*ff018026SKenneth E. Jansen   PetscCall(DMGetLocalSection(*dm, &section));
77*ff018026SKenneth E. Jansen   PetscCall(PetscSectionSetFieldName(section, 0, ""));
78*ff018026SKenneth E. Jansen   PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure"));
79*ff018026SKenneth E. Jansen   PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX"));
80*ff018026SKenneth E. Jansen   PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY"));
81*ff018026SKenneth E. Jansen   PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ"));
82*ff018026SKenneth E. Jansen   PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature"));
83*ff018026SKenneth E. Jansen   PetscFunctionReturn(PETSC_SUCCESS);
84*ff018026SKenneth E. Jansen }
85*ff018026SKenneth E. Jansen 
86*ff018026SKenneth E. Jansen int main(int argc, char **argv)
87*ff018026SKenneth E. Jansen {
88*ff018026SKenneth E. Jansen   AppCtx      user;
89*ff018026SKenneth E. Jansen   MPI_Comm    comm;
90*ff018026SKenneth E. Jansen   const char *infilename;
91*ff018026SKenneth E. Jansen 
92*ff018026SKenneth E. Jansen   PetscFunctionBeginUser;
93*ff018026SKenneth E. Jansen   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
94*ff018026SKenneth E. Jansen   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
95*ff018026SKenneth E. Jansen   infilename = user.infile;
96*ff018026SKenneth E. Jansen 
97*ff018026SKenneth E. Jansen   DM          dm;
98*ff018026SKenneth E. Jansen   Vec         V;
99*ff018026SKenneth E. Jansen   PetscViewer viewer;
100*ff018026SKenneth E. Jansen   const char *name;
101*ff018026SKenneth E. Jansen   PetscReal   time;
102*ff018026SKenneth E. Jansen   PetscBool   set;
103*ff018026SKenneth E. Jansen   comm = PETSC_COMM_WORLD;
104*ff018026SKenneth E. Jansen 
105*ff018026SKenneth E. Jansen   // Load DM from CGNS file
106*ff018026SKenneth E. Jansen   PetscCall(ReadCGNSDM(comm, infilename, &dm));
107*ff018026SKenneth E. Jansen   PetscCall(DMSetOptionsPrefix(dm, "loaded_"));
108*ff018026SKenneth E. Jansen   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
109*ff018026SKenneth E. Jansen 
110*ff018026SKenneth E. Jansen   // Load solution from CGNS file
111*ff018026SKenneth E. Jansen   PetscCall(PetscViewerCGNSOpen(comm, infilename, FILE_MODE_READ, &viewer));
112*ff018026SKenneth E. Jansen   PetscCall(DMGetGlobalVector(dm, &V));
113*ff018026SKenneth E. Jansen   PetscCall(PetscViewerCGNSSetSolutionIndex(viewer, 1));
114*ff018026SKenneth E. Jansen   PetscCall(PetscViewerCGNSGetSolutionName(viewer, &name));
115*ff018026SKenneth E. Jansen   PetscCall(PetscViewerCGNSGetSolutionTime(viewer, &time, &set));
116*ff018026SKenneth E. Jansen   PetscCall(PetscPrintf(comm, "Solution Name: %s, and time %g\n", name, (double)time));
117*ff018026SKenneth E. Jansen   PetscCall(VecLoad(V, viewer));
118*ff018026SKenneth E. Jansen   PetscCall(PetscViewerDestroy(&viewer));
119*ff018026SKenneth E. Jansen 
120*ff018026SKenneth E. Jansen   // Write loaded solution (e.g. example in TEST below is to CGNS file)
121*ff018026SKenneth E. Jansen   PetscCall(VecViewFromOptions(V, NULL, "-vec_view"));
122*ff018026SKenneth E. Jansen 
123*ff018026SKenneth E. Jansen   PetscCall(DMRestoreGlobalVector(dm, &V));
124*ff018026SKenneth E. Jansen   PetscCall(DMDestroy(&dm));
125*ff018026SKenneth E. Jansen 
126*ff018026SKenneth E. Jansen   PetscCall(PetscFinalize());
127*ff018026SKenneth E. Jansen   return 0;
128*ff018026SKenneth E. Jansen }
129*ff018026SKenneth E. Jansen 
130*ff018026SKenneth E. Jansen /*TEST
131*ff018026SKenneth E. Jansen   build:
132*ff018026SKenneth E. Jansen     requires: cgns
133*ff018026SKenneth E. Jansen   testset:
134*ff018026SKenneth E. Jansen     suffix: cgns
135*ff018026SKenneth E. Jansen     requires: !complex
136*ff018026SKenneth E. Jansen     nsize: 4
137*ff018026SKenneth E. Jansen     args: -infile ${wPETSC_DIR}/share/petsc/datafiles/meshes/2x2x2_Q3_wave.cgns
138*ff018026SKenneth E. Jansen     args: -dm_plex_cgns_parallel -loaded_dm_view
139*ff018026SKenneth E. Jansen     test:
140*ff018026SKenneth E. Jansen       suffix: simple
141*ff018026SKenneth E. Jansen       args: -vec_view cgns:2x2x2_Q3Vecview.cgns
142*ff018026SKenneth E. Jansen TEST*/
143