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