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