1 // Copyright (c) 2017-2022, 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 problems using the Newtonian Qfunction 10 11 #include "../qfunctions/newtonian.h" 12 13 #include <ceed.h> 14 #include <petscdm.h> 15 16 #include "../navierstokes.h" 17 #include "../qfunctions/setupgeo.h" 18 19 // For use with PetscOptionsEnum 20 static const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "StateVariable", "STATEVAR_", NULL}; 21 22 // Compute relative error |a - b|/|s| 23 static PetscErrorCode CheckPrimitiveWithTolerance(StatePrimitive sY, StatePrimitive aY, StatePrimitive bY, const char *name, PetscReal rtol_pressure, 24 PetscReal rtol_velocity, PetscReal rtol_temperature) { 25 StatePrimitive eY; // relative error 26 27 PetscFunctionBeginUser; 28 eY.pressure = (aY.pressure - bY.pressure) / sY.pressure; 29 PetscScalar u = sqrt(Square(sY.velocity[0]) + Square(sY.velocity[1]) + Square(sY.velocity[2])); 30 for (int j = 0; j < 3; j++) eY.velocity[j] = (aY.velocity[j] - bY.velocity[j]) / u; 31 eY.temperature = (aY.temperature - bY.temperature) / sY.temperature; 32 if (fabs(eY.pressure) > rtol_pressure) printf("%s: pressure error %g\n", name, eY.pressure); 33 for (int j = 0; j < 3; j++) { 34 if (fabs(eY.velocity[j]) > rtol_velocity) printf("%s: velocity[%d] error %g\n", name, j, eY.velocity[j]); 35 } 36 if (fabs(eY.temperature) > rtol_temperature) printf("%s: temperature error %g\n", name, eY.temperature); 37 PetscFunctionReturn(PETSC_SUCCESS); 38 } 39 40 static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) { 41 Units units = user->units; 42 const CeedScalar eps = 1e-6; 43 const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, Pascal = units->Pascal; 44 PetscFunctionBeginUser; 45 const CeedScalar rho = 1.2 * kg / (m * m * m), u = 40 * m / sec; 46 CeedScalar U[5] = {rho, rho * u, rho * u * 1.1, rho * u * 1.2, 250e3 * Pascal + .5 * rho * u * u}; 47 State s = StateFromU(gas, U); 48 for (int i = 0; i < 8; i++) { 49 CeedScalar dU[5] = {0}; 50 if (i < 5) dU[i] = U[i]; 51 State ds = StateFromU_fwd(gas, s, dU); 52 for (int j = 0; j < 5; j++) dU[j] = (1 + eps * (i == j)) * U[j]; 53 State t = StateFromU(gas, dU); 54 StatePrimitive dY; 55 dY.pressure = (t.Y.pressure - s.Y.pressure) / eps; 56 for (int j = 0; j < 3; j++) dY.velocity[j] = (t.Y.velocity[j] - s.Y.velocity[j]) / eps; 57 dY.temperature = (t.Y.temperature - s.Y.temperature) / eps; 58 char buf[128]; 59 snprintf(buf, sizeof buf, "StateFromU_fwd i=%d", i); 60 PetscCall(CheckPrimitiveWithTolerance(dY, ds.Y, dY, buf, 5e-6, 1e-6, 1e-6)); 61 } 62 PetscFunctionReturn(PETSC_SUCCESS); 63 } 64 65 PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) { 66 SetupContext setup_context; 67 User user = *(User *)ctx; 68 CeedInt degree = user->app_ctx->degree; 69 StabilizationType stab; 70 StateVariable state_var; 71 MPI_Comm comm = user->comm; 72 Ceed ceed = user->ceed; 73 PetscBool implicit; 74 PetscBool unit_tests; 75 NewtonianIdealGasContext newtonian_ig_ctx; 76 CeedQFunctionContext newtonian_ig_context; 77 78 PetscFunctionBeginUser; 79 PetscCall(PetscCalloc1(1, &setup_context)); 80 PetscCall(PetscCalloc1(1, &newtonian_ig_ctx)); 81 82 // ------------------------------------------------------ 83 // Setup Generic Newtonian IG Problem 84 // ------------------------------------------------------ 85 problem->dim = 3; 86 problem->q_data_size_vol = 10; 87 problem->q_data_size_sur = 10; 88 problem->jac_data_size_sur = 11; 89 problem->setup_vol.qfunction = Setup; 90 problem->setup_vol.qfunction_loc = Setup_loc; 91 problem->setup_sur.qfunction = SetupBoundary; 92 problem->setup_sur.qfunction_loc = SetupBoundary_loc; 93 problem->non_zero_time = PETSC_FALSE; 94 problem->print_info = PRINT_NEWTONIAN; 95 problem->uses_newtonian = PETSC_TRUE; 96 97 // ------------------------------------------------------ 98 // Create the libCEED context 99 // ------------------------------------------------------ 100 CeedScalar cv = 717.; // J/(kg K) 101 CeedScalar cp = 1004.; // J/(kg K) 102 CeedScalar g[3] = {0, 0, 0}; // m/s^2 103 CeedScalar lambda = -2. / 3.; // - 104 CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 105 CeedScalar k = 0.02638; // W/(m K) 106 CeedScalar c_tau = 0.5 / degree; // - 107 CeedScalar Ctau_t = 1.0; // - 108 CeedScalar Cv_func[3] = {36, 60, 128}; 109 CeedScalar Ctau_v = Cv_func[(CeedInt)Min(3, degree) - 1]; 110 CeedScalar Ctau_C = 0.25 / degree; 111 CeedScalar Ctau_M = 0.25 / degree; 112 CeedScalar Ctau_E = 0.125; 113 PetscReal domain_min[3], domain_max[3], domain_size[3]; 114 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 115 for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 116 117 StatePrimitive reference = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15}; 118 CeedScalar idl_decay_time = -1, idl_start = 0, idl_length = 0; 119 PetscBool idl_enable = PETSC_FALSE; 120 121 // ------------------------------------------------------ 122 // Create the PETSc context 123 // ------------------------------------------------------ 124 PetscScalar meter = 1; // 1 meter in scaled length units 125 PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 126 PetscScalar second = 1; // 1 second in scaled time units 127 PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 128 PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 129 130 // ------------------------------------------------------ 131 // Command line Options 132 // ------------------------------------------------------ 133 PetscBool given_option = PETSC_FALSE; 134 PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL); 135 // -- Conservative vs Primitive variables 136 PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE), 137 (PetscEnum *)&state_var, NULL)); 138 139 switch (state_var) { 140 case STATEVAR_CONSERVATIVE: 141 problem->ics.qfunction = ICsNewtonianIG_Conserv; 142 problem->ics.qfunction_loc = ICsNewtonianIG_Conserv_loc; 143 problem->apply_vol_rhs.qfunction = RHSFunction_Newtonian; 144 problem->apply_vol_rhs.qfunction_loc = RHSFunction_Newtonian_loc; 145 problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Conserv; 146 problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Conserv_loc; 147 problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Conserv; 148 problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Conserv_loc; 149 problem->apply_inflow.qfunction = BoundaryIntegral_Conserv; 150 problem->apply_inflow.qfunction_loc = BoundaryIntegral_Conserv_loc; 151 problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Conserv; 152 problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc; 153 break; 154 155 case STATEVAR_PRIMITIVE: 156 problem->ics.qfunction = ICsNewtonianIG_Prim; 157 problem->ics.qfunction_loc = ICsNewtonianIG_Prim_loc; 158 problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Prim; 159 problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Prim_loc; 160 problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Prim; 161 problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Prim_loc; 162 problem->apply_inflow.qfunction = BoundaryIntegral_Prim; 163 problem->apply_inflow.qfunction_loc = BoundaryIntegral_Prim_loc; 164 problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Prim; 165 problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc; 166 break; 167 } 168 169 // -- Physics 170 PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL)); 171 PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL)); 172 PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL)); 173 PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL)); 174 PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL)); 175 176 PetscInt dim = problem->dim; 177 PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL)); 178 PetscCall(PetscOptionsRealArray("-g", "Gravitational acceleration vector", NULL, g, &dim, &given_option)); 179 dim = problem->dim; 180 PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option)); 181 if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim); 182 183 PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 184 PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 185 PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 186 PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL)); 187 PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL)); 188 PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL)); 189 PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL)); 190 PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 191 PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL)); 192 193 dim = 3; 194 PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL)); 195 PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL)); 196 PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL)); 197 198 // -- Units 199 PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 200 meter = fabs(meter); 201 PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 202 kilogram = fabs(kilogram); 203 PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 204 second = fabs(second); 205 PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL)); 206 Kelvin = fabs(Kelvin); 207 208 // -- Warnings 209 if (stab == STAB_SUPG && !implicit) { 210 PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n")); 211 } 212 PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP, 213 "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n"); 214 215 PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point", 216 NULL, idl_decay_time, &idl_decay_time, &idl_enable)); 217 PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero."); 218 if (idl_decay_time < 0) idl_enable = PETSC_FALSE; 219 PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL)); 220 PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL)); 221 PetscOptionsEnd(); 222 223 // ------------------------------------------------------ 224 // Set up the PETSc context 225 // ------------------------------------------------------ 226 // -- Define derived units 227 Pascal = kilogram / (meter * PetscSqr(second)); 228 J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 229 m_per_squared_s = meter / PetscSqr(second); 230 W_per_m_K = kilogram * meter / (pow(second, 3) * Kelvin); 231 232 user->units->meter = meter; 233 user->units->kilogram = kilogram; 234 user->units->second = second; 235 user->units->Kelvin = Kelvin; 236 user->units->Pascal = Pascal; 237 user->units->J_per_kg_K = J_per_kg_K; 238 user->units->m_per_squared_s = m_per_squared_s; 239 user->units->W_per_m_K = W_per_m_K; 240 241 // ------------------------------------------------------ 242 // Set up the libCEED context 243 // ------------------------------------------------------ 244 // -- Scale variables to desired units 245 cv *= J_per_kg_K; 246 cp *= J_per_kg_K; 247 mu *= Pascal * second; 248 k *= W_per_m_K; 249 for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter; 250 for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s; 251 reference.pressure *= Pascal; 252 for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second; 253 reference.temperature *= Kelvin; 254 problem->dm_scale = meter; 255 256 // -- Solver Settings 257 user->phys->implicit = implicit; 258 user->phys->state_var = state_var; 259 260 // -- QFunction Context 261 newtonian_ig_ctx->lambda = lambda; 262 newtonian_ig_ctx->mu = mu; 263 newtonian_ig_ctx->k = k; 264 newtonian_ig_ctx->cv = cv; 265 newtonian_ig_ctx->cp = cp; 266 newtonian_ig_ctx->c_tau = c_tau; 267 newtonian_ig_ctx->Ctau_t = Ctau_t; 268 newtonian_ig_ctx->Ctau_v = Ctau_v; 269 newtonian_ig_ctx->Ctau_C = Ctau_C; 270 newtonian_ig_ctx->Ctau_M = Ctau_M; 271 newtonian_ig_ctx->Ctau_E = Ctau_E; 272 newtonian_ig_ctx->P0 = reference.pressure; 273 newtonian_ig_ctx->stabilization = stab; 274 newtonian_ig_ctx->P0 = reference.pressure; 275 newtonian_ig_ctx->is_implicit = implicit; 276 newtonian_ig_ctx->state_var = state_var; 277 newtonian_ig_ctx->idl_enable = idl_enable; 278 newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * second); 279 newtonian_ig_ctx->idl_start = idl_start * meter; 280 newtonian_ig_ctx->idl_length = idl_length * meter; 281 PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3)); 282 283 // -- Setup Context 284 setup_context->reference = reference; 285 setup_context->gas = *newtonian_ig_ctx; 286 setup_context->lx = domain_size[0]; 287 setup_context->ly = domain_size[1]; 288 setup_context->lz = domain_size[2]; 289 setup_context->time = 0; 290 291 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context)); 292 PetscCallCeed(ceed, 293 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 294 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc)); 295 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, "evaluation time", offsetof(struct SetupContext_, time), 1, 296 "Time of evaluation")); 297 298 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context)); 299 PetscCallCeed(ceed, 300 CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx)); 301 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, FreeContextPetsc)); 302 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1, 303 "Size of timestep, delta t")); 304 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift", 305 offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1, 306 "Shift for mass matrix in IJacobian")); 307 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1, 308 "Current solution time")); 309 310 problem->apply_vol_rhs.qfunction_context = newtonian_ig_context; 311 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ifunction.qfunction_context)); 312 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ijacobian.qfunction_context)); 313 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow.qfunction_context)); 314 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow_jacobian.qfunction_context)); 315 316 if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 317 if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 318 if (bc->num_slip > 0) PetscCall(SlipBCSetup(problem, dm, ctx, newtonian_ig_context)); 319 320 if (unit_tests) { 321 PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx)); 322 } 323 PetscFunctionReturn(PETSC_SUCCESS); 324 } 325 326 PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData *problem, AppCtx app_ctx) { 327 MPI_Comm comm = user->comm; 328 Ceed ceed = user->ceed; 329 NewtonianIdealGasContext newtonian_ctx; 330 331 PetscFunctionBeginUser; 332 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ctx)); 333 PetscCall(PetscPrintf(comm, 334 " Problem:\n" 335 " Problem Name : %s\n" 336 " Stabilization : %s\n", 337 app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization])); 338 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ctx)); 339 PetscFunctionReturn(PETSC_SUCCESS); 340 } 341