xref: /libCEED/examples/fluids/navierstokes.h (revision f179a3328604d643d4b2fd6d7256d83a63ecfb4c)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 #ifndef libceed_fluids_examples_navier_stokes_h
9 #define libceed_fluids_examples_navier_stokes_h
10 
11 #include <ceed.h>
12 #include <petscdm.h>
13 #include <petscdmplex.h>
14 #include <petscsys.h>
15 #include <petscts.h>
16 #include <stdbool.h>
17 
18 #include "qfunctions/newtonian_types.h"
19 #include "qfunctions/stabilization_types.h"
20 
21 // -----------------------------------------------------------------------------
22 // PETSc Version
23 // -----------------------------------------------------------------------------
24 #if PETSC_VERSION_LT(3, 17, 0)
25 #error "PETSc v3.17 or later is required"
26 #endif
27 
28 // -----------------------------------------------------------------------------
29 // Enums
30 // -----------------------------------------------------------------------------
31 // Translate PetscMemType to CeedMemType
32 static inline CeedMemType MemTypeP2C(PetscMemType mem_type) { return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; }
33 
34 // Advection - Wind Options
35 typedef enum {
36   WIND_ROTATION    = 0,
37   WIND_TRANSLATION = 1,
38 } WindType;
39 static const char *const WindTypes[] = {"rotation", "translation", "WindType", "WIND_", NULL};
40 
41 // Advection - Bubble Types
42 typedef enum {
43   BUBBLE_SPHERE   = 0,  // dim=3
44   BUBBLE_CYLINDER = 1,  // dim=2
45 } BubbleType;
46 static const char *const BubbleTypes[] = {"sphere", "cylinder", "BubbleType", "BUBBLE_", NULL};
47 
48 // Advection - Bubble Continuity Types
49 typedef enum {
50   BUBBLE_CONTINUITY_SMOOTH     = 0,  // Original continuous, smooth shape
51   BUBBLE_CONTINUITY_BACK_SHARP = 1,  // Discontinuous, sharp back half shape
52   BUBBLE_CONTINUITY_THICK      = 2,  // Define a finite thickness
53 } BubbleContinuityType;
54 static const char *const BubbleContinuityTypes[] = {"smooth", "back_sharp", "thick", "BubbleContinuityType", "BUBBLE_CONTINUITY_", NULL};
55 
56 // Euler - test cases
57 typedef enum {
58   EULER_TEST_ISENTROPIC_VORTEX = 0,
59   EULER_TEST_1                 = 1,
60   EULER_TEST_2                 = 2,
61   EULER_TEST_3                 = 3,
62   EULER_TEST_4                 = 4,
63   EULER_TEST_5                 = 5,
64 } EulerTestType;
65 static const char *const EulerTestTypes[] = {"isentropic_vortex", "test_1",      "test_2", "test_3", "test_4", "test_5",
66                                              "EulerTestType",     "EULER_TEST_", NULL};
67 
68 // Stabilization methods
69 static const char *const StabilizationTypes[] = {"none", "SU", "SUPG", "StabilizationType", "STAB_", NULL};
70 
71 // -----------------------------------------------------------------------------
72 // Structs
73 // -----------------------------------------------------------------------------
74 // Structs declarations
75 typedef struct AppCtx_private   *AppCtx;
76 typedef struct CeedData_private *CeedData;
77 typedef struct User_private     *User;
78 typedef struct Units_private    *Units;
79 typedef struct SimpleBC_private *SimpleBC;
80 typedef struct Physics_private  *Physics;
81 
82 // Application context from user command line options
83 struct AppCtx_private {
84   // libCEED arguments
85   char     ceed_resource[PETSC_MAX_PATH_LEN];  // libCEED backend
86   PetscInt degree;
87   PetscInt q_extra;
88   // Solver arguments
89   MatType   amat_type;
90   PetscBool pmat_pbdiagonal;
91   // Post-processing arguments
92   PetscInt  checkpoint_interval;
93   PetscInt  viz_refine;
94   PetscInt  cont_steps;
95   PetscReal cont_time;
96   char      cont_file[PETSC_MAX_PATH_LEN];
97   char      cont_time_file[PETSC_MAX_PATH_LEN];
98   char      output_dir[PETSC_MAX_PATH_LEN];
99   PetscBool add_stepnum2bin;
100   PetscBool checkpoint_vtk;
101   // Problem type arguments
102   PetscFunctionList problems;
103   char              problem_name[PETSC_MAX_PATH_LEN];
104   // Test mode arguments
105   PetscBool   test_mode;
106   PetscScalar test_tol;
107   char        file_path[PETSC_MAX_PATH_LEN];
108 };
109 
110 // libCEED data struct
111 struct CeedData_private {
112   CeedVector          x_coord, q_data;
113   CeedBasis           basis_x, basis_xc, basis_q, basis_x_sur, basis_q_sur, basis_xc_sur;
114   CeedElemRestriction elem_restr_x, elem_restr_q, elem_restr_qd_i;
115   CeedOperator        op_setup_vol, op_ics;
116   CeedQFunction       qf_setup_vol, qf_ics, qf_rhs_vol, qf_ifunction_vol, qf_setup_sur, qf_apply_inflow, qf_apply_inflow_jacobian, qf_apply_outflow,
117       qf_apply_outflow_jacobian, qf_apply_freestream, qf_apply_freestream_jacobian;
118   struct {
119     CeedElemRestriction elem_restr_parent_x, elem_restr_parent_stats, elem_restr_parent_qd, elem_restr_parent_colloc, elem_restr_child_colloc;
120     CeedBasis           basis_x, basis_stats;
121     CeedVector          x_coord, q_data;
122   } spanstats;
123 };
124 
125 typedef struct {
126   DM           dm;
127   PetscSF      sf;  // For communicating child data to parents
128   CeedOperator op_stats_collect, op_stats_proj;
129   CeedVector   stats_ceed, child_stats, parent_stats;
130   PetscInt     num_comp_stats;
131 } Span_Stats;
132 
133 // PETSc user data
134 struct User_private {
135   MPI_Comm     comm;
136   DM           dm;
137   DM           dm_viz;
138   Mat          interp_viz;
139   Ceed         ceed;
140   Units        units;
141   Vec          M, Q_loc, Q_dot_loc;
142   Physics      phys;
143   AppCtx       app_ctx;
144   CeedVector   q_ceed, q_dot_ceed, g_ceed, coo_values_amat, coo_values_pmat, x_ceed;
145   CeedOperator op_rhs_vol, op_rhs, op_ifunction_vol, op_ifunction, op_ijacobian, op_dirichlet;
146   bool         matrices_set_up;
147   CeedScalar   time, dt;
148   Span_Stats   spanstats;
149 };
150 
151 // Units
152 struct Units_private {
153   // fundamental units
154   PetscScalar meter;
155   PetscScalar kilogram;
156   PetscScalar second;
157   PetscScalar Kelvin;
158   // derived units
159   PetscScalar Pascal;
160   PetscScalar J_per_kg_K;
161   PetscScalar m_per_squared_s;
162   PetscScalar W_per_m_K;
163   PetscScalar Joule;
164 };
165 
166 // Boundary conditions
167 struct SimpleBC_private {
168   PetscInt num_wall,  // Number of faces with wall BCs
169       wall_comps[5],  // An array of constrained component numbers
170       num_comps,
171       num_slip[3],  // Number of faces with slip BCs
172       num_inflow, num_outflow, num_freestream;
173   PetscInt  walls[16], slips[3][16], inflows[16], outflows[16], freestreams[16];
174   PetscBool user_bc;
175 };
176 
177 // Struct that contains all enums and structs used for the physics of all problems
178 struct Physics_private {
179   WindType              wind_type;
180   BubbleType            bubble_type;
181   BubbleContinuityType  bubble_continuity_type;
182   EulerTestType         euler_test;
183   StabilizationType     stab;
184   PetscBool             implicit;
185   StateVariable         state_var;
186   PetscBool             has_curr_time;
187   PetscBool             has_neumann;
188   CeedContextFieldLabel solution_time_label;
189   CeedContextFieldLabel stg_solution_time_label;
190   CeedContextFieldLabel timestep_size_label;
191   CeedContextFieldLabel ics_time_label;
192   CeedContextFieldLabel ijacobian_time_shift_label;
193 };
194 
195 typedef struct {
196   CeedQFunctionUser    qfunction;
197   const char          *qfunction_loc;
198   CeedQFunctionContext qfunction_context;
199 } ProblemQFunctionSpec;
200 
201 // Problem specific data
202 typedef struct ProblemData_private ProblemData;
203 struct ProblemData_private {
204   CeedInt              dim, q_data_size_vol, q_data_size_sur, jac_data_size_sur;
205   CeedScalar           dm_scale;
206   ProblemQFunctionSpec setup_vol, setup_sur, ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow,
207       apply_freestream, apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian;
208   bool non_zero_time;
209   PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
210   void     *bc_ctx;
211   PetscBool bc_from_ics, use_dirichlet_ceed;
212   PetscErrorCode (*print_info)(ProblemData *, AppCtx);
213 };
214 
215 extern int FreeContextPetsc(void *);
216 
217 // -----------------------------------------------------------------------------
218 // Set up problems
219 // -----------------------------------------------------------------------------
220 // Set up function for each problem
221 extern PetscErrorCode NS_NEWTONIAN_WAVE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
222 extern PetscErrorCode NS_CHANNEL(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
223 extern PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
224 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
225 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
226 extern PetscErrorCode NS_EULER_VORTEX(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
227 extern PetscErrorCode NS_SHOCKTUBE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
228 extern PetscErrorCode NS_ADVECTION(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
229 extern PetscErrorCode NS_ADVECTION2D(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
230 
231 // Print function for each problem
232 extern PetscErrorCode PRINT_NEWTONIAN(ProblemData *problem, AppCtx app_ctx);
233 
234 extern PetscErrorCode PRINT_EULER_VORTEX(ProblemData *problem, AppCtx app_ctx);
235 
236 extern PetscErrorCode PRINT_SHOCKTUBE(ProblemData *problem, AppCtx app_ctx);
237 
238 extern PetscErrorCode PRINT_ADVECTION(ProblemData *problem, AppCtx app_ctx);
239 
240 extern PetscErrorCode PRINT_ADVECTION2D(ProblemData *problem, AppCtx app_ctx);
241 
242 // -----------------------------------------------------------------------------
243 // libCEED functions
244 // -----------------------------------------------------------------------------
245 // Utility function - essential BC dofs are encoded in closure indices as -(i+1).
246 PetscInt Involute(PetscInt i);
247 
248 // Utility function to create local CEED restriction
249 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr);
250 
251 // Utility function to get Ceed Restriction for each domain
252 PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, PetscInt value, CeedInt Q, CeedInt q_data_size,
253                                        CeedElemRestriction *elem_restr_q, CeedElemRestriction *elem_restr_x, CeedElemRestriction *elem_restr_qd_i);
254 
255 // Utility function to create CEED Composite Operator for the entire domain
256 PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol,
257                                        CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur,
258                                        CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian);
259 
260 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc);
261 
262 // -----------------------------------------------------------------------------
263 // Time-stepping functions
264 // -----------------------------------------------------------------------------
265 // Compute mass matrix for explicit scheme
266 PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M);
267 
268 // RHS (Explicit time-stepper) function setup
269 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data);
270 
271 // Implicit time-stepper function setup
272 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data);
273 
274 // User provided TS Monitor
275 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx);
276 
277 // TS: Create, setup, and solve
278 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts);
279 
280 // -----------------------------------------------------------------------------
281 // Setup DM
282 // -----------------------------------------------------------------------------
283 // Create mesh
284 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType, VecType, DM *dm);
285 
286 // Set up DM
287 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, SimpleBC bc, Physics phys);
288 
289 // Refine DM for high-order viz
290 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys);
291 
292 // -----------------------------------------------------------------------------
293 // Process command line options
294 // -----------------------------------------------------------------------------
295 // Register problems to be available on the command line
296 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx);
297 
298 // Process general command line options
299 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc);
300 
301 // -----------------------------------------------------------------------------
302 // Miscellaneous utility functions
303 // -----------------------------------------------------------------------------
304 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time);
305 
306 PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
307                                              Vec grad_FVM);
308 
309 // Compare reference solution values with current test run for CI
310 PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q);
311 
312 // Get error for problems with exact solutions
313 PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time);
314 
315 // Post-processing
316 PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time);
317 
318 // -- Gather initial Q values in case of continuation of simulation
319 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q);
320 
321 // Record boundary values from initial condition
322 PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc);
323 
324 // Versioning token for binary checkpoints
325 extern const PetscInt FLUIDS_FILE_TOKEN;
326 
327 // Create appropriate mass qfunction based on number of components N
328 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf);
329 
330 PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree, SimpleBC bc);
331 
332 PetscErrorCode SetupStatsCollection(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem);
333 
334 // -----------------------------------------------------------------------------
335 // Boundary Condition Related Functions
336 // -----------------------------------------------------------------------------
337 
338 // Setup StrongBCs that use QFunctions
339 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc, CeedInt Q_sur,
340                                   CeedInt q_data_size_sur);
341 
342 PetscErrorCode FreestreamBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
343 PetscErrorCode OutflowBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
344 
345 #endif  // libceed_fluids_examples_navier_stokes_h
346