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