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