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