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(User user, CeedOperator *op_mass) { 18 Ceed ceed = user->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(user->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(user->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, user->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] user `User` object 82 @param[in] ceed_data `CeedData` object 83 @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 84 @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 85 **/ 86 static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_AdvDif(User user, CeedData ceed_data, DivDiffFluxProjectionData diff_flux_proj, 87 CeedOperator *op_rhs) { 88 Ceed ceed = user->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(ceed_data->basis_q, &num_comp_q)); 99 100 { // Get newtonian QF context 101 CeedOperator *sub_ops; 102 PetscInt sub_op_index = 0; // will be 0 for the volume op 103 104 if (user->op_ifunction) PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_ifunction, &sub_ops)); 105 else PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->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, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, 119 &elem_restr_qd, &q_data, &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, user->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", ceed_data->elem_restr_q, ceed_data->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(user->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, user->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q)); 199 PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, user->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(QDataBoundaryGradientGet(ceed, user->dm, face_orientation_label, orientation, ceed_data->x_coord, &elem_restr_qdata, &q_data, 204 &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 PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 233 AdvDifWindType wind_type; 234 AdvDifICType advectionic_type; 235 AdvDifBubbleContinuityType bubble_continuity_type = -1; 236 StabilizationType stab; 237 StabilizationTauType stab_tau; 238 SetupContextAdv setup_context; 239 User user = *(User *)ctx; 240 MPI_Comm comm = user->comm; 241 Ceed ceed = user->ceed; 242 PetscBool implicit; 243 AdvectionContext advection_ctx; 244 CeedQFunctionContext advection_qfctx; 245 PetscInt dim; 246 247 PetscFunctionBeginUser; 248 PetscCall(PetscCalloc1(1, &setup_context)); 249 PetscCall(PetscCalloc1(1, &advection_ctx)); 250 PetscCall(DMGetDimension(dm, &dim)); 251 252 // ------------------------------------------------------ 253 // SET UP ADVECTION 254 // ------------------------------------------------------ 255 problem->print_info = PRINT_ADVECTION; 256 problem->jac_data_size_vol = 0; 257 problem->jac_data_size_sur = 0; 258 switch (dim) { 259 case 2: 260 problem->ics.qf_func_ptr = ICsAdvection2d; 261 problem->ics.qf_loc = ICsAdvection2d_loc; 262 problem->apply_vol_rhs.qf_func_ptr = RHS_Advection2d; 263 problem->apply_vol_rhs.qf_loc = RHS_Advection2d_loc; 264 problem->apply_vol_ifunction.qf_func_ptr = IFunction_Advection2d; 265 problem->apply_vol_ifunction.qf_loc = IFunction_Advection2d_loc; 266 problem->apply_inflow.qf_func_ptr = Advection2d_InOutFlow; 267 problem->apply_inflow.qf_loc = Advection2d_InOutFlow_loc; 268 problem->compute_exact_solution_error = PETSC_TRUE; 269 break; 270 case 3: 271 problem->ics.qf_func_ptr = ICsAdvection; 272 problem->ics.qf_loc = ICsAdvection_loc; 273 problem->apply_vol_rhs.qf_func_ptr = RHS_Advection; 274 problem->apply_vol_rhs.qf_loc = RHS_Advection_loc; 275 problem->apply_vol_ifunction.qf_func_ptr = IFunction_Advection; 276 problem->apply_vol_ifunction.qf_loc = IFunction_Advection_loc; 277 problem->apply_inflow.qf_func_ptr = Advection_InOutFlow; 278 problem->apply_inflow.qf_loc = Advection_InOutFlow_loc; 279 problem->compute_exact_solution_error = PETSC_FALSE; 280 break; 281 } 282 283 PetscCall(DivDiffFluxProjectionCreate(user, 1, &user->diff_flux_proj)); 284 if (user->diff_flux_proj) { 285 DivDiffFluxProjectionData diff_flux_proj = user->diff_flux_proj; 286 NodalProjectionData projection = diff_flux_proj->projection; 287 288 diff_flux_proj->CreateRHSOperator_Direct = DivDiffFluxProjectionCreateRHS_Direct_AdvDif; 289 diff_flux_proj->CreateRHSOperator_Indirect = NULL; 290 291 switch (user->diff_flux_proj->method) { 292 case DIV_DIFF_FLUX_PROJ_DIRECT: { 293 PetscSection section; 294 295 PetscCall(DMGetLocalSection(projection->dm, §ion)); 296 PetscCall(PetscSectionSetFieldName(section, 0, "")); 297 PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_Scalar")); 298 } break; 299 case DIV_DIFF_FLUX_PROJ_INDIRECT: { 300 PetscSection section; 301 302 PetscCall(DMGetLocalSection(projection->dm, §ion)); 303 PetscCall(PetscSectionSetFieldName(section, 0, "")); 304 PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_ScalarX")); 305 PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_ScalarY")); 306 PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_ScalarZ")); 307 } break; 308 case DIV_DIFF_FLUX_PROJ_NONE: 309 SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 310 DivDiffFluxProjectionMethods[user->app_ctx->divFdiffproj_method]); 311 break; 312 } 313 } 314 315 // ------------------------------------------------------ 316 // Create the libCEED context 317 // ------------------------------------------------------ 318 CeedScalar rc = 1000.; // m (Radius of bubble) 319 CeedScalar CtauS = 0.; // dimensionless 320 PetscBool strong_form = PETSC_FALSE; 321 CeedScalar E_wind = 1.e6; // J 322 CeedScalar Ctau_a = PetscPowScalarInt(user->app_ctx->degree, 2); 323 CeedScalar Ctau_d = PetscPowScalarInt(user->app_ctx->degree, 4); 324 CeedScalar Ctau_t = 0.; 325 PetscReal wind[3] = {1., 0, 0}; // m/s 326 CeedScalar diffusion_coeff = 0.; 327 CeedScalar wave_frequency = 2 * M_PI; 328 CeedScalar wave_phase = 0; 329 AdvDifWaveType wave_type = -1; 330 PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.}; 331 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 332 for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 333 334 // ------------------------------------------------------ 335 // Create the PETSc context 336 // ------------------------------------------------------ 337 PetscScalar meter = 1e-2; // 1 meter in scaled length units 338 PetscScalar kilogram = 1e-6; // 1 kilogram in scaled mass units 339 PetscScalar second = 1e-2; // 1 second in scaled time units 340 PetscScalar Joule; 341 342 // ------------------------------------------------------ 343 // Command line Options 344 // ------------------------------------------------------ 345 PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL); 346 // -- Physics 347 PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL)); 348 PetscBool translation; 349 PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, AdvDifWindTypes, (PetscEnum)(wind_type = ADVDIF_WIND_ROTATION), 350 (PetscEnum *)&wind_type, &translation)); 351 PetscInt n = dim; 352 PetscBool user_wind; 353 PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind)); 354 PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL)); 355 PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL)); 356 PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form, 357 NULL)); 358 PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL)); 359 PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvDifICTypes, 360 (PetscEnum)(advectionic_type = ADVDIF_IC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL)); 361 PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 362 PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU), 363 (PetscEnum *)&stab_tau, NULL)); 364 PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 365 PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization, advection component", NULL, Ctau_a, &Ctau_a, NULL)); 366 PetscCall(PetscOptionsScalar("-Ctau_d", "Coefficient for the stabilization, diffusion component", NULL, Ctau_d, &Ctau_d, NULL)); 367 PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 368 369 if (advectionic_type == ADVDIF_IC_WAVE) { 370 PetscCall(PetscOptionsEnum("-wave_type", "Type of wave", NULL, AdvDifWaveTypes, (PetscEnum)(wave_type = ADVDIF_WAVE_SINE), 371 (PetscEnum *)&wave_type, NULL)); 372 PetscCall(PetscOptionsScalar("-wave_frequency", "Frequency of sine wave", NULL, wave_frequency, &wave_frequency, NULL)); 373 PetscCall(PetscOptionsScalar("-wave_phase", "Length correction", NULL, wave_phase, &wave_phase, NULL)); 374 } 375 376 if (advectionic_type == ADVDIF_IC_BUBBLE_CYLINDER || advectionic_type == ADVDIF_IC_BUBBLE_SPHERE) { 377 bubble_continuity_type = dim == 3 ? ADVDIF_BUBBLE_CONTINUITY_SMOOTH : ADVDIF_BUBBLE_CONTINUITY_COSINE; 378 PetscCall(PetscOptionsEnum("-bubble_continuity", "Smooth, back_sharp, or thick", NULL, AdvDifBubbleContinuityTypes, 379 (PetscEnum)bubble_continuity_type, (PetscEnum *)&bubble_continuity_type, NULL)); 380 } 381 382 // -- Units 383 PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 384 meter = fabs(meter); 385 PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 386 kilogram = fabs(kilogram); 387 PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 388 second = fabs(second); 389 390 // -- Warnings 391 if (wind_type == ADVDIF_WIND_ROTATION && user_wind) { 392 PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n")); 393 } 394 if (wind_type == ADVDIF_WIND_TRANSLATION && advectionic_type == ADVDIF_IC_BUBBLE_CYLINDER && wind[2] != 0.) { 395 wind[2] = 0; 396 PetscCall( 397 PetscPrintf(comm, "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n")); 398 } 399 if (stab == STAB_NONE && CtauS != 0) { 400 PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n")); 401 } 402 PetscOptionsEnd(); 403 404 if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized; 405 406 // ------------------------------------------------------ 407 // Set up the PETSc context 408 // ------------------------------------------------------ 409 // -- Define derived units 410 Joule = kilogram * PetscSqr(meter) / PetscSqr(second); 411 412 user->units->meter = meter; 413 user->units->kilogram = kilogram; 414 user->units->second = second; 415 user->units->Joule = Joule; 416 417 // ------------------------------------------------------ 418 // Set up the libCEED context 419 // ------------------------------------------------------ 420 // -- Scale variables to desired units 421 E_wind *= Joule; 422 rc = fabs(rc) * meter; 423 for (PetscInt i = 0; i < dim; i++) { 424 wind[i] *= (meter / second); 425 domain_size[i] *= meter; 426 } 427 428 // -- Setup Context 429 setup_context->rc = rc; 430 setup_context->lx = domain_size[0]; 431 setup_context->ly = domain_size[1]; 432 setup_context->lz = dim == 3 ? domain_size[2] : 0.; 433 setup_context->wind[0] = wind[0]; 434 setup_context->wind[1] = wind[1]; 435 setup_context->wind[2] = dim == 3 ? wind[2] : 0.; 436 setup_context->wind_type = wind_type; 437 setup_context->initial_condition_type = advectionic_type; 438 setup_context->bubble_continuity_type = bubble_continuity_type; 439 setup_context->time = 0; 440 setup_context->wave_frequency = wave_frequency; 441 setup_context->wave_phase = wave_phase; 442 setup_context->wave_type = wave_type; 443 444 // -- QFunction Context 445 user->phys->implicit = implicit; 446 advection_ctx->CtauS = CtauS; 447 advection_ctx->E_wind = E_wind; 448 advection_ctx->implicit = implicit; 449 advection_ctx->strong_form = strong_form; 450 advection_ctx->stabilization = stab; 451 advection_ctx->stabilization_tau = stab_tau; 452 advection_ctx->Ctau_a = Ctau_a; 453 advection_ctx->Ctau_d = Ctau_d; 454 advection_ctx->Ctau_t = Ctau_t; 455 advection_ctx->diffusion_coeff = diffusion_coeff; 456 advection_ctx->divFdiff_method = user->app_ctx->divFdiffproj_method; 457 458 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfctx)); 459 PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 460 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc)); 461 462 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &advection_qfctx)); 463 PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx)); 464 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 465 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_qfctx, "timestep size", offsetof(struct AdvectionContext_, dt), 1, 466 "Size of timestep, delta t")); 467 problem->apply_vol_rhs.qfctx = advection_qfctx; 468 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_vol_ifunction.qfctx)); 469 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_inflow.qfctx)); 470 PetscFunctionReturn(PETSC_SUCCESS); 471 } 472 473 PetscErrorCode PRINT_ADVECTION(User user, ProblemData problem, AppCtx app_ctx) { 474 MPI_Comm comm = user->comm; 475 Ceed ceed = user->ceed; 476 SetupContextAdv setup_ctx; 477 AdvectionContext advection_ctx; 478 PetscInt dim; 479 480 PetscFunctionBeginUser; 481 PetscCall(DMGetDimension(user->dm, &dim)); 482 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &setup_ctx)); 483 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &advection_ctx)); 484 PetscCall(PetscPrintf(comm, 485 " Problem:\n" 486 " Problem Name : %s\n" 487 " Stabilization : %s\n" 488 " Stabilization Tau : %s\n" 489 " Wind Type : %s\n", 490 app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization], 491 StabilizationTauTypes[advection_ctx->stabilization_tau], AdvDifWindTypes[setup_ctx->wind_type])); 492 493 if (setup_ctx->wind_type == ADVDIF_WIND_TRANSLATION) { 494 CeedScalar *wind = setup_ctx->wind; 495 switch (dim) { 496 case 2: 497 PetscCall(PetscPrintf(comm, " Background Wind : %f,%f\n", wind[0], wind[1])); 498 break; 499 case 3: 500 PetscCall(PetscPrintf(comm, " Background Wind : %f,%f,%f\n", wind[0], wind[1], wind[2])); 501 break; 502 } 503 } 504 505 PetscCall(PetscPrintf(comm, " Initial Condition Type : %s\n", AdvDifICTypes[setup_ctx->initial_condition_type])); 506 switch (setup_ctx->initial_condition_type) { 507 case ADVDIF_IC_SKEW: 508 case ADVDIF_IC_COSINE_HILL: 509 break; 510 case ADVDIF_IC_BUBBLE_SPHERE: 511 case ADVDIF_IC_BUBBLE_CYLINDER: 512 PetscCall(PetscPrintf(comm, " Bubble Continuity : %s\n", AdvDifBubbleContinuityTypes[setup_ctx->bubble_continuity_type])); 513 break; 514 case ADVDIF_IC_WAVE: 515 PetscCall(PetscPrintf(comm, " Wave Type : %s\n", AdvDifWaveTypes[setup_ctx->wave_type])); 516 break; 517 } 518 519 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &setup_ctx)); 520 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &advection_ctx)); 521 PetscFunctionReturn(PETSC_SUCCESS); 522 } 523