xref: /libCEED/examples/fluids/navierstokes.h (revision 8fb33541c0006446646a5f49fd7454feff98577f)
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   struct {
134     CeedElemRestriction elem_restr_parent_x, elem_restr_parent_stats, elem_restr_parent_qd, elem_restr_parent_colloc, elem_restr_child_colloc;
135     CeedBasis           basis_x, basis_stats;
136     CeedVector          x_coord, q_data;
137     CeedQFunction       qf_stats_collect, qf_stats_proj;
138   } spanstats;
139 };
140 
141 typedef struct {
142   DM                dm;
143   PetscSF           sf;  // For communicating child data to parents
144   CeedOperator      op_stats_collect, op_stats_proj;
145   PetscInt          num_comp_stats;
146   CeedVector        child_inst_stats, child_stats, parent_stats;  // collocated statistics data
147   CeedVector        rhs_ceed, x_ceed, y_ceed;
148   Vec               M_inv;  // Lumped Mass matrix inverse
149   MatopApplyContext M_ctx, mms_error_ctx;
150   KSP               ksp;         // For the L^2 projection solve
151   CeedScalar        span_width;  // spanwise width of the child domain
152   PetscScalar       prev_time;
153   PetscBool         do_mms_test;
154 } Span_Stats;
155 
156 // PETSc user data
157 struct User_private {
158   MPI_Comm     comm;
159   DM           dm;
160   DM           dm_viz;
161   Mat          interp_viz;
162   Ceed         ceed;
163   Units        units;
164   Vec          M, Q_loc, Q_dot_loc;
165   Physics      phys;
166   AppCtx       app_ctx;
167   CeedVector   q_ceed, q_dot_ceed, g_ceed, coo_values_amat, coo_values_pmat, x_ceed;
168   CeedOperator op_rhs_vol, op_rhs, op_ifunction_vol, op_ifunction, op_ijacobian, op_dirichlet;
169   bool         matrices_set_up;
170   CeedScalar   time, dt, time_bc_set;
171   Span_Stats   spanstats;
172 };
173 
174 // Units
175 struct Units_private {
176   // fundamental units
177   PetscScalar meter;
178   PetscScalar kilogram;
179   PetscScalar second;
180   PetscScalar Kelvin;
181   // derived units
182   PetscScalar Pascal;
183   PetscScalar J_per_kg_K;
184   PetscScalar m_per_squared_s;
185   PetscScalar W_per_m_K;
186   PetscScalar Joule;
187 };
188 
189 // Boundary conditions
190 struct SimpleBC_private {
191   PetscInt num_wall,  // Number of faces with wall BCs
192       wall_comps[5],  // An array of constrained component numbers
193       num_comps,
194       num_slip[3],  // Number of faces with slip BCs
195       num_inflow, num_outflow, num_freestream;
196   PetscInt  walls[16], slips[3][16], inflows[16], outflows[16], freestreams[16];
197   PetscBool user_bc;
198 };
199 
200 // Struct that contains all enums and structs used for the physics of all problems
201 struct Physics_private {
202   WindType              wind_type;
203   BubbleType            bubble_type;
204   BubbleContinuityType  bubble_continuity_type;
205   EulerTestType         euler_test;
206   StabilizationType     stab;
207   PetscBool             implicit;
208   StateVariable         state_var;
209   PetscBool             has_curr_time;
210   PetscBool             has_neumann;
211   CeedContextFieldLabel solution_time_label;
212   CeedContextFieldLabel stg_solution_time_label;
213   CeedContextFieldLabel timestep_size_label;
214   CeedContextFieldLabel ics_time_label;
215   CeedContextFieldLabel ijacobian_time_shift_label;
216 };
217 
218 typedef struct {
219   CeedQFunctionUser    qfunction;
220   const char          *qfunction_loc;
221   CeedQFunctionContext qfunction_context;
222 } ProblemQFunctionSpec;
223 
224 // Problem specific data
225 typedef struct ProblemData_private ProblemData;
226 struct ProblemData_private {
227   CeedInt              dim, q_data_size_vol, q_data_size_sur, jac_data_size_sur;
228   CeedScalar           dm_scale;
229   ProblemQFunctionSpec setup_vol, setup_sur, ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow,
230       apply_freestream, apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian;
231   bool non_zero_time;
232   PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
233   void     *bc_ctx;
234   PetscBool bc_from_ics, use_dirichlet_ceed;
235   PetscErrorCode (*print_info)(ProblemData *, AppCtx);
236 };
237 
238 extern int FreeContextPetsc(void *);
239 
240 // -----------------------------------------------------------------------------
241 // Set up problems
242 // -----------------------------------------------------------------------------
243 // Set up function for each problem
244 extern PetscErrorCode NS_NEWTONIAN_WAVE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
245 extern PetscErrorCode NS_CHANNEL(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
246 extern PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
247 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
248 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
249 extern PetscErrorCode NS_EULER_VORTEX(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
250 extern PetscErrorCode NS_SHOCKTUBE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
251 extern PetscErrorCode NS_ADVECTION(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
252 extern PetscErrorCode NS_ADVECTION2D(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
253 
254 // Print function for each problem
255 extern PetscErrorCode PRINT_NEWTONIAN(ProblemData *problem, AppCtx app_ctx);
256 
257 extern PetscErrorCode PRINT_EULER_VORTEX(ProblemData *problem, AppCtx app_ctx);
258 
259 extern PetscErrorCode PRINT_SHOCKTUBE(ProblemData *problem, AppCtx app_ctx);
260 
261 extern PetscErrorCode PRINT_ADVECTION(ProblemData *problem, AppCtx app_ctx);
262 
263 extern PetscErrorCode PRINT_ADVECTION2D(ProblemData *problem, AppCtx app_ctx);
264 
265 // -----------------------------------------------------------------------------
266 // libCEED functions
267 // -----------------------------------------------------------------------------
268 // Utility function - essential BC dofs are encoded in closure indices as -(i+1).
269 PetscInt Involute(PetscInt i);
270 
271 // Utility function to create local CEED restriction
272 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr);
273 
274 // Utility function to get Ceed Restriction for each domain
275 PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, PetscInt value, CeedInt Q, CeedInt q_data_size,
276                                        CeedElemRestriction *elem_restr_q, CeedElemRestriction *elem_restr_x, CeedElemRestriction *elem_restr_qd_i);
277 
278 // Utility function to create CEED Composite Operator for the entire domain
279 PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol,
280                                        CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur,
281                                        CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian);
282 
283 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc);
284 
285 // -----------------------------------------------------------------------------
286 // Time-stepping functions
287 // -----------------------------------------------------------------------------
288 // Compute mass matrix for explicit scheme
289 PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M);
290 
291 // RHS (Explicit time-stepper) function setup
292 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data);
293 
294 // Implicit time-stepper function setup
295 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data);
296 
297 // User provided TS Monitor
298 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx);
299 
300 // TS: Create, setup, and solve
301 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts);
302 
303 // Update Boundary Values when time has changed
304 PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t);
305 
306 // -----------------------------------------------------------------------------
307 // Setup DM
308 // -----------------------------------------------------------------------------
309 // Create mesh
310 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType, VecType, DM *dm);
311 
312 // Set up DM
313 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, SimpleBC bc, Physics phys);
314 
315 // Refine DM for high-order viz
316 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys);
317 
318 // -----------------------------------------------------------------------------
319 // Process command line options
320 // -----------------------------------------------------------------------------
321 // Register problems to be available on the command line
322 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx);
323 
324 // Process general command line options
325 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc);
326 
327 // -----------------------------------------------------------------------------
328 // Miscellaneous utility functions
329 // -----------------------------------------------------------------------------
330 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time);
331 
332 PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
333                                              Vec grad_FVM);
334 
335 // Compare reference solution values with current test run for CI
336 PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q);
337 
338 // Get error for problems with exact solutions
339 PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time);
340 
341 // Post-processing
342 PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time);
343 
344 // -- Gather initial Q values in case of continuation of simulation
345 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q);
346 
347 // Record boundary values from initial condition
348 PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc);
349 
350 // Versioning token for binary checkpoints
351 extern const PetscInt FLUIDS_FILE_TOKEN;
352 
353 // Create appropriate mass qfunction based on number of components N
354 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf);
355 
356 PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree, SimpleBC bc);
357 
358 PetscErrorCode SetupStatsCollection(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem);
359 
360 PetscErrorCode TSMonitor_Statistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
361 
362 PetscErrorCode CleanupStats(User user, CeedData ceed_data);
363 
364 // -----------------------------------------------------------------------------
365 // Boundary Condition Related Functions
366 // -----------------------------------------------------------------------------
367 
368 // Setup StrongBCs that use QFunctions
369 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc, CeedInt Q_sur,
370                                   CeedInt q_data_size_sur);
371 
372 PetscErrorCode FreestreamBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
373 PetscErrorCode OutflowBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
374 
375 #endif  // libceed_fluids_examples_navier_stokes_h
376