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