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