1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 /// @file 9 /// Utility functions for setting up ADVECTION 10 11 #include "../qfunctions/advection.h" 12 13 #include <ceed.h> 14 #include <petscdm.h> 15 16 #include "../navierstokes.h" 17 #include "../qfunctions/setupgeo.h" 18 #include "../qfunctions/setupgeo2d.h" 19 20 // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping 21 // 22 // Only used for SUPG stabilization 23 PetscErrorCode CreateKSPMassOperator_AdvectionStabilized(User user, CeedOperator *op_mass) { 24 Ceed ceed = user->ceed; 25 CeedInt num_comp_q, q_data_size; 26 CeedQFunction qf_mass; 27 CeedElemRestriction elem_restr_q, elem_restr_qd_i; 28 CeedBasis basis_q; 29 CeedVector q_data; 30 CeedQFunctionContext qf_ctx = NULL; 31 PetscInt dim; 32 33 PetscFunctionBeginUser; 34 PetscCall(DMGetDimension(user->dm, &dim)); 35 { // Get restriction and basis from the RHS function 36 CeedOperator *sub_ops; 37 CeedOperatorField field; 38 PetscInt sub_op_index = 0; // will be 0 for the volume op 39 40 PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops)); 41 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field)); 42 PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q)); 43 PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q)); 44 45 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field)); 46 PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i)); 47 PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data)); 48 49 PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx)); 50 } 51 52 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 53 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size)); 54 55 switch (dim) { 56 case 2: 57 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection2D, MassFunction_Advection2D_loc, &qf_mass)); 58 break; 59 case 3: 60 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection, MassFunction_Advection_loc, &qf_mass)); 61 break; 62 } 63 64 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx)); 65 PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0)); 66 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP)); 67 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP)); 68 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE)); 69 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP)); 70 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD)); 71 72 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 73 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 74 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed)); 75 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data)); 76 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 77 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 78 79 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 80 PetscFunctionReturn(PETSC_SUCCESS); 81 } 82 83 PetscErrorCode NS_ADVECTION(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) { 84 WindType wind_type; 85 AdvectionICType advectionic_type; 86 BubbleContinuityType bubble_continuity_type; 87 StabilizationType stab; 88 StabilizationTauType stab_tau; 89 SetupContextAdv setup_context; 90 User user = *(User *)ctx; 91 MPI_Comm comm = user->comm; 92 Ceed ceed = user->ceed; 93 PetscBool implicit; 94 AdvectionContext advection_ctx; 95 CeedQFunctionContext advection_context; 96 PetscInt dim; 97 98 PetscFunctionBeginUser; 99 PetscCall(PetscCalloc1(1, &setup_context)); 100 PetscCall(PetscCalloc1(1, &advection_ctx)); 101 PetscCall(DMGetDimension(dm, &dim)); 102 103 // ------------------------------------------------------ 104 // SET UP ADVECTION 105 // ------------------------------------------------------ 106 switch (dim) { 107 case 2: 108 problem->dim = 2; 109 problem->q_data_size_vol = 5; 110 problem->q_data_size_sur = 3; 111 problem->setup_vol.qfunction = Setup2d; 112 problem->setup_vol.qfunction_loc = Setup2d_loc; 113 problem->setup_sur.qfunction = SetupBoundary2d; 114 problem->setup_sur.qfunction_loc = SetupBoundary2d_loc; 115 problem->ics.qfunction = ICsAdvection2d; 116 problem->ics.qfunction_loc = ICsAdvection2d_loc; 117 problem->apply_vol_rhs.qfunction = RHS_Advection2d; 118 problem->apply_vol_rhs.qfunction_loc = RHS_Advection2d_loc; 119 problem->apply_vol_ifunction.qfunction = IFunction_Advection2d; 120 problem->apply_vol_ifunction.qfunction_loc = IFunction_Advection2d_loc; 121 problem->apply_inflow.qfunction = Advection2d_InOutFlow; 122 problem->apply_inflow.qfunction_loc = Advection2d_InOutFlow_loc; 123 problem->non_zero_time = PETSC_TRUE; 124 problem->print_info = PRINT_ADVECTION; 125 break; 126 case 3: 127 problem->dim = 3; 128 problem->q_data_size_vol = 10; 129 problem->q_data_size_sur = 10; 130 problem->setup_vol.qfunction = Setup; 131 problem->setup_vol.qfunction_loc = Setup_loc; 132 problem->setup_sur.qfunction = SetupBoundary; 133 problem->setup_sur.qfunction_loc = SetupBoundary_loc; 134 problem->ics.qfunction = ICsAdvection; 135 problem->ics.qfunction_loc = ICsAdvection_loc; 136 problem->apply_vol_rhs.qfunction = RHS_Advection; 137 problem->apply_vol_rhs.qfunction_loc = RHS_Advection_loc; 138 problem->apply_vol_ifunction.qfunction = IFunction_Advection; 139 problem->apply_vol_ifunction.qfunction_loc = IFunction_Advection_loc; 140 problem->apply_inflow.qfunction = Advection_InOutFlow; 141 problem->apply_inflow.qfunction_loc = Advection_InOutFlow_loc; 142 problem->non_zero_time = PETSC_FALSE; 143 problem->print_info = PRINT_ADVECTION; 144 break; 145 } 146 147 // ------------------------------------------------------ 148 // Create the libCEED context 149 // ------------------------------------------------------ 150 CeedScalar rc = 1000.; // m (Radius of bubble) 151 CeedScalar CtauS = 0.; // dimensionless 152 PetscBool strong_form = PETSC_FALSE; 153 CeedScalar E_wind = 1.e6; // J 154 CeedScalar Ctau_a = PetscPowScalarInt(user->app_ctx->degree, 2); 155 CeedScalar Ctau_t = 0.; 156 PetscReal wind[3] = {1., 0, 0}; // m/s 157 PetscReal domain_min[3], domain_max[3], domain_size[3]; 158 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 159 for (PetscInt i = 0; i < problem->dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 160 161 // ------------------------------------------------------ 162 // Create the PETSc context 163 // ------------------------------------------------------ 164 PetscScalar meter = 1e-2; // 1 meter in scaled length units 165 PetscScalar kilogram = 1e-6; // 1 kilogram in scaled mass units 166 PetscScalar second = 1e-2; // 1 second in scaled time units 167 PetscScalar Joule; 168 169 // ------------------------------------------------------ 170 // Command line Options 171 // ------------------------------------------------------ 172 PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL); 173 // -- Physics 174 PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL)); 175 PetscBool translation; 176 PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, WindTypes, (PetscEnum)(wind_type = WIND_ROTATION), (PetscEnum *)&wind_type, 177 &translation)); 178 PetscInt n = problem->dim; 179 PetscBool user_wind; 180 PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind)); 181 PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL)); 182 PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form, 183 NULL)); 184 PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL)); 185 PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvectionICTypes, 186 (PetscEnum)(advectionic_type = ADVECTIONIC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL)); 187 bubble_continuity_type = problem->dim == 3 ? BUBBLE_CONTINUITY_SMOOTH : BUBBLE_CONTINUITY_COSINE; 188 PetscCall(PetscOptionsEnum("-bubble_continuity", "Smooth, back_sharp, or thick", NULL, BubbleContinuityTypes, (PetscEnum)bubble_continuity_type, 189 (PetscEnum *)&bubble_continuity_type, NULL)); 190 PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 191 PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU), 192 (PetscEnum *)&stab_tau, NULL)); 193 PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 194 PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization ", NULL, Ctau_a, &Ctau_a, NULL)); 195 PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 196 197 // -- Units 198 PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 199 meter = fabs(meter); 200 PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 201 kilogram = fabs(kilogram); 202 PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 203 second = fabs(second); 204 205 // -- Warnings 206 if (wind_type == WIND_ROTATION && user_wind) { 207 PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n")); 208 } 209 if (wind_type == WIND_TRANSLATION && advectionic_type == ADVECTIONIC_BUBBLE_CYLINDER && wind[2] != 0.) { 210 wind[2] = 0; 211 PetscCall( 212 PetscPrintf(comm, "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n")); 213 } 214 if (stab == STAB_NONE && CtauS != 0) { 215 PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n")); 216 } 217 PetscOptionsEnd(); 218 219 if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized; 220 221 // ------------------------------------------------------ 222 // Set up the PETSc context 223 // ------------------------------------------------------ 224 // -- Define derived units 225 Joule = kilogram * PetscSqr(meter) / PetscSqr(second); 226 227 user->units->meter = meter; 228 user->units->kilogram = kilogram; 229 user->units->second = second; 230 user->units->Joule = Joule; 231 232 // ------------------------------------------------------ 233 // Set up the libCEED context 234 // ------------------------------------------------------ 235 // -- Scale variables to desired units 236 E_wind *= Joule; 237 rc = fabs(rc) * meter; 238 for (PetscInt i = 0; i < problem->dim; i++) { 239 wind[i] *= (meter / second); 240 domain_size[i] *= meter; 241 } 242 problem->dm_scale = meter; 243 244 // -- Setup Context 245 setup_context->rc = rc; 246 setup_context->lx = domain_size[0]; 247 setup_context->ly = domain_size[1]; 248 setup_context->lz = problem->dim == 3 ? domain_size[2] : 0.; 249 setup_context->wind[0] = wind[0]; 250 setup_context->wind[1] = wind[1]; 251 setup_context->wind[2] = problem->dim == 3 ? wind[2] : 0.; 252 setup_context->wind_type = wind_type; 253 setup_context->initial_condition_type = advectionic_type; 254 setup_context->bubble_continuity_type = bubble_continuity_type; 255 setup_context->time = 0; 256 257 // -- QFunction Context 258 user->phys->implicit = implicit; 259 advection_ctx->CtauS = CtauS; 260 advection_ctx->E_wind = E_wind; 261 advection_ctx->implicit = implicit; 262 advection_ctx->strong_form = strong_form; 263 advection_ctx->stabilization = stab; 264 advection_ctx->stabilization_tau = stab_tau; 265 advection_ctx->Ctau_a = Ctau_a; 266 advection_ctx->Ctau_t = Ctau_t; 267 268 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context)); 269 PetscCallCeed(ceed, 270 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 271 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc)); 272 273 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &advection_context)); 274 PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx)); 275 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_context, CEED_MEM_HOST, FreeContextPetsc)); 276 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_context, "timestep size", offsetof(struct AdvectionContext_, dt), 1, 277 "Size of timestep, delta t")); 278 problem->apply_vol_rhs.qfunction_context = advection_context; 279 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_vol_ifunction.qfunction_context)); 280 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_inflow.qfunction_context)); 281 PetscFunctionReturn(PETSC_SUCCESS); 282 } 283 284 PetscErrorCode PRINT_ADVECTION(User user, ProblemData *problem, AppCtx app_ctx) { 285 MPI_Comm comm = user->comm; 286 Ceed ceed = user->ceed; 287 SetupContextAdv setup_ctx; 288 AdvectionContext advection_ctx; 289 290 PetscFunctionBeginUser; 291 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfunction_context, CEED_MEM_HOST, &setup_ctx)); 292 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &advection_ctx)); 293 PetscCall(PetscPrintf(comm, 294 " Problem:\n" 295 " Problem Name : %s\n" 296 " Stabilization : %s\n" 297 " Initial Condition Type : %s\n" 298 " Bubble Continuity : %s\n" 299 " Wind Type : %s\n", 300 app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization], AdvectionICTypes[setup_ctx->initial_condition_type], 301 BubbleContinuityTypes[setup_ctx->bubble_continuity_type], WindTypes[setup_ctx->wind_type])); 302 303 if (setup_ctx->wind_type == WIND_TRANSLATION) { 304 switch (problem->dim) { 305 case 2: 306 PetscCall(PetscPrintf(comm, " Background Wind : %f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1])); 307 break; 308 case 3: 309 PetscCall( 310 PetscPrintf(comm, " Background Wind : %f,%f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1], setup_ctx->wind[2])); 311 break; 312 } 313 } 314 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfunction_context, &setup_ctx)); 315 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &advection_ctx)); 316 PetscFunctionReturn(PETSC_SUCCESS); 317 } 318