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