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 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 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 **/ 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 CeedInt num_comp_q; 150 PetscInt dim, label_value = 0; 151 DMLabel domain_label = NULL; 152 CeedQFunctionContext advection_qfctx = NULL; 153 154 PetscFunctionBeginUser; 155 // -- Get Pre-requisite things 156 PetscCall(DMGetDimension(projection->dm, &dim)); 157 PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q)); 158 159 { // Get advection-diffusion QF context 160 CeedOperator *sub_ops; 161 PetscInt sub_op_index = 0; // will be 0 for the volume op 162 163 if (honee->op_ifunction) PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_ifunction, &sub_ops)); 164 else PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 165 PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx)); 166 } 167 PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, op_rhs)); 168 { // Add the volume integral CeedOperator 169 CeedQFunction qf_rhs_volume = NULL; 170 CeedOperator op_rhs_volume; 171 CeedVector q_data; 172 CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL; 173 CeedBasis basis_diff_flux = NULL; 174 CeedInt q_data_size; 175 176 PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL)); 177 PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data, 178 &q_data_size)); 179 switch (dim) { 180 case 2: 181 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_2D, DivDiffusiveFluxVolumeRHS_AdvDif_2D_loc, 182 &qf_rhs_volume)); 183 break; 184 case 3: 185 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_3D, DivDiffusiveFluxVolumeRHS_AdvDif_3D_loc, 186 &qf_rhs_volume)); 187 break; 188 } 189 PetscCheck(qf_rhs_volume, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim); 190 191 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, advection_qfctx)); 192 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 193 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE)); 194 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD)); 195 196 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume)); 197 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 198 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 199 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE)); 200 201 PetscCallCeed(ceed, CeedOperatorCompositeAddSub(*op_rhs, op_rhs_volume)); 202 203 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 204 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 205 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume)); 206 PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 207 PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume)); 208 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume)); 209 } 210 211 { // Add the boundary integral CeedOperator 212 CeedQFunction qf_rhs_boundary; 213 DMLabel face_sets_label; 214 PetscInt num_face_set_values, *face_set_values; 215 CeedInt q_data_size; 216 217 // -- Build RHS operator 218 switch (dim) { 219 case 2: 220 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_2D, DivDiffusiveFluxBoundaryRHS_AdvDif_2D_loc, 221 &qf_rhs_boundary)); 222 break; 223 case 3: 224 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_3D, DivDiffusiveFluxBoundaryRHS_AdvDif_3D_loc, 225 &qf_rhs_boundary)); 226 break; 227 } 228 229 PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size)); 230 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, advection_qfctx)); 231 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 232 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE)); 233 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP)); 234 235 PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label)); 236 PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values)); 237 for (PetscInt f = 0; f < num_face_set_values; f++) { 238 DMLabel face_orientation_label; 239 PetscInt num_orientations_values, *orientation_values; 240 241 { 242 char *face_orientation_label_name; 243 244 PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name)); 245 PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label)); 246 PetscCall(PetscFree(face_orientation_label_name)); 247 } 248 PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values)); 249 for (PetscInt o = 0; o < num_orientations_values; o++) { 250 CeedOperator op_rhs_boundary; 251 CeedBasis basis_q, basis_diff_flux_boundary; 252 CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary; 253 CeedVector q_data; 254 CeedInt q_data_size; 255 PetscInt orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1; 256 257 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q)); 258 PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q)); 259 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, 260 &elem_restr_diff_flux_boundary)); 261 PetscCall(CreateBasisFromPlex(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary)); 262 PetscCall(QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, honee->x_coord, &elem_restr_qdata, &q_data, 263 &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 **/ 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; 302 CeedElemRestriction elem_restr_diff_flux, elem_restr_qd; 303 CeedVector q_data; 304 CeedInt num_comp_q, q_data_size; 305 PetscInt dim; 306 PetscInt label_value = 0, height = 0, dm_field = 0; 307 DMLabel domain_label = NULL; 308 CeedQFunction qf_rhs = NULL; 309 CeedQFunctionContext advection_qfctx = NULL; 310 311 PetscFunctionBeginUser; 312 PetscCall(DMGetDimension(projection->dm, &dim)); 313 PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q)); 314 315 { // Get advection-diffusion QF context 316 CeedOperator *sub_ops; 317 PetscInt sub_op_index = 0; // will be 0 for the volume op 318 319 if (honee->op_ifunction) PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_ifunction, &sub_ops)); 320 else PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 321 PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx)); 322 } 323 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux)); 324 PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux)); 325 PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data, 326 &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", honee->elem_restr_q, honee->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, CeedBasisDestroy(&basis_diff_flux)); 351 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 352 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 353 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 354 PetscFunctionReturn(PETSC_SUCCESS); 355 } 356 357 static PetscErrorCode AdvectionInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) { 358 HoneeBCStruct honee_bc; 359 DM dm; 360 PetscInt dim; 361 362 PetscFunctionBeginUser; 363 PetscCall(BCDefinitionGetContext(bc_def, &honee_bc)); 364 PetscCall(BCDefinitionGetDM(bc_def, &dm)); 365 PetscCall(DMGetDimension(dm, &dim)); 366 switch (dim) { 367 case 2: 368 PetscCall(HoneeBCCreateIFunctionQF(bc_def, Advection2d_InOutFlow, Advection2d_InOutFlow_loc, honee_bc->qfctx, qf)); 369 break; 370 case 3: 371 PetscCall(HoneeBCCreateIFunctionQF(bc_def, Advection_InOutFlow, Advection_InOutFlow_loc, honee_bc->qfctx, qf)); 372 break; 373 } 374 PetscFunctionReturn(PETSC_SUCCESS); 375 } 376 377 static const char *const component_names[] = {"Density", "MomentumX", "MomentumY", "MomentumZ", "TotalEnergy"}; 378 379 PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx) { 380 AdvDifWindType wind_type; 381 AdvDifICType advectionic_type; 382 AdvDifBubbleContinuityType bubble_continuity_type = -1; 383 StabilizationType stab; 384 StabilizationTauType stab_tau; 385 SetupContextAdv setup_context; 386 Honee honee = *(Honee *)ctx; 387 MPI_Comm comm = honee->comm; 388 Ceed ceed = honee->ceed; 389 PetscBool implicit; 390 AdvectionContext advection_ctx; 391 CeedQFunctionContext advection_qfctx, ics_qfctx; 392 PetscInt dim; 393 394 PetscFunctionBeginUser; 395 PetscCall(PetscNew(&setup_context)); 396 PetscCall(PetscNew(&advection_ctx)); 397 PetscCall(DMGetDimension(dm, &dim)); 398 399 // ------------------------------------------------------ 400 // Create the QFunction context 401 // ------------------------------------------------------ 402 CeedScalar rc = 1000.; // m (Radius of bubble) 403 CeedScalar CtauS = 0.; // dimensionless 404 PetscBool strong_form = PETSC_FALSE; 405 CeedScalar E_wind = 1.e6; // J 406 CeedScalar Ctau_a = PetscPowScalarInt(honee->app_ctx->degree, 2); 407 CeedScalar Ctau_d = PetscPowScalarInt(honee->app_ctx->degree, 4); 408 CeedScalar Ctau_t = 0.; 409 PetscReal wind[3] = {1., 0, 0}; // m/s 410 CeedScalar diffusion_coeff = 0.; 411 CeedScalar wave_frequency = 2 * M_PI; 412 CeedScalar wave_phase = 0; 413 AdvDifWaveType wave_type = -1; 414 PetscScalar bl_height_factor = 1.; 415 PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.}; 416 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 417 for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 418 419 // ------------------------------------------------------ 420 // Command line Options 421 // ------------------------------------------------------ 422 PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL); 423 // -- Physics 424 PetscBool translation; 425 PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, AdvDifWindTypes, (PetscEnum)(wind_type = ADVDIF_WIND_ROTATION), 426 (PetscEnum *)&wind_type, &translation)); 427 PetscInt n = dim; 428 PetscBool user_wind; 429 PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind)); 430 PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL)); 431 PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL)); 432 PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form, 433 NULL)); 434 PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL)); 435 PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 436 PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU), 437 (PetscEnum *)&stab_tau, NULL)); 438 PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 439 PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization, advection component", NULL, Ctau_a, &Ctau_a, NULL)); 440 PetscCall(PetscOptionsScalar("-Ctau_d", "Coefficient for the stabilization, diffusion component", NULL, Ctau_d, &Ctau_d, NULL)); 441 PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 442 PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvDifICTypes, 443 (PetscEnum)(advectionic_type = ADVDIF_IC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL)); 444 // IC-specific options 445 switch (advectionic_type) { 446 case ADVDIF_IC_WAVE: 447 PetscCall(PetscOptionsDeprecated("-wave_type", "-advection_ic_wave_type", "HONEE 0.0", NULL)); 448 PetscCall(PetscOptionsDeprecated("-wave_frequency", "-advection_ic_wave_frequency", "HONEE 0.0", NULL)); 449 PetscCall(PetscOptionsDeprecated("-wave_phase", "-advection_ic_wave_phase", "HONEE 0.0", NULL)); 450 PetscCall(PetscOptionsEnum("-advection_ic_wave_type", "Type of wave", NULL, AdvDifWaveTypes, (PetscEnum)(wave_type = ADVDIF_WAVE_SINE), 451 (PetscEnum *)&wave_type, NULL)); 452 PetscCall(PetscOptionsScalar("-advection_ic_wave_frequency", "Frequency of sine wave", NULL, wave_frequency, &wave_frequency, NULL)); 453 PetscCall(PetscOptionsScalar("-advection_ic_wave_phase", "Length correction", NULL, wave_phase, &wave_phase, NULL)); 454 break; 455 case ADVDIF_IC_BOUNDARY_LAYER: 456 PetscCall(PetscOptionsScalar("-advection_ic_bl_height_factor", "Height of boundary layer in IC", NULL, bl_height_factor, &bl_height_factor, 457 NULL)); 458 break; 459 case ADVDIF_IC_BUBBLE_CYLINDER: 460 case ADVDIF_IC_BUBBLE_SPHERE: 461 PetscCall(PetscOptionsDeprecated("-rc", "-advection_ic_bubble_rc", "HONEE 0.0", NULL)); 462 PetscCall(PetscOptionsDeprecated("-bubble_continuity", "-advection_ic_bubble_continuity", "HONEE 0.0", NULL)); 463 PetscCall(PetscOptionsScalar("-advection_ic_bubble_rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL)); 464 bubble_continuity_type = dim == 3 ? ADVDIF_BUBBLE_CONTINUITY_SMOOTH : ADVDIF_BUBBLE_CONTINUITY_COSINE; 465 PetscCall(PetscOptionsEnum("-advection_ic_bubble_continuity", "Smooth, back_sharp, or thick", NULL, AdvDifBubbleContinuityTypes, 466 (PetscEnum)bubble_continuity_type, (PetscEnum *)&bubble_continuity_type, NULL)); 467 break; 468 case ADVDIF_IC_SKEW: 469 case ADVDIF_IC_COSINE_HILL: 470 break; 471 } 472 473 // -- Warnings 474 if (wind_type == ADVDIF_WIND_ROTATION && user_wind) { 475 PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n")); 476 } 477 if (wind_type == ADVDIF_WIND_TRANSLATION && advectionic_type == ADVDIF_IC_BUBBLE_CYLINDER && wind[2] != 0.) { 478 wind[2] = 0; 479 PetscCall(PetscPrintf(comm, 480 "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n")); 481 } 482 if (stab == STAB_NONE && CtauS != 0) { 483 PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n")); 484 } 485 PetscOptionsEnd(); 486 487 if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized; 488 489 // ------------------------------------------------------ 490 // Set up the QFunction contexts 491 // ------------------------------------------------------ 492 // -- Scale variables to desired units 493 Units units = honee->units; 494 E_wind *= units->Joule; 495 rc = fabs(rc) * units->meter; 496 for (PetscInt i = 0; i < dim; i++) { 497 wind[i] *= (units->meter / units->second); 498 domain_size[i] *= units->meter; 499 } 500 501 // -- Setup Context 502 setup_context->rc = rc; 503 setup_context->lx = domain_size[0]; 504 setup_context->ly = domain_size[1]; 505 setup_context->lz = dim == 3 ? domain_size[2] : 0.; 506 setup_context->wind[0] = wind[0]; 507 setup_context->wind[1] = wind[1]; 508 setup_context->wind[2] = dim == 3 ? wind[2] : 0.; 509 setup_context->wind_type = wind_type; 510 setup_context->initial_condition_type = advectionic_type; 511 setup_context->bubble_continuity_type = bubble_continuity_type; 512 setup_context->time = 0; 513 setup_context->wave_frequency = wave_frequency; 514 setup_context->wave_phase = wave_phase; 515 setup_context->wave_type = wave_type; 516 setup_context->bl_height_factor = bl_height_factor; 517 518 // -- QFunction Context 519 honee->phys->implicit = implicit; 520 advection_ctx->CtauS = CtauS; 521 advection_ctx->E_wind = E_wind; 522 advection_ctx->implicit = implicit; 523 advection_ctx->strong_form = strong_form; 524 advection_ctx->stabilization = stab; 525 advection_ctx->stabilization_tau = stab_tau; 526 advection_ctx->Ctau_a = Ctau_a; 527 advection_ctx->Ctau_d = Ctau_d; 528 advection_ctx->Ctau_t = Ctau_t; 529 advection_ctx->diffusion_coeff = diffusion_coeff; 530 advection_ctx->divFdiff_method = honee->app_ctx->divFdiffproj_method; 531 532 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &ics_qfctx)); 533 PetscCallCeed(ceed, CeedQFunctionContextSetData(ics_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 534 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(ics_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 535 536 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &advection_qfctx)); 537 PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx)); 538 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 539 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_qfctx, "timestep size", offsetof(struct AdvectionContext_, dt), 1, 540 "Size of timestep, delta t")); 541 542 // ------------------------------------------------------ 543 // SET UP ADVECTION 544 // ------------------------------------------------------ 545 problem->print_info = PRINT_ADVECTION; 546 problem->num_comps_jac_data = 0; 547 switch (dim) { 548 case 2: 549 problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsAdvection2d, .qf_loc = ICsAdvection2d_loc}; 550 problem->apply_vol_rhs = (HoneeQFSpec){.qf_func_ptr = RHS_Advection2d, .qf_loc = RHS_Advection2d_loc}; 551 problem->apply_vol_ifunction = (HoneeQFSpec){.qf_func_ptr = IFunction_Advection2d, .qf_loc = IFunction_Advection2d_loc}; 552 problem->compute_exact_solution_error = PETSC_TRUE; 553 break; 554 case 3: 555 problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsAdvection, .qf_loc = ICsAdvection_loc}; 556 problem->apply_vol_rhs = (HoneeQFSpec){.qf_func_ptr = RHS_Advection, .qf_loc = RHS_Advection_loc}; 557 problem->apply_vol_ifunction = (HoneeQFSpec){.qf_func_ptr = IFunction_Advection, .qf_loc = IFunction_Advection_loc}; 558 problem->compute_exact_solution_error = PETSC_FALSE; 559 break; 560 } 561 problem->ics.qfctx = ics_qfctx; 562 problem->apply_vol_rhs.qfctx = advection_qfctx; 563 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_vol_ifunction.qfctx)); 564 565 problem->num_components = 5; 566 PetscCall(PetscMalloc1(problem->num_components, &problem->component_names)); 567 for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(component_names[i], &problem->component_names[i])); 568 569 PetscCall(DivDiffFluxProjectionCreate(honee, honee->app_ctx->divFdiffproj_method, 1, &honee->diff_flux_proj)); 570 if (honee->diff_flux_proj) { 571 DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj; 572 NodalProjectionData projection = diff_flux_proj->projection; 573 PetscSection section; 574 575 diff_flux_proj->CreateRHSOperator_Direct = DivDiffFluxProjectionCreateRHS_Direct_AdvDif; 576 diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_AdvDif; 577 PetscCall(DMGetLocalSection(projection->dm, §ion)); 578 switch (honee->diff_flux_proj->method) { 579 case DIV_DIFF_FLUX_PROJ_DIRECT: { 580 PetscCall(PetscSectionSetFieldName(section, 0, "")); 581 PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_Scalar")); 582 } break; 583 case DIV_DIFF_FLUX_PROJ_INDIRECT: { 584 PetscCall(PetscSectionSetFieldName(section, 0, "")); 585 PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_ScalarX")); 586 PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_ScalarY")); 587 if (dim >= 3) PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_ScalarZ")); 588 } break; 589 case DIV_DIFF_FLUX_PROJ_NONE: 590 SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 591 DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]); 592 break; 593 } 594 } 595 596 for (PetscCount b = 0; b < problem->num_bc_defs; b++) { 597 BCDefinition bc_def = problem->bc_defs[b]; 598 const char *name; 599 600 PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL)); 601 if (!strcmp(name, "inflow")) { 602 HoneeBCStruct honee_bc; 603 604 PetscCall(PetscNew(&honee_bc)); 605 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &honee_bc->qfctx)); 606 honee_bc->honee = honee; 607 honee_bc->num_comps_jac_data = 0; 608 PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc)); 609 610 PetscCall(BCDefinitionSetIFunction(bc_def, AdvectionInflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp)); 611 PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL)); 612 } 613 } 614 PetscFunctionReturn(PETSC_SUCCESS); 615 } 616