xref: /libCEED/examples/fluids/navierstokes.h (revision 3737f8325a1e09ea5db4e0ce05b04230e5d3648e)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 #ifndef libceed_fluids_examples_navier_stokes_h
9 #define libceed_fluids_examples_navier_stokes_h
10 
11 #include <ceed.h>
12 #include <petscts.h>
13 #include <stdbool.h>
14 
15 #include "./include/petsc_ops.h"
16 #include "qfunctions/newtonian_types.h"
17 #include "qfunctions/stabilization_types.h"
18 
19 #if PETSC_VERSION_LT(3, 20, 0)
20 #error "PETSc v3.20 or later is required"
21 #endif
22 
23 #if PETSC_VERSION_LT(3, 21, 0)
24 #define DMSetCoordinateDisc(a, b, c) DMProjectCoordinates(a, b)
25 #endif
26 
27 #define PetscCallCeed(ceed, ...)                                    \
28   do {                                                              \
29     int ierr = __VA_ARGS__;                                         \
30     if (ierr != CEED_ERROR_SUCCESS) {                               \
31       const char *error_message;                                    \
32       CeedGetErrorMessage(ceed, &error_message);                    \
33       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "%s", error_message); \
34     }                                                               \
35   } while (0)
36 
37 // -----------------------------------------------------------------------------
38 // Enums
39 // -----------------------------------------------------------------------------
40 // Translate PetscMemType to CeedMemType
41 static inline CeedMemType MemTypeP2C(PetscMemType mem_type) { return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; }
42 
43 // Euler - test cases
44 typedef enum {
45   EULER_TEST_ISENTROPIC_VORTEX = 0,
46   EULER_TEST_1                 = 1,
47   EULER_TEST_2                 = 2,
48   EULER_TEST_3                 = 3,
49   EULER_TEST_4                 = 4,
50   EULER_TEST_5                 = 5,
51 } EulerTestType;
52 static const char *const EulerTestTypes[] = {"isentropic_vortex", "test_1",      "test_2", "test_3", "test_4", "test_5",
53                                              "EulerTestType",     "EULER_TEST_", NULL};
54 
55 // Advection - Wind types
56 static const char *const WindTypes[] = {"rotation", "translation", "WindType", "WIND_", NULL};
57 
58 // Advection - Initial Condition Types
59 static const char *const AdvectionICTypes[] = {"sphere", "cylinder", "cosine_hill", "skew", "AdvectionICType", "ADVECTIONIC_", NULL};
60 
61 // Advection - Bubble Continuity Types
62 static const char *const BubbleContinuityTypes[] = {"smooth", "back_sharp", "thick", "cosine", "BubbleContinuityType", "BUBBLE_CONTINUITY_", NULL};
63 
64 // Stabilization methods
65 static const char *const StabilizationTypes[] = {"none", "SU", "SUPG", "StabilizationType", "STAB_", NULL};
66 
67 // Test mode type
68 typedef enum {
69   TESTTYPE_NONE           = 0,
70   TESTTYPE_SOLVER         = 1,
71   TESTTYPE_TURB_SPANSTATS = 2,
72   TESTTYPE_DIFF_FILTER    = 3,
73 } TestType;
74 static const char *const TestTypes[] = {"none", "solver", "turb_spanstats", "diff_filter", "TestType", "TESTTYPE_", NULL};
75 
76 // Subgrid-Stress mode type
77 typedef enum {
78   SGS_MODEL_NONE        = 0,
79   SGS_MODEL_DATA_DRIVEN = 1,
80 } SGSModelType;
81 static const char *const SGSModelTypes[] = {"none", "data_driven", "SGSModelType", "SGS_MODEL_", NULL};
82 
83 // Mesh transformation type
84 typedef enum {
85   MESH_TRANSFORM_NONE      = 0,
86   MESH_TRANSFORM_PLATEMESH = 1,
87 } MeshTransformType;
88 static const char *const MeshTransformTypes[] = {"none", "platemesh", "MeshTransformType", "MESH_TRANSFORM_", NULL};
89 
90 static const char *const DifferentialFilterDampingFunctions[] = {
91     "none", "van_driest", "mms", "DifferentialFilterDampingFunction", "DIFF_FILTER_DAMP_", NULL};
92 
93 // -----------------------------------------------------------------------------
94 // Log Events
95 // -----------------------------------------------------------------------------
96 extern PetscLogEvent FLUIDS_CeedOperatorApply;
97 extern PetscLogEvent FLUIDS_CeedOperatorAssemble;
98 extern PetscLogEvent FLUIDS_CeedOperatorAssembleDiagonal;
99 extern PetscLogEvent FLUIDS_CeedOperatorAssemblePointBlockDiagonal;
100 extern PetscLogEvent FLUIDS_SmartRedis_Init;
101 extern PetscLogEvent FLUIDS_SmartRedis_Meta;
102 extern PetscLogEvent FLUIDS_SmartRedis_Train;
103 extern PetscLogEvent FLUIDS_TrainDataCompute;
104 extern PetscLogEvent FLUIDS_DifferentialFilter;
105 extern PetscLogEvent FLUIDS_VelocityGradientProjection;
106 PetscErrorCode       RegisterLogEvents();
107 
108 // -----------------------------------------------------------------------------
109 // Structs
110 // -----------------------------------------------------------------------------
111 // Structs declarations
112 typedef struct AppCtx_private   *AppCtx;
113 typedef struct CeedData_private *CeedData;
114 typedef struct User_private     *User;
115 typedef struct Units_private    *Units;
116 typedef struct SimpleBC_private *SimpleBC;
117 typedef struct Physics_private  *Physics;
118 
119 // Application context from user command line options
120 struct AppCtx_private {
121   // libCEED arguments
122   char     ceed_resource[PETSC_MAX_PATH_LEN];  // libCEED backend
123   PetscInt degree;
124   PetscInt q_extra;
125   // Solver arguments
126   MatType   amat_type;
127   PetscBool pmat_pbdiagonal;
128   // Post-processing arguments
129   PetscInt  checkpoint_interval;
130   PetscInt  viz_refine;
131   PetscInt  cont_steps;
132   PetscReal cont_time;
133   char      cont_file[PETSC_MAX_PATH_LEN];
134   char      cont_time_file[PETSC_MAX_PATH_LEN];
135   char      output_dir[PETSC_MAX_PATH_LEN];
136   PetscBool add_stepnum2bin;
137   PetscBool checkpoint_vtk;
138   // Problem type arguments
139   PetscFunctionList problems;
140   char              problem_name[PETSC_MAX_PATH_LEN];
141   // Test mode arguments
142   TestType    test_type;
143   PetscScalar test_tol;
144   char        test_file_path[PETSC_MAX_PATH_LEN];
145   // Turbulent spanwise statistics
146   PetscBool         turb_spanstats_enable;
147   PetscInt          turb_spanstats_collect_interval;
148   PetscInt          turb_spanstats_viewer_interval;
149   PetscViewer       turb_spanstats_viewer;
150   PetscViewerFormat turb_spanstats_viewer_format;
151   // Wall forces
152   struct {
153     PetscInt          num_wall;
154     PetscInt         *walls;
155     PetscViewer       viewer;
156     PetscViewerFormat viewer_format;
157     PetscBool         header_written;
158   } wall_forces;
159   // Subgrid Stress Model
160   SGSModelType sgs_model_type;
161   PetscBool    sgs_train_enable;
162   // Differential Filtering
163   PetscBool         diff_filter_monitor;
164   MeshTransformType mesh_transform_type;
165 };
166 
167 // libCEED data struct
168 struct CeedData_private {
169   CeedVector           x_coord, q_data;
170   CeedBasis            basis_x, basis_xc, basis_q, basis_x_sur, basis_q_sur, basis_xc_sur;
171   CeedElemRestriction  elem_restr_x, elem_restr_q, elem_restr_qd_i;
172   CeedOperator         op_setup_vol;
173   OperatorApplyContext op_ics_ctx;
174   CeedQFunction        qf_setup_vol, qf_ics, qf_rhs_vol, qf_ifunction_vol, qf_setup_sur, qf_apply_inflow, qf_apply_inflow_jacobian, qf_apply_outflow,
175       qf_apply_outflow_jacobian, qf_apply_freestream, qf_apply_freestream_jacobian;
176 };
177 
178 typedef struct {
179   DM                    dm;
180   PetscSF               sf;  // For communicating child data to parents
181   OperatorApplyContext  op_stats_collect_ctx, op_proj_rhs_ctx;
182   PetscInt              num_comp_stats;
183   Vec                   Child_Stats_loc, Parent_Stats_loc;
184   KSP                   ksp;         // For the L^2 projection solve
185   CeedScalar            span_width;  // spanwise width of the child domain
186   PetscBool             do_mms_test;
187   OperatorApplyContext  mms_error_ctx;
188   CeedContextFieldLabel solution_time_label, previous_time_label;
189 } SpanStatsData;
190 
191 typedef struct {
192   DM                   dm;
193   PetscInt             num_comp;
194   OperatorApplyContext l2_rhs_ctx;
195   KSP                  ksp;
196 } *NodalProjectionData;
197 
198 typedef struct {
199   DM                   dm_sgs;
200   PetscInt             num_comp_sgs;
201   OperatorApplyContext op_nodal_evaluation_ctx, op_sgs_apply_ctx;
202   CeedVector           sgs_nodal_ceed;
203 } *SgsDDData;
204 
205 typedef struct {
206   DM                   dm_dd_training;
207   PetscInt             num_comp_dd_inputs, write_data_interval;
208   OperatorApplyContext op_training_data_calc_ctx;
209   NodalProjectionData  filtered_grad_velo_proj;
210   size_t               training_data_array_dims[2];
211   PetscBool            overwrite_training_data;
212 } *SGS_DD_TrainingData;
213 
214 typedef struct {
215   DM                   dm_filter;
216   PetscInt             num_filtered_fields;
217   CeedInt             *num_field_components;
218   PetscInt             field_prim_state, field_velo_prod;
219   OperatorApplyContext op_rhs_ctx;
220   KSP                  ksp;
221   PetscBool            do_mms_test;
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                  M_inv, Q_loc, Q_dot_loc;
239   Physics              phys;
240   AppCtx               app_ctx;
241   CeedVector           q_ceed, q_dot_ceed, g_ceed, coo_values_amat, coo_values_pmat, x_ceed;
242   CeedOperator         op_rhs_vol, op_ifunction_vol, op_ifunction, op_ijacobian;
243   OperatorApplyContext op_rhs_ctx, op_strong_bc_ctx;
244   bool                 matrices_set_up;
245   CeedScalar           time_bc_set;
246   SpanStatsData        spanstats;
247   NodalProjectionData  grad_velo_proj;
248   SgsDDData            sgs_dd_data;
249   DiffFilterData       diff_filter;
250   SmartSimData         smartsim;
251   SGS_DD_TrainingData  sgs_dd_train;
252 };
253 
254 // Units
255 struct Units_private {
256   // fundamental units
257   PetscScalar meter;
258   PetscScalar kilogram;
259   PetscScalar second;
260   PetscScalar Kelvin;
261   // derived units
262   PetscScalar Pascal;
263   PetscScalar J_per_kg_K;
264   PetscScalar m_per_squared_s;
265   PetscScalar W_per_m_K;
266   PetscScalar Joule;
267 };
268 
269 // Boundary conditions
270 struct SimpleBC_private {
271   PetscInt num_wall,  // Number of faces with wall BCs
272       wall_comps[5],  // An array of constrained component numbers
273       num_comps,
274       num_slip[3],  // Number of faces with slip BCs
275       num_inflow, num_outflow, num_freestream;
276   PetscInt  walls[16], slips[3][16], inflows[16], outflows[16], freestreams[16];
277   PetscBool user_bc;
278 };
279 
280 // Struct that contains all enums and structs used for the physics of all problems
281 struct Physics_private {
282   PetscBool             implicit;
283   StateVariable         state_var;
284   CeedContextFieldLabel solution_time_label;
285   CeedContextFieldLabel stg_solution_time_label;
286   CeedContextFieldLabel timestep_size_label;
287   CeedContextFieldLabel ics_time_label;
288   CeedContextFieldLabel ijacobian_time_shift_label;
289 };
290 
291 typedef struct {
292   CeedQFunctionUser    qfunction;
293   const char          *qfunction_loc;
294   CeedQFunctionContext qfunction_context;
295 } ProblemQFunctionSpec;
296 
297 // Problem specific data
298 typedef struct ProblemData_private ProblemData;
299 struct ProblemData_private {
300   CeedInt              dim, q_data_size_vol, q_data_size_sur, jac_data_size_sur;
301   CeedScalar           dm_scale;
302   ProblemQFunctionSpec setup_vol, setup_sur, ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow,
303       apply_freestream, apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian;
304   bool      non_zero_time;
305   PetscBool bc_from_ics, use_strong_bc_ceed, uses_newtonian;
306   PetscErrorCode (*print_info)(User, ProblemData *, AppCtx);
307 };
308 
309 extern int FreeContextPetsc(void *);
310 
311 // -----------------------------------------------------------------------------
312 // Set up problems
313 // -----------------------------------------------------------------------------
314 // Set up function for each problem
315 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
316 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
317 extern PetscErrorCode NS_CHANNEL(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
318 extern PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
319 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
320 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
321 extern PetscErrorCode NS_EULER_VORTEX(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
322 extern PetscErrorCode NS_SHOCKTUBE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
323 extern PetscErrorCode NS_ADVECTION(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
324 extern PetscErrorCode NS_ADVECTION2D(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
325 
326 // Print function for each problem
327 extern PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData *problem, AppCtx app_ctx);
328 
329 extern PetscErrorCode PRINT_EULER_VORTEX(User user, ProblemData *problem, AppCtx app_ctx);
330 
331 extern PetscErrorCode PRINT_SHOCKTUBE(User user, ProblemData *problem, AppCtx app_ctx);
332 
333 extern PetscErrorCode PRINT_ADVECTION(User user, ProblemData *problem, AppCtx app_ctx);
334 
335 extern PetscErrorCode PRINT_ADVECTION2D(User user, ProblemData *problem, AppCtx app_ctx);
336 
337 PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData *problem, MPI_Comm comm);
338 
339 // -----------------------------------------------------------------------------
340 // libCEED functions
341 // -----------------------------------------------------------------------------
342 // Utility function to create local CEED restriction
343 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt label_value, PetscInt dm_field,
344                                          CeedElemRestriction *elem_restr);
345 
346 PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
347                                                CeedElemRestriction *restriction);
348 PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
349                                                          CeedElemRestriction *restriction);
350 PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
351                                                     PetscInt q_data_size, CeedElemRestriction *restriction);
352 PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
353                                                          PetscInt q_data_size, CeedElemRestriction *restriction);
354 
355 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis);
356 
357 // Utility function to create CEED Composite Operator for the entire domain
358 PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol,
359                                        CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur,
360                                        CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian);
361 
362 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc);
363 
364 // -----------------------------------------------------------------------------
365 // Time-stepping functions
366 // -----------------------------------------------------------------------------
367 // Compute mass matrix for explicit scheme
368 PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M);
369 
370 // RHS (Explicit time-stepper) function setup
371 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data);
372 
373 // Implicit time-stepper function setup
374 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data);
375 
376 // User provided TS Monitor
377 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx);
378 
379 // TS: Create, setup, and solve
380 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts);
381 
382 // Update Boundary Values when time has changed
383 PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t);
384 
385 // -----------------------------------------------------------------------------
386 // Setup DM
387 // -----------------------------------------------------------------------------
388 // Create mesh
389 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType, VecType, DM *dm);
390 
391 // Set up DM
392 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys);
393 PetscErrorCode DMSetupByOrderBegin_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 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm);
396 PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
397                                   PetscInt num_fields, const PetscInt *field_sizes, DM dm);
398 
399 // Refine DM for high-order viz
400 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys);
401 
402 // -----------------------------------------------------------------------------
403 // Process command line options
404 // -----------------------------------------------------------------------------
405 // Register problems to be available on the command line
406 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx);
407 
408 // Process general command line options
409 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc);
410 
411 // -----------------------------------------------------------------------------
412 // Miscellaneous utility functions
413 // -----------------------------------------------------------------------------
414 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time);
415 
416 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
417                                                   Vec grad_FVM);
418 
419 // Compare reference solution values with current test run for CI
420 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q);
421 
422 // Get error for problems with exact solutions
423 PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time);
424 
425 // Post-processing
426 PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time);
427 
428 // -- Gather initial Q values in case of continuation of simulation
429 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q);
430 
431 // Record boundary values from initial condition
432 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc);
433 
434 // Versioning token for binary checkpoints
435 extern const PetscInt32 FLUIDS_FILE_TOKEN;  // for backwards compatibility
436 extern const PetscInt32 FLUIDS_FILE_TOKEN_32;
437 extern const PetscInt32 FLUIDS_FILE_TOKEN_64;
438 
439 // Create appropriate mass qfunction based on number of components N
440 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf);
441 
442 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context);
443 
444 PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
445                                  FILE **fp);
446 
447 PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows);
448 
449 PetscErrorCode PhastaDatFileReadToArrayReal(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]);
450 
451 PetscErrorCode IntArrayC2P(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc);
452 PetscErrorCode IntArrayP2C(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed);
453 
454 // -----------------------------------------------------------------------------
455 // Turbulence Statistics Collection Functions
456 // -----------------------------------------------------------------------------
457 
458 PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem);
459 PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
460 PetscErrorCode TurbulenceStatisticsDestroy(User user, CeedData ceed_data);
461 
462 // -----------------------------------------------------------------------------
463 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions
464 // -----------------------------------------------------------------------------
465 
466 PetscErrorCode SgsDDModelSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem);
467 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data);
468 PetscErrorCode SgsDDModelApplyIFunction(User user, const Vec Q_loc, Vec G_loc);
469 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem, StateVariable state_var_input,
470                                                CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj);
471 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient);
472 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso,
473                                                         CeedVector *grid_aniso_vector);
474 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso,
475                                                              CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso);
476 
477 // -----------------------------------------------------------------------------
478 // Boundary Condition Related Functions
479 // -----------------------------------------------------------------------------
480 
481 // Setup StrongBCs that use QFunctions
482 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, ProblemData *problem, SimpleBC bc);
483 
484 PetscErrorCode FreestreamBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
485 PetscErrorCode OutflowBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
486 
487 // -----------------------------------------------------------------------------
488 // Differential Filtering Functions
489 // -----------------------------------------------------------------------------
490 
491 PetscErrorCode DifferentialFilterSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem);
492 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData diff_filter);
493 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
494 PetscErrorCode DifferentialFilterApply(User user, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution);
495 PetscErrorCode DifferentialFilterMmsICSetup(ProblemData *problem);
496 
497 // -----------------------------------------------------------------------------
498 // SGS Data-Driven Training via SmartSim
499 // -----------------------------------------------------------------------------
500 PetscErrorCode SmartSimSetup(User user);
501 PetscErrorCode SmartSimDataDestroy(SmartSimData smartsim);
502 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem);
503 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx);
504 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts);
505 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train);
506 
507 #endif  // libceed_fluids_examples_navier_stokes_h
508