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