xref: /libCEED/examples/fluids/navierstokes.h (revision 999ff5c7ca91b692a440ecf149a8ab458a4dc7a6)
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 // -----------------------------------------------------------------------------
20 // PETSc Version
21 // -----------------------------------------------------------------------------
22 #if PETSC_VERSION_LT(3, 19, 0)
23 #error "PETSc v3.19 or later is required"
24 #endif
25 
26 // -----------------------------------------------------------------------------
27 // Enums
28 // -----------------------------------------------------------------------------
29 // Translate PetscMemType to CeedMemType
30 static inline CeedMemType MemTypeP2C(PetscMemType mem_type) { return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; }
31 
32 // Advection - Wind Options
33 typedef enum {
34   WIND_ROTATION    = 0,
35   WIND_TRANSLATION = 1,
36 } WindType;
37 static const char *const WindTypes[] = {"rotation", "translation", "WindType", "WIND_", NULL};
38 
39 // Advection - Bubble Types
40 typedef enum {
41   BUBBLE_SPHERE   = 0,  // dim=3
42   BUBBLE_CYLINDER = 1,  // dim=2
43 } BubbleType;
44 static const char *const BubbleTypes[] = {"sphere", "cylinder", "BubbleType", "BUBBLE_", NULL};
45 
46 // Advection - Bubble Continuity Types
47 typedef enum {
48   BUBBLE_CONTINUITY_SMOOTH     = 0,  // Original continuous, smooth shape
49   BUBBLE_CONTINUITY_BACK_SHARP = 1,  // Discontinuous, sharp back half shape
50   BUBBLE_CONTINUITY_THICK      = 2,  // Define a finite thickness
51 } BubbleContinuityType;
52 static const char *const BubbleContinuityTypes[] = {"smooth", "back_sharp", "thick", "BubbleContinuityType", "BUBBLE_CONTINUITY_", NULL};
53 
54 // Euler - test cases
55 typedef enum {
56   EULER_TEST_ISENTROPIC_VORTEX = 0,
57   EULER_TEST_1                 = 1,
58   EULER_TEST_2                 = 2,
59   EULER_TEST_3                 = 3,
60   EULER_TEST_4                 = 4,
61   EULER_TEST_5                 = 5,
62 } EulerTestType;
63 static const char *const EulerTestTypes[] = {"isentropic_vortex", "test_1",      "test_2", "test_3", "test_4", "test_5",
64                                              "EulerTestType",     "EULER_TEST_", NULL};
65 
66 // Stabilization methods
67 static const char *const StabilizationTypes[] = {"none", "SU", "SUPG", "StabilizationType", "STAB_", NULL};
68 
69 // Euler - test cases
70 typedef enum {
71   TESTTYPE_NONE           = 0,
72   TESTTYPE_SOLVER         = 1,
73   TESTTYPE_TURB_SPANSTATS = 2,
74 } TestType;
75 static const char *const TestTypes[] = {"none", "solver", "turb_spanstats", "TestType", "TESTTYPE_", NULL};
76 
77 // -----------------------------------------------------------------------------
78 // Structs
79 // -----------------------------------------------------------------------------
80 // Structs declarations
81 typedef struct AppCtx_private   *AppCtx;
82 typedef struct CeedData_private *CeedData;
83 typedef struct User_private     *User;
84 typedef struct Units_private    *Units;
85 typedef struct SimpleBC_private *SimpleBC;
86 typedef struct Physics_private  *Physics;
87 
88 // Application context from user command line options
89 struct AppCtx_private {
90   // libCEED arguments
91   char     ceed_resource[PETSC_MAX_PATH_LEN];  // libCEED backend
92   PetscInt degree;
93   PetscInt q_extra;
94   // Solver arguments
95   MatType   amat_type;
96   PetscBool pmat_pbdiagonal;
97   // Post-processing arguments
98   PetscInt  checkpoint_interval;
99   PetscInt  viz_refine;
100   PetscInt  cont_steps;
101   PetscReal cont_time;
102   char      cont_file[PETSC_MAX_PATH_LEN];
103   char      cont_time_file[PETSC_MAX_PATH_LEN];
104   char      output_dir[PETSC_MAX_PATH_LEN];
105   PetscBool add_stepnum2bin;
106   PetscBool checkpoint_vtk;
107   // Problem type arguments
108   PetscFunctionList problems;
109   char              problem_name[PETSC_MAX_PATH_LEN];
110   // Test mode arguments
111   TestType    test_type;
112   PetscScalar test_tol;
113   char        test_file_path[PETSC_MAX_PATH_LEN];
114   // Turbulent spanwise statistics
115   PetscBool         turb_spanstats_enable;
116   PetscInt          turb_spanstats_collect_interval;
117   PetscInt          turb_spanstats_viewer_interval;
118   PetscViewer       turb_spanstats_viewer;
119   PetscViewerFormat turb_spanstats_viewer_format;
120   // Wall forces
121   struct {
122     PetscInt          num_wall;
123     PetscInt         *walls;
124     PetscViewer       viewer;
125     PetscViewerFormat viewer_format;
126     PetscBool         header_written;
127   } wall_forces;
128 };
129 
130 // libCEED data struct
131 struct CeedData_private {
132   CeedVector          x_coord, q_data;
133   CeedBasis           basis_x, basis_xc, basis_q, basis_x_sur, basis_q_sur, basis_xc_sur;
134   CeedElemRestriction elem_restr_x, elem_restr_q, elem_restr_qd_i;
135   CeedOperator        op_setup_vol, op_ics;
136   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,
137       qf_apply_outflow_jacobian, qf_apply_freestream, qf_apply_freestream_jacobian;
138 };
139 
140 typedef struct {
141   DM                    dm;
142   PetscSF               sf;  // For communicating child data to parents
143   CeedOperator          op_stats_collect, op_stats_proj;
144   PetscInt              num_comp_stats;
145   CeedVector            child_stats, parent_stats;  // collocated statistics data
146   CeedVector            rhs_ceed;
147   KSP                   ksp;         // For the L^2 projection solve
148   CeedScalar            span_width;  // spanwise width of the child domain
149   PetscBool             do_mms_test;
150   OperatorApplyContext  mms_error_ctx;
151   CeedContextFieldLabel solution_time_label, previous_time_label;
152 } Span_Stats;
153 
154 typedef struct {
155   DM                   dm;
156   PetscInt             num_comp;
157   OperatorApplyContext l2_rhs_ctx;
158   KSP                  ksp;
159 } *NodalProjectionData;
160 
161 // PETSc user data
162 struct User_private {
163   MPI_Comm            comm;
164   DM                  dm;
165   DM                  dm_viz;
166   Mat                 interp_viz;
167   Ceed                ceed;
168   Units               units;
169   Vec                 M, Q_loc, Q_dot_loc;
170   Physics             phys;
171   AppCtx              app_ctx;
172   CeedVector          q_ceed, q_dot_ceed, g_ceed, coo_values_amat, coo_values_pmat, x_ceed;
173   CeedOperator        op_rhs_vol, op_rhs, op_ifunction_vol, op_ifunction, op_ijacobian, op_dirichlet;
174   bool                matrices_set_up;
175   CeedScalar          time_bc_set;
176   Span_Stats          spanstats;
177   NodalProjectionData grad_velo_proj;
178 };
179 
180 // Units
181 struct Units_private {
182   // fundamental units
183   PetscScalar meter;
184   PetscScalar kilogram;
185   PetscScalar second;
186   PetscScalar Kelvin;
187   // derived units
188   PetscScalar Pascal;
189   PetscScalar J_per_kg_K;
190   PetscScalar m_per_squared_s;
191   PetscScalar W_per_m_K;
192   PetscScalar Joule;
193 };
194 
195 // Boundary conditions
196 struct SimpleBC_private {
197   PetscInt num_wall,  // Number of faces with wall BCs
198       wall_comps[5],  // An array of constrained component numbers
199       num_comps,
200       num_slip[3],  // Number of faces with slip BCs
201       num_inflow, num_outflow, num_freestream;
202   PetscInt  walls[16], slips[3][16], inflows[16], outflows[16], freestreams[16];
203   PetscBool user_bc;
204 };
205 
206 // Struct that contains all enums and structs used for the physics of all problems
207 struct Physics_private {
208   WindType              wind_type;
209   BubbleType            bubble_type;
210   BubbleContinuityType  bubble_continuity_type;
211   EulerTestType         euler_test;
212   StabilizationType     stab;
213   PetscBool             implicit;
214   StateVariable         state_var;
215   PetscBool             has_curr_time;
216   PetscBool             has_neumann;
217   CeedContextFieldLabel solution_time_label;
218   CeedContextFieldLabel stg_solution_time_label;
219   CeedContextFieldLabel timestep_size_label;
220   CeedContextFieldLabel ics_time_label;
221   CeedContextFieldLabel ijacobian_time_shift_label;
222 };
223 
224 typedef struct {
225   CeedQFunctionUser    qfunction;
226   const char          *qfunction_loc;
227   CeedQFunctionContext qfunction_context;
228 } ProblemQFunctionSpec;
229 
230 // Problem specific data
231 typedef struct ProblemData_private ProblemData;
232 struct ProblemData_private {
233   CeedInt              dim, q_data_size_vol, q_data_size_sur, jac_data_size_sur;
234   CeedScalar           dm_scale;
235   ProblemQFunctionSpec setup_vol, setup_sur, ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow,
236       apply_freestream, apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian;
237   bool non_zero_time;
238   PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
239   void     *bc_ctx;
240   PetscBool bc_from_ics, use_dirichlet_ceed;
241   PetscErrorCode (*print_info)(ProblemData *, AppCtx);
242 };
243 
244 extern int FreeContextPetsc(void *);
245 
246 // -----------------------------------------------------------------------------
247 // Set up problems
248 // -----------------------------------------------------------------------------
249 // Set up function for each problem
250 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
251 extern PetscErrorCode NS_CHANNEL(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
252 extern PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
253 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
254 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
255 extern PetscErrorCode NS_EULER_VORTEX(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
256 extern PetscErrorCode NS_SHOCKTUBE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
257 extern PetscErrorCode NS_ADVECTION(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
258 extern PetscErrorCode NS_ADVECTION2D(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
259 
260 // Print function for each problem
261 extern PetscErrorCode PRINT_NEWTONIAN(ProblemData *problem, AppCtx app_ctx);
262 
263 extern PetscErrorCode PRINT_EULER_VORTEX(ProblemData *problem, AppCtx app_ctx);
264 
265 extern PetscErrorCode PRINT_SHOCKTUBE(ProblemData *problem, AppCtx app_ctx);
266 
267 extern PetscErrorCode PRINT_ADVECTION(ProblemData *problem, AppCtx app_ctx);
268 
269 extern PetscErrorCode PRINT_ADVECTION2D(ProblemData *problem, AppCtx app_ctx);
270 
271 // -----------------------------------------------------------------------------
272 // libCEED functions
273 // -----------------------------------------------------------------------------
274 // Utility function - essential BC dofs are encoded in closure indices as -(i+1).
275 PetscInt Involute(PetscInt i);
276 
277 // Utility function to create local CEED restriction
278 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr);
279 
280 // Utility function to get Ceed Restriction for each domain
281 PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, PetscInt value, CeedInt Q, CeedInt q_data_size,
282                                        CeedElemRestriction *elem_restr_q, CeedElemRestriction *elem_restr_x, CeedElemRestriction *elem_restr_qd_i);
283 
284 // Utility function to create CEED Composite Operator for the entire domain
285 PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol,
286                                        CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur,
287                                        CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian);
288 
289 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc);
290 
291 // -----------------------------------------------------------------------------
292 // Time-stepping functions
293 // -----------------------------------------------------------------------------
294 // Compute mass matrix for explicit scheme
295 PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M);
296 
297 // RHS (Explicit time-stepper) function setup
298 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data);
299 
300 // Implicit time-stepper function setup
301 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data);
302 
303 // User provided TS Monitor
304 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx);
305 
306 // TS: Create, setup, and solve
307 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts);
308 
309 // Update Boundary Values when time has changed
310 PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t);
311 
312 // -----------------------------------------------------------------------------
313 // Setup DM
314 // -----------------------------------------------------------------------------
315 // Create mesh
316 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType, VecType, DM *dm);
317 
318 // Set up DM
319 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, SimpleBC bc, Physics phys);
320 
321 // Refine DM for high-order viz
322 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys);
323 
324 // -----------------------------------------------------------------------------
325 // Process command line options
326 // -----------------------------------------------------------------------------
327 // Register problems to be available on the command line
328 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx);
329 
330 // Process general command line options
331 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc);
332 
333 // -----------------------------------------------------------------------------
334 // Miscellaneous utility functions
335 // -----------------------------------------------------------------------------
336 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time);
337 
338 PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
339                                              Vec grad_FVM);
340 
341 // Compare reference solution values with current test run for CI
342 PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q);
343 
344 // Get error for problems with exact solutions
345 PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time);
346 
347 // Post-processing
348 PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time);
349 
350 // -- Gather initial Q values in case of continuation of simulation
351 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q);
352 
353 // Record boundary values from initial condition
354 PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc);
355 
356 // Versioning token for binary checkpoints
357 extern const PetscInt FLUIDS_FILE_TOKEN;
358 
359 // Create appropriate mass qfunction based on number of components N
360 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf);
361 
362 PetscErrorCode ComputeL2Projection(Vec source_vec, Vec target_vec, OperatorApplyContext rhs_matop_ctx, KSP ksp);
363 
364 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context);
365 
366 PetscErrorCode PHASTADatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
367                                  FILE **fp);
368 
369 PetscErrorCode PHASTADatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows);
370 
371 PetscErrorCode PHASTADatFileReadToArrayReal(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]);
372 
373 // -----------------------------------------------------------------------------
374 // Turbulence Statistics Collection Functions
375 // -----------------------------------------------------------------------------
376 
377 PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree, SimpleBC bc);
378 PetscErrorCode SetupStatsCollection(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem);
379 PetscErrorCode TSMonitor_Statistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
380 PetscErrorCode DestroyStats(User user, CeedData ceed_data);
381 
382 // -----------------------------------------------------------------------------
383 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions
384 // -----------------------------------------------------------------------------
385 
386 PetscErrorCode SGS_DD_ModelSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem);
387 
388 // -----------------------------------------------------------------------------
389 // Boundary Condition Related Functions
390 // -----------------------------------------------------------------------------
391 
392 // Setup StrongBCs that use QFunctions
393 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc, CeedInt Q_sur,
394                                   CeedInt q_data_size_sur);
395 
396 PetscErrorCode FreestreamBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
397 PetscErrorCode OutflowBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
398 
399 #endif  // libceed_fluids_examples_navier_stokes_h
400