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