xref: /libCEED/examples/fluids/navierstokes.h (revision 990fdeb6bb8fc9af2da4472bdc0d1f57da5da0e5)
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 <petscdm.h>
13 #include <petscdmplex.h>
14 #include <petscsys.h>
15 #include <petscts.h>
16 #include <stdbool.h>
17 #include "qfunctions/stabilization_types.h"
18 
19 // -----------------------------------------------------------------------------
20 // PETSc Version
21 // -----------------------------------------------------------------------------
22 #if PETSC_VERSION_LT(3,17,0)
23 #error "PETSc v3.17 or later is required"
24 #endif
25 
26 // -----------------------------------------------------------------------------
27 // Enums
28 // -----------------------------------------------------------------------------
29 // Translate PetscMemType to CeedMemType
30 static inline CeedMemType MemTypeP2C(PetscMemType mem_type) {
31   return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST;
32 }
33 
34 // Advection - Wind Options
35 typedef enum {
36   WIND_ROTATION    = 0,
37   WIND_TRANSLATION = 1,
38 } WindType;
39 static const char *const WindTypes[] = {
40   "rotation",
41   "translation",
42   "WindType", "WIND_", NULL
43 };
44 
45 // Advection - Bubble Types
46 typedef enum {
47   BUBBLE_SPHERE   = 0, // dim=3
48   BUBBLE_CYLINDER = 1, // dim=2
49 } BubbleType;
50 static const char *const BubbleTypes[] = {
51   "sphere",
52   "cylinder",
53   "BubbleType", "BUBBLE_", NULL
54 };
55 
56 // Advection - Bubble Continuity Types
57 typedef enum {
58   BUBBLE_CONTINUITY_SMOOTH     = 0,  // Original continuous, smooth shape
59   BUBBLE_CONTINUITY_BACK_SHARP = 1,  // Discontinuous, sharp back half shape
60   BUBBLE_CONTINUITY_THICK      = 2,  // Define a finite thickness
61 } BubbleContinuityType;
62 static const char *const BubbleContinuityTypes[] = {
63   "smooth",
64   "back_sharp",
65   "thick",
66   "BubbleContinuityType", "BUBBLE_CONTINUITY_", NULL
67 };
68 
69 // Euler - test cases
70 typedef enum {
71   EULER_TEST_ISENTROPIC_VORTEX = 0,
72   EULER_TEST_1 = 1,
73   EULER_TEST_2 = 2,
74   EULER_TEST_3 = 3,
75   EULER_TEST_4 = 4,
76   EULER_TEST_5 = 5,
77 } EulerTestType;
78 static const char *const EulerTestTypes[] = {
79   "isentropic_vortex",
80   "test_1",
81   "test_2",
82   "test_3",
83   "test_4",
84   "test_5",
85   "EulerTestType", "EULER_TEST_", NULL
86 };
87 
88 // Stabilization methods
89 static const char *const StabilizationTypes[] = {
90   "none",
91   "SU",
92   "SUPG",
93   "StabilizationType", "STAB_", NULL
94 };
95 
96 // -----------------------------------------------------------------------------
97 // Structs
98 // -----------------------------------------------------------------------------
99 // Structs declarations
100 typedef struct AppCtx_private    *AppCtx;
101 typedef struct CeedData_private  *CeedData;
102 typedef struct User_private      *User;
103 typedef struct Units_private     *Units;
104 typedef struct SimpleBC_private  *SimpleBC;
105 typedef struct Physics_private   *Physics;
106 
107 // Application context from user command line options
108 struct AppCtx_private {
109   // libCEED arguments
110   char              ceed_resource[PETSC_MAX_PATH_LEN]; // libCEED backend
111   PetscInt          degree;
112   PetscInt          q_extra;
113   // Post-processing arguments
114   PetscInt          output_freq;
115   PetscInt          viz_refine;
116   PetscInt          cont_steps;
117   char              output_dir[PETSC_MAX_PATH_LEN];
118   // Problem type arguments
119   PetscFunctionList problems;
120   char              problem_name[PETSC_MAX_PATH_LEN];
121   // Test mode arguments
122   PetscBool         test_mode;
123   PetscScalar       test_tol;
124   char              file_path[PETSC_MAX_PATH_LEN];
125 };
126 
127 // libCEED data struct
128 struct CeedData_private {
129   CeedVector           x_coord, q_data;
130   CeedQFunction        qf_setup_vol, qf_ics, qf_rhs_vol, qf_ifunction_vol,
131                        qf_setup_sur,
132                        qf_apply_inflow, qf_apply_inflow_jacobian,
133                        qf_apply_outflow, qf_apply_outflow_jacobian;
134   CeedBasis            basis_x, basis_xc, basis_q, basis_x_sur, basis_q_sur;
135   CeedElemRestriction  elem_restr_x, elem_restr_q, elem_restr_qd_i;
136   CeedOperator         op_setup_vol, op_ics;
137 };
138 
139 // PETSc user data
140 struct User_private {
141   MPI_Comm     comm;
142   DM           dm;
143   DM           dm_viz;
144   Mat          interp_viz;
145   Ceed         ceed;
146   Units        units;
147   Vec          M, Q_loc, Q_dot_loc;
148   Physics      phys;
149   AppCtx       app_ctx;
150   CeedVector   q_ceed, q_dot_ceed, g_ceed, coo_values;
151   CeedOperator op_rhs_vol, op_rhs, op_ifunction_vol, op_ifunction, op_ijacobian;
152   bool matrices_set_up;
153   CeedScalar time, dt;
154 };
155 
156 // Units
157 struct Units_private {
158   // fundamental units
159   PetscScalar meter;
160   PetscScalar kilogram;
161   PetscScalar second;
162   PetscScalar Kelvin;
163   // derived units
164   PetscScalar Pascal;
165   PetscScalar J_per_kg_K;
166   PetscScalar m_per_squared_s;
167   PetscScalar W_per_m_K;
168   PetscScalar Joule;
169 };
170 
171 // Boundary conditions
172 struct SimpleBC_private {
173   PetscInt  num_wall,    // Number of faces with wall BCs
174             wall_comps[5], // An array of constrained component numbers
175             num_comps,
176             num_slip[3], // Number of faces with slip BCs
177             num_inflow,
178             num_outflow;
179   PetscInt  walls[16], slips[3][16], inflows[16], outflows[16];
180   PetscBool user_bc;
181 };
182 
183 // Struct that contains all enums and structs used for the physics of all problems
184 struct Physics_private {
185   WindType                 wind_type;
186   BubbleType               bubble_type;
187   BubbleContinuityType     bubble_continuity_type;
188   EulerTestType            euler_test;
189   StabilizationType        stab;
190   PetscBool                implicit;
191   PetscBool                has_curr_time;
192   PetscBool                has_neumann;
193   CeedContextFieldLabel    solution_time_label;
194   CeedContextFieldLabel    timestep_size_label;
195   CeedContextFieldLabel    ics_time_label;
196   CeedContextFieldLabel    ijacobian_time_shift_label;
197 };
198 
199 typedef struct {
200   CeedQFunctionUser    qfunction;
201   const char           *qfunction_loc;
202   CeedQFunctionContext qfunction_context;
203 } ProblemQFunctionSpec;
204 
205 // Problem specific data
206 // *INDENT-OFF*
207 typedef struct ProblemData_private ProblemData;
208 struct ProblemData_private {
209   CeedInt           dim, q_data_size_vol, q_data_size_sur, jac_data_size_sur;
210   CeedScalar        dm_scale;
211   ProblemQFunctionSpec setup_vol, setup_sur, ics, apply_vol_rhs, apply_vol_ifunction,
212     apply_vol_ijacobian, apply_inflow, apply_outflow,
213     apply_inflow_jacobian, apply_outflow_jacobian;
214   bool              non_zero_time;
215   PetscErrorCode    (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt,
216                           PetscScalar[], void *);
217   void *bc_ctx;
218   PetscBool bc_from_ics;
219   PetscErrorCode    (*print_info)(ProblemData*, AppCtx);
220 };
221 // *INDENT-ON*
222 
223 extern int FreeContextPetsc(void *);
224 
225 // -----------------------------------------------------------------------------
226 // Set up problems
227 // -----------------------------------------------------------------------------
228 // Set up function for each problem
229 extern PetscErrorCode NS_CHANNEL(ProblemData *problem, DM dm,
230                                  void *ctx);
231 extern PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm,
232                                  void *ctx);
233 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm,
234                                       void *ctx);
235 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm,
236     void *ctx);
237 
238 extern PetscErrorCode NS_EULER_VORTEX(ProblemData *problem, DM dm,
239                                       void *ctx);
240 extern PetscErrorCode NS_SHOCKTUBE(ProblemData *problem, DM dm,
241                                    void *ctx);
242 extern PetscErrorCode NS_ADVECTION(ProblemData *problem, DM dm,
243                                    void *ctx);
244 extern PetscErrorCode NS_ADVECTION2D(ProblemData *problem, DM dm,
245                                      void *ctx);
246 
247 // Print function for each problem
248 extern PetscErrorCode PRINT_DENSITY_CURRENT(ProblemData *problem,
249     AppCtx app_ctx);
250 
251 extern PetscErrorCode PRINT_EULER_VORTEX(ProblemData *problem,
252     AppCtx app_ctx);
253 
254 extern PetscErrorCode PRINT_SHOCKTUBE(ProblemData *problem,
255                                       AppCtx app_ctx);
256 
257 extern PetscErrorCode PRINT_ADVECTION(ProblemData *problem,
258                                       AppCtx app_ctx);
259 
260 extern PetscErrorCode PRINT_ADVECTION2D(ProblemData *problem,
261                                         AppCtx app_ctx);
262 
263 // -----------------------------------------------------------------------------
264 // libCEED functions
265 // -----------------------------------------------------------------------------
266 // Utility function - essential BC dofs are encoded in closure indices as -(i+1).
267 PetscInt Involute(PetscInt i);
268 
269 // Utility function to create local CEED restriction
270 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height,
271     DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr);
272 
273 // Utility function to get Ceed Restriction for each domain
274 PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height,
275                                        DMLabel domain_label, PetscInt value,
276                                        CeedInt Q, CeedInt q_data_size,
277                                        CeedElemRestriction *elem_restr_q,
278                                        CeedElemRestriction *elem_restr_x,
279                                        CeedElemRestriction *elem_restr_qd_i);
280 
281 // Utility function to create CEED Composite Operator for the entire domain
282 PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc,
283                                        CeedData ceed_data, Physics phys,
284                                        CeedOperator op_apply_vol,
285                                        CeedOperator op_apply_ijacobian_vol,
286                                        CeedInt height,
287                                        CeedInt P_sur, CeedInt Q_sur,
288                                        CeedInt q_data_size_sur, CeedInt jac_data_size_sur,
289                                        CeedOperator *op_apply, CeedOperator *op_apply_ijacobian);
290 
291 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user,
292                             AppCtx app_ctx, ProblemData *problem, SimpleBC bc);
293 
294 // -----------------------------------------------------------------------------
295 // Time-stepping functions
296 // -----------------------------------------------------------------------------
297 // Compute mass matrix for explicit scheme
298 PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data,
299                                        Vec M);
300 
301 // RHS (Explicit time-stepper) function setup
302 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data);
303 
304 // Implicit time-stepper function setup
305 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G,
306                             void *user_data);
307 
308 // User provided TS Monitor
309 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q,
310                             void *ctx);
311 
312 // TS: Create, setup, and solve
313 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys,
314                           Vec *Q, PetscScalar *f_time, TS *ts);
315 
316 // -----------------------------------------------------------------------------
317 // Setup DM
318 // -----------------------------------------------------------------------------
319 // Create mesh
320 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem,
321                         MatType, VecType, DM *dm);
322 
323 // Set up DM
324 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree,
325                        SimpleBC bc, Physics phys);
326 
327 // Refine DM for high-order viz
328 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem,
329                            SimpleBC bc, Physics phys);
330 
331 // -----------------------------------------------------------------------------
332 // Process command line options
333 // -----------------------------------------------------------------------------
334 // Register problems to be available on the command line
335 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx);
336 
337 // Process general command line options
338 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx,
339     SimpleBC bc);
340 
341 // -----------------------------------------------------------------------------
342 // Miscellaneous utility functions
343 // -----------------------------------------------------------------------------
344 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user,
345                                    Vec Q_loc, Vec Q,
346                                    CeedScalar time);
347 
348 PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm,
349     PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM,
350     Vec cell_geom_FVM, Vec grad_FVM);
351 
352 // Compare reference solution values with current test run for CI
353 PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q);
354 
355 // Get error for problems with exact solutions
356 PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q,
357                            PetscScalar final_time);
358 
359 // Post-processing
360 PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm,
361                               ProblemData *problem, User user,
362                               Vec Q, PetscScalar final_time);
363 
364 // -- Gather initial Q values in case of continuation of simulation
365 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q);
366 
367 // Record boundary values from initial condition
368 PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc);
369 
370 // -----------------------------------------------------------------------------
371 
372 #endif // libceed_fluids_examples_navier_stokes_h
373