xref: /libCEED/examples/fluids/navierstokes.h (revision 3b38d1dfe7a8d0ddfedd40f2e4c35e536721e848)
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 PetscErrorCode (*SgsDDNodalStressEval)(User user, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc);
199 typedef PetscErrorCode (*SgsDDNodalStressInference)(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx);
200 typedef struct {
201   DM                        dm_sgs, dm_dd_inputs, dm_dd_outputs;
202   PetscInt                  num_comp_sgs;
203   CeedInt                   num_comp_inputs, num_comp_outputs;
204   OperatorApplyContext      op_nodal_evaluation_ctx, op_nodal_dd_inputs_ctx, op_nodal_dd_outputs_ctx, op_sgs_apply_ctx;
205   CeedVector                sgs_nodal_ceed, grad_velo_ceed;
206   SgsDDNodalStressEval      sgs_nodal_eval;
207   SgsDDNodalStressInference sgs_nodal_inference;
208   void                     *sgs_nodal_inference_ctx;
209   PetscErrorCode (*sgs_nodal_inference_ctx_destroy)(void *ctx);
210 } *SgsDDData;
211 
212 typedef struct {
213   DM                   dm_dd_training;
214   PetscInt             num_comp_dd_inputs, write_data_interval;
215   OperatorApplyContext op_training_data_calc_ctx;
216   NodalProjectionData  filtered_grad_velo_proj;
217   size_t               training_data_array_dims[2];
218   PetscBool            overwrite_training_data;
219 } *SGS_DD_TrainingData;
220 
221 typedef struct {
222   DM                   dm_filter;
223   PetscInt             num_filtered_fields;
224   CeedInt             *num_field_components;
225   PetscInt             field_prim_state, field_velo_prod;
226   OperatorApplyContext op_rhs_ctx;
227   KSP                  ksp;
228   PetscBool            do_mms_test;
229 } *DiffFilterData;
230 
231 typedef struct {
232   void    *client;
233   char     rank_id_name[16];
234   PetscInt collocated_database_num_ranks;
235 } *SmartSimData;
236 
237 // PETSc user data
238 struct User_private {
239   MPI_Comm             comm;
240   DM                   dm;
241   DM                   dm_viz;
242   Mat                  interp_viz;
243   Ceed                 ceed;
244   Units                units;
245   Vec                  Q_loc, Q_dot_loc;
246   Physics              phys;
247   AppCtx               app_ctx;
248   CeedVector           q_ceed, q_dot_ceed, g_ceed, coo_values_amat, coo_values_pmat, x_ceed;
249   CeedOperator         op_rhs_vol, op_ifunction_vol, op_ifunction, op_ijacobian;
250   KSP                  mass_ksp;
251   OperatorApplyContext op_rhs_ctx, op_strong_bc_ctx;
252   bool                 matrices_set_up;
253   CeedScalar           time_bc_set;
254   SpanStatsData        spanstats;
255   NodalProjectionData  grad_velo_proj;
256   SgsDDData            sgs_dd_data;
257   DiffFilterData       diff_filter;
258   SmartSimData         smartsim;
259   SGS_DD_TrainingData  sgs_dd_train;
260 };
261 
262 // Units
263 struct Units_private {
264   // fundamental units
265   PetscScalar meter;
266   PetscScalar kilogram;
267   PetscScalar second;
268   PetscScalar Kelvin;
269   // derived units
270   PetscScalar Pascal;
271   PetscScalar J_per_kg_K;
272   PetscScalar m_per_squared_s;
273   PetscScalar W_per_m_K;
274   PetscScalar Joule;
275 };
276 
277 // Boundary conditions
278 struct SimpleBC_private {
279   PetscInt num_wall,  // Number of faces with wall BCs
280       wall_comps[5],  // An array of constrained component numbers
281       num_comps,
282       num_slip[3],  // Number of faces with slip BCs
283       num_inflow, num_outflow, num_freestream;
284   PetscInt  walls[16], slips[3][16], inflows[16], outflows[16], freestreams[16];
285   PetscBool user_bc;
286 };
287 
288 // Struct that contains all enums and structs used for the physics of all problems
289 struct Physics_private {
290   PetscBool             implicit;
291   StateVariable         state_var;
292   CeedContextFieldLabel solution_time_label;
293   CeedContextFieldLabel stg_solution_time_label;
294   CeedContextFieldLabel timestep_size_label;
295   CeedContextFieldLabel ics_time_label;
296   CeedContextFieldLabel ijacobian_time_shift_label;
297 };
298 
299 typedef struct {
300   CeedQFunctionUser    qfunction;
301   const char          *qfunction_loc;
302   CeedQFunctionContext qfunction_context;
303 } ProblemQFunctionSpec;
304 
305 // Problem specific data
306 typedef struct ProblemData_private ProblemData;
307 struct ProblemData_private {
308   CeedInt              dim, q_data_size_vol, q_data_size_sur, jac_data_size_sur;
309   CeedScalar           dm_scale;
310   ProblemQFunctionSpec setup_vol, setup_sur, ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow,
311       apply_freestream, apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian;
312   bool      non_zero_time;
313   PetscBool bc_from_ics, use_strong_bc_ceed, uses_newtonian;
314   PetscErrorCode (*print_info)(User, ProblemData *, AppCtx);
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(User user, ProblemData *problem, AppCtx app_ctx);
336 
337 extern PetscErrorCode PRINT_EULER_VORTEX(User user, ProblemData *problem, AppCtx app_ctx);
338 
339 extern PetscErrorCode PRINT_SHOCKTUBE(User user, ProblemData *problem, AppCtx app_ctx);
340 
341 extern PetscErrorCode PRINT_ADVECTION(User user, ProblemData *problem, AppCtx app_ctx);
342 
343 extern PetscErrorCode PRINT_ADVECTION2D(User user, ProblemData *problem, AppCtx app_ctx);
344 
345 PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData *problem, MPI_Comm comm);
346 
347 // -----------------------------------------------------------------------------
348 // libCEED functions
349 // -----------------------------------------------------------------------------
350 // Utility function to create local CEED restriction
351 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt label_value, PetscInt dm_field,
352                                          CeedElemRestriction *elem_restr);
353 
354 PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
355                                                CeedElemRestriction *restriction);
356 PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
357                                                          CeedElemRestriction *restriction);
358 PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
359                                                     PetscInt q_data_size, CeedElemRestriction *restriction);
360 PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
361                                                          PetscInt q_data_size, CeedElemRestriction *restriction);
362 
363 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis);
364 
365 // Utility function to create CEED Composite Operator for the entire domain
366 PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol,
367                                        CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur,
368                                        CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian);
369 
370 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc);
371 
372 // -----------------------------------------------------------------------------
373 // Time-stepping functions
374 // -----------------------------------------------------------------------------
375 // Create KSP to solve the inverse mass operator for explicit time stepping schemes
376 PetscErrorCode CreateKSPMassOperator(User user, CeedData ceed_data);
377 
378 // RHS (Explicit time-stepper) function setup
379 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data);
380 
381 // Implicit time-stepper function setup
382 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data);
383 
384 // User provided TS Monitor
385 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx);
386 
387 // TS: Create, setup, and solve
388 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts);
389 
390 // Update Boundary Values when time has changed
391 PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t);
392 
393 // -----------------------------------------------------------------------------
394 // Setup DM
395 // -----------------------------------------------------------------------------
396 // Create mesh
397 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType, VecType, DM *dm);
398 
399 // Set up DM
400 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys);
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 // Refine DM for high-order viz
408 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys);
409 
410 // -----------------------------------------------------------------------------
411 // Process command line options
412 // -----------------------------------------------------------------------------
413 // Register problems to be available on the command line
414 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx);
415 
416 // Process general command line options
417 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc);
418 
419 // -----------------------------------------------------------------------------
420 // Miscellaneous utility functions
421 // -----------------------------------------------------------------------------
422 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time);
423 
424 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
425                                                   Vec grad_FVM);
426 
427 // Compare reference solution values with current test run for CI
428 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q);
429 
430 // Get error for problems with exact solutions
431 PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time);
432 
433 // Post-processing
434 PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time);
435 
436 // -- Gather initial Q values in case of continuation of simulation
437 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q);
438 
439 // Record boundary values from initial condition
440 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc);
441 
442 // Versioning token for binary checkpoints
443 extern const PetscInt32 FLUIDS_FILE_TOKEN;  // for backwards compatibility
444 extern const PetscInt32 FLUIDS_FILE_TOKEN_32;
445 extern const PetscInt32 FLUIDS_FILE_TOKEN_64;
446 
447 // Create appropriate mass qfunction based on number of components N
448 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf);
449 
450 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context);
451 
452 PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
453                                  FILE **fp);
454 
455 PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows);
456 
457 PetscErrorCode PhastaDatFileReadToArrayReal(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]);
458 
459 PetscErrorCode IntArrayC2P(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc);
460 PetscErrorCode IntArrayP2C(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed);
461 
462 // -----------------------------------------------------------------------------
463 // Turbulence Statistics Collection Functions
464 // -----------------------------------------------------------------------------
465 
466 PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem);
467 PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
468 PetscErrorCode TurbulenceStatisticsDestroy(User user, CeedData ceed_data);
469 
470 // -----------------------------------------------------------------------------
471 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions
472 // -----------------------------------------------------------------------------
473 
474 PetscErrorCode SgsDDSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem);
475 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data);
476 PetscErrorCode SgsDDApplyIFunction(User user, const Vec Q_loc, Vec G_loc);
477 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem, StateVariable state_var_input,
478                                                CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj);
479 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient);
480 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso,
481                                                         CeedVector *grid_aniso_vector);
482 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso,
483                                                              CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso);
484 
485 // -----------------------------------------------------------------------------
486 // Boundary Condition Related Functions
487 // -----------------------------------------------------------------------------
488 
489 // Setup StrongBCs that use QFunctions
490 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, ProblemData *problem, SimpleBC bc);
491 
492 PetscErrorCode FreestreamBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
493 PetscErrorCode OutflowBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
494 
495 // -----------------------------------------------------------------------------
496 // Differential Filtering Functions
497 // -----------------------------------------------------------------------------
498 
499 PetscErrorCode DifferentialFilterSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem);
500 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData diff_filter);
501 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
502 PetscErrorCode DifferentialFilterApply(User user, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution);
503 PetscErrorCode DifferentialFilterMmsICSetup(ProblemData *problem);
504 
505 // -----------------------------------------------------------------------------
506 // SGS Data-Driven Training via SmartSim
507 // -----------------------------------------------------------------------------
508 PetscErrorCode SmartSimSetup(User user);
509 PetscErrorCode SmartSimDataDestroy(SmartSimData smartsim);
510 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem);
511 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx);
512 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts);
513 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train);
514 
515 #endif  // libceed_fluids_examples_navier_stokes_h
516