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