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