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