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