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 qf_ctx = 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], &qf_ctx)); 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, qf_ctx)); 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(&qf_ctx)); 74 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 75 PetscFunctionReturn(PETSC_SUCCESS); 76 } 77 78 PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 79 WindType wind_type; 80 AdvectionICType advectionic_type; 81 BubbleContinuityType bubble_continuity_type; 82 StabilizationType stab; 83 StabilizationTauType stab_tau; 84 SetupContextAdv setup_context; 85 User user = *(User *)ctx; 86 MPI_Comm comm = user->comm; 87 Ceed ceed = user->ceed; 88 PetscBool implicit; 89 AdvectionContext advection_ctx; 90 CeedQFunctionContext advection_qfctx; 91 PetscInt dim; 92 93 PetscFunctionBeginUser; 94 PetscCall(PetscCalloc1(1, &setup_context)); 95 PetscCall(PetscCalloc1(1, &advection_ctx)); 96 PetscCall(DMGetDimension(dm, &dim)); 97 98 // ------------------------------------------------------ 99 // SET UP ADVECTION 100 // ------------------------------------------------------ 101 switch (dim) { 102 case 2: 103 problem->ics.qf_func_ptr = ICsAdvection2d; 104 problem->ics.qf_loc = ICsAdvection2d_loc; 105 problem->apply_vol_rhs.qf_func_ptr = RHS_Advection2d; 106 problem->apply_vol_rhs.qf_loc = RHS_Advection2d_loc; 107 problem->apply_vol_ifunction.qf_func_ptr = IFunction_Advection2d; 108 problem->apply_vol_ifunction.qf_loc = IFunction_Advection2d_loc; 109 problem->apply_inflow.qf_func_ptr = Advection2d_InOutFlow; 110 problem->apply_inflow.qf_loc = Advection2d_InOutFlow_loc; 111 problem->compute_exact_solution_error = PETSC_TRUE; 112 problem->print_info = PRINT_ADVECTION; 113 break; 114 case 3: 115 problem->ics.qf_func_ptr = ICsAdvection; 116 problem->ics.qf_loc = ICsAdvection_loc; 117 problem->apply_vol_rhs.qf_func_ptr = RHS_Advection; 118 problem->apply_vol_rhs.qf_loc = RHS_Advection_loc; 119 problem->apply_vol_ifunction.qf_func_ptr = IFunction_Advection; 120 problem->apply_vol_ifunction.qf_loc = IFunction_Advection_loc; 121 problem->apply_inflow.qf_func_ptr = Advection_InOutFlow; 122 problem->apply_inflow.qf_loc = Advection_InOutFlow_loc; 123 problem->compute_exact_solution_error = PETSC_FALSE; 124 problem->print_info = PRINT_ADVECTION; 125 break; 126 } 127 128 // ------------------------------------------------------ 129 // Create the libCEED context 130 // ------------------------------------------------------ 131 CeedScalar rc = 1000.; // m (Radius of bubble) 132 CeedScalar CtauS = 0.; // dimensionless 133 PetscBool strong_form = PETSC_FALSE; 134 CeedScalar E_wind = 1.e6; // J 135 CeedScalar Ctau_a = PetscPowScalarInt(user->app_ctx->degree, 2); 136 CeedScalar Ctau_t = 0.; 137 PetscReal wind[3] = {1., 0, 0}; // m/s 138 CeedScalar diffusion_coeff = 0.; 139 PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.}; 140 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 141 for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 142 143 // ------------------------------------------------------ 144 // Create the PETSc context 145 // ------------------------------------------------------ 146 PetscScalar meter = 1e-2; // 1 meter in scaled length units 147 PetscScalar kilogram = 1e-6; // 1 kilogram in scaled mass units 148 PetscScalar second = 1e-2; // 1 second in scaled time units 149 PetscScalar Joule; 150 151 // ------------------------------------------------------ 152 // Command line Options 153 // ------------------------------------------------------ 154 PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL); 155 // -- Physics 156 PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL)); 157 PetscBool translation; 158 PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, WindTypes, (PetscEnum)(wind_type = WIND_ROTATION), (PetscEnum *)&wind_type, 159 &translation)); 160 PetscInt n = dim; 161 PetscBool user_wind; 162 PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind)); 163 PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL)); 164 PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL)); 165 PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form, 166 NULL)); 167 PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL)); 168 PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvectionICTypes, 169 (PetscEnum)(advectionic_type = ADVECTIONIC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL)); 170 bubble_continuity_type = dim == 3 ? BUBBLE_CONTINUITY_SMOOTH : BUBBLE_CONTINUITY_COSINE; 171 PetscCall(PetscOptionsEnum("-bubble_continuity", "Smooth, back_sharp, or thick", NULL, BubbleContinuityTypes, (PetscEnum)bubble_continuity_type, 172 (PetscEnum *)&bubble_continuity_type, NULL)); 173 PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 174 PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU), 175 (PetscEnum *)&stab_tau, NULL)); 176 PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 177 PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization ", NULL, Ctau_a, &Ctau_a, NULL)); 178 PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 179 180 // -- Units 181 PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 182 meter = fabs(meter); 183 PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 184 kilogram = fabs(kilogram); 185 PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 186 second = fabs(second); 187 188 // -- Warnings 189 if (wind_type == WIND_ROTATION && user_wind) { 190 PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n")); 191 } 192 if (wind_type == WIND_TRANSLATION && advectionic_type == ADVECTIONIC_BUBBLE_CYLINDER && wind[2] != 0.) { 193 wind[2] = 0; 194 PetscCall( 195 PetscPrintf(comm, "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n")); 196 } 197 if (stab == STAB_NONE && CtauS != 0) { 198 PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n")); 199 } 200 PetscOptionsEnd(); 201 202 if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized; 203 204 // ------------------------------------------------------ 205 // Set up the PETSc context 206 // ------------------------------------------------------ 207 // -- Define derived units 208 Joule = kilogram * PetscSqr(meter) / PetscSqr(second); 209 210 user->units->meter = meter; 211 user->units->kilogram = kilogram; 212 user->units->second = second; 213 user->units->Joule = Joule; 214 215 // ------------------------------------------------------ 216 // Set up the libCEED context 217 // ------------------------------------------------------ 218 // -- Scale variables to desired units 219 E_wind *= Joule; 220 rc = fabs(rc) * meter; 221 for (PetscInt i = 0; i < dim; i++) { 222 wind[i] *= (meter / second); 223 domain_size[i] *= meter; 224 } 225 226 // -- Setup Context 227 setup_context->rc = rc; 228 setup_context->lx = domain_size[0]; 229 setup_context->ly = domain_size[1]; 230 setup_context->lz = dim == 3 ? domain_size[2] : 0.; 231 setup_context->wind[0] = wind[0]; 232 setup_context->wind[1] = wind[1]; 233 setup_context->wind[2] = dim == 3 ? wind[2] : 0.; 234 setup_context->wind_type = wind_type; 235 setup_context->initial_condition_type = advectionic_type; 236 setup_context->bubble_continuity_type = bubble_continuity_type; 237 setup_context->time = 0; 238 239 // -- QFunction Context 240 user->phys->implicit = implicit; 241 advection_ctx->CtauS = CtauS; 242 advection_ctx->E_wind = E_wind; 243 advection_ctx->implicit = implicit; 244 advection_ctx->strong_form = strong_form; 245 advection_ctx->stabilization = stab; 246 advection_ctx->stabilization_tau = stab_tau; 247 advection_ctx->Ctau_a = Ctau_a; 248 advection_ctx->Ctau_t = Ctau_t; 249 advection_ctx->diffusion_coeff = diffusion_coeff; 250 251 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfctx)); 252 PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 253 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc)); 254 255 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &advection_qfctx)); 256 PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx)); 257 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 258 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_qfctx, "timestep size", offsetof(struct AdvectionContext_, dt), 1, 259 "Size of timestep, delta t")); 260 problem->apply_vol_rhs.qfctx = advection_qfctx; 261 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_vol_ifunction.qfctx)); 262 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_inflow.qfctx)); 263 PetscFunctionReturn(PETSC_SUCCESS); 264 } 265 266 PetscErrorCode PRINT_ADVECTION(User user, ProblemData problem, AppCtx app_ctx) { 267 MPI_Comm comm = user->comm; 268 Ceed ceed = user->ceed; 269 SetupContextAdv setup_ctx; 270 AdvectionContext advection_ctx; 271 PetscInt dim; 272 273 PetscFunctionBeginUser; 274 PetscCall(DMGetDimension(user->dm, &dim)); 275 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &setup_ctx)); 276 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &advection_ctx)); 277 PetscCall(PetscPrintf(comm, 278 " Problem:\n" 279 " Problem Name : %s\n" 280 " Stabilization : %s\n" 281 " Initial Condition Type : %s\n" 282 " Bubble Continuity : %s\n" 283 " Wind Type : %s\n", 284 app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization], AdvectionICTypes[setup_ctx->initial_condition_type], 285 BubbleContinuityTypes[setup_ctx->bubble_continuity_type], WindTypes[setup_ctx->wind_type])); 286 287 if (setup_ctx->wind_type == WIND_TRANSLATION) { 288 switch (dim) { 289 case 2: 290 PetscCall(PetscPrintf(comm, " Background Wind : %f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1])); 291 break; 292 case 3: 293 PetscCall( 294 PetscPrintf(comm, " Background Wind : %f,%f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1], setup_ctx->wind[2])); 295 break; 296 } 297 } 298 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &setup_ctx)); 299 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &advection_ctx)); 300 PetscFunctionReturn(PETSC_SUCCESS); 301 } 302