xref: /libCEED/examples/fluids/navierstokes.h (revision 9a3a46e1d4300a3ad2fa9b1d35445254cd604d16)
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;
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 };
154 
155 // Units
156 struct Units_private {
157   // fundamental units
158   PetscScalar meter;
159   PetscScalar kilogram;
160   PetscScalar second;
161   PetscScalar Kelvin;
162   // derived units
163   PetscScalar Pascal;
164   PetscScalar J_per_kg_K;
165   PetscScalar m_per_squared_s;
166   PetscScalar W_per_m_K;
167   PetscScalar Joule;
168 };
169 
170 // Boundary conditions
171 struct SimpleBC_private {
172   PetscInt  num_wall,    // Number of faces with wall BCs
173             wall_comps[5], // An array of constrained component numbers
174             num_comps,
175             num_slip[3], // Number of faces with slip BCs
176             num_inflow,
177             num_outflow;
178   PetscInt  walls[16], slips[3][16], inflows[16], outflows[16];
179   PetscBool user_bc;
180 };
181 
182 // Struct that contains all enums and structs used for the physics of all problems
183 struct Physics_private {
184   WindType                 wind_type;
185   BubbleType               bubble_type;
186   BubbleContinuityType     bubble_continuity_type;
187   EulerTestType            euler_test;
188   StabilizationType        stab;
189   PetscBool                implicit;
190   PetscBool                has_curr_time;
191   PetscBool                has_neumann;
192   CeedContextFieldLabel    solution_time_label;
193   CeedContextFieldLabel    timestep_size_label;
194   CeedContextFieldLabel    ics_time_label;
195   CeedContextFieldLabel    ijacobian_time_shift_label;
196 };
197 
198 typedef struct {
199   CeedQFunctionUser    qfunction;
200   const char           *qfunction_loc;
201   CeedQFunctionContext qfunction_context;
202 } ProblemQFunctionSpec;
203 
204 // Problem specific data
205 // *INDENT-OFF*
206 typedef struct ProblemData_private ProblemData;
207 struct ProblemData_private {
208   CeedInt           dim, q_data_size_vol, q_data_size_sur, jac_data_size_sur;
209   CeedScalar        dm_scale;
210   ProblemQFunctionSpec setup_vol, setup_sur, ics, apply_vol_rhs, apply_vol_ifunction,
211     apply_vol_ijacobian, apply_inflow, apply_outflow,
212     apply_inflow_jacobian, apply_outflow_jacobian;
213   bool              non_zero_time;
214   PetscErrorCode    (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt,
215                           PetscScalar[], void *);
216   void *bc_ctx;
217   PetscBool bc_from_ics;
218   PetscErrorCode    (*print_info)(ProblemData*, AppCtx);
219 };
220 // *INDENT-ON*
221 
222 extern int FreeContextPetsc(void *);
223 
224 // -----------------------------------------------------------------------------
225 // Set up problems
226 // -----------------------------------------------------------------------------
227 // Set up function for each problem
228 extern PetscErrorCode NS_CHANNEL(ProblemData *problem, DM dm,
229                                  void *ctx);
230 extern PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm,
231                                  void *ctx);
232 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm,
233                                       void *ctx);
234 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm,
235     void *ctx);
236 
237 extern PetscErrorCode NS_EULER_VORTEX(ProblemData *problem, DM dm,
238                                       void *ctx);
239 extern PetscErrorCode NS_SHOCKTUBE(ProblemData *problem, DM dm,
240                                    void *ctx);
241 extern PetscErrorCode NS_ADVECTION(ProblemData *problem, DM dm,
242                                    void *ctx);
243 extern PetscErrorCode NS_ADVECTION2D(ProblemData *problem, DM dm,
244                                      void *ctx);
245 
246 // Print function for each problem
247 extern PetscErrorCode PRINT_DENSITY_CURRENT(ProblemData *problem,
248     AppCtx app_ctx);
249 
250 extern PetscErrorCode PRINT_EULER_VORTEX(ProblemData *problem,
251     AppCtx app_ctx);
252 
253 extern PetscErrorCode PRINT_SHOCKTUBE(ProblemData *problem,
254                                       AppCtx app_ctx);
255 
256 extern PetscErrorCode PRINT_ADVECTION(ProblemData *problem,
257                                       AppCtx app_ctx);
258 
259 extern PetscErrorCode PRINT_ADVECTION2D(ProblemData *problem,
260                                         AppCtx app_ctx);
261 
262 // -----------------------------------------------------------------------------
263 // libCEED functions
264 // -----------------------------------------------------------------------------
265 // Utility function - essential BC dofs are encoded in closure indices as -(i+1).
266 PetscInt Involute(PetscInt i);
267 
268 // Utility function to create local CEED restriction
269 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height,
270     DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr);
271 
272 // Utility function to get Ceed Restriction for each domain
273 PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height,
274                                        DMLabel domain_label, PetscInt value,
275                                        CeedInt Q, CeedInt q_data_size,
276                                        CeedElemRestriction *elem_restr_q,
277                                        CeedElemRestriction *elem_restr_x,
278                                        CeedElemRestriction *elem_restr_qd_i);
279 
280 // Utility function to create CEED Composite Operator for the entire domain
281 PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc,
282                                        CeedData ceed_data, Physics phys,
283                                        CeedOperator op_apply_vol,
284                                        CeedOperator op_apply_ijacobian_vol,
285                                        CeedInt height,
286                                        CeedInt P_sur, CeedInt Q_sur,
287                                        CeedInt q_data_size_sur, CeedInt jac_data_size_sur,
288                                        CeedOperator *op_apply, CeedOperator *op_apply_ijacobian);
289 
290 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user,
291                             AppCtx app_ctx, ProblemData *problem, SimpleBC bc);
292 
293 // -----------------------------------------------------------------------------
294 // Time-stepping functions
295 // -----------------------------------------------------------------------------
296 // Compute mass matrix for explicit scheme
297 PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data,
298                                        Vec M);
299 
300 // RHS (Explicit time-stepper) function setup
301 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data);
302 
303 // Implicit time-stepper function setup
304 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G,
305                             void *user_data);
306 
307 // User provided TS Monitor
308 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q,
309                             void *ctx);
310 
311 // TS: Create, setup, and solve
312 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys,
313                           Vec *Q, PetscScalar *f_time, TS *ts);
314 
315 // -----------------------------------------------------------------------------
316 // Setup DM
317 // -----------------------------------------------------------------------------
318 // Create mesh
319 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem,
320                         MatType, VecType, DM *dm);
321 
322 // Set up DM
323 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree,
324                        SimpleBC bc, Physics phys);
325 
326 // Refine DM for high-order viz
327 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem,
328                            SimpleBC bc, Physics phys);
329 
330 // -----------------------------------------------------------------------------
331 // Process command line options
332 // -----------------------------------------------------------------------------
333 // Register problems to be available on the command line
334 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx);
335 
336 // Process general command line options
337 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx,
338     SimpleBC bc);
339 
340 // -----------------------------------------------------------------------------
341 // Miscellaneous utility functions
342 // -----------------------------------------------------------------------------
343 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user,
344                                    Vec Q_loc, Vec Q,
345                                    CeedScalar time);
346 
347 PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm,
348     PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM,
349     Vec cell_geom_FVM, Vec grad_FVM);
350 
351 // Compare reference solution values with current test run for CI
352 PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q);
353 
354 // Get error for problems with exact solutions
355 PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q,
356                            PetscScalar final_time);
357 
358 // Post-processing
359 PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm,
360                               ProblemData *problem, User user,
361                               Vec Q, PetscScalar final_time);
362 
363 // -- Gather initial Q values in case of continuation of simulation
364 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q);
365 
366 // Record boundary values from initial condition
367 PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc);
368 
369 // -----------------------------------------------------------------------------
370 
371 #endif // libceed_fluids_examples_navier_stokes_h
372