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