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 const char *const AdvDifWindTypes[] = {"ROTATION", "TRANSLATION", "BOUNDARY_LAYER", "WindType", "ADVDIF_WIND_", NULL};
15 const char *const AdvDifICTypes[] = {"SPHERE", "CYLINDER", "COSINE_HILL", "SKEW", "WAVE", "BOUNDARY_LAYER", "AdvectionICType", "ADVDIF_IC_", NULL};
16 const char *const AdvDifWaveTypes[] = {"SINE", "SQUARE", "AdvDifWaveType", "ADVDIF_WAVE_", NULL};
17 const char *const AdvDifBubbleContinuityTypes[] = {"SMOOTH", "BACK_SHARP", "THICK", "COSINE", "BubbleContinuityType", "ADVDIF_BUBBLE_CONTINUITY_",
18 NULL};
19 const char *const StabilizationTauTypes[] = {"CTAU", "ADVDIFF_SHAKIB", "StabilizationTauType", "STAB_TAU_", NULL};
20
PRINT_ADVECTION(Honee honee,ProblemData problem,AppCtx app_ctx)21 static PetscErrorCode PRINT_ADVECTION(Honee honee, ProblemData problem, AppCtx app_ctx) {
22 MPI_Comm comm = honee->comm;
23 Ceed ceed = honee->ceed;
24 SetupContextAdv setup_ctx;
25 AdvectionContext advection_ctx;
26 PetscInt dim;
27
28 PetscFunctionBeginUser;
29 PetscCall(DMGetDimension(honee->dm, &dim));
30 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &setup_ctx));
31 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &advection_ctx));
32 PetscCall(PetscPrintf(comm,
33 " Problem:\n"
34 " Problem Name : %s\n"
35 " Stabilization : %s\n"
36 " Stabilization Tau : %s\n"
37 " Wind Type : %s\n",
38 app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization],
39 StabilizationTauTypes[advection_ctx->stabilization_tau], AdvDifWindTypes[setup_ctx->wind_type]));
40
41 if (setup_ctx->wind_type == ADVDIF_WIND_TRANSLATION) {
42 CeedScalar *wind = setup_ctx->wind;
43 switch (dim) {
44 case 2:
45 PetscCall(PetscPrintf(comm, " Background Wind : %f,%f\n", wind[0], wind[1]));
46 break;
47 case 3:
48 PetscCall(PetscPrintf(comm, " Background Wind : %f,%f,%f\n", wind[0], wind[1], wind[2]));
49 break;
50 }
51 }
52
53 PetscCall(PetscPrintf(comm, " Initial Condition Type : %s\n", AdvDifICTypes[setup_ctx->initial_condition_type]));
54 switch (setup_ctx->initial_condition_type) {
55 case ADVDIF_IC_SKEW:
56 case ADVDIF_IC_COSINE_HILL:
57 case ADVDIF_IC_BOUNDARY_LAYER:
58 break;
59 case ADVDIF_IC_BUBBLE_SPHERE:
60 case ADVDIF_IC_BUBBLE_CYLINDER:
61 PetscCall(PetscPrintf(comm, " Bubble Continuity : %s\n", AdvDifBubbleContinuityTypes[setup_ctx->bubble_continuity_type]));
62 break;
63 case ADVDIF_IC_WAVE:
64 PetscCall(PetscPrintf(comm, " Wave Type : %s\n", AdvDifWaveTypes[setup_ctx->wave_type]));
65 break;
66 }
67
68 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &setup_ctx));
69 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &advection_ctx));
70 PetscFunctionReturn(PETSC_SUCCESS);
71 }
72
73 // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping
74 //
75 // Only used for SUPG stabilization
CreateKSPMassOperator_AdvectionStabilized(Honee honee,CeedOperator * op_mass)76 PetscErrorCode CreateKSPMassOperator_AdvectionStabilized(Honee honee, CeedOperator *op_mass) {
77 Ceed ceed = honee->ceed;
78 CeedInt num_comp_q, q_data_size;
79 CeedQFunction qf_mass = NULL;
80 CeedElemRestriction elem_restr_q, elem_restr_qd;
81 CeedBasis basis_q;
82 CeedVector q_data;
83 CeedQFunctionContext qfctx = NULL;
84 PetscInt dim;
85
86 PetscFunctionBeginUser;
87 PetscCall(DMGetDimension(honee->dm, &dim));
88 { // Get restriction and basis from the RHS function
89 CeedOperator *sub_ops;
90 CeedOperatorField op_field;
91 PetscInt sub_op_index = 0; // will be 0 for the volume op
92
93 PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops));
94 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field));
95 PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL));
96 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &op_field));
97 PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_qd, NULL, &q_data));
98
99 PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qfctx));
100 }
101
102 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
103 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd, &q_data_size));
104
105 switch (dim) {
106 case 2:
107 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection2D, MassFunction_Advection2D_loc, &qf_mass));
108 break;
109 case 3:
110 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection, MassFunction_Advection_loc, &qf_mass));
111 break;
112 }
113
114 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qfctx));
115 PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
116 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
117 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
118 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
119 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
120 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
121
122 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
123 PetscCallCeed(ceed, CeedOperatorSetName(*op_mass, "RHS Mass Operator, Advection-Diffusion Stabilized"));
124 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
125 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, honee->q_ceed));
126 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
127 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
128 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
129
130 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
131 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
132 PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
133 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
134 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qfctx));
135 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
136 PetscFunctionReturn(PETSC_SUCCESS);
137 }
138
139 /**
140 @brief Create RHS CeedOperator for direct projection of divergence of diffusive flux
141
142 @param[in] honee `Honee` context
143 @param[in] diff_flux_proj `DivDiffFluxProjectionData` object
144 @param[out] op_rhs Operator to calculate the RHS of the L^2 projection
145 **/
DivDiffFluxProjectionCreateRHS_Direct_AdvDif(Honee honee,DivDiffFluxProjectionData diff_flux_proj,CeedOperator * op_rhs)146 static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_AdvDif(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) {
147 Ceed ceed = honee->ceed;
148 NodalProjectionData projection = diff_flux_proj->projection;
149 PetscInt dim, num_comp_q;
150 CeedQFunctionContext advection_qfctx = NULL;
151
152 PetscFunctionBeginUser;
153 // -- Get Pre-requisite things
154 PetscCall(DMGetDimension(projection->dm, &dim));
155 PetscCall(DMGetFieldNumComps(honee->dm, 0, &num_comp_q));
156
157 { // Get advection-diffusion QF context
158 CeedOperator *sub_ops;
159 PetscInt sub_op_index = 0; // will be 0 for the volume op
160
161 if (honee->op_ifunction) PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_ifunction, &sub_ops));
162 else PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops));
163 PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx));
164 }
165 PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, op_rhs));
166 { // Add the volume integral CeedOperator
167 CeedQFunction qf_rhs_volume = NULL;
168 CeedOperator op_rhs_volume;
169 CeedVector q_data;
170 CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL, elem_restr_q;
171 CeedBasis basis_diff_flux = NULL, basis_q;
172 CeedInt q_data_size;
173
174 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q));
175 PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q));
176 PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL));
177 PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, &elem_restr_qd, &q_data, &q_data_size));
178 switch (dim) {
179 case 2:
180 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_2D, DivDiffusiveFluxVolumeRHS_AdvDif_2D_loc,
181 &qf_rhs_volume));
182 break;
183 case 3:
184 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_3D, DivDiffusiveFluxVolumeRHS_AdvDif_3D_loc,
185 &qf_rhs_volume));
186 break;
187 }
188 PetscCheck(qf_rhs_volume, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
189
190 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, advection_qfctx));
191 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
192 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE));
193 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD));
194
195 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume));
196 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
197 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
198 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE));
199
200 PetscCallCeed(ceed, CeedOperatorCompositeAddSub(*op_rhs, op_rhs_volume));
201
202 PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
203 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
204 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
205 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
206 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume));
207 PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
208 PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume));
209 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume));
210 }
211
212 { // Add the boundary integral CeedOperator
213 CeedQFunction qf_rhs_boundary;
214 DMLabel face_sets_label;
215 PetscInt num_face_set_values, *face_set_values;
216 CeedInt q_data_size;
217
218 // -- Build RHS operator
219 switch (dim) {
220 case 2:
221 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_2D, DivDiffusiveFluxBoundaryRHS_AdvDif_2D_loc,
222 &qf_rhs_boundary));
223 break;
224 case 3:
225 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_3D, DivDiffusiveFluxBoundaryRHS_AdvDif_3D_loc,
226 &qf_rhs_boundary));
227 break;
228 }
229
230 PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size));
231 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, advection_qfctx));
232 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
233 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE));
234 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP));
235
236 PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label));
237 PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values));
238 for (PetscInt f = 0; f < num_face_set_values; f++) {
239 DMLabel face_orientation_label;
240 PetscInt num_orientations_values, *orientation_values;
241
242 {
243 char *face_orientation_label_name;
244
245 PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name));
246 PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label));
247 PetscCall(PetscFree(face_orientation_label_name));
248 }
249 PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values));
250 for (PetscInt o = 0; o < num_orientations_values; o++) {
251 CeedOperator op_rhs_boundary;
252 CeedBasis basis_q, basis_diff_flux_boundary;
253 CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary;
254 CeedVector q_data;
255 CeedInt q_data_size;
256 PetscInt orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1;
257
258 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q));
259 PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q));
260 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0,
261 &elem_restr_diff_flux_boundary));
262 PetscCall(DMPlexCeedBasisCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary));
263 PetscCall(QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, &elem_restr_qdata, &q_data, &q_data_size));
264
265 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary));
266 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
267 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data));
268 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary,
269 CEED_VECTOR_ACTIVE));
270
271 PetscCallCeed(ceed, CeedOperatorCompositeAddSub(*op_rhs, op_rhs_boundary));
272
273 PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary));
274 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata));
275 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
276 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary));
277 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
278 PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary));
279 PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
280 }
281 PetscCall(PetscFree(orientation_values));
282 }
283 PetscCall(PetscFree(face_set_values));
284 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary));
285 }
286
287 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx));
288 PetscFunctionReturn(PETSC_SUCCESS);
289 }
290
291 /**
292 @brief Create RHS CeedOperator for indirect projection of divergence of diffusive flux
293
294 @param[in] honee `Honee` context
295 @param[in] diff_flux_proj `DivDiffFluxProjectionData` object
296 @param[out] op_rhs Operator to calculate the RHS of the L^2 projection
297 **/
DivDiffFluxProjectionCreateRHS_Indirect_AdvDif(Honee honee,DivDiffFluxProjectionData diff_flux_proj,CeedOperator * op_rhs)298 static PetscErrorCode DivDiffFluxProjectionCreateRHS_Indirect_AdvDif(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) {
299 Ceed ceed = honee->ceed;
300 NodalProjectionData projection = diff_flux_proj->projection;
301 CeedBasis basis_diff_flux, basis_q;
302 CeedElemRestriction elem_restr_diff_flux, elem_restr_qd, elem_restr_q;
303 CeedVector q_data;
304 CeedInt q_data_size;
305 PetscInt dim, num_comp_q;
306 PetscInt height = 0, dm_field = 0;
307 CeedQFunction qf_rhs = NULL;
308 CeedQFunctionContext advection_qfctx = NULL;
309
310 PetscFunctionBeginUser;
311 PetscCall(DMGetDimension(projection->dm, &dim));
312 PetscCall(DMGetFieldNumComps(honee->dm, 0, &num_comp_q));
313
314 { // Get advection-diffusion QF context
315 CeedOperator *sub_ops;
316 PetscInt sub_op_index = 0; // will be 0 for the volume op
317
318 if (honee->op_ifunction) PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_ifunction, &sub_ops));
319 else PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops));
320 PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx));
321 }
322 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q));
323 PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q));
324 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &elem_restr_diff_flux));
325 PetscCall(DMPlexCeedBasisCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &basis_diff_flux));
326 PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, &elem_restr_qd, &q_data, &q_data_size));
327
328 switch (dim) {
329 case 2:
330 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_2D, DiffusiveFluxRHS_AdvDif_2D_loc, &qf_rhs));
331 break;
332 case 3:
333 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_3D, DiffusiveFluxRHS_AdvDif_3D_loc, &qf_rhs));
334 break;
335 }
336 PetscCheck(qf_rhs, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
337
338 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, advection_qfctx));
339 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
340 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE));
341 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP));
342
343 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, op_rhs));
344 PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
345 PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
346 PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
347
348 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs));
349 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx));
350 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
351 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
352 PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
353 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
354 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux));
355 PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
356 PetscFunctionReturn(PETSC_SUCCESS);
357 }
358
AdvectionInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def,CeedQFunction * qf)359 static PetscErrorCode AdvectionInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) {
360 HoneeBCStruct honee_bc;
361 DM dm;
362 PetscInt dim;
363
364 PetscFunctionBeginUser;
365 PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
366 PetscCall(BCDefinitionGetDM(bc_def, &dm));
367 PetscCall(DMGetDimension(dm, &dim));
368 switch (dim) {
369 case 2:
370 PetscCall(HoneeBCCreateIFunctionQF(bc_def, Advection2d_InOutFlow, Advection2d_InOutFlow_loc, honee_bc->qfctx, qf));
371 break;
372 case 3:
373 PetscCall(HoneeBCCreateIFunctionQF(bc_def, Advection_InOutFlow, Advection_InOutFlow_loc, honee_bc->qfctx, qf));
374 break;
375 }
376 PetscFunctionReturn(PETSC_SUCCESS);
377 }
378
379 static const char *const component_names[] = {"Density", "MomentumX", "MomentumY", "MomentumZ", "TotalEnergy"};
380
NS_ADVECTION(ProblemData problem,DM dm,void * ctx)381 PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx) {
382 AdvDifWindType wind_type;
383 AdvDifICType advectionic_type;
384 AdvDifBubbleContinuityType bubble_continuity_type = -1;
385 StabilizationType stab;
386 StabilizationTauType stab_tau;
387 SetupContextAdv setup_context;
388 Honee honee = *(Honee *)ctx;
389 MPI_Comm comm = honee->comm;
390 Ceed ceed = honee->ceed;
391 PetscBool implicit;
392 AdvectionContext advection_ctx;
393 CeedQFunctionContext advection_qfctx, ics_qfctx;
394 PetscInt dim;
395
396 PetscFunctionBeginUser;
397 PetscCall(PetscNew(&setup_context));
398 PetscCall(PetscNew(&advection_ctx));
399 PetscCall(DMGetDimension(dm, &dim));
400
401 // ------------------------------------------------------
402 // Create the QFunction context
403 // ------------------------------------------------------
404 CeedScalar rc = 1000.; // m (Radius of bubble)
405 CeedScalar CtauS = 0.; // dimensionless
406 PetscBool strong_form = PETSC_FALSE;
407 CeedScalar E_wind = 1.e6; // J
408 CeedScalar Ctau_a = PetscPowScalarInt(honee->app_ctx->degree, 2);
409 CeedScalar Ctau_d = PetscPowScalarInt(honee->app_ctx->degree, 4);
410 CeedScalar Ctau_t = 0.;
411 PetscReal wind[3] = {1., 0, 0}; // m/s
412 CeedScalar diffusion_coeff = 0.;
413 CeedScalar wave_frequency = 2 * M_PI;
414 CeedScalar wave_phase = 0;
415 AdvDifWaveType wave_type = -1;
416 PetscScalar bl_height_factor = 1.;
417 PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.};
418 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
419 for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i];
420
421 // ------------------------------------------------------
422 // Command line Options
423 // ------------------------------------------------------
424 PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL);
425 // -- Physics
426 PetscBool translation;
427 PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, AdvDifWindTypes, (PetscEnum)(wind_type = ADVDIF_WIND_ROTATION),
428 (PetscEnum *)&wind_type, &translation));
429 PetscInt n = dim;
430 PetscBool user_wind;
431 PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind));
432 PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL));
433 PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL));
434 PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form,
435 NULL));
436 PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL));
437 PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
438 PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU),
439 (PetscEnum *)&stab_tau, NULL));
440 PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
441 PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization, advection component", NULL, Ctau_a, &Ctau_a, NULL));
442 PetscCall(PetscOptionsScalar("-Ctau_d", "Coefficient for the stabilization, diffusion component", NULL, Ctau_d, &Ctau_d, NULL));
443 PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
444 PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvDifICTypes,
445 (PetscEnum)(advectionic_type = ADVDIF_IC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL));
446 // IC-specific options
447 switch (advectionic_type) {
448 case ADVDIF_IC_WAVE:
449 PetscCall(PetscOptionsDeprecated("-wave_type", "-advection_ic_wave_type", "HONEE 0.0", NULL));
450 PetscCall(PetscOptionsDeprecated("-wave_frequency", "-advection_ic_wave_frequency", "HONEE 0.0", NULL));
451 PetscCall(PetscOptionsDeprecated("-wave_phase", "-advection_ic_wave_phase", "HONEE 0.0", NULL));
452 PetscCall(PetscOptionsEnum("-advection_ic_wave_type", "Type of wave", NULL, AdvDifWaveTypes, (PetscEnum)(wave_type = ADVDIF_WAVE_SINE),
453 (PetscEnum *)&wave_type, NULL));
454 PetscCall(PetscOptionsScalar("-advection_ic_wave_frequency", "Frequency of sine wave", NULL, wave_frequency, &wave_frequency, NULL));
455 PetscCall(PetscOptionsScalar("-advection_ic_wave_phase", "Length correction", NULL, wave_phase, &wave_phase, NULL));
456 break;
457 case ADVDIF_IC_BOUNDARY_LAYER:
458 PetscCall(PetscOptionsScalar("-advection_ic_bl_height_factor", "Height of boundary layer in IC", NULL, bl_height_factor, &bl_height_factor,
459 NULL));
460 break;
461 case ADVDIF_IC_BUBBLE_CYLINDER:
462 case ADVDIF_IC_BUBBLE_SPHERE:
463 PetscCall(PetscOptionsDeprecated("-rc", "-advection_ic_bubble_rc", "HONEE 0.0", NULL));
464 PetscCall(PetscOptionsDeprecated("-bubble_continuity", "-advection_ic_bubble_continuity", "HONEE 0.0", NULL));
465 PetscCall(PetscOptionsScalar("-advection_ic_bubble_rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL));
466 bubble_continuity_type = dim == 3 ? ADVDIF_BUBBLE_CONTINUITY_SMOOTH : ADVDIF_BUBBLE_CONTINUITY_COSINE;
467 PetscCall(PetscOptionsEnum("-advection_ic_bubble_continuity", "Smooth, back_sharp, or thick", NULL, AdvDifBubbleContinuityTypes,
468 (PetscEnum)bubble_continuity_type, (PetscEnum *)&bubble_continuity_type, NULL));
469 break;
470 case ADVDIF_IC_SKEW:
471 case ADVDIF_IC_COSINE_HILL:
472 break;
473 }
474
475 // -- Warnings
476 if (wind_type == ADVDIF_WIND_ROTATION && user_wind) {
477 PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n"));
478 }
479 if (wind_type == ADVDIF_WIND_TRANSLATION && advectionic_type == ADVDIF_IC_BUBBLE_CYLINDER && wind[2] != 0.) {
480 wind[2] = 0;
481 PetscCall(PetscPrintf(comm,
482 "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n"));
483 }
484 if (stab == STAB_NONE && CtauS != 0) {
485 PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n"));
486 }
487 PetscOptionsEnd();
488
489 if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized;
490
491 // ------------------------------------------------------
492 // Set up the QFunction contexts
493 // ------------------------------------------------------
494 // -- Scale variables to desired units
495 Units units = honee->units;
496 E_wind *= units->Joule;
497 rc = fabs(rc) * units->meter;
498 for (PetscInt i = 0; i < dim; i++) {
499 wind[i] *= (units->meter / units->second);
500 }
501
502 // -- Setup Context
503 setup_context->rc = rc;
504 setup_context->lx = domain_size[0];
505 setup_context->ly = domain_size[1];
506 setup_context->lz = dim == 3 ? domain_size[2] : 0.;
507 setup_context->wind[0] = wind[0];
508 setup_context->wind[1] = wind[1];
509 setup_context->wind[2] = dim == 3 ? wind[2] : 0.;
510 setup_context->wind_type = wind_type;
511 setup_context->initial_condition_type = advectionic_type;
512 setup_context->bubble_continuity_type = bubble_continuity_type;
513 setup_context->time = 0;
514 setup_context->wave_frequency = wave_frequency;
515 setup_context->wave_phase = wave_phase;
516 setup_context->wave_type = wave_type;
517 setup_context->bl_height_factor = bl_height_factor;
518
519 // -- QFunction Context
520 honee->phys->implicit = implicit;
521 advection_ctx->CtauS = CtauS;
522 advection_ctx->E_wind = E_wind;
523 advection_ctx->implicit = implicit;
524 advection_ctx->strong_form = strong_form;
525 advection_ctx->stabilization = stab;
526 advection_ctx->stabilization_tau = stab_tau;
527 advection_ctx->Ctau_a = Ctau_a;
528 advection_ctx->Ctau_d = Ctau_d;
529 advection_ctx->Ctau_t = Ctau_t;
530 advection_ctx->diffusion_coeff = diffusion_coeff;
531 advection_ctx->divFdiff_method = honee->app_ctx->divFdiffproj_method;
532
533 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &ics_qfctx));
534 PetscCallCeed(ceed, CeedQFunctionContextSetData(ics_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
535 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(ics_qfctx, CEED_MEM_HOST, FreeContextPetsc));
536
537 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &advection_qfctx));
538 PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx));
539 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_qfctx, CEED_MEM_HOST, FreeContextPetsc));
540 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_qfctx, "timestep size", offsetof(struct AdvectionContext_, dt), 1,
541 "Size of timestep, delta t"));
542
543 // ------------------------------------------------------
544 // SET UP ADVECTION
545 // ------------------------------------------------------
546 problem->print_info = PRINT_ADVECTION;
547 problem->num_comps_jac_data = 0;
548 switch (dim) {
549 case 2:
550 problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsAdvection2d, .qf_loc = ICsAdvection2d_loc};
551 problem->apply_vol_rhs = (HoneeQFSpec){.qf_func_ptr = RHS_Advection2d, .qf_loc = RHS_Advection2d_loc};
552 problem->apply_vol_ifunction = (HoneeQFSpec){.qf_func_ptr = IFunction_Advection2d, .qf_loc = IFunction_Advection2d_loc};
553 problem->compute_exact_solution_error = PETSC_TRUE;
554 break;
555 case 3:
556 problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsAdvection, .qf_loc = ICsAdvection_loc};
557 problem->apply_vol_rhs = (HoneeQFSpec){.qf_func_ptr = RHS_Advection, .qf_loc = RHS_Advection_loc};
558 problem->apply_vol_ifunction = (HoneeQFSpec){.qf_func_ptr = IFunction_Advection, .qf_loc = IFunction_Advection_loc};
559 problem->compute_exact_solution_error = PETSC_FALSE;
560 break;
561 }
562 problem->ics.qfctx = ics_qfctx;
563 problem->apply_vol_rhs.qfctx = advection_qfctx;
564 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_vol_ifunction.qfctx));
565
566 problem->num_components = 5;
567 PetscCall(PetscMalloc1(problem->num_components, &problem->component_names));
568 for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(component_names[i], &problem->component_names[i]));
569
570 PetscCall(DivDiffFluxProjectionCreate(honee, honee->app_ctx->divFdiffproj_method, 1, &honee->diff_flux_proj));
571 if (honee->diff_flux_proj) {
572 DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj;
573 NodalProjectionData projection = diff_flux_proj->projection;
574 PetscSection section;
575
576 diff_flux_proj->CreateRHSOperator_Direct = DivDiffFluxProjectionCreateRHS_Direct_AdvDif;
577 diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_AdvDif;
578 PetscCall(DMGetLocalSection(projection->dm, §ion));
579 switch (honee->diff_flux_proj->method) {
580 case DIV_DIFF_FLUX_PROJ_DIRECT: {
581 PetscCall(PetscSectionSetFieldName(section, 0, ""));
582 PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_Scalar"));
583 } break;
584 case DIV_DIFF_FLUX_PROJ_INDIRECT: {
585 PetscCall(PetscSectionSetFieldName(section, 0, ""));
586 PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_ScalarX"));
587 PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_ScalarY"));
588 if (dim >= 3) PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_ScalarZ"));
589 } break;
590 case DIV_DIFF_FLUX_PROJ_NONE:
591 SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
592 DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]);
593 break;
594 }
595 }
596
597 for (PetscCount b = 0; b < problem->num_bc_defs; b++) {
598 BCDefinition bc_def = problem->bc_defs[b];
599 const char *name;
600
601 PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL));
602 if (!strcmp(name, "inflow")) {
603 HoneeBCStruct honee_bc;
604
605 PetscCall(PetscNew(&honee_bc));
606 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &honee_bc->qfctx));
607 honee_bc->honee = honee;
608 honee_bc->num_comps_jac_data = 0;
609 PetscCall(BCDefinitionSetContext(bc_def, (PetscCtxDestroyFn *)HoneeBCDestroy, honee_bc));
610
611 PetscCall(BCDefinitionSetIFunction(bc_def, AdvectionInflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
612 PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL));
613 }
614 }
615 PetscFunctionReturn(PETSC_SUCCESS);
616 }
617