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 problems using the Newtonian Qfunction 6 7 #include "../qfunctions/newtonian.h" 8 9 #include <ceed.h> 10 #include <petscdm.h> 11 12 #include <navierstokes.h> 13 14 // For use with PetscOptionsEnum 15 static const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "ENTROPY", "StateVariable", "STATEVAR_", NULL}; 16 17 static PetscErrorCode CheckQWithTolerance(const CeedScalar Q_s[5], const CeedScalar Q_a[5], const CeedScalar Q_b[5], const char *name, 18 PetscReal rtol_0, PetscReal rtol_u, PetscReal rtol_4) { 19 CeedScalar relative_error[5]; // relative error 20 CeedScalar divisor_threshold = 10 * CEED_EPSILON; 21 22 PetscFunctionBeginUser; 23 relative_error[0] = (Q_a[0] - Q_b[0]) / (fabs(Q_s[0]) > divisor_threshold ? Q_s[0] : 1); 24 relative_error[4] = (Q_a[4] - Q_b[4]) / (fabs(Q_s[4]) > divisor_threshold ? Q_s[4] : 1); 25 26 CeedScalar u_magnitude = sqrt(Square(Q_s[1]) + Square(Q_s[2]) + Square(Q_s[3])); 27 CeedScalar u_divisor = u_magnitude > divisor_threshold ? u_magnitude : 1; 28 for (int i = 1; i < 4; i++) { 29 relative_error[i] = (Q_a[i] - Q_b[i]) / u_divisor; 30 } 31 32 if (fabs(relative_error[0]) >= rtol_0) { 33 printf("%s[0] error %g (expected %.10e, got %.10e)\n", name, relative_error[0], Q_s[0], Q_a[0]); 34 } 35 for (int i = 1; i < 4; i++) { 36 if (fabs(relative_error[i]) >= rtol_u) { 37 printf("%s[%d] error %g (expected %.10e, got %.10e)\n", name, i, relative_error[i], Q_s[i], Q_a[i]); 38 } 39 } 40 if (fabs(relative_error[4]) >= rtol_4) { 41 printf("%s[4] error %g (expected %.10e, got %.10e)\n", name, relative_error[4], Q_s[4], Q_a[4]); 42 } 43 PetscFunctionReturn(PETSC_SUCCESS); 44 } 45 46 // @brief Verify `StateFromQ` by converting A0 -> B0 -> A0_test, where A0 should equal A0_test 47 static PetscErrorCode TestState(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5], 48 CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) { 49 CeedScalar B0[5], A0_test[5]; 50 char buf[128]; 51 const char *const StateVariables_Initial[] = {"U", "Y", "V"}; 52 53 PetscFunctionBeginUser; 54 const char *A_initial = StateVariables_Initial[state_var_A]; 55 const char *B_initial = StateVariables_Initial[state_var_B]; 56 57 State state_A0 = StateFromQ(gas, A0, state_var_A); 58 StateToQ(gas, state_A0, B0, state_var_B); 59 State state_B0 = StateFromQ(gas, B0, state_var_B); 60 StateToQ(gas, state_B0, A0_test, state_var_A); 61 62 snprintf(buf, sizeof buf, "%s->%s->%s: %s", A_initial, B_initial, A_initial, A_initial); 63 PetscCall(CheckQWithTolerance(A0, A0_test, A0, buf, rtol_0, rtol_u, rtol_4)); 64 PetscFunctionReturn(PETSC_SUCCESS); 65 } 66 67 // @brief Verify `StateFromQ_fwd` via a finite difference approximation 68 static PetscErrorCode TestState_fwd(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5], 69 CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) { 70 CeedScalar eps = 4e-7; // Finite difference step 71 char buf[128]; 72 const char *const StateVariables_Initial[] = {"U", "Y", "V"}; 73 74 PetscFunctionBeginUser; 75 const char *A_initial = StateVariables_Initial[state_var_A]; 76 const char *B_initial = StateVariables_Initial[state_var_B]; 77 State state_0 = StateFromQ(gas, A0, state_var_A); 78 79 for (int i = 0; i < 5; i++) { 80 CeedScalar dB[5] = {0.}, dB_fd[5] = {0.}; 81 { // Calculate dB using State functions 82 CeedScalar dA[5] = {0}; 83 84 dA[i] = A0[i]; 85 State dstate_0 = StateFromQ_fwd(gas, state_0, dA, state_var_A); 86 StateToQ_fwd(gas, state_0, dstate_0, dB, state_var_B); 87 } 88 89 { // Calculate dB_fd via finite difference approximation 90 CeedScalar A1[5], B0[5], B1[5]; 91 92 for (int j = 0; j < 5; j++) A1[j] = (1 + eps * (i == j)) * A0[j]; 93 State state_1 = StateFromQ(gas, A1, state_var_A); 94 StateToQ(gas, state_0, B0, state_var_B); 95 StateToQ(gas, state_1, B1, state_var_B); 96 for (int j = 0; j < 5; j++) dB_fd[j] = (B1[j] - B0[j]) / eps; 97 } 98 99 snprintf(buf, sizeof buf, "d%s->d%s: StateFrom%s_fwd i=%d: d%s", A_initial, B_initial, A_initial, i, B_initial); 100 PetscCall(CheckQWithTolerance(dB_fd, dB, dB_fd, buf, rtol_0, rtol_u, rtol_4)); 101 } 102 PetscFunctionReturn(PETSC_SUCCESS); 103 } 104 105 // @brief Test the Newtonian State transformation functions, `StateFrom*` 106 static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) { 107 Units units = user->units; 108 const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, K = units->Kelvin; 109 110 PetscFunctionBeginUser; 111 const CeedScalar T = 200 * K; 112 const CeedScalar rho = 1.2 * kg / Cube(m); 113 const CeedScalar P = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T; 114 const CeedScalar u_base = 40 * m / sec; 115 const CeedScalar u[3] = {u_base, u_base * 1.1, u_base * 1.2}; 116 const CeedScalar e_kinetic = 0.5 * Dot3(u, u); 117 const CeedScalar e_internal = gas->cv * T; 118 const CeedScalar e_total = e_kinetic + e_internal; 119 const CeedScalar gamma = HeatCapacityRatio(gas); 120 const CeedScalar entropy = log(P) - gamma * log(rho); 121 const CeedScalar rho_div_p = rho / P; 122 const CeedScalar Y0[5] = {P, u[0], u[1], u[2], T}; 123 const CeedScalar U0[5] = {rho, rho * u[0], rho * u[1], rho * u[2], rho * e_total}; 124 const CeedScalar V0[5] = {(gamma - entropy) / (gamma - 1) - rho_div_p * (e_kinetic), rho_div_p * u[0], rho_div_p * u[1], rho_div_p * u[2], 125 -rho_div_p}; 126 127 { 128 CeedScalar rtol = 20 * CEED_EPSILON; 129 130 PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol)); 131 PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol)); 132 PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol)); 133 PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, rtol, rtol, rtol)); 134 PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, rtol, rtol, rtol)); 135 PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, rtol, rtol, rtol)); 136 } 137 138 { 139 CeedScalar rtol = 5e-6; 140 141 PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol)); 142 PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol)); 143 PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol)); 144 PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, 10 * rtol, rtol, rtol)); 145 PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, 5 * rtol, rtol, rtol)); 146 PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, 5 * rtol, 5 * rtol, 5 * rtol)); 147 } 148 PetscFunctionReturn(PETSC_SUCCESS); 149 } 150 151 // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping 152 // 153 // Only used for SUPG stabilization 154 PetscErrorCode CreateKSPMassOperator_NewtonianStabilized(User user, CeedOperator *op_mass) { 155 Ceed ceed = user->ceed; 156 CeedInt num_comp_q, q_data_size; 157 CeedQFunction qf_mass; 158 CeedElemRestriction elem_restr_q, elem_restr_qd_i; 159 CeedBasis basis_q; 160 CeedVector q_data; 161 CeedQFunctionContext qf_ctx = NULL; 162 PetscInt dim = 3; 163 164 PetscFunctionBeginUser; 165 { // Get restriction and basis from the RHS function 166 CeedOperator *sub_ops; 167 CeedOperatorField field; 168 PetscInt sub_op_index = 0; // will be 0 for the volume op 169 170 PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops)); 171 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field)); 172 PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q)); 173 PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q)); 174 175 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field)); 176 PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i)); 177 PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data)); 178 179 PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx)); 180 } 181 182 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 183 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size)); 184 185 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Newtonian_Conserv, MassFunction_Newtonian_Conserv_loc, &qf_mass)); 186 187 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx)); 188 PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0)); 189 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP)); 190 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP)); 191 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE)); 192 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP)); 193 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD)); 194 195 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 196 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 197 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed)); 198 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data)); 199 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 200 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 201 202 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qf_ctx)); 203 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 204 PetscFunctionReturn(PETSC_SUCCESS); 205 } 206 207 PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 208 SetupContext setup_context; 209 User user = *(User *)ctx; 210 CeedInt degree = user->app_ctx->degree; 211 StabilizationType stab; 212 StateVariable state_var; 213 MPI_Comm comm = user->comm; 214 Ceed ceed = user->ceed; 215 PetscBool implicit; 216 PetscBool unit_tests; 217 NewtonianIdealGasContext newtonian_ig_ctx; 218 CeedQFunctionContext newtonian_ig_qfctx; 219 220 PetscFunctionBeginUser; 221 PetscCall(PetscCalloc1(1, &setup_context)); 222 PetscCall(PetscCalloc1(1, &newtonian_ig_ctx)); 223 224 // ------------------------------------------------------ 225 // Setup Generic Newtonian IG Problem 226 // ------------------------------------------------------ 227 problem->dim = 3; 228 problem->jac_data_size_sur = 11; 229 problem->compute_exact_solution_error = PETSC_FALSE; 230 problem->print_info = PRINT_NEWTONIAN; 231 problem->uses_newtonian = PETSC_TRUE; 232 233 // ------------------------------------------------------ 234 // Create the libCEED context 235 // ------------------------------------------------------ 236 CeedScalar cv = 717.; // J/(kg K) 237 CeedScalar cp = 1004.; // J/(kg K) 238 CeedScalar g[3] = {0, 0, 0}; // m/s^2 239 CeedScalar lambda = -2. / 3.; // - 240 CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 241 CeedScalar k = 0.02638; // W/(m K) 242 CeedScalar c_tau = 0.5 / degree; // - 243 CeedScalar Ctau_t = 1.0; // - 244 CeedScalar Cv_func[3] = {36, 60, 128}; 245 CeedScalar Ctau_v = Cv_func[(CeedInt)Min(3, degree) - 1]; 246 CeedScalar Ctau_C = 0.25 / degree; 247 CeedScalar Ctau_M = 0.25 / degree; 248 CeedScalar Ctau_E = 0.125; 249 PetscReal domain_min[3], domain_max[3], domain_size[3]; 250 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 251 for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 252 253 StatePrimitive reference = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15}; 254 CeedScalar idl_decay_time = -1, idl_start = 0, idl_length = 0, idl_pressure = reference.pressure; 255 PetscBool idl_enable = PETSC_FALSE; 256 257 // ------------------------------------------------------ 258 // Create the PETSc context 259 // ------------------------------------------------------ 260 PetscScalar meter = 1; // 1 meter in scaled length units 261 PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 262 PetscScalar second = 1; // 1 second in scaled time units 263 PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 264 PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 265 266 // ------------------------------------------------------ 267 // Command line Options 268 // ------------------------------------------------------ 269 PetscBool given_option = PETSC_FALSE; 270 PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL); 271 // -- Conservative vs Primitive variables 272 PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE), 273 (PetscEnum *)&state_var, NULL)); 274 275 switch (state_var) { 276 case STATEVAR_CONSERVATIVE: 277 problem->ics.qf_func_ptr = ICsNewtonianIG_Conserv; 278 problem->ics.qf_loc = ICsNewtonianIG_Conserv_loc; 279 problem->apply_vol_rhs.qf_func_ptr = RHSFunction_Newtonian; 280 problem->apply_vol_rhs.qf_loc = RHSFunction_Newtonian_loc; 281 problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Conserv; 282 problem->apply_vol_ifunction.qf_loc = IFunction_Newtonian_Conserv_loc; 283 problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Conserv; 284 problem->apply_vol_ijacobian.qf_loc = IJacobian_Newtonian_Conserv_loc; 285 problem->apply_inflow.qf_func_ptr = BoundaryIntegral_Conserv; 286 problem->apply_inflow.qf_loc = BoundaryIntegral_Conserv_loc; 287 problem->apply_inflow_jacobian.qf_func_ptr = BoundaryIntegral_Jacobian_Conserv; 288 problem->apply_inflow_jacobian.qf_loc = BoundaryIntegral_Jacobian_Conserv_loc; 289 break; 290 case STATEVAR_PRIMITIVE: 291 problem->ics.qf_func_ptr = ICsNewtonianIG_Prim; 292 problem->ics.qf_loc = ICsNewtonianIG_Prim_loc; 293 problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Prim; 294 problem->apply_vol_ifunction.qf_loc = IFunction_Newtonian_Prim_loc; 295 problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Prim; 296 problem->apply_vol_ijacobian.qf_loc = IJacobian_Newtonian_Prim_loc; 297 problem->apply_inflow.qf_func_ptr = BoundaryIntegral_Prim; 298 problem->apply_inflow.qf_loc = BoundaryIntegral_Prim_loc; 299 problem->apply_inflow_jacobian.qf_func_ptr = BoundaryIntegral_Jacobian_Prim; 300 problem->apply_inflow_jacobian.qf_loc = BoundaryIntegral_Jacobian_Prim_loc; 301 break; 302 case STATEVAR_ENTROPY: 303 problem->ics.qf_func_ptr = ICsNewtonianIG_Entropy; 304 problem->ics.qf_loc = ICsNewtonianIG_Entropy_loc; 305 problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Entropy; 306 problem->apply_vol_ifunction.qf_loc = IFunction_Newtonian_Entropy_loc; 307 problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Entropy; 308 problem->apply_vol_ijacobian.qf_loc = IJacobian_Newtonian_Entropy_loc; 309 problem->apply_inflow.qf_func_ptr = BoundaryIntegral_Entropy; 310 problem->apply_inflow.qf_loc = BoundaryIntegral_Entropy_loc; 311 problem->apply_inflow_jacobian.qf_func_ptr = BoundaryIntegral_Jacobian_Entropy; 312 problem->apply_inflow_jacobian.qf_loc = BoundaryIntegral_Jacobian_Entropy_loc; 313 break; 314 } 315 316 // -- Physics 317 PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL)); 318 PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL)); 319 PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL)); 320 PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL)); 321 PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL)); 322 323 PetscInt dim = problem->dim; 324 PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL)); 325 PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option)); 326 if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim); 327 328 PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 329 PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 330 PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 331 PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL)); 332 PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL)); 333 PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL)); 334 PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL)); 335 PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 336 PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL)); 337 338 dim = 3; 339 PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL)); 340 PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL)); 341 PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL)); 342 343 // -- Units 344 PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 345 meter = fabs(meter); 346 PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 347 kilogram = fabs(kilogram); 348 PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 349 second = fabs(second); 350 PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL)); 351 Kelvin = fabs(Kelvin); 352 353 // -- Warnings 354 PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP, 355 "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n"); 356 357 PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point", 358 NULL, idl_decay_time, &idl_decay_time, &idl_enable)); 359 PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero."); 360 if (idl_decay_time < 0) idl_enable = PETSC_FALSE; 361 PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL)); 362 PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL)); 363 idl_pressure = reference.pressure; 364 PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-reference_pressure`)", NULL, idl_pressure, 365 &idl_pressure, NULL)); 366 PetscOptionsEnd(); 367 368 if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized; 369 370 // ------------------------------------------------------ 371 // Set up the PETSc context 372 // ------------------------------------------------------ 373 // -- Define derived units 374 Pascal = kilogram / (meter * PetscSqr(second)); 375 J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 376 m_per_squared_s = meter / PetscSqr(second); 377 W_per_m_K = kilogram * meter / (pow(second, 3) * Kelvin); 378 379 user->units->meter = meter; 380 user->units->kilogram = kilogram; 381 user->units->second = second; 382 user->units->Kelvin = Kelvin; 383 user->units->Pascal = Pascal; 384 user->units->J_per_kg_K = J_per_kg_K; 385 user->units->m_per_squared_s = m_per_squared_s; 386 user->units->W_per_m_K = W_per_m_K; 387 388 // ------------------------------------------------------ 389 // Set up the libCEED context 390 // ------------------------------------------------------ 391 // -- Scale variables to desired units 392 cv *= J_per_kg_K; 393 cp *= J_per_kg_K; 394 mu *= Pascal * second; 395 k *= W_per_m_K; 396 for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter; 397 for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s; 398 reference.pressure *= Pascal; 399 for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second; 400 reference.temperature *= Kelvin; 401 problem->dm_scale = meter; 402 403 // -- Solver Settings 404 user->phys->implicit = implicit; 405 user->phys->state_var = state_var; 406 407 // -- QFunction Context 408 newtonian_ig_ctx->lambda = lambda; 409 newtonian_ig_ctx->mu = mu; 410 newtonian_ig_ctx->k = k; 411 newtonian_ig_ctx->cv = cv; 412 newtonian_ig_ctx->cp = cp; 413 newtonian_ig_ctx->c_tau = c_tau; 414 newtonian_ig_ctx->Ctau_t = Ctau_t; 415 newtonian_ig_ctx->Ctau_v = Ctau_v; 416 newtonian_ig_ctx->Ctau_C = Ctau_C; 417 newtonian_ig_ctx->Ctau_M = Ctau_M; 418 newtonian_ig_ctx->Ctau_E = Ctau_E; 419 newtonian_ig_ctx->stabilization = stab; 420 newtonian_ig_ctx->is_implicit = implicit; 421 newtonian_ig_ctx->state_var = state_var; 422 newtonian_ig_ctx->idl_enable = idl_enable; 423 newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * second); 424 newtonian_ig_ctx->idl_start = idl_start * meter; 425 newtonian_ig_ctx->idl_length = idl_length * meter; 426 newtonian_ig_ctx->idl_pressure = idl_pressure; 427 PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3)); 428 429 // -- Setup Context 430 setup_context->reference = reference; 431 setup_context->gas = *newtonian_ig_ctx; 432 setup_context->lx = domain_size[0]; 433 setup_context->ly = domain_size[1]; 434 setup_context->lz = domain_size[2]; 435 setup_context->time = 0; 436 437 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfctx)); 438 PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 439 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc)); 440 PetscCallCeed( 441 ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfctx, "evaluation time", offsetof(struct SetupContext_, time), 1, "Time of evaluation")); 442 443 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &newtonian_ig_qfctx)); 444 PetscCallCeed(ceed, CeedQFunctionContextSetData(newtonian_ig_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx)); 445 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 446 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1, 447 "Size of timestep, delta t")); 448 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "ijacobian time shift", 449 offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1, 450 "Shift for mass matrix in IJacobian")); 451 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1, 452 "Current solution time")); 453 454 problem->apply_vol_rhs.qfctx = newtonian_ig_qfctx; 455 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ifunction.qfctx)); 456 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ijacobian.qfctx)); 457 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_inflow.qfctx)); 458 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_inflow_jacobian.qfctx)); 459 460 if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 461 if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 462 if (bc->num_slip > 0) PetscCall(SlipBCSetup(problem, dm, ctx, newtonian_ig_qfctx)); 463 464 if (unit_tests) { 465 PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx)); 466 } 467 PetscFunctionReturn(PETSC_SUCCESS); 468 } 469 470 PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx) { 471 MPI_Comm comm = user->comm; 472 Ceed ceed = user->ceed; 473 NewtonianIdealGasContext newtonian_ctx; 474 475 PetscFunctionBeginUser; 476 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ctx)); 477 PetscCall(PetscPrintf(comm, 478 " Problem:\n" 479 " Problem Name : %s\n" 480 " Stabilization : %s\n", 481 app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization])); 482 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ctx)); 483 PetscFunctionReturn(PETSC_SUCCESS); 484 } 485