xref: /honee/include/navierstokes.h (revision 07126c9a638dba5b0bb1e3c5c6d0702032cbfb5e)
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 // -----------------------------------------------------------------------------
68 // Structs
69 // -----------------------------------------------------------------------------
70 // Structs declarations
71 typedef struct AppCtx_private      *AppCtx;
72 typedef struct Units_private       *Units;
73 typedef struct SimpleBC_private    *SimpleBC;
74 typedef struct Physics_private     *Physics;
75 typedef struct ProblemData_private *ProblemData;
76 
77 // Application context from user command line options
78 struct AppCtx_private {
79   // libCEED arguments
80   char     ceed_resource[PETSC_MAX_PATH_LEN];  // libCEED backend
81   PetscInt degree;
82   PetscInt q_extra;
83   // Solver arguments
84   MatType amat_type;
85   // Post-processing arguments
86   PetscInt  checkpoint_interval;
87   PetscInt  viz_refine;
88   PetscBool use_continue_file;
89   PetscInt  cont_steps;
90   PetscReal cont_time;
91   char      cont_file[PETSC_MAX_PATH_LEN];
92   char      output_dir[PETSC_MAX_PATH_LEN];
93   PetscBool add_stepnum2bin;
94   PetscBool checkpoint_vtk;
95   // Problem type arguments
96   PetscFunctionList problems;
97   char              problem_name[PETSC_MAX_PATH_LEN];
98   // Test mode arguments
99   TestType    test_type;
100   PetscScalar test_tol;
101   char        test_file_path[PETSC_MAX_PATH_LEN];
102   // Wall forces
103   struct {
104     PetscInt          num_wall;
105     PetscInt         *walls;
106     PetscViewer       viewer;
107     PetscViewerFormat viewer_format;
108     PetscBool         header_written;
109   } wall_forces;
110   // Subgrid Stress Model
111   SGSModelType sgs_model_type;
112   PetscBool    sgs_train_enable;
113 
114   // Divergence of Diffusive Flux Projection
115   DivDiffFluxProjectionMethod divFdiffproj_method;
116 
117   PetscInt check_step_interval;
118 };
119 
120 typedef struct DivDiffFluxProjectionData_ *DivDiffFluxProjectionData;
121 struct DivDiffFluxProjectionData_ {
122   PetscInt                    num_diff_flux_comps;
123   DivDiffFluxProjectionMethod method;
124   NodalProjectionData         projection;
125 
126   // CeedOperator Objects
127   CeedElemRestriction elem_restr_div_diff_flux;
128   CeedBasis           basis_div_diff_flux;
129   CeedEvalMode        eval_mode_div_diff_flux;
130   CeedVector          div_diff_flux_ceed;
131 
132   // Problem specific setup functions
133   PetscErrorCode (*CreateRHSOperator_Direct)(Honee, DivDiffFluxProjectionData, CeedOperator *);
134   PetscErrorCode (*CreateRHSOperator_Indirect)(Honee, DivDiffFluxProjectionData, CeedOperator *);
135 
136   // Only used for direct method:
137   Vec          DivDiffFlux_loc;
138   PetscMemType DivDiffFlux_memtype;
139   PetscBool    ceed_vec_has_array;
140 
141   // Only used for indirect method:
142   OperatorApplyContext calc_div_diff_flux;
143 };
144 
145 typedef struct _HoneeOps *HoneeOps;
146 struct _HoneeOps {};
147 
148 PetscErrorCode HoneeInit(MPI_Comm comm, Honee *honee);
149 PetscErrorCode HoneeDestroy(Honee *honee);
150 
151 // PETSc user data
152 struct Honee_private {
153   PETSCHEADER(struct _HoneeOps);
154   MPI_Comm                  comm;
155   DM                        dm;
156   DM                        dm_viz;
157   Mat                       interp_viz;
158   Ceed                      ceed;
159   Units                     units;
160   Vec                       Q_loc, Q_dot_loc;
161   Physics                   phys;
162   AppCtx                    app_ctx;
163   CeedVector                q_ceed, q_dot_ceed, g_ceed, x_ceed;
164   CeedOperator              op_ifunction;
165   Mat                       mat_ijacobian;
166   KSP                       mass_ksp;
167   OperatorApplyContext      op_rhs_ctx, op_strong_bc_ctx;
168   CeedScalar                time_bc_set;
169   DivDiffFluxProjectionData diff_flux_proj;
170 
171   ProblemData problem_data;
172 
173   // old CeedData
174   CeedVector           x_coord;
175   CeedBasis            basis_x;
176   CeedElemRestriction  elem_restr_x;
177   OperatorApplyContext op_ics_ctx;
178 
179   PetscBool set_poststep;
180   time_t    start_time;
181   time_t    max_wall_time;
182   PetscInt  max_wall_time_interval;
183 };
184 
185 // Units
186 struct Units_private {
187   // fundamental units
188   PetscScalar meter;
189   PetscScalar kilogram;
190   PetscScalar second;
191   PetscScalar Kelvin;
192   // derived units
193   PetscScalar Pascal;
194   PetscScalar J_per_kg_K;
195   PetscScalar m_per_squared_s;
196   PetscScalar W_per_m_K;
197   PetscScalar Joule;
198 };
199 
200 // Struct that contains all enums and structs used for the physics of all problems
201 struct Physics_private {
202   PetscBool             implicit;
203   StateVariable         state_var;
204   CeedContextFieldLabel solution_time_label;
205   CeedContextFieldLabel stg_solution_time_label;
206   CeedContextFieldLabel timestep_size_label;
207   CeedContextFieldLabel ics_time_label;
208 };
209 
210 typedef struct HoneeBCStruct_ *HoneeBCStruct;
211 struct HoneeBCStruct_ {
212   Honee                honee;
213   CeedInt              num_comps_jac_data;
214   CeedQFunctionContext qfctx;
215   void                *ctx;
216   PetscCtxDestroyFn   *DestroyCtx;
217 };
218 
219 PetscErrorCode BoundaryConditionSetUp(Honee honee, ProblemData problem, AppCtx app_ctx);
220 PetscErrorCode HoneeBCDestroy(void **ctx);
221 PetscErrorCode HoneeBCCreateIFunctionQF(BCDefinition bc_def, CeedQFunctionUser qf_func_ptr, const char *qf_loc, CeedQFunctionContext qfctx,
222                                         CeedQFunction *qf_ifunc);
223 PetscErrorCode HoneeBCCreateIJacobianQF(BCDefinition bc_def, CeedQFunctionUser qf_func_ptr, const char *qf_loc, CeedQFunctionContext qfctx,
224                                         CeedQFunction *qf_ijac);
225 PetscErrorCode HoneeBCAddIFunctionOp(BCDefinition bc_def, DMLabel domain_label, PetscInt label_value, CeedQFunction qf_ifunc, CeedOperator op_ifunc,
226                                      CeedOperator *sub_op_ifunc);
227 PetscErrorCode HoneeBCAddIJacobianOp(BCDefinition bc_def, CeedOperator sub_op_ifunc, DMLabel domain_label, PetscInt label_value,
228                                      CeedQFunction qf_ijac, CeedOperator op_ijac);
229 
230 typedef struct {
231   CeedQFunctionUser    qf_func_ptr;  // !< QFunction function pointer
232   const char          *qf_loc;       // !< Absolute path to QFunction source file
233   CeedQFunctionContext qfctx;        // !< QFunctionContext to attach to QFunction
234 } HoneeQFSpec;
235 
236 // Problem specific data
237 struct ProblemData_private {
238   // DM Field Settings
239   PetscInt num_components;
240   char   **component_names;
241 
242   CeedInt     num_comps_jac_data;
243   HoneeQFSpec ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian;
244   bool        compute_exact_solution_error;
245   PetscBool   set_bc_from_ics, use_strong_bc_ceed;
246 
247   // BC Definitions
248   PetscCount    num_bc_defs;
249   BCDefinition *bc_defs;
250 
251   PetscErrorCode (*print_info)(Honee, ProblemData, AppCtx);
252   PetscErrorCode (*create_mass_operator)(Honee, CeedOperator *);
253 };
254 
255 extern int FreeContextPetsc(void *);
256 
257 // -----------------------------------------------------------------------------
258 // Set up problems
259 // -----------------------------------------------------------------------------
260 // Set up function for each problem
261 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData problem, DM dm, void *ctx);
262 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx);
263 extern PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx);
264 extern PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx);
265 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx);
266 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx);
267 extern PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx);
268 extern PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx);
269 extern PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx);
270 
271 PetscErrorCode PrintRunInfo(Honee honee, Physics phys_ctx, ProblemData problem, TS ts);
272 
273 // -----------------------------------------------------------------------------
274 // libCEED functions
275 // -----------------------------------------------------------------------------
276 PetscErrorCode SetupLibceed(Ceed ceed, DM dm, Honee honee, AppCtx app_ctx, ProblemData problem);
277 
278 PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction *elem_restr_qd, CeedVector *q_data,
279                         CeedInt *q_data_size);
280 PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size);
281 PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction *elem_restr_qd, CeedVector *q_data,
282                                 CeedInt *q_data_size);
283 PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size);
284 PetscErrorCode QDataBoundaryGradientGetNumComponents(DM dm, CeedInt *q_data_size);
285 PetscErrorCode QDataBoundaryGradientGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedVector x_coord,
286                                         CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
287 PetscErrorCode QDataClearStoredData();
288 // -----------------------------------------------------------------------------
289 // Time-stepping functions
290 // -----------------------------------------------------------------------------
291 PetscErrorCode TSSolve_NS(DM dm, Honee honee, AppCtx app_ctx, Physics phys, ProblemData problem, Vec Q, PetscScalar *f_time, TS *ts);
292 PetscErrorCode UpdateBoundaryValues(Honee honee, Vec Q_loc, PetscReal t);
293 
294 // -----------------------------------------------------------------------------
295 // Setup DM
296 // -----------------------------------------------------------------------------
297 PetscErrorCode CreateDM(Honee honee, ProblemData problem, MatType, VecType, DM *dm);
298 PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, Physics phys);
299 PetscErrorCode VizRefineDM(DM dm, Honee honee, ProblemData problem, Physics phys);
300 
301 // -----------------------------------------------------------------------------
302 // Process command line options
303 // -----------------------------------------------------------------------------
304 PetscErrorCode ProcessCommandLineOptions(Honee honee);
305 PetscErrorCode HoneeOptionsSetValueDefault(PetscOptions options, const char name[], const char value[]);
306 
307 // -----------------------------------------------------------------------------
308 // Miscellaneous utility functions
309 // -----------------------------------------------------------------------------
310 PetscErrorCode GetInverseMultiplicity(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
311                                       PetscBool get_global_multiplicity, CeedElemRestriction *elem_restr_inv_multiplicity,
312                                       CeedVector *inv_multiplicity);
313 PetscErrorCode ICs_FixMultiplicity(DM dm, Honee honee, Vec Q_loc, Vec Q, CeedScalar time);
314 
315 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
316                                                   Vec grad_FVM);
317 
318 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q);
319 PetscErrorCode PrintError(DM dm, Honee honee, Vec Q, PetscScalar final_time);
320 PetscErrorCode PostProcess(TS ts, DM dm, ProblemData problem, Honee honee, Vec Q, PetscScalar final_time);
321 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc);
322 PetscErrorCode HoneeMassQFunctionCreate(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf);
323 PetscErrorCode HoneeCalculateDomainSize(Honee honee, PetscScalar *volume);
324 
325 // -----------------------------------------------------------------------------
326 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions
327 // -----------------------------------------------------------------------------
328 PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, ProblemData problem);
329 PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc);
330 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, Honee honee, ProblemData problem, StateVariable state_var_input,
331                                                CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj);
332 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient);
333 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso,
334                                                         CeedVector *grid_aniso_vector);
335 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso,
336                                                              CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso);
337 
338 // -----------------------------------------------------------------------------
339 // Boundary Condition Related Functions
340 // -----------------------------------------------------------------------------
341 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, DM dm, Honee honee, ProblemData problem);
342 PetscErrorCode FreestreamBCSetup(BCDefinition bc_def, ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx,
343                                  const StatePrimitive *reference);
344 PetscErrorCode OutflowBCSetup(BCDefinition bc_def, ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx,
345                               const StatePrimitive *reference);
346 PetscErrorCode SlipBCSetup(BCDefinition bc_def, ProblemData problem, DM dm, void *ctx, CeedQFunctionContext newtonian_ig_qfctx);
347 
348 // -----------------------------------------------------------------------------
349 // Divergence of Diffusive Flux Projection
350 // -----------------------------------------------------------------------------
351 PetscErrorCode DivDiffFluxProjectionCreate(Honee honee, DivDiffFluxProjectionMethod divFdiffproj_method, PetscInt num_diff_flux_comps,
352                                            DivDiffFluxProjectionData *pdiff_flux_proj);
353 PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis,
354                                                          CeedVector *vector, CeedEvalMode *eval_mode);
355 PetscErrorCode DivDiffFluxProjectionSetup(Honee honee, DivDiffFluxProjectionData diff_flux_proj);
356 PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc);
357 PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj);
358 
359 PetscErrorCode SetupMontiorTotalKineticEnergy(TS ts, PetscViewerAndFormat *ctx);
360 PetscErrorCode TSMonitor_TotalKineticEnergy(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx);
361 
362 PetscErrorCode SetupMontiorCfl(TS ts, PetscViewerAndFormat *ctx);
363 PetscErrorCode TSMonitor_Cfl(TS ts, PetscInt step, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx);
364 
365 PetscErrorCode KSPPostSolve_Honee(KSP ksp, Vec rhs, Vec x, void *ctx);
366