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