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