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