xref: /honee/include/navierstokes.h (revision 224fc8c867bfdf53209e49be2150ec3eb3d75900)
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, num_outflow, num_freestream, num_slip;
283   PetscInt inflows[16], outflows[16], freestreams[16], slips[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 PetscErrorCode BoundaryConditionSetUp(Honee honee, ProblemData problem, AppCtx app_ctx, SimpleBC bc);
297 
298 typedef struct {
299   CeedQFunctionUser    qf_func_ptr;  // !< QFunction function pointer
300   const char          *qf_loc;       // !< Absolute path to QFunction source file
301   CeedQFunctionContext qfctx;        // !< QFunctionContext to attach to QFunction
302 } ProblemQFunctionSpec;
303 
304 // Problem specific data
305 struct ProblemData_private {
306   CeedInt              jac_data_size_vol, jac_data_size_sur;
307   ProblemQFunctionSpec ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow, apply_freestream, apply_slip,
308       apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian, apply_slip_jacobian;
309   bool          compute_exact_solution_error;
310   PetscBool     set_bc_from_ics, use_strong_bc_ceed;
311   PetscCount    num_bc_defs;
312   BCDefinition *bc_defs;
313   PetscErrorCode (*print_info)(Honee, ProblemData, AppCtx);
314   PetscErrorCode (*create_mass_operator)(Honee, CeedOperator *);
315 };
316 
317 extern int FreeContextPetsc(void *);
318 
319 // -----------------------------------------------------------------------------
320 // Set up problems
321 // -----------------------------------------------------------------------------
322 // Set up function for each problem
323 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
324 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
325 extern PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
326 extern PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
327 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
328 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
329 extern PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
330 extern PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
331 extern PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
332 extern PetscErrorCode NS_ADVECTION2D(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
333 
334 // Print function for each problem
335 extern PetscErrorCode PRINT_NEWTONIAN(Honee honee, ProblemData problem, AppCtx app_ctx);
336 extern PetscErrorCode PRINT_EULER_VORTEX(Honee honee, ProblemData problem, AppCtx app_ctx);
337 extern PetscErrorCode PRINT_SHOCKTUBE(Honee honee, ProblemData problem, AppCtx app_ctx);
338 extern PetscErrorCode PRINT_ADVECTION(Honee honee, ProblemData problem, AppCtx app_ctx);
339 extern PetscErrorCode PRINT_ADVECTION2D(Honee honee, ProblemData problem, AppCtx app_ctx);
340 
341 PetscErrorCode PrintRunInfo(Honee honee, Physics phys_ctx, ProblemData problem, TS ts);
342 
343 // -----------------------------------------------------------------------------
344 // libCEED functions
345 // -----------------------------------------------------------------------------
346 PetscErrorCode SetupLibceed(Ceed ceed, DM dm, Honee honee, AppCtx app_ctx, ProblemData problem, SimpleBC bc);
347 
348 PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x,
349                         CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
350 PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size);
351 PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x,
352                                 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
353 PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size);
354 PetscErrorCode QDataBoundaryGradientGetNumComponents(DM dm, CeedInt *q_data_size);
355 PetscErrorCode QDataBoundaryGradientGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedVector x_coord,
356                                         CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
357 PetscErrorCode QDataClearStoredData();
358 // -----------------------------------------------------------------------------
359 // Time-stepping functions
360 // -----------------------------------------------------------------------------
361 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data);
362 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data);
363 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx);
364 PetscErrorCode TSSolve_NS(DM dm, Honee honee, AppCtx app_ctx, Physics phys, ProblemData problem, Vec Q, PetscScalar *f_time, TS *ts);
365 PetscErrorCode UpdateBoundaryValues(Honee honee, Vec Q_loc, PetscReal t);
366 
367 // -----------------------------------------------------------------------------
368 // Setup DM
369 // -----------------------------------------------------------------------------
370 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType, VecType, DM *dm);
371 PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys);
372 PetscErrorCode VizRefineDM(DM dm, Honee honee, ProblemData problem, SimpleBC bc, Physics phys);
373 
374 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
375                                        PetscInt num_fields, const PetscInt *field_sizes, DM dm);
376 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm);
377 PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
378                                   PetscInt num_fields, const PetscInt *field_sizes, DM dm);
379 
380 // -----------------------------------------------------------------------------
381 // Process command line options
382 // -----------------------------------------------------------------------------
383 PetscErrorCode ProcessCommandLineOptions(Honee honee, SimpleBC bc);
384 PetscErrorCode HoneeOptionsSetValueDefault(PetscOptions options, const char name[], const char value[]);
385 
386 // -----------------------------------------------------------------------------
387 // Miscellaneous utility functions
388 // -----------------------------------------------------------------------------
389 PetscErrorCode GetInverseMultiplicity(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
390                                       PetscBool get_global_multiplicity, CeedElemRestriction *elem_restr_inv_multiplicity,
391                                       CeedVector *inv_multiplicity);
392 PetscErrorCode ICs_FixMultiplicity(DM dm, Honee honee, Vec Q_loc, Vec Q, CeedScalar time);
393 
394 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
395                                                   Vec grad_FVM);
396 
397 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q);
398 PetscErrorCode PrintError(DM dm, Honee honee, Vec Q, PetscScalar final_time);
399 PetscErrorCode PostProcess(TS ts, DM dm, ProblemData problem, Honee honee, Vec Q, PetscScalar final_time);
400 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc);
401 PetscErrorCode HoneeMassQFunctionCreate(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf);
402 PetscErrorCode HoneeCalculateDomainSize(Honee honee, PetscScalar *volume);
403 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context);
404 
405 // -----------------------------------------------------------------------------
406 // Turbulence Statistics Collection Functions
407 // -----------------------------------------------------------------------------
408 PetscErrorCode SpanwiseStatisticsSetup_Turbulence(TS ts, PetscViewerAndFormat *ctx);
409 PetscErrorCode TSMonitor_TurbulenceSpanwiseStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx);
410 
411 // -----------------------------------------------------------------------------
412 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions
413 // -----------------------------------------------------------------------------
414 PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, ProblemData problem);
415 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data);
416 PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc);
417 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, Honee honee, ProblemData problem, StateVariable state_var_input,
418                                                CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj);
419 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient);
420 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso,
421                                                         CeedVector *grid_aniso_vector);
422 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso,
423                                                              CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso);
424 
425 // -----------------------------------------------------------------------------
426 // Boundary Condition Related Functions
427 // -----------------------------------------------------------------------------
428 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, DM dm, Honee honee, ProblemData problem, SimpleBC bc);
429 PetscErrorCode FreestreamBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
430 PetscErrorCode OutflowBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
431 PetscErrorCode SlipBCSetup(ProblemData problem, DM dm, void *ctx, CeedQFunctionContext newtonian_ig_qfctx);
432 
433 // -----------------------------------------------------------------------------
434 // Differential Filtering Functions
435 // -----------------------------------------------------------------------------
436 PetscErrorCode DifferentialFilterSetup(Ceed ceed, Honee honee, ProblemData problem);
437 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData diff_filter);
438 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
439 PetscErrorCode DifferentialFilterApply(Honee honee, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution);
440 PetscErrorCode DifferentialFilterMmsICSetup(ProblemData problem);
441 
442 // -----------------------------------------------------------------------------
443 // SGS Data-Driven Training via SmartSim
444 // -----------------------------------------------------------------------------
445 PetscErrorCode SmartSimSetup(Honee honee);
446 PetscErrorCode SmartSimDataDestroy(SmartSimData smartsim);
447 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, Honee honee, ProblemData problem);
448 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx);
449 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts);
450 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train);
451 
452 // -----------------------------------------------------------------------------
453 // Divergence of Diffusive Flux Projection
454 // -----------------------------------------------------------------------------
455 PetscErrorCode DivDiffFluxProjectionCreate(Honee honee, PetscInt num_diff_flux_comps, DivDiffFluxProjectionData *pdiff_flux_proj);
456 PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis,
457                                                          CeedVector *vector, CeedEvalMode *eval_mode);
458 PetscErrorCode DivDiffFluxProjectionSetup(Honee honee, DivDiffFluxProjectionData diff_flux_proj);
459 PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc);
460 PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj);
461 
462 PetscErrorCode SetupMontiorTotalKineticEnergy(TS ts, PetscViewerAndFormat *ctx);
463 PetscErrorCode TSMonitor_TotalKineticEnergy(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx);
464 
465 PetscErrorCode SetupMontiorCfl(TS ts, PetscViewerAndFormat *ctx);
466 PetscErrorCode TSMonitor_Cfl(TS ts, PetscInt step, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx);
467 
468 PetscErrorCode KSPPostSolve_Honee(KSP ksp, Vec rhs, Vec x, void *ctx);
469