xref: /libCEED/examples/fluids/problems/advection.c (revision 2e88d319acd75658c85aee6ae8a14c6e519e582d)
1 // Copyright (c) 2017-2024, 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 /// @file
9 /// Utility functions for setting up ADVECTION
10 
11 #include "../qfunctions/advection.h"
12 
13 #include <ceed.h>
14 #include <petscdm.h>
15 
16 #include "../navierstokes.h"
17 
18 // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping
19 //
20 // Only used for SUPG stabilization
21 PetscErrorCode CreateKSPMassOperator_AdvectionStabilized(User user, CeedOperator *op_mass) {
22   Ceed                 ceed = user->ceed;
23   CeedInt              num_comp_q, q_data_size;
24   CeedQFunction        qf_mass;
25   CeedElemRestriction  elem_restr_q, elem_restr_qd_i;
26   CeedBasis            basis_q;
27   CeedVector           q_data;
28   CeedQFunctionContext qf_ctx = NULL;
29   PetscInt             dim;
30 
31   PetscFunctionBeginUser;
32   PetscCall(DMGetDimension(user->dm, &dim));
33   {  // Get restriction and basis from the RHS function
34     CeedOperator     *sub_ops;
35     CeedOperatorField field;
36     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
37 
38     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops));
39     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
40     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
41     PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q));
42 
43     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field));
44     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i));
45     PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data));
46 
47     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx));
48   }
49 
50   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
51   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size));
52 
53   switch (dim) {
54     case 2:
55       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection2D, MassFunction_Advection2D_loc, &qf_mass));
56       break;
57     case 3:
58       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection, MassFunction_Advection_loc, &qf_mass));
59       break;
60   }
61 
62   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx));
63   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
64   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
65   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
66   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
67   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
68   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
69 
70   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
71   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
72   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed));
73   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data));
74   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
75   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
76 
77   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qf_ctx));
78   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
79   PetscFunctionReturn(PETSC_SUCCESS);
80 }
81 
82 PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
83   WindType             wind_type;
84   AdvectionICType      advectionic_type;
85   BubbleContinuityType bubble_continuity_type;
86   StabilizationType    stab;
87   StabilizationTauType stab_tau;
88   SetupContextAdv      setup_context;
89   User                 user = *(User *)ctx;
90   MPI_Comm             comm = user->comm;
91   Ceed                 ceed = user->ceed;
92   PetscBool            implicit;
93   AdvectionContext     advection_ctx;
94   CeedQFunctionContext advection_context;
95   PetscInt             dim;
96 
97   PetscFunctionBeginUser;
98   PetscCall(PetscCalloc1(1, &setup_context));
99   PetscCall(PetscCalloc1(1, &advection_ctx));
100   PetscCall(DMGetDimension(dm, &dim));
101 
102   // ------------------------------------------------------
103   //               SET UP ADVECTION
104   // ------------------------------------------------------
105   switch (dim) {
106     case 2:
107       problem->dim                               = 2;
108       problem->ics.qfunction                     = ICsAdvection2d;
109       problem->ics.qfunction_loc                 = ICsAdvection2d_loc;
110       problem->apply_vol_rhs.qfunction           = RHS_Advection2d;
111       problem->apply_vol_rhs.qfunction_loc       = RHS_Advection2d_loc;
112       problem->apply_vol_ifunction.qfunction     = IFunction_Advection2d;
113       problem->apply_vol_ifunction.qfunction_loc = IFunction_Advection2d_loc;
114       problem->apply_inflow.qfunction            = Advection2d_InOutFlow;
115       problem->apply_inflow.qfunction_loc        = Advection2d_InOutFlow_loc;
116       problem->compute_exact_solution_error      = PETSC_TRUE;
117       problem->print_info                        = PRINT_ADVECTION;
118       break;
119     case 3:
120       problem->dim                               = 3;
121       problem->ics.qfunction                     = ICsAdvection;
122       problem->ics.qfunction_loc                 = ICsAdvection_loc;
123       problem->apply_vol_rhs.qfunction           = RHS_Advection;
124       problem->apply_vol_rhs.qfunction_loc       = RHS_Advection_loc;
125       problem->apply_vol_ifunction.qfunction     = IFunction_Advection;
126       problem->apply_vol_ifunction.qfunction_loc = IFunction_Advection_loc;
127       problem->apply_inflow.qfunction            = Advection_InOutFlow;
128       problem->apply_inflow.qfunction_loc        = Advection_InOutFlow_loc;
129       problem->compute_exact_solution_error      = PETSC_FALSE;
130       problem->print_info                        = PRINT_ADVECTION;
131       break;
132   }
133 
134   // ------------------------------------------------------
135   //             Create the libCEED context
136   // ------------------------------------------------------
137   CeedScalar rc              = 1000.;  // m (Radius of bubble)
138   CeedScalar CtauS           = 0.;     // dimensionless
139   PetscBool  strong_form     = PETSC_FALSE;
140   CeedScalar E_wind          = 1.e6;  // J
141   CeedScalar Ctau_a          = PetscPowScalarInt(user->app_ctx->degree, 2);
142   CeedScalar Ctau_t          = 0.;
143   PetscReal  wind[3]         = {1., 0, 0};  // m/s
144   CeedScalar diffusion_coeff = 0.;
145   PetscReal  domain_min[3], domain_max[3], domain_size[3];
146   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
147   for (PetscInt i = 0; i < problem->dim; i++) domain_size[i] = domain_max[i] - domain_min[i];
148 
149   // ------------------------------------------------------
150   //             Create the PETSc context
151   // ------------------------------------------------------
152   PetscScalar meter    = 1e-2;  // 1 meter in scaled length units
153   PetscScalar kilogram = 1e-6;  // 1 kilogram in scaled mass units
154   PetscScalar second   = 1e-2;  // 1 second in scaled time units
155   PetscScalar Joule;
156 
157   // ------------------------------------------------------
158   //              Command line Options
159   // ------------------------------------------------------
160   PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL);
161   // -- Physics
162   PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL));
163   PetscBool translation;
164   PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, WindTypes, (PetscEnum)(wind_type = WIND_ROTATION), (PetscEnum *)&wind_type,
165                              &translation));
166   PetscInt  n = problem->dim;
167   PetscBool user_wind;
168   PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind));
169   PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL));
170   PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL));
171   PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form,
172                              NULL));
173   PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL));
174   PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvectionICTypes,
175                              (PetscEnum)(advectionic_type = ADVECTIONIC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL));
176   bubble_continuity_type = problem->dim == 3 ? BUBBLE_CONTINUITY_SMOOTH : BUBBLE_CONTINUITY_COSINE;
177   PetscCall(PetscOptionsEnum("-bubble_continuity", "Smooth, back_sharp, or thick", NULL, BubbleContinuityTypes, (PetscEnum)bubble_continuity_type,
178                              (PetscEnum *)&bubble_continuity_type, NULL));
179   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
180   PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU),
181                              (PetscEnum *)&stab_tau, NULL));
182   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
183   PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization ", NULL, Ctau_a, &Ctau_a, NULL));
184   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
185 
186   // -- Units
187   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
188   meter = fabs(meter);
189   PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL));
190   kilogram = fabs(kilogram);
191   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
192   second = fabs(second);
193 
194   // -- Warnings
195   if (wind_type == WIND_ROTATION && user_wind) {
196     PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n"));
197   }
198   if (wind_type == WIND_TRANSLATION && advectionic_type == ADVECTIONIC_BUBBLE_CYLINDER && wind[2] != 0.) {
199     wind[2] = 0;
200     PetscCall(
201         PetscPrintf(comm, "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n"));
202   }
203   if (stab == STAB_NONE && CtauS != 0) {
204     PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n"));
205   }
206   PetscOptionsEnd();
207 
208   if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized;
209 
210   // ------------------------------------------------------
211   //           Set up the PETSc context
212   // ------------------------------------------------------
213   // -- Define derived units
214   Joule = kilogram * PetscSqr(meter) / PetscSqr(second);
215 
216   user->units->meter    = meter;
217   user->units->kilogram = kilogram;
218   user->units->second   = second;
219   user->units->Joule    = Joule;
220 
221   // ------------------------------------------------------
222   //           Set up the libCEED context
223   // ------------------------------------------------------
224   // -- Scale variables to desired units
225   E_wind *= Joule;
226   rc = fabs(rc) * meter;
227   for (PetscInt i = 0; i < problem->dim; i++) {
228     wind[i] *= (meter / second);
229     domain_size[i] *= meter;
230   }
231   problem->dm_scale = meter;
232 
233   // -- Setup Context
234   setup_context->rc                     = rc;
235   setup_context->lx                     = domain_size[0];
236   setup_context->ly                     = domain_size[1];
237   setup_context->lz                     = problem->dim == 3 ? domain_size[2] : 0.;
238   setup_context->wind[0]                = wind[0];
239   setup_context->wind[1]                = wind[1];
240   setup_context->wind[2]                = problem->dim == 3 ? wind[2] : 0.;
241   setup_context->wind_type              = wind_type;
242   setup_context->initial_condition_type = advectionic_type;
243   setup_context->bubble_continuity_type = bubble_continuity_type;
244   setup_context->time                   = 0;
245 
246   // -- QFunction Context
247   user->phys->implicit             = implicit;
248   advection_ctx->CtauS             = CtauS;
249   advection_ctx->E_wind            = E_wind;
250   advection_ctx->implicit          = implicit;
251   advection_ctx->strong_form       = strong_form;
252   advection_ctx->stabilization     = stab;
253   advection_ctx->stabilization_tau = stab_tau;
254   advection_ctx->Ctau_a            = Ctau_a;
255   advection_ctx->Ctau_t            = Ctau_t;
256   advection_ctx->diffusion_coeff   = diffusion_coeff;
257 
258   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context));
259   PetscCallCeed(ceed,
260                 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
261   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc));
262 
263   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &advection_context));
264   PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx));
265   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_context, CEED_MEM_HOST, FreeContextPetsc));
266   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_context, "timestep size", offsetof(struct AdvectionContext_, dt), 1,
267                                                          "Size of timestep, delta t"));
268   problem->apply_vol_rhs.qfunction_context = advection_context;
269   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_vol_ifunction.qfunction_context));
270   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_inflow.qfunction_context));
271   PetscFunctionReturn(PETSC_SUCCESS);
272 }
273 
274 PetscErrorCode PRINT_ADVECTION(User user, ProblemData problem, AppCtx app_ctx) {
275   MPI_Comm         comm = user->comm;
276   Ceed             ceed = user->ceed;
277   SetupContextAdv  setup_ctx;
278   AdvectionContext advection_ctx;
279 
280   PetscFunctionBeginUser;
281   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfunction_context, CEED_MEM_HOST, &setup_ctx));
282   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &advection_ctx));
283   PetscCall(PetscPrintf(comm,
284                         "  Problem:\n"
285                         "    Problem Name                       : %s\n"
286                         "    Stabilization                      : %s\n"
287                         "    Initial Condition Type             : %s\n"
288                         "    Bubble Continuity                  : %s\n"
289                         "    Wind Type                          : %s\n",
290                         app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization], AdvectionICTypes[setup_ctx->initial_condition_type],
291                         BubbleContinuityTypes[setup_ctx->bubble_continuity_type], WindTypes[setup_ctx->wind_type]));
292 
293   if (setup_ctx->wind_type == WIND_TRANSLATION) {
294     switch (problem->dim) {
295       case 2:
296         PetscCall(PetscPrintf(comm, "    Background Wind                    : %f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1]));
297         break;
298       case 3:
299         PetscCall(
300             PetscPrintf(comm, "    Background Wind                    : %f,%f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1], setup_ctx->wind[2]));
301         break;
302     }
303   }
304   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfunction_context, &setup_ctx));
305   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &advection_ctx));
306   PetscFunctionReturn(PETSC_SUCCESS);
307 }
308