xref: /honee/include/navierstokes.h (revision fcefbc99808dbd582b56f5c769be3dc8c9d159d6)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 #pragma once
4 
5 #include <bc_definition.h>
6 #include <ceed.h>
7 #include <dm-utils.h>
8 #include <honee-file.h>
9 #include <honee.h>
10 #include <log_events.h>
11 #include <mat-ceed.h>
12 #include <petsc-ceed-utils.h>
13 #include <petscts.h>
14 #include <stdbool.h>
15 #include <time.h>
16 
17 #include <nodal_projection.h>
18 
19 #include <petsc_ops.h>
20 #include "../qfunctions/newtonian_types.h"
21 
22 #if PETSC_VERSION_LT(3, 23, 0)
23 #error "PETSc v3.23 or later is required"
24 #endif
25 
26 // -----------------------------------------------------------------------------
27 // Enums
28 // -----------------------------------------------------------------------------
29 
30 // Euler - test cases
31 typedef enum {
32   EULER_TEST_ISENTROPIC_VORTEX = 0,
33   EULER_TEST_1                 = 1,
34   EULER_TEST_2                 = 2,
35   EULER_TEST_3                 = 3,
36   EULER_TEST_4                 = 4,
37   EULER_TEST_5                 = 5,
38 } EulerTestType;
39 static const char *const EulerTestTypes[] = {"ISENTROPIC_VORTEX", "1", "2", "3", "4", "5", "EulerTestType", "EULER_TEST_", NULL};
40 
41 // Test mode type
42 typedef enum {
43   TESTTYPE_NONE           = 0,
44   TESTTYPE_SOLVER         = 1,
45   TESTTYPE_TURB_SPANSTATS = 2,
46   TESTTYPE_DIFF_FILTER    = 3,
47   TESTTYPE_SPANSTATS      = 4,
48 } TestType;
49 static const char *const TestTypes[] = {"NONE", "SOLVER", "TURB_SPANSTATS", "DIFF_FILTER", "SPANSTATS", "TestType", "TESTTYPE_", NULL};
50 
51 // Subgrid-Stress mode type
52 typedef enum {
53   SGS_MODEL_NONE        = 0,
54   SGS_MODEL_DATA_DRIVEN = 1,
55 } SGSModelType;
56 static const char *const SGSModelTypes[] = {"NONE", "DATA_DRIVEN", "SGSModelType", "SGS_MODEL_", NULL};
57 
58 // Subgrid-Stress mode type
59 typedef enum {
60   SGS_MODEL_DD_FUSED           = 0,
61   SGS_MODEL_DD_SEQENTIAL_CEED  = 1,
62   SGS_MODEL_DD_SEQENTIAL_TORCH = 2,
63 } SGSModelDDImplementation;
64 static const char *const SGSModelDDImplementations[] = {"FUSED", "SEQUENTIAL_CEED", "SEQUENTIAL_TORCH", "SGSModelDDImplementation", "SGS_MODEL_DD_",
65                                                         NULL};
66 
67 // Mesh transformation type
68 typedef enum {
69   MESH_TRANSFORM_NONE      = 0,
70   MESH_TRANSFORM_PLATEMESH = 1,
71 } MeshTransformType;
72 static const char *const MeshTransformTypes[] = {"NONE", "PLATEMESH", "MeshTransformType", "MESH_TRANSFORM_", NULL};
73 
74 // -----------------------------------------------------------------------------
75 // Structs
76 // -----------------------------------------------------------------------------
77 // Structs declarations
78 typedef struct AppCtx_private      *AppCtx;
79 typedef struct Units_private       *Units;
80 typedef struct SimpleBC_private    *SimpleBC;
81 typedef struct Physics_private     *Physics;
82 typedef struct ProblemData_private *ProblemData;
83 
84 // Application context from user command line options
85 struct AppCtx_private {
86   // libCEED arguments
87   char     ceed_resource[PETSC_MAX_PATH_LEN];  // libCEED backend
88   PetscInt degree;
89   PetscInt q_extra;
90   // Solver arguments
91   MatType amat_type;
92   // Post-processing arguments
93   PetscInt  checkpoint_interval;
94   PetscInt  viz_refine;
95   PetscBool use_continue_file;
96   PetscInt  cont_steps;
97   PetscReal cont_time;
98   char      cont_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   TestType    test_type;
107   PetscScalar test_tol;
108   char        test_file_path[PETSC_MAX_PATH_LEN];
109   // Wall forces
110   struct {
111     PetscInt          num_wall;
112     PetscInt         *walls;
113     PetscViewer       viewer;
114     PetscViewerFormat viewer_format;
115     PetscBool         header_written;
116   } wall_forces;
117   // Subgrid Stress Model
118   SGSModelType sgs_model_type;
119   PetscBool    sgs_train_enable;
120 
121   MeshTransformType mesh_transform_type;
122   // Divergence of Diffusive Flux Projection
123   DivDiffFluxProjectionMethod divFdiffproj_method;
124 
125   PetscInt check_step_interval;
126 };
127 
128 typedef struct DivDiffFluxProjectionData_ *DivDiffFluxProjectionData;
129 struct DivDiffFluxProjectionData_ {
130   PetscInt                    num_diff_flux_comps;
131   DivDiffFluxProjectionMethod method;
132   NodalProjectionData         projection;
133 
134   // CeedOperator Objects
135   CeedElemRestriction elem_restr_div_diff_flux;
136   CeedBasis           basis_div_diff_flux;
137   CeedEvalMode        eval_mode_div_diff_flux;
138   CeedVector          div_diff_flux_ceed;
139 
140   // Problem specific setup functions
141   PetscErrorCode (*CreateRHSOperator_Direct)(Honee, DivDiffFluxProjectionData, CeedOperator *);
142   PetscErrorCode (*CreateRHSOperator_Indirect)(Honee, DivDiffFluxProjectionData, CeedOperator *);
143 
144   // Only used for direct method:
145   Vec          DivDiffFlux_loc;
146   PetscMemType DivDiffFlux_memtype;
147   PetscBool    ceed_vec_has_array;
148 
149   // Only used for indirect method:
150   OperatorApplyContext calc_div_diff_flux;
151 };
152 
153 typedef struct _HoneeOps *HoneeOps;
154 struct _HoneeOps {};
155 
156 PetscErrorCode HoneeInit(MPI_Comm comm, Honee *honee);
157 PetscErrorCode HoneeDestroy(Honee *honee);
158 
159 // PETSc user data
160 struct Honee_private {
161   PETSCHEADER(struct _HoneeOps);
162   MPI_Comm                  comm;
163   DM                        dm;
164   DM                        dm_viz;
165   Mat                       interp_viz;
166   Ceed                      ceed;
167   Units                     units;
168   Vec                       Q_loc, Q_dot_loc;
169   Physics                   phys;
170   AppCtx                    app_ctx;
171   CeedVector                q_ceed, q_dot_ceed, g_ceed, x_ceed;
172   CeedOperator              op_ifunction;
173   Mat                       mat_ijacobian;
174   KSP                       mass_ksp;
175   OperatorApplyContext      op_rhs_ctx, op_strong_bc_ctx;
176   CeedScalar                time_bc_set;
177   DivDiffFluxProjectionData diff_flux_proj;
178 
179   ProblemData problem_data;
180 
181   // old CeedData
182   CeedVector           x_coord;
183   CeedBasis            basis_x, basis_q;
184   CeedElemRestriction  elem_restr_x, elem_restr_q;
185   OperatorApplyContext op_ics_ctx;
186 
187   PetscBool set_poststep;
188   time_t    start_time;
189   time_t    max_wall_time;
190   PetscInt  max_wall_time_interval;
191 };
192 
193 // Units
194 struct Units_private {
195   // fundamental units
196   PetscScalar meter;
197   PetscScalar kilogram;
198   PetscScalar second;
199   PetscScalar Kelvin;
200   // derived units
201   PetscScalar Pascal;
202   PetscScalar J_per_kg_K;
203   PetscScalar m_per_squared_s;
204   PetscScalar W_per_m_K;
205   PetscScalar Joule;
206 };
207 
208 // Struct that contains all enums and structs used for the physics of all problems
209 struct Physics_private {
210   PetscBool             implicit;
211   StateVariable         state_var;
212   CeedContextFieldLabel solution_time_label;
213   CeedContextFieldLabel stg_solution_time_label;
214   CeedContextFieldLabel timestep_size_label;
215   CeedContextFieldLabel ics_time_label;
216 };
217 
218 typedef struct HoneeBCStruct_ *HoneeBCStruct;
219 struct HoneeBCStruct_ {
220   Honee                honee;
221   CeedInt              num_comps_jac_data;
222   CeedQFunctionContext qfctx;
223   void                *ctx;
224   PetscCtxDestroyFn   *DestroyCtx;
225 };
226 
227 PetscErrorCode BoundaryConditionSetUp(Honee honee, ProblemData problem, AppCtx app_ctx);
228 PetscErrorCode HoneeBCDestroy(void **ctx);
229 PetscErrorCode HoneeBCCreateIFunctionQF(BCDefinition bc_def, CeedQFunctionUser qf_func_ptr, const char *qf_loc, CeedQFunctionContext qfctx,
230                                         CeedQFunction *qf_ifunc);
231 PetscErrorCode HoneeBCCreateIJacobianQF(BCDefinition bc_def, CeedQFunctionUser qf_func_ptr, const char *qf_loc, CeedQFunctionContext qfctx,
232                                         CeedQFunction *qf_ijac);
233 PetscErrorCode HoneeBCAddIFunctionOp(BCDefinition bc_def, DMLabel domain_label, PetscInt label_value, CeedQFunction qf_ifunc, CeedOperator op_ifunc,
234                                      CeedOperator *sub_op_ifunc);
235 PetscErrorCode HoneeBCAddIJacobianOp(BCDefinition bc_def, CeedOperator sub_op_ifunc, DMLabel domain_label, PetscInt label_value,
236                                      CeedQFunction qf_ijac, CeedOperator op_ijac);
237 
238 typedef struct {
239   CeedQFunctionUser    qf_func_ptr;  // !< QFunction function pointer
240   const char          *qf_loc;       // !< Absolute path to QFunction source file
241   CeedQFunctionContext qfctx;        // !< QFunctionContext to attach to QFunction
242 } HoneeQFSpec;
243 
244 // Problem specific data
245 struct ProblemData_private {
246   // DM Field Settings
247   PetscInt num_components;
248   char   **component_names;
249 
250   CeedInt     num_comps_jac_data;
251   HoneeQFSpec ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian;
252   bool        compute_exact_solution_error;
253   PetscBool   set_bc_from_ics, use_strong_bc_ceed;
254 
255   // BC Definitions
256   PetscCount    num_bc_defs;
257   BCDefinition *bc_defs;
258 
259   PetscErrorCode (*print_info)(Honee, ProblemData, AppCtx);
260   PetscErrorCode (*create_mass_operator)(Honee, CeedOperator *);
261 };
262 
263 extern int FreeContextPetsc(void *);
264 
265 // -----------------------------------------------------------------------------
266 // Set up problems
267 // -----------------------------------------------------------------------------
268 // Set up function for each problem
269 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData problem, DM dm, void *ctx);
270 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx);
271 extern PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx);
272 extern PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx);
273 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx);
274 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx);
275 extern PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx);
276 extern PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx);
277 extern PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx);
278 
279 PetscErrorCode PrintRunInfo(Honee honee, Physics phys_ctx, ProblemData problem, TS ts);
280 
281 // -----------------------------------------------------------------------------
282 // libCEED functions
283 // -----------------------------------------------------------------------------
284 PetscErrorCode SetupLibceed(Ceed ceed, DM dm, Honee honee, AppCtx app_ctx, ProblemData problem);
285 
286 PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x,
287                         CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
288 PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size);
289 PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x,
290                                 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
291 PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size);
292 PetscErrorCode QDataBoundaryGradientGetNumComponents(DM dm, CeedInt *q_data_size);
293 PetscErrorCode QDataBoundaryGradientGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedVector x_coord,
294                                         CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
295 PetscErrorCode QDataClearStoredData();
296 // -----------------------------------------------------------------------------
297 // Time-stepping functions
298 // -----------------------------------------------------------------------------
299 PetscErrorCode TSSolve_NS(DM dm, Honee honee, AppCtx app_ctx, Physics phys, ProblemData problem, Vec Q, PetscScalar *f_time, TS *ts);
300 PetscErrorCode UpdateBoundaryValues(Honee honee, Vec Q_loc, PetscReal t);
301 
302 // -----------------------------------------------------------------------------
303 // Setup DM
304 // -----------------------------------------------------------------------------
305 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType, VecType, DM *dm);
306 PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, Physics phys);
307 PetscErrorCode VizRefineDM(DM dm, Honee honee, ProblemData problem, Physics phys);
308 
309 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
310                                        PetscInt num_fields, const PetscInt *field_sizes, DM dm);
311 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm);
312 PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
313                                   PetscInt num_fields, const PetscInt *field_sizes, DM dm);
314 
315 // -----------------------------------------------------------------------------
316 // Process command line options
317 // -----------------------------------------------------------------------------
318 PetscErrorCode ProcessCommandLineOptions(Honee honee);
319 PetscErrorCode HoneeOptionsSetValueDefault(PetscOptions options, const char name[], const char value[]);
320 
321 // -----------------------------------------------------------------------------
322 // Miscellaneous utility functions
323 // -----------------------------------------------------------------------------
324 PetscErrorCode GetInverseMultiplicity(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
325                                       PetscBool get_global_multiplicity, CeedElemRestriction *elem_restr_inv_multiplicity,
326                                       CeedVector *inv_multiplicity);
327 PetscErrorCode ICs_FixMultiplicity(DM dm, Honee honee, Vec Q_loc, Vec Q, CeedScalar time);
328 
329 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
330                                                   Vec grad_FVM);
331 
332 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q);
333 PetscErrorCode PrintError(DM dm, Honee honee, Vec Q, PetscScalar final_time);
334 PetscErrorCode PostProcess(TS ts, DM dm, ProblemData problem, Honee honee, Vec Q, PetscScalar final_time);
335 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc);
336 PetscErrorCode HoneeMassQFunctionCreate(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf);
337 PetscErrorCode HoneeCalculateDomainSize(Honee honee, PetscScalar *volume);
338 
339 // -----------------------------------------------------------------------------
340 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions
341 // -----------------------------------------------------------------------------
342 PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, ProblemData problem);
343 PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc);
344 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, Honee honee, ProblemData problem, StateVariable state_var_input,
345                                                CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj);
346 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient);
347 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso,
348                                                         CeedVector *grid_aniso_vector);
349 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso,
350                                                              CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso);
351 
352 // -----------------------------------------------------------------------------
353 // Boundary Condition Related Functions
354 // -----------------------------------------------------------------------------
355 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, DM dm, Honee honee, ProblemData problem);
356 PetscErrorCode FreestreamBCSetup(BCDefinition bc_def, ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx,
357                                  const StatePrimitive *reference);
358 PetscErrorCode OutflowBCSetup(BCDefinition bc_def, ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx,
359                               const StatePrimitive *reference);
360 PetscErrorCode SlipBCSetup(BCDefinition bc_def, ProblemData problem, DM dm, void *ctx, CeedQFunctionContext newtonian_ig_qfctx);
361 
362 // -----------------------------------------------------------------------------
363 // Divergence of Diffusive Flux Projection
364 // -----------------------------------------------------------------------------
365 PetscErrorCode DivDiffFluxProjectionCreate(Honee honee, DivDiffFluxProjectionMethod divFdiffproj_method, PetscInt num_diff_flux_comps,
366                                            DivDiffFluxProjectionData *pdiff_flux_proj);
367 PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis,
368                                                          CeedVector *vector, CeedEvalMode *eval_mode);
369 PetscErrorCode DivDiffFluxProjectionSetup(Honee honee, DivDiffFluxProjectionData diff_flux_proj);
370 PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc);
371 PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj);
372 
373 PetscErrorCode SetupMontiorTotalKineticEnergy(TS ts, PetscViewerAndFormat *ctx);
374 PetscErrorCode TSMonitor_TotalKineticEnergy(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx);
375 
376 PetscErrorCode SetupMontiorCfl(TS ts, PetscViewerAndFormat *ctx);
377 PetscErrorCode TSMonitor_Cfl(TS ts, PetscInt step, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx);
378 
379 PetscErrorCode KSPPostSolve_Honee(KSP ksp, Vec rhs, Vec x, void *ctx);
380