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