xref: /libCEED/examples/fluids/navierstokes.h (revision c0b5abf0f23b15c4f0ada76f8abe9f8d2b6fa247)
1 // Copyright (c) 2017-2024, 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 #pragma once
8 
9 #include <ceed.h>
10 #include <bc_definition.h>
11 #include <log_events.h>
12 #include <mat-ceed.h>
13 #include <petsc-ceed-utils.h>
14 #include <petscts.h>
15 #include <stdbool.h>
16 
17 #include "./include/petsc_ops.h"
18 #include "qfunctions/newtonian_types.h"
19 
20 #if PETSC_VERSION_LT(3, 22, 0)
21 #error "PETSc v3.22 or later is required"
22 #endif
23 
24 // -----------------------------------------------------------------------------
25 // Enums
26 // -----------------------------------------------------------------------------
27 
28 // Euler - test cases
29 typedef enum {
30   EULER_TEST_ISENTROPIC_VORTEX = 0,
31   EULER_TEST_1                 = 1,
32   EULER_TEST_2                 = 2,
33   EULER_TEST_3                 = 3,
34   EULER_TEST_4                 = 4,
35   EULER_TEST_5                 = 5,
36 } EulerTestType;
37 static const char *const EulerTestTypes[] = {"ISENTROPIC_VORTEX", "1", "2", "3", "4", "5", "EulerTestType", "EULER_TEST_", NULL};
38 
39 // Advection - Wind types
40 static const char *const WindTypes[] = {"ROTATION", "TRANSLATION", "WindType", "WIND_", NULL};
41 
42 // Advection - Initial Condition Types
43 static const char *const AdvectionICTypes[] = {"SPHERE", "CYLINDER", "COSINE_HILL", "SKEW", "AdvectionICType", "ADVECTIONIC_", NULL};
44 
45 // Advection - Bubble Continuity Types
46 static const char *const BubbleContinuityTypes[] = {"SMOOTH", "BACK_SHARP", "THICK", "COSINE", "BubbleContinuityType", "BUBBLE_CONTINUITY_", NULL};
47 
48 // Stabilization methods
49 static const char *const StabilizationTypes[] = {"NONE", "SU", "SUPG", "StabilizationType", "STAB_", NULL};
50 
51 // Stabilization tau constants
52 static const char *const StabilizationTauTypes[] = {"CTAU", "ADVDIFF_SHAKIB", "ADVDIFF_SHAKIB_P", "StabilizationTauType", "STAB_TAU_", NULL};
53 
54 // Test mode type
55 typedef enum {
56   TESTTYPE_NONE           = 0,
57   TESTTYPE_SOLVER         = 1,
58   TESTTYPE_TURB_SPANSTATS = 2,
59   TESTTYPE_DIFF_FILTER    = 3,
60 } TestType;
61 static const char *const TestTypes[] = {"NONE", "SOLVER", "TURB_SPANSTATS", "DIFF_FILTER", "TestType", "TESTTYPE_", NULL};
62 
63 // Mesh transformation type
64 typedef enum {
65   MESH_TRANSFORM_NONE      = 0,
66   MESH_TRANSFORM_PLATEMESH = 1,
67 } MeshTransformType;
68 static const char *const MeshTransformTypes[] = {"NONE", "PLATEMESH", "MeshTransformType", "MESH_TRANSFORM_", NULL};
69 
70 static const char *const DifferentialFilterDampingFunctions[] = {
71     "NONE", "VAN_DRIEST", "MMS", "DifferentialFilterDampingFunction", "DIFF_FILTER_DAMP_", NULL};
72 
73 // -----------------------------------------------------------------------------
74 // Structs
75 // -----------------------------------------------------------------------------
76 // Structs declarations
77 typedef struct AppCtx_private      *AppCtx;
78 typedef struct CeedData_private    *CeedData;
79 typedef struct User_private        *User;
80 typedef struct Units_private       *Units;
81 typedef struct SimpleBC_private    *SimpleBC;
82 typedef struct Physics_private     *Physics;
83 typedef struct ProblemData_private *ProblemData;
84 
85 // Application context from user command line options
86 struct AppCtx_private {
87   // libCEED arguments
88   char     ceed_resource[PETSC_MAX_PATH_LEN];  // libCEED backend
89   PetscInt degree;
90   PetscInt q_extra;
91   // Solver arguments
92   MatType amat_type;
93   // Post-processing arguments
94   PetscInt  checkpoint_interval;
95   PetscInt  viz_refine;
96   PetscInt  cont_steps;
97   PetscReal cont_time;
98   char      cont_file[PETSC_MAX_PATH_LEN];
99   char      cont_time_file[PETSC_MAX_PATH_LEN];
100   char      output_dir[PETSC_MAX_PATH_LEN];
101   PetscBool add_stepnum2bin;
102   PetscBool checkpoint_vtk;
103   // Problem type arguments
104   PetscFunctionList problems;
105   char              problem_name[PETSC_MAX_PATH_LEN];
106   // Test mode arguments
107   TestType    test_type;
108   PetscScalar test_tol;
109   char        test_file_path[PETSC_MAX_PATH_LEN];
110   // Turbulent spanwise statistics
111   PetscBool         turb_spanstats_enable;
112   PetscInt          turb_spanstats_collect_interval;
113   PetscInt          turb_spanstats_viewer_interval;
114   PetscViewer       turb_spanstats_viewer;
115   PetscViewerFormat turb_spanstats_viewer_format;
116   // Wall forces
117   struct {
118     PetscInt          num_wall;
119     PetscInt         *walls;
120     PetscViewer       viewer;
121     PetscViewerFormat viewer_format;
122     PetscBool         header_written;
123   } wall_forces;
124   // Differential Filtering
125   PetscBool         diff_filter_monitor;
126   MeshTransformType mesh_transform_type;
127 };
128 
129 // libCEED data struct
130 struct CeedData_private {
131   CeedVector           x_coord, q_data;
132   CeedBasis            basis_x, basis_q;
133   CeedElemRestriction  elem_restr_x, elem_restr_q, elem_restr_qd_i;
134   OperatorApplyContext op_ics_ctx;
135 };
136 
137 typedef struct {
138   DM                    dm;
139   PetscSF               sf;  // For communicating child data to parents
140   OperatorApplyContext  op_stats_collect_ctx, op_proj_rhs_ctx;
141   PetscInt              num_comp_stats;
142   Vec                   Child_Stats_loc, Parent_Stats_loc;
143   KSP                   ksp;         // For the L^2 projection solve
144   CeedScalar            span_width;  // spanwise width of the child domain
145   PetscBool             do_mms_test;
146   OperatorApplyContext  mms_error_ctx;
147   CeedContextFieldLabel solution_time_label, previous_time_label;
148 } SpanStatsData;
149 
150 typedef struct {
151   DM                   dm;
152   PetscInt             num_comp;
153   OperatorApplyContext l2_rhs_ctx;
154   KSP                  ksp;
155 } *NodalProjectionData;
156 
157 typedef struct {
158   DM                    dm_filter;
159   PetscInt              num_filtered_fields;
160   CeedInt              *num_field_components;
161   PetscInt              field_prim_state, field_velo_prod;
162   OperatorApplyContext  op_rhs_ctx;
163   KSP                   ksp;
164   PetscObjectState      X_loc_state;
165   PetscBool             do_mms_test;
166   CeedContextFieldLabel filter_width_scaling_label;
167 } *DiffFilterData;
168 
169 // PETSc user data
170 struct User_private {
171   MPI_Comm             comm;
172   DM                   dm;
173   DM                   dm_viz;
174   Mat                  interp_viz;
175   Ceed                 ceed;
176   Units                units;
177   Vec                  Q_loc, Q_dot_loc;
178   Physics              phys;
179   AppCtx               app_ctx;
180   CeedVector           q_ceed, q_dot_ceed, g_ceed, x_ceed;
181   CeedOperator         op_ifunction;
182   Mat                  mat_ijacobian;
183   KSP                  mass_ksp;
184   OperatorApplyContext op_rhs_ctx, op_strong_bc_ctx;
185   CeedScalar           time_bc_set;
186   SpanStatsData        spanstats;
187   NodalProjectionData  grad_velo_proj;
188   DiffFilterData       diff_filter;
189 };
190 
191 // Units
192 struct Units_private {
193   // fundamental units
194   PetscScalar meter;
195   PetscScalar kilogram;
196   PetscScalar second;
197   PetscScalar Kelvin;
198   // derived units
199   PetscScalar Pascal;
200   PetscScalar J_per_kg_K;
201   PetscScalar m_per_squared_s;
202   PetscScalar W_per_m_K;
203   PetscScalar Joule;
204 };
205 
206 // Boundary conditions
207 struct SimpleBC_private {
208   PetscInt num_inflow, num_outflow, num_freestream, num_slip;
209   PetscInt inflows[16], outflows[16], freestreams[16], slips[16];
210 };
211 
212 // Struct that contains all enums and structs used for the physics of all problems
213 struct Physics_private {
214   PetscBool             implicit;
215   StateVariable         state_var;
216   CeedContextFieldLabel solution_time_label;
217   CeedContextFieldLabel stg_solution_time_label;
218   CeedContextFieldLabel timestep_size_label;
219   CeedContextFieldLabel ics_time_label;
220 };
221 
222 PetscErrorCode BoundaryConditionSetUp(User user, ProblemData problem, AppCtx app_ctx, SimpleBC bc);
223 
224 typedef struct {
225   CeedQFunctionUser    qfunction;
226   const char          *qfunction_loc;
227   CeedQFunctionContext qfunction_context;
228 } ProblemQFunctionSpec;
229 
230 // Problem specific data
231 struct ProblemData_private {
232   CeedInt              dim, q_data_size_vol, q_data_size_sur, jac_data_size_sur;
233   CeedScalar           dm_scale;
234   ProblemQFunctionSpec ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow, apply_freestream, apply_slip,
235       apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian, apply_slip_jacobian;
236   bool          compute_exact_solution_error;
237   PetscBool     set_bc_from_ics, use_strong_bc_ceed, uses_newtonian;
238   PetscCount    num_bc_defs;
239   BCDefinition *bc_defs;
240   PetscErrorCode (*print_info)(User, ProblemData, AppCtx);
241   PetscErrorCode (*create_mass_operator)(User, CeedOperator *);
242 };
243 
244 extern int FreeContextPetsc(void *);
245 
246 // -----------------------------------------------------------------------------
247 // Set up problems
248 // -----------------------------------------------------------------------------
249 // Set up function for each problem
250 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
251 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
252 extern PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
253 extern PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
254 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
255 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
256 extern PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
257 extern PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
258 extern PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
259 extern PetscErrorCode NS_ADVECTION2D(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
260 
261 // Print function for each problem
262 extern PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx);
263 
264 extern PetscErrorCode PRINT_EULER_VORTEX(User user, ProblemData problem, AppCtx app_ctx);
265 
266 extern PetscErrorCode PRINT_SHOCKTUBE(User user, ProblemData problem, AppCtx app_ctx);
267 
268 extern PetscErrorCode PRINT_ADVECTION(User user, ProblemData problem, AppCtx app_ctx);
269 
270 extern PetscErrorCode PRINT_ADVECTION2D(User user, ProblemData problem, AppCtx app_ctx);
271 
272 PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData problem, TS ts);
273 
274 // -----------------------------------------------------------------------------
275 // libCEED functions
276 // -----------------------------------------------------------------------------
277 // Utility function to create local CEED restriction
278 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt label_value, PetscInt dm_field,
279                                          CeedElemRestriction *elem_restr);
280 
281 PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
282                                                CeedElemRestriction *restriction);
283 PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
284                                                          CeedElemRestriction *restriction);
285 PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
286                                                     PetscInt q_data_size, CeedElemRestriction *restriction);
287 PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
288                                                          PetscInt q_data_size, CeedElemRestriction *restriction);
289 
290 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis);
291 
292 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData problem, SimpleBC bc);
293 
294 PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x,
295                         CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
296 PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size);
297 PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x,
298                                 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
299 PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size);
300 // -----------------------------------------------------------------------------
301 // Time-stepping functions
302 // -----------------------------------------------------------------------------
303 // RHS (Explicit time-stepper) function setup
304 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data);
305 
306 // Implicit time-stepper function setup
307 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data);
308 
309 // User provided TS Monitor
310 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx);
311 
312 // TS: Create, setup, and solve
313 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, ProblemData problem, Vec *Q, PetscScalar *f_time, TS *ts);
314 
315 // Update Boundary Values when time has changed
316 PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t);
317 
318 // -----------------------------------------------------------------------------
319 // Setup DM
320 // -----------------------------------------------------------------------------
321 // Create mesh
322 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType, VecType, DM *dm);
323 
324 // Set up DM
325 PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys);
326 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
327                                        PetscInt num_fields, const PetscInt *field_sizes, DM dm);
328 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm);
329 PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
330                                   PetscInt num_fields, const PetscInt *field_sizes, DM dm);
331 
332 // Refine DM for high-order viz
333 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData problem, SimpleBC bc, Physics phys);
334 
335 // -----------------------------------------------------------------------------
336 // Process command line options
337 // -----------------------------------------------------------------------------
338 // Register problems to be available on the command line
339 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx);
340 
341 // Process general command line options
342 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc);
343 
344 // -----------------------------------------------------------------------------
345 // Miscellaneous utility functions
346 // -----------------------------------------------------------------------------
347 PetscErrorCode GetInverseMultiplicity(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
348                                       PetscBool get_global_multiplicity, CeedElemRestriction *elem_restr_inv_multiplicity,
349                                       CeedVector *inv_multiplicity);
350 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time);
351 
352 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
353                                                   Vec grad_FVM);
354 
355 // Compare reference solution values with current test run for CI
356 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q);
357 
358 // Get error for problems with exact solutions
359 PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time);
360 
361 // Post-processing
362 PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData problem, User user, Vec Q, PetscScalar final_time);
363 
364 // -- Gather initial Q values in case of continuation of simulation
365 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q);
366 
367 // Record boundary values from initial condition
368 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc);
369 
370 // Versioning token for binary checkpoints
371 extern const PetscInt32 FLUIDS_FILE_TOKEN;  // for backwards compatibility
372 extern const PetscInt32 FLUIDS_FILE_TOKEN_32;
373 extern const PetscInt32 FLUIDS_FILE_TOKEN_64;
374 
375 // Create appropriate mass qfunction based on number of components N
376 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf);
377 
378 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context);
379 
380 PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
381                                  FILE **fp);
382 
383 PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows);
384 
385 PetscErrorCode PhastaDatFileReadToArrayReal(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]);
386 
387 // -----------------------------------------------------------------------------
388 // Turbulence Statistics Collection Functions
389 // -----------------------------------------------------------------------------
390 
391 PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem);
392 PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
393 PetscErrorCode TurbulenceStatisticsDestroy(User user, CeedData ceed_data);
394 
395 // -----------------------------------------------------------------------------
396 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions
397 // -----------------------------------------------------------------------------
398 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem, StateVariable state_var_input,
399                                                CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj);
400 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient);
401 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso,
402                                                         CeedVector *grid_aniso_vector);
403 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso,
404                                                              CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso);
405 
406 // -----------------------------------------------------------------------------
407 // Boundary Condition Related Functions
408 // -----------------------------------------------------------------------------
409 
410 // Setup StrongBCs that use QFunctions
411 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, ProblemData problem, SimpleBC bc);
412 
413 PetscErrorCode FreestreamBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
414 PetscErrorCode OutflowBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
415 PetscErrorCode SlipBCSetup(ProblemData problem, DM dm, void *ctx, CeedQFunctionContext newtonian_ig_qfctx);
416 
417 // -----------------------------------------------------------------------------
418 // Differential Filtering Functions
419 // -----------------------------------------------------------------------------
420 
421 PetscErrorCode DifferentialFilterSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem);
422 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData diff_filter);
423 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
424 PetscErrorCode DifferentialFilterApply(User user, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution);
425 PetscErrorCode DifferentialFilterMmsICSetup(ProblemData problem);
426