xref: /honee/include/navierstokes.h (revision f978755d57fb6574cfe1ff974b5424124ae3c75e)
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 <ceed.h>
6 #include <bc_definition.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 <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 AdvDifWindTypes[] = {"ROTATION", "TRANSLATION", "BOUNDARY_LAYER", "WindType", "ADVDIF_WIND_", NULL};
41 
42 // Advection - Initial Condition Types
43 static const char *const AdvDifICTypes[] = {"SPHERE",         "CYLINDER",        "COSINE_HILL", "SKEW", "WAVE",
44                                             "BOUNDARY_LAYER", "AdvectionICType", "ADVDIF_IC_",  NULL};
45 
46 // Advection - Wave Types
47 static const char *const AdvDifWaveTypes[] = {"SINE", "SQUARE", "AdvDifWaveType", "ADVDIF_WAVE_", NULL};
48 
49 // Advection - Bubble Continuity Types
50 static const char *const AdvDifBubbleContinuityTypes[] = {
51     "SMOOTH", "BACK_SHARP", "THICK", "COSINE", "BubbleContinuityType", "ADVDIF_BUBBLE_CONTINUITY_", NULL};
52 
53 // Stabilization methods
54 static const char *const StabilizationTypes[] = {"NONE", "SU", "SUPG", "StabilizationType", "STAB_", NULL};
55 
56 // Stabilization tau constants
57 static const char *const StabilizationTauTypes[] = {"CTAU", "ADVDIFF_SHAKIB", "StabilizationTauType", "STAB_TAU_", NULL};
58 
59 // Test mode type
60 typedef enum {
61   TESTTYPE_NONE           = 0,
62   TESTTYPE_SOLVER         = 1,
63   TESTTYPE_TURB_SPANSTATS = 2,
64   TESTTYPE_DIFF_FILTER    = 3,
65 } TestType;
66 static const char *const TestTypes[] = {"NONE", "SOLVER", "TURB_SPANSTATS", "DIFF_FILTER", "TestType", "TESTTYPE_", NULL};
67 
68 // Subgrid-Stress mode type
69 typedef enum {
70   SGS_MODEL_NONE        = 0,
71   SGS_MODEL_DATA_DRIVEN = 1,
72 } SGSModelType;
73 static const char *const SGSModelTypes[] = {"NONE", "DATA_DRIVEN", "SGSModelType", "SGS_MODEL_", NULL};
74 
75 // Subgrid-Stress mode type
76 typedef enum {
77   SGS_MODEL_DD_FUSED           = 0,
78   SGS_MODEL_DD_SEQENTIAL_CEED  = 1,
79   SGS_MODEL_DD_SEQENTIAL_TORCH = 2,
80 } SGSModelDDImplementation;
81 static const char *const SGSModelDDImplementations[] = {"FUSED", "SEQUENTIAL_CEED", "SEQUENTIAL_TORCH", "SGSModelDDImplementation", "SGS_MODEL_DD_",
82                                                         NULL};
83 
84 // Mesh transformation type
85 typedef enum {
86   MESH_TRANSFORM_NONE      = 0,
87   MESH_TRANSFORM_PLATEMESH = 1,
88 } MeshTransformType;
89 static const char *const MeshTransformTypes[] = {"NONE", "PLATEMESH", "MeshTransformType", "MESH_TRANSFORM_", NULL};
90 
91 static const char *const DifferentialFilterDampingFunctions[] = {
92     "NONE", "VAN_DRIEST", "MMS", "DifferentialFilterDampingFunction", "DIFF_FILTER_DAMP_", NULL};
93 
94 static const char *const DivDiffFluxProjectionMethods[] = {"NONE", "DIRECT", "INDIRECT", "DivDiffFluxProjectionMethod", "DIV_DIFF_FLUX_PROJ_", NULL};
95 
96 static const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "ENTROPY", "StateVariable", "STATEVAR_", NULL};
97 
98 // -----------------------------------------------------------------------------
99 // Structs
100 // -----------------------------------------------------------------------------
101 // Structs declarations
102 typedef struct AppCtx_private      *AppCtx;
103 typedef struct Units_private       *Units;
104 typedef struct SimpleBC_private    *SimpleBC;
105 typedef struct Physics_private     *Physics;
106 typedef struct ProblemData_private *ProblemData;
107 
108 // Application context from user command line options
109 struct AppCtx_private {
110   // libCEED arguments
111   char     ceed_resource[PETSC_MAX_PATH_LEN];  // libCEED backend
112   PetscInt degree;
113   PetscInt q_extra;
114   // Solver arguments
115   MatType amat_type;
116   // Post-processing arguments
117   PetscInt  checkpoint_interval;
118   PetscInt  viz_refine;
119   PetscBool use_continue_file;
120   PetscInt  cont_steps;
121   PetscReal cont_time;
122   char      cont_file[PETSC_MAX_PATH_LEN];
123   char      output_dir[PETSC_MAX_PATH_LEN];
124   PetscBool add_stepnum2bin;
125   PetscBool checkpoint_vtk;
126   // Problem type arguments
127   PetscFunctionList problems;
128   char              problem_name[PETSC_MAX_PATH_LEN];
129   // Test mode arguments
130   TestType    test_type;
131   PetscScalar test_tol;
132   char        test_file_path[PETSC_MAX_PATH_LEN];
133   // Wall forces
134   struct {
135     PetscInt          num_wall;
136     PetscInt         *walls;
137     PetscViewer       viewer;
138     PetscViewerFormat viewer_format;
139     PetscBool         header_written;
140   } wall_forces;
141   // Subgrid Stress Model
142   SGSModelType sgs_model_type;
143   PetscBool    sgs_train_enable;
144   // Differential Filtering
145   PetscBool         diff_filter_monitor;
146   MeshTransformType mesh_transform_type;
147   // Divergence of Diffusive Flux Projection
148   DivDiffFluxProjectionMethod divFdiffproj_method;
149 
150   PetscInt check_step_interval;
151 };
152 
153 typedef struct {
154   DM                   dm;
155   PetscInt             num_comp;
156   OperatorApplyContext l2_rhs_ctx;
157   KSP                  ksp;
158 } *NodalProjectionData;
159 
160 typedef struct DivDiffFluxProjectionData_private *DivDiffFluxProjectionData;
161 
162 struct DivDiffFluxProjectionData_private {
163   PetscInt                    num_diff_flux_comps;
164   DivDiffFluxProjectionMethod method;
165   NodalProjectionData         projection;
166 
167   // CeedOperator Objects
168   CeedElemRestriction elem_restr_div_diff_flux;
169   CeedBasis           basis_div_diff_flux;
170   CeedEvalMode        eval_mode_div_diff_flux;
171   CeedVector          div_diff_flux_ceed;
172 
173   // Problem specific setup functions
174   PetscErrorCode (*CreateRHSOperator_Direct)(Honee, DivDiffFluxProjectionData, CeedOperator *);
175   PetscErrorCode (*CreateRHSOperator_Indirect)(Honee, DivDiffFluxProjectionData, CeedOperator *);
176 
177   // Only used for direct method:
178   Vec          DivDiffFlux_loc;
179   PetscMemType DivDiffFlux_memtype;
180   PetscBool    ceed_vec_has_array;
181 
182   // Only used for indirect method:
183   OperatorApplyContext calc_div_diff_flux;
184 };
185 
186 typedef PetscErrorCode (*SgsDDNodalStressEval)(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc);
187 typedef PetscErrorCode (*SgsDDNodalStressInference)(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx);
188 typedef struct {
189   DM                        dm_sgs, dm_dd_inputs, dm_dd_outputs;
190   PetscInt                  num_comp_sgs, num_comp_inputs, num_comp_outputs;
191   OperatorApplyContext      op_nodal_evaluation_ctx, op_nodal_dd_inputs_ctx, op_nodal_dd_outputs_ctx, op_sgs_apply_ctx;
192   CeedVector                sgs_nodal_ceed, grad_velo_ceed;
193   SgsDDNodalStressEval      sgs_nodal_eval;
194   SgsDDNodalStressInference sgs_nodal_inference;
195   void                     *sgs_nodal_inference_ctx;
196   PetscErrorCode (*sgs_nodal_inference_ctx_destroy)(void *ctx);
197 } *SgsDDData;
198 
199 typedef struct {
200   DM                   dm_dd_training;
201   PetscInt             num_comp_dd_inputs, write_data_interval, num_filter_widths;
202   PetscScalar          filter_widths[16];
203   OperatorApplyContext op_training_data_calc_ctx;
204   NodalProjectionData  filtered_grad_velo_proj;
205   size_t               training_data_array_dims[2];
206   PetscBool            overwrite_training_data;
207 } *SGS_DD_TrainingData;
208 
209 typedef struct {
210   DM                    dm_filter;
211   PetscInt              num_filtered_fields;
212   CeedInt              *num_field_components;
213   PetscInt              field_prim_state, field_velo_prod;
214   OperatorApplyContext  op_rhs_ctx;
215   KSP                   ksp;
216   PetscObjectState      X_loc_state;
217   PetscBool             do_mms_test;
218   CeedContextFieldLabel filter_width_scaling_label;
219 } *DiffFilterData;
220 
221 typedef struct {
222   void    *client;
223   char     rank_id_name[16];
224   PetscInt collocated_database_num_ranks;
225 } *SmartSimData;
226 
227 // PETSc user data
228 struct Honee_private {
229   MPI_Comm                  comm;
230   DM                        dm;
231   DM                        dm_viz;
232   Mat                       interp_viz;
233   Ceed                      ceed;
234   Units                     units;
235   Vec                       Q_loc, Q_dot_loc;
236   Physics                   phys;
237   AppCtx                    app_ctx;
238   CeedVector                q_ceed, q_dot_ceed, g_ceed, x_ceed;
239   CeedOperator              op_ifunction;
240   Mat                       mat_ijacobian;
241   KSP                       mass_ksp;
242   OperatorApplyContext      op_rhs_ctx, op_strong_bc_ctx;
243   CeedScalar                time_bc_set;
244   NodalProjectionData       grad_velo_proj;
245   SgsDDData                 sgs_dd_data;
246   DiffFilterData            diff_filter;
247   SmartSimData              smartsim;
248   SGS_DD_TrainingData       sgs_dd_train;
249   DivDiffFluxProjectionData diff_flux_proj;
250 
251   ProblemData problem_data;
252 
253   // old CeedData
254   CeedVector           x_coord;
255   CeedBasis            basis_x, basis_q;
256   CeedElemRestriction  elem_restr_x, elem_restr_q;
257   OperatorApplyContext op_ics_ctx;
258 
259   PetscBool set_poststep;
260   time_t    start_time;
261   time_t    max_wall_time;
262   PetscInt  max_wall_time_interval;
263 };
264 
265 // Units
266 struct Units_private {
267   // fundamental units
268   PetscScalar meter;
269   PetscScalar kilogram;
270   PetscScalar second;
271   PetscScalar Kelvin;
272   // derived units
273   PetscScalar Pascal;
274   PetscScalar J_per_kg_K;
275   PetscScalar m_per_squared_s;
276   PetscScalar W_per_m_K;
277   PetscScalar Joule;
278 };
279 
280 // Boundary conditions
281 struct SimpleBC_private {
282   PetscInt num_inflow;
283   PetscInt inflows[16];
284 };
285 
286 // Struct that contains all enums and structs used for the physics of all problems
287 struct Physics_private {
288   PetscBool             implicit;
289   StateVariable         state_var;
290   CeedContextFieldLabel solution_time_label;
291   CeedContextFieldLabel stg_solution_time_label;
292   CeedContextFieldLabel timestep_size_label;
293   CeedContextFieldLabel ics_time_label;
294 };
295 
296 typedef struct {
297   Honee                honee;
298   CeedInt              jac_data_size_sur;
299   CeedQFunctionContext qfctx;
300   void                *ctx;
301   PetscCtxDestroyFn   *DestroyCtx;
302 } *HoneeBCStruct;
303 
304 PetscErrorCode BoundaryConditionSetUp(Honee honee, ProblemData problem, AppCtx app_ctx, SimpleBC bc);
305 PetscErrorCode HoneeBCDestroy(void **ctx);
306 PetscErrorCode HoneeBCCreateIFunctionQF(BCDefinition bc_def, CeedQFunctionUser qf_func_ptr, const char *qf_loc, CeedQFunctionContext qfctx,
307                                         CeedQFunction *qf_ifunc);
308 PetscErrorCode HoneeBCCreateIJacobianQF(BCDefinition bc_def, CeedQFunctionUser qf_func_ptr, const char *qf_loc, CeedQFunctionContext qfctx,
309                                         CeedQFunction *qf_ijac);
310 PetscErrorCode HoneeBCAddIFunctionOp(BCDefinition bc_def, DMLabel domain_label, PetscInt label_value, CeedQFunction qf_ifunc, CeedOperator op_ifunc,
311                                      CeedOperator *sub_op_ifunc);
312 PetscErrorCode HoneeBCAddIJacobianOp(BCDefinition bc_def, CeedOperator sub_op_ifunc, DMLabel domain_label, PetscInt label_value,
313                                      CeedQFunction qf_ijac, CeedOperator op_ijac);
314 
315 typedef struct {
316   CeedQFunctionUser    qf_func_ptr;  // !< QFunction function pointer
317   const char          *qf_loc;       // !< Absolute path to QFunction source file
318   CeedQFunctionContext qfctx;        // !< QFunctionContext to attach to QFunction
319 } ProblemQFunctionSpec;
320 
321 // Problem specific data
322 struct ProblemData_private {
323   CeedInt              jac_data_size_vol, jac_data_size_sur;
324   ProblemQFunctionSpec ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_inflow_jacobian;
325   bool                 compute_exact_solution_error;
326   PetscBool            set_bc_from_ics, use_strong_bc_ceed;
327   PetscCount           num_bc_defs;
328   BCDefinition        *bc_defs;
329   PetscErrorCode (*print_info)(Honee, ProblemData, AppCtx);
330   PetscErrorCode (*create_mass_operator)(Honee, CeedOperator *);
331 };
332 
333 extern int FreeContextPetsc(void *);
334 
335 // -----------------------------------------------------------------------------
336 // Set up problems
337 // -----------------------------------------------------------------------------
338 // Set up function for each problem
339 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
340 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
341 extern PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
342 extern PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
343 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
344 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
345 extern PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
346 extern PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
347 extern PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
348 extern PetscErrorCode NS_ADVECTION2D(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
349 
350 // Print function for each problem
351 extern PetscErrorCode PRINT_NEWTONIAN(Honee honee, ProblemData problem, AppCtx app_ctx);
352 extern PetscErrorCode PRINT_EULER_VORTEX(Honee honee, ProblemData problem, AppCtx app_ctx);
353 extern PetscErrorCode PRINT_SHOCKTUBE(Honee honee, ProblemData problem, AppCtx app_ctx);
354 extern PetscErrorCode PRINT_ADVECTION(Honee honee, ProblemData problem, AppCtx app_ctx);
355 extern PetscErrorCode PRINT_ADVECTION2D(Honee honee, ProblemData problem, AppCtx app_ctx);
356 
357 PetscErrorCode PrintRunInfo(Honee honee, Physics phys_ctx, ProblemData problem, TS ts);
358 
359 // -----------------------------------------------------------------------------
360 // libCEED functions
361 // -----------------------------------------------------------------------------
362 PetscErrorCode SetupLibceed(Ceed ceed, DM dm, Honee honee, AppCtx app_ctx, ProblemData problem, SimpleBC bc);
363 
364 PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x,
365                         CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
366 PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size);
367 PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x,
368                                 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
369 PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size);
370 PetscErrorCode QDataBoundaryGradientGetNumComponents(DM dm, CeedInt *q_data_size);
371 PetscErrorCode QDataBoundaryGradientGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedVector x_coord,
372                                         CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
373 PetscErrorCode QDataClearStoredData();
374 // -----------------------------------------------------------------------------
375 // Time-stepping functions
376 // -----------------------------------------------------------------------------
377 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data);
378 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data);
379 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx);
380 PetscErrorCode TSSolve_NS(DM dm, Honee honee, AppCtx app_ctx, Physics phys, ProblemData problem, Vec Q, PetscScalar *f_time, TS *ts);
381 PetscErrorCode UpdateBoundaryValues(Honee honee, Vec Q_loc, PetscReal t);
382 
383 // -----------------------------------------------------------------------------
384 // Setup DM
385 // -----------------------------------------------------------------------------
386 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType, VecType, DM *dm);
387 PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys);
388 PetscErrorCode VizRefineDM(DM dm, Honee honee, ProblemData problem, SimpleBC bc, Physics phys);
389 
390 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
391                                        PetscInt num_fields, const PetscInt *field_sizes, DM dm);
392 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm);
393 PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
394                                   PetscInt num_fields, const PetscInt *field_sizes, DM dm);
395 
396 // -----------------------------------------------------------------------------
397 // Process command line options
398 // -----------------------------------------------------------------------------
399 PetscErrorCode ProcessCommandLineOptions(Honee honee, SimpleBC bc);
400 PetscErrorCode HoneeOptionsSetValueDefault(PetscOptions options, const char name[], const char value[]);
401 
402 // -----------------------------------------------------------------------------
403 // Miscellaneous utility functions
404 // -----------------------------------------------------------------------------
405 PetscErrorCode GetInverseMultiplicity(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
406                                       PetscBool get_global_multiplicity, CeedElemRestriction *elem_restr_inv_multiplicity,
407                                       CeedVector *inv_multiplicity);
408 PetscErrorCode ICs_FixMultiplicity(DM dm, Honee honee, Vec Q_loc, Vec Q, CeedScalar time);
409 
410 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
411                                                   Vec grad_FVM);
412 
413 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q);
414 PetscErrorCode PrintError(DM dm, Honee honee, Vec Q, PetscScalar final_time);
415 PetscErrorCode PostProcess(TS ts, DM dm, ProblemData problem, Honee honee, Vec Q, PetscScalar final_time);
416 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc);
417 PetscErrorCode HoneeMassQFunctionCreate(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf);
418 PetscErrorCode HoneeCalculateDomainSize(Honee honee, PetscScalar *volume);
419 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context);
420 
421 // -----------------------------------------------------------------------------
422 // Turbulence Statistics Collection Functions
423 // -----------------------------------------------------------------------------
424 PetscErrorCode SpanwiseStatisticsSetup_Turbulence(TS ts, PetscViewerAndFormat *ctx);
425 PetscErrorCode TSMonitor_TurbulenceSpanwiseStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx);
426 
427 // -----------------------------------------------------------------------------
428 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions
429 // -----------------------------------------------------------------------------
430 PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, ProblemData problem);
431 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data);
432 PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc);
433 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, Honee honee, ProblemData problem, StateVariable state_var_input,
434                                                CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj);
435 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient);
436 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso,
437                                                         CeedVector *grid_aniso_vector);
438 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso,
439                                                              CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso);
440 
441 // -----------------------------------------------------------------------------
442 // Boundary Condition Related Functions
443 // -----------------------------------------------------------------------------
444 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, DM dm, Honee honee, ProblemData problem, SimpleBC bc);
445 PetscErrorCode FreestreamBCSetup(BCDefinition bc_def, ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx,
446                                  const StatePrimitive *reference);
447 PetscErrorCode OutflowBCSetup(BCDefinition bc_def, ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx,
448                               const StatePrimitive *reference);
449 PetscErrorCode SlipBCSetup(BCDefinition bc_def, ProblemData problem, DM dm, void *ctx, CeedQFunctionContext newtonian_ig_qfctx);
450 
451 // -----------------------------------------------------------------------------
452 // Differential Filtering Functions
453 // -----------------------------------------------------------------------------
454 PetscErrorCode DifferentialFilterSetup(Ceed ceed, Honee honee, ProblemData problem);
455 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData diff_filter);
456 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
457 PetscErrorCode DifferentialFilterApply(Honee honee, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution);
458 PetscErrorCode DifferentialFilterMmsICSetup(ProblemData problem);
459 
460 // -----------------------------------------------------------------------------
461 // SGS Data-Driven Training via SmartSim
462 // -----------------------------------------------------------------------------
463 PetscErrorCode SmartSimSetup(Honee honee);
464 PetscErrorCode SmartSimDataDestroy(SmartSimData smartsim);
465 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, Honee honee, ProblemData problem);
466 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx);
467 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts);
468 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train);
469 
470 // -----------------------------------------------------------------------------
471 // Divergence of Diffusive Flux Projection
472 // -----------------------------------------------------------------------------
473 PetscErrorCode DivDiffFluxProjectionCreate(Honee honee, PetscInt num_diff_flux_comps, DivDiffFluxProjectionData *pdiff_flux_proj);
474 PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis,
475                                                          CeedVector *vector, CeedEvalMode *eval_mode);
476 PetscErrorCode DivDiffFluxProjectionSetup(Honee honee, DivDiffFluxProjectionData diff_flux_proj);
477 PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc);
478 PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj);
479 
480 PetscErrorCode SetupMontiorTotalKineticEnergy(TS ts, PetscViewerAndFormat *ctx);
481 PetscErrorCode TSMonitor_TotalKineticEnergy(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx);
482 
483 PetscErrorCode SetupMontiorCfl(TS ts, PetscViewerAndFormat *ctx);
484 PetscErrorCode TSMonitor_Cfl(TS ts, PetscInt step, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx);
485 
486 PetscErrorCode KSPPostSolve_Honee(KSP ksp, Vec rhs, Vec x, void *ctx);
487