xref: /libCEED/examples/fluids/navierstokes.h (revision 7650ae9a66c1de2783569eed4c328204687c633e)
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  output_freq;
93   PetscInt  viz_refine;
94   PetscInt  cont_steps;
95   char      cont_file[PETSC_MAX_PATH_LEN];
96   char      cont_time_file[PETSC_MAX_PATH_LEN];
97   char      output_dir[PETSC_MAX_PATH_LEN];
98   PetscBool add_stepnum2bin;
99   // Problem type arguments
100   PetscFunctionList problems;
101   char              problem_name[PETSC_MAX_PATH_LEN];
102   // Test mode arguments
103   PetscBool   test_mode;
104   PetscScalar test_tol;
105   char        file_path[PETSC_MAX_PATH_LEN];
106 };
107 
108 // libCEED data struct
109 struct CeedData_private {
110   CeedVector    x_coord, q_data;
111   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,
112       qf_apply_outflow_jacobian, qf_apply_freestream, qf_apply_freestream_jacobian;
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 };
117 
118 // PETSc user data
119 struct User_private {
120   MPI_Comm     comm;
121   DM           dm;
122   DM           dm_viz;
123   Mat          interp_viz;
124   Ceed         ceed;
125   Units        units;
126   Vec          M, Q_loc, Q_dot_loc;
127   Physics      phys;
128   AppCtx       app_ctx;
129   CeedVector   q_ceed, q_dot_ceed, g_ceed, coo_values_amat, coo_values_pmat, x_ceed;
130   CeedOperator op_rhs_vol, op_rhs, op_ifunction_vol, op_ifunction, op_ijacobian, op_dirichlet;
131   bool         matrices_set_up;
132   CeedScalar   time, dt;
133 };
134 
135 // Units
136 struct Units_private {
137   // fundamental units
138   PetscScalar meter;
139   PetscScalar kilogram;
140   PetscScalar second;
141   PetscScalar Kelvin;
142   // derived units
143   PetscScalar Pascal;
144   PetscScalar J_per_kg_K;
145   PetscScalar m_per_squared_s;
146   PetscScalar W_per_m_K;
147   PetscScalar Joule;
148 };
149 
150 // Boundary conditions
151 struct SimpleBC_private {
152   PetscInt num_wall,  // Number of faces with wall BCs
153       wall_comps[5],  // An array of constrained component numbers
154       num_comps,
155       num_slip[3],  // Number of faces with slip BCs
156       num_inflow, num_outflow, num_freestream;
157   PetscInt  walls[16], slips[3][16], inflows[16], outflows[16], freestreams[16];
158   PetscBool user_bc;
159 };
160 
161 // Struct that contains all enums and structs used for the physics of all problems
162 struct Physics_private {
163   WindType              wind_type;
164   BubbleType            bubble_type;
165   BubbleContinuityType  bubble_continuity_type;
166   EulerTestType         euler_test;
167   StabilizationType     stab;
168   PetscBool             implicit;
169   StateVariable         state_var;
170   PetscBool             has_curr_time;
171   PetscBool             has_neumann;
172   CeedContextFieldLabel solution_time_label;
173   CeedContextFieldLabel stg_solution_time_label;
174   CeedContextFieldLabel timestep_size_label;
175   CeedContextFieldLabel ics_time_label;
176   CeedContextFieldLabel ijacobian_time_shift_label;
177 };
178 
179 typedef struct {
180   CeedQFunctionUser    qfunction;
181   const char          *qfunction_loc;
182   CeedQFunctionContext qfunction_context;
183 } ProblemQFunctionSpec;
184 
185 // Problem specific data
186 typedef struct ProblemData_private ProblemData;
187 struct ProblemData_private {
188   CeedInt              dim, q_data_size_vol, q_data_size_sur, jac_data_size_sur;
189   CeedScalar           dm_scale;
190   ProblemQFunctionSpec setup_vol, setup_sur, ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow,
191       apply_freestream, apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian;
192   bool non_zero_time;
193   PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
194   void     *bc_ctx;
195   PetscBool bc_from_ics, use_dirichlet_ceed;
196   PetscErrorCode (*print_info)(ProblemData *, AppCtx);
197 };
198 
199 extern int FreeContextPetsc(void *);
200 
201 // -----------------------------------------------------------------------------
202 // Set up problems
203 // -----------------------------------------------------------------------------
204 // Set up function for each problem
205 extern PetscErrorCode NS_NEWTONIAN_WAVE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
206 extern PetscErrorCode NS_CHANNEL(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
207 extern PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
208 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
209 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
210 
211 extern PetscErrorCode NS_EULER_VORTEX(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
212 extern PetscErrorCode NS_SHOCKTUBE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
213 extern PetscErrorCode NS_ADVECTION(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
214 extern PetscErrorCode NS_ADVECTION2D(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
215 
216 // Print function for each problem
217 extern PetscErrorCode PRINT_NEWTONIAN(ProblemData *problem, AppCtx app_ctx);
218 
219 extern PetscErrorCode PRINT_EULER_VORTEX(ProblemData *problem, AppCtx app_ctx);
220 
221 extern PetscErrorCode PRINT_SHOCKTUBE(ProblemData *problem, AppCtx app_ctx);
222 
223 extern PetscErrorCode PRINT_ADVECTION(ProblemData *problem, AppCtx app_ctx);
224 
225 extern PetscErrorCode PRINT_ADVECTION2D(ProblemData *problem, AppCtx app_ctx);
226 
227 // -----------------------------------------------------------------------------
228 // libCEED functions
229 // -----------------------------------------------------------------------------
230 // Utility function - essential BC dofs are encoded in closure indices as -(i+1).
231 PetscInt Involute(PetscInt i);
232 
233 // Utility function to create local CEED restriction
234 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr);
235 
236 // Utility function to get Ceed Restriction for each domain
237 PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, PetscInt value, CeedInt Q, CeedInt q_data_size,
238                                        CeedElemRestriction *elem_restr_q, CeedElemRestriction *elem_restr_x, CeedElemRestriction *elem_restr_qd_i);
239 
240 // Utility function to create CEED Composite Operator for the entire domain
241 PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol,
242                                        CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur,
243                                        CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian);
244 
245 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc);
246 
247 // -----------------------------------------------------------------------------
248 // Time-stepping functions
249 // -----------------------------------------------------------------------------
250 // Compute mass matrix for explicit scheme
251 PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M);
252 
253 // RHS (Explicit time-stepper) function setup
254 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data);
255 
256 // Implicit time-stepper function setup
257 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data);
258 
259 // User provided TS Monitor
260 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx);
261 
262 // TS: Create, setup, and solve
263 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts);
264 
265 // -----------------------------------------------------------------------------
266 // Setup DM
267 // -----------------------------------------------------------------------------
268 // Create mesh
269 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType, VecType, DM *dm);
270 
271 // Set up DM
272 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, SimpleBC bc, Physics phys);
273 
274 // Refine DM for high-order viz
275 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys);
276 
277 // -----------------------------------------------------------------------------
278 // Process command line options
279 // -----------------------------------------------------------------------------
280 // Register problems to be available on the command line
281 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx);
282 
283 // Process general command line options
284 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc);
285 
286 // -----------------------------------------------------------------------------
287 // Miscellaneous utility functions
288 // -----------------------------------------------------------------------------
289 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time);
290 
291 PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
292                                              Vec grad_FVM);
293 
294 // Compare reference solution values with current test run for CI
295 PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q);
296 
297 // Get error for problems with exact solutions
298 PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time);
299 
300 // Post-processing
301 PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time);
302 
303 // -- Gather initial Q values in case of continuation of simulation
304 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q);
305 
306 // Record boundary values from initial condition
307 PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc);
308 
309 // -----------------------------------------------------------------------------
310 // Boundary Condition Related Functions
311 // -----------------------------------------------------------------------------
312 
313 // Setup StrongBCs that use QFunctions
314 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc, CeedInt Q_sur,
315                                   CeedInt q_data_size_sur);
316 
317 PetscErrorCode FreestreamBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx);
318 
319 #endif  // libceed_fluids_examples_navier_stokes_h
320