1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 33a8779fbSJames Wright // 4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 53a8779fbSJames Wright // 6727da7e7SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 73a8779fbSJames Wright 83a8779fbSJames Wright /// @file 93a8779fbSJames Wright /// Utility functions for setting up problems using the Newtonian Qfunction 103a8779fbSJames Wright 113a8779fbSJames Wright #include "../qfunctions/newtonian.h" 123a8779fbSJames Wright 13e419654dSJeremy L Thompson #include <ceed.h> 14e419654dSJeremy L Thompson #include <petscdm.h> 15e419654dSJeremy L Thompson 162b916ea7SJeremy L Thompson #include "../navierstokes.h" 176dd99bedSLeila Ghaffari 18f7534ca8SJed Brown // For use with PetscOptionsEnum 19f7534ca8SJed Brown static const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "StateVariable", "STATEVAR_", NULL}; 20f7534ca8SJed Brown 212b916ea7SJeremy L Thompson // Compute relative error |a - b|/|s| 222b916ea7SJeremy L Thompson static PetscErrorCode CheckPrimitiveWithTolerance(StatePrimitive sY, StatePrimitive aY, StatePrimitive bY, const char *name, PetscReal rtol_pressure, 232b916ea7SJeremy L Thompson PetscReal rtol_velocity, PetscReal rtol_temperature) { 24f0b65372SJed Brown StatePrimitive eY; // relative error 2506f41313SJames Wright 2606f41313SJames Wright PetscFunctionBeginUser; 27f0b65372SJed Brown eY.pressure = (aY.pressure - bY.pressure) / sY.pressure; 282b916ea7SJeremy L Thompson PetscScalar u = sqrt(Square(sY.velocity[0]) + Square(sY.velocity[1]) + Square(sY.velocity[2])); 29f0b65372SJed Brown for (int j = 0; j < 3; j++) eY.velocity[j] = (aY.velocity[j] - bY.velocity[j]) / u; 30f0b65372SJed Brown eY.temperature = (aY.temperature - bY.temperature) / sY.temperature; 312b916ea7SJeremy L Thompson if (fabs(eY.pressure) > rtol_pressure) printf("%s: pressure error %g\n", name, eY.pressure); 322b916ea7SJeremy L Thompson for (int j = 0; j < 3; j++) { 332b916ea7SJeremy L Thompson if (fabs(eY.velocity[j]) > rtol_velocity) printf("%s: velocity[%d] error %g\n", name, j, eY.velocity[j]); 342b916ea7SJeremy L Thompson } 352b916ea7SJeremy L Thompson if (fabs(eY.temperature) > rtol_temperature) printf("%s: temperature error %g\n", name, eY.temperature); 36d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 37f0b65372SJed Brown } 38f0b65372SJed Brown 392b916ea7SJeremy L Thompson static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) { 40f0b65372SJed Brown Units units = user->units; 41f0b65372SJed Brown const CeedScalar eps = 1e-6; 422b916ea7SJeremy L Thompson const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, Pascal = units->Pascal; 43f0b65372SJed Brown PetscFunctionBeginUser; 442b916ea7SJeremy L Thompson const CeedScalar rho = 1.2 * kg / (m * m * m), u = 40 * m / sec; 45f0b65372SJed Brown CeedScalar U[5] = {rho, rho * u, rho * u * 1.1, rho * u * 1.2, 250e3 * Pascal + .5 * rho * u * u}; 46edcfef1bSKenneth E. Jansen State s = StateFromU(gas, U); 47f0b65372SJed Brown for (int i = 0; i < 8; i++) { 48c4f7305dSKenneth E. Jansen CeedScalar dU[5] = {0}; 49f0b65372SJed Brown if (i < 5) dU[i] = U[i]; 50edcfef1bSKenneth E. Jansen State ds = StateFromU_fwd(gas, s, dU); 51f0b65372SJed Brown for (int j = 0; j < 5; j++) dU[j] = (1 + eps * (i == j)) * U[j]; 52edcfef1bSKenneth E. Jansen State t = StateFromU(gas, dU); 53f0b65372SJed Brown StatePrimitive dY; 54f0b65372SJed Brown dY.pressure = (t.Y.pressure - s.Y.pressure) / eps; 552b916ea7SJeremy L Thompson for (int j = 0; j < 3; j++) dY.velocity[j] = (t.Y.velocity[j] - s.Y.velocity[j]) / eps; 56f0b65372SJed Brown dY.temperature = (t.Y.temperature - s.Y.temperature) / eps; 57f0b65372SJed Brown char buf[128]; 58f0b65372SJed Brown snprintf(buf, sizeof buf, "StateFromU_fwd i=%d", i); 59f0b65372SJed Brown PetscCall(CheckPrimitiveWithTolerance(dY, ds.Y, dY, buf, 5e-6, 1e-6, 1e-6)); 60f0b65372SJed Brown } 61d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 62f0b65372SJed Brown } 63f0b65372SJed Brown 6465dee3d2SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping 6565dee3d2SJames Wright // 6665dee3d2SJames Wright // Only used for SUPG stabilization 6765dee3d2SJames Wright PetscErrorCode CreateKSPMassOperator_NewtonianStabilized(User user, CeedOperator *op_mass) { 6865dee3d2SJames Wright Ceed ceed = user->ceed; 6965dee3d2SJames Wright CeedInt num_comp_q, q_data_size; 7065dee3d2SJames Wright CeedQFunction qf_mass; 7165dee3d2SJames Wright CeedElemRestriction elem_restr_q, elem_restr_qd_i; 7265dee3d2SJames Wright CeedBasis basis_q; 7365dee3d2SJames Wright CeedVector q_data; 7465dee3d2SJames Wright CeedQFunctionContext qf_ctx = NULL; 7565dee3d2SJames Wright PetscInt dim = 3; 7665dee3d2SJames Wright 7765dee3d2SJames Wright PetscFunctionBeginUser; 7865dee3d2SJames Wright { // Get restriction and basis from the RHS function 7965dee3d2SJames Wright CeedOperator *sub_ops; 8065dee3d2SJames Wright CeedOperatorField field; 8165dee3d2SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 8265dee3d2SJames Wright 8365dee3d2SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops)); 8465dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field)); 8565dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q)); 8665dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q)); 8765dee3d2SJames Wright 8865dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field)); 8965dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i)); 9065dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data)); 9165dee3d2SJames Wright 9265dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx)); 9365dee3d2SJames Wright } 9465dee3d2SJames Wright 9565dee3d2SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 9665dee3d2SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size)); 9765dee3d2SJames Wright 9865dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Newtonian_Conserv, MassFunction_Newtonian_Conserv_loc, &qf_mass)); 9965dee3d2SJames Wright 10065dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx)); 10165dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0)); 10265dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP)); 10365dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP)); 10465dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE)); 10565dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP)); 10665dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD)); 10765dee3d2SJames Wright 10865dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 10965dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 11065dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed)); 11165dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data)); 11265dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 11365dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 11465dee3d2SJames Wright 11565dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 11665dee3d2SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 11765dee3d2SJames Wright } 118*4c07ec22SJames Wright 119991aef52SJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 120b7f03f12SJed Brown SetupContext setup_context; 1213a8779fbSJames Wright User user = *(User *)ctx; 1224351d27fSLeila Ghaffari CeedInt degree = user->app_ctx->degree; 1233a8779fbSJames Wright StabilizationType stab; 1243636f6a4SJames Wright StateVariable state_var; 125b4c37c5cSJames Wright MPI_Comm comm = user->comm; 126b4c37c5cSJames Wright Ceed ceed = user->ceed; 1273a8779fbSJames Wright PetscBool implicit; 1280d8cd818SJames Wright PetscBool unit_tests; 12915a3537eSJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 13015a3537eSJed Brown CeedQFunctionContext newtonian_ig_context; 1313a8779fbSJames Wright 13215a3537eSJed Brown PetscFunctionBeginUser; 1332b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &setup_context)); 1342b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &newtonian_ig_ctx)); 1353a8779fbSJames Wright 1363a8779fbSJames Wright // ------------------------------------------------------ 1373a8779fbSJames Wright // Setup Generic Newtonian IG Problem 1383a8779fbSJames Wright // ------------------------------------------------------ 1393a8779fbSJames Wright problem->dim = 3; 140b01ba163SJames Wright problem->jac_data_size_sur = 11; 1413a8779fbSJames Wright problem->non_zero_time = PETSC_FALSE; 142cbe60e31SLeila Ghaffari problem->print_info = PRINT_NEWTONIAN; 1437d8a615bSJames Wright problem->uses_newtonian = PETSC_TRUE; 1443a8779fbSJames Wright 1453a8779fbSJames Wright // ------------------------------------------------------ 1463a8779fbSJames Wright // Create the libCEED context 1473a8779fbSJames Wright // ------------------------------------------------------ 1483a8779fbSJames Wright CeedScalar cv = 717.; // J/(kg K) 1493a8779fbSJames Wright CeedScalar cp = 1004.; // J/(kg K) 150d9bb1cdbSJames Wright CeedScalar g[3] = {0, 0, 0}; // m/s^2 1513a8779fbSJames Wright CeedScalar lambda = -2. / 3.; // - 152bb8a0c61SJames Wright CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 1533a8779fbSJames Wright CeedScalar k = 0.02638; // W/(m K) 1544351d27fSLeila Ghaffari CeedScalar c_tau = 0.5 / degree; // - 155bb8a0c61SJames Wright CeedScalar Ctau_t = 1.0; // - 156b5786772SLeila Ghaffari CeedScalar Cv_func[3] = {36, 60, 128}; 157c176988bSLeila Ghaffari CeedScalar Ctau_v = Cv_func[(CeedInt)Min(3, degree) - 1]; 1582a0eee6eSLeila Ghaffari CeedScalar Ctau_C = 0.25 / degree; 1592a0eee6eSLeila Ghaffari CeedScalar Ctau_M = 0.25 / degree; 1602a0eee6eSLeila Ghaffari CeedScalar Ctau_E = 0.125; 1613a8779fbSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 1622b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 163493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 1643a8779fbSJames Wright 165b8fb7609SAdeleke O. Bankole StatePrimitive reference = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15}; 166a213b8aaSJames Wright CeedScalar idl_decay_time = -1, idl_start = 0, idl_length = 0, idl_pressure = reference.pressure; 167e7754af5SKenneth E. Jansen PetscBool idl_enable = PETSC_FALSE; 168b8fb7609SAdeleke O. Bankole 1693a8779fbSJames Wright // ------------------------------------------------------ 1703a8779fbSJames Wright // Create the PETSc context 1713a8779fbSJames Wright // ------------------------------------------------------ 172bb8a0c61SJames Wright PetscScalar meter = 1; // 1 meter in scaled length units 173bb8a0c61SJames Wright PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 174bb8a0c61SJames Wright PetscScalar second = 1; // 1 second in scaled time units 1753a8779fbSJames Wright PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 1763a8779fbSJames Wright PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 1773a8779fbSJames Wright 1783a8779fbSJames Wright // ------------------------------------------------------ 1793a8779fbSJames Wright // Command line Options 1803a8779fbSJames Wright // ------------------------------------------------------ 181d9bb1cdbSJames Wright PetscBool given_option = PETSC_FALSE; 1822b916ea7SJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL); 183cbe60e31SLeila Ghaffari // -- Conservative vs Primitive variables 1842b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE), 1852b916ea7SJeremy L Thompson (PetscEnum *)&state_var, NULL)); 1863636f6a4SJames Wright 1873636f6a4SJames Wright switch (state_var) { 1883636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 189b8fb7609SAdeleke O. Bankole problem->ics.qfunction = ICsNewtonianIG_Conserv; 190b8fb7609SAdeleke O. Bankole problem->ics.qfunction_loc = ICsNewtonianIG_Conserv_loc; 191cbe60e31SLeila Ghaffari problem->apply_vol_rhs.qfunction = RHSFunction_Newtonian; 192cbe60e31SLeila Ghaffari problem->apply_vol_rhs.qfunction_loc = RHSFunction_Newtonian_loc; 19376555becSJames Wright problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Conserv; 19476555becSJames Wright problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Conserv_loc; 19576555becSJames Wright problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Conserv; 19676555becSJames Wright problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Conserv_loc; 197d4559bbeSJames Wright problem->apply_inflow.qfunction = BoundaryIntegral_Conserv; 198d4559bbeSJames Wright problem->apply_inflow.qfunction_loc = BoundaryIntegral_Conserv_loc; 199d4559bbeSJames Wright problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Conserv; 200d4559bbeSJames Wright problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc; 2013636f6a4SJames Wright break; 2023636f6a4SJames Wright 2033636f6a4SJames Wright case STATEVAR_PRIMITIVE: 2043636f6a4SJames Wright problem->ics.qfunction = ICsNewtonianIG_Prim; 2053636f6a4SJames Wright problem->ics.qfunction_loc = ICsNewtonianIG_Prim_loc; 2063636f6a4SJames Wright problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Prim; 2073636f6a4SJames Wright problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Prim_loc; 2083636f6a4SJames Wright problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Prim; 2093636f6a4SJames Wright problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Prim_loc; 2103636f6a4SJames Wright problem->apply_inflow.qfunction = BoundaryIntegral_Prim; 2113636f6a4SJames Wright problem->apply_inflow.qfunction_loc = BoundaryIntegral_Prim_loc; 2123636f6a4SJames Wright problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Prim; 2133636f6a4SJames Wright problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc; 2143636f6a4SJames Wright break; 215cbe60e31SLeila Ghaffari } 2161485969bSJeremy L Thompson 2173a8779fbSJames Wright // -- Physics 2182b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL)); 2192b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL)); 2202b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL)); 2212b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL)); 2222b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL)); 2233a8779fbSJames Wright 224bb8a0c61SJames Wright PetscInt dim = problem->dim; 225d9bb1cdbSJames Wright PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL)); 226d9bb1cdbSJames Wright PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option)); 227d9bb1cdbSJames Wright if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim); 228d9bb1cdbSJames Wright 2292b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 2302b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 2312b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 2322b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL)); 2332b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL)); 2342b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL)); 2352b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL)); 2362b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 2372b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL)); 2383a8779fbSJames Wright 239db95602eSJames Wright dim = 3; 240b8fb7609SAdeleke O. Bankole PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL)); 241b8fb7609SAdeleke O. Bankole PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL)); 242b8fb7609SAdeleke O. Bankole PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL)); 243b8fb7609SAdeleke O. Bankole 2443a8779fbSJames Wright // -- Units 2452b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 2463a8779fbSJames Wright meter = fabs(meter); 2472b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 2483a8779fbSJames Wright kilogram = fabs(kilogram); 2492b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 2503a8779fbSJames Wright second = fabs(second); 2512b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL)); 2523a8779fbSJames Wright Kelvin = fabs(Kelvin); 2533a8779fbSJames Wright 2543a8779fbSJames Wright // -- Warnings 2555d28dccaSJames Wright PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP, 2565d28dccaSJames Wright "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n"); 257e7754af5SKenneth E. Jansen 258e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point", 259e7754af5SKenneth E. Jansen NULL, idl_decay_time, &idl_decay_time, &idl_enable)); 2609cbdf780SJames Wright PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero."); 2619cbdf780SJames Wright if (idl_decay_time < 0) idl_enable = PETSC_FALSE; 262e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL)); 263e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL)); 264a213b8aaSJames Wright idl_pressure = reference.pressure; 265a213b8aaSJames Wright PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-reference_pressure`)", NULL, idl_pressure, 266a213b8aaSJames Wright &idl_pressure, NULL)); 2671485969bSJeremy L Thompson PetscOptionsEnd(); 2683a8779fbSJames Wright 26965dee3d2SJames Wright if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized; 27065dee3d2SJames Wright 2713a8779fbSJames Wright // ------------------------------------------------------ 2723a8779fbSJames Wright // Set up the PETSc context 2733a8779fbSJames Wright // ------------------------------------------------------ 2743a8779fbSJames Wright // -- Define derived units 2753a8779fbSJames Wright Pascal = kilogram / (meter * PetscSqr(second)); 2763a8779fbSJames Wright J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 2773a8779fbSJames Wright m_per_squared_s = meter / PetscSqr(second); 2783a8779fbSJames Wright W_per_m_K = kilogram * meter / (pow(second, 3) * Kelvin); 2793a8779fbSJames Wright 2803a8779fbSJames Wright user->units->meter = meter; 2813a8779fbSJames Wright user->units->kilogram = kilogram; 2823a8779fbSJames Wright user->units->second = second; 2833a8779fbSJames Wright user->units->Kelvin = Kelvin; 2843a8779fbSJames Wright user->units->Pascal = Pascal; 2853a8779fbSJames Wright user->units->J_per_kg_K = J_per_kg_K; 2863a8779fbSJames Wright user->units->m_per_squared_s = m_per_squared_s; 2873a8779fbSJames Wright user->units->W_per_m_K = W_per_m_K; 2883a8779fbSJames Wright 2893a8779fbSJames Wright // ------------------------------------------------------ 2903a8779fbSJames Wright // Set up the libCEED context 2913a8779fbSJames Wright // ------------------------------------------------------ 2923a8779fbSJames Wright // -- Scale variables to desired units 2933a8779fbSJames Wright cv *= J_per_kg_K; 2943a8779fbSJames Wright cp *= J_per_kg_K; 2953a8779fbSJames Wright mu *= Pascal * second; 2963a8779fbSJames Wright k *= W_per_m_K; 297493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter; 298493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s; 299b8fb7609SAdeleke O. Bankole reference.pressure *= Pascal; 300b8fb7609SAdeleke O. Bankole for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second; 301b8fb7609SAdeleke O. Bankole reference.temperature *= Kelvin; 3023a8779fbSJames Wright problem->dm_scale = meter; 3033a8779fbSJames Wright 3043a8779fbSJames Wright // -- Solver Settings 3053a8779fbSJames Wright user->phys->implicit = implicit; 3063636f6a4SJames Wright user->phys->state_var = state_var; 3073a8779fbSJames Wright 3083a8779fbSJames Wright // -- QFunction Context 30915a3537eSJed Brown newtonian_ig_ctx->lambda = lambda; 31015a3537eSJed Brown newtonian_ig_ctx->mu = mu; 31115a3537eSJed Brown newtonian_ig_ctx->k = k; 31215a3537eSJed Brown newtonian_ig_ctx->cv = cv; 31315a3537eSJed Brown newtonian_ig_ctx->cp = cp; 31415a3537eSJed Brown newtonian_ig_ctx->c_tau = c_tau; 31515a3537eSJed Brown newtonian_ig_ctx->Ctau_t = Ctau_t; 31615a3537eSJed Brown newtonian_ig_ctx->Ctau_v = Ctau_v; 31715a3537eSJed Brown newtonian_ig_ctx->Ctau_C = Ctau_C; 31815a3537eSJed Brown newtonian_ig_ctx->Ctau_M = Ctau_M; 31915a3537eSJed Brown newtonian_ig_ctx->Ctau_E = Ctau_E; 32015a3537eSJed Brown newtonian_ig_ctx->stabilization = stab; 3218085925cSJames Wright newtonian_ig_ctx->is_implicit = implicit; 3223636f6a4SJames Wright newtonian_ig_ctx->state_var = state_var; 323e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_enable = idl_enable; 324e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * second); 325e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_start = idl_start * meter; 326e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_length = idl_length * meter; 327a213b8aaSJames Wright newtonian_ig_ctx->idl_pressure = idl_pressure; 3282b916ea7SJeremy L Thompson PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3)); 3293a8779fbSJames Wright 330b8fb7609SAdeleke O. Bankole // -- Setup Context 331b8fb7609SAdeleke O. Bankole setup_context->reference = reference; 332b8fb7609SAdeleke O. Bankole setup_context->gas = *newtonian_ig_ctx; 333b8fb7609SAdeleke O. Bankole setup_context->lx = domain_size[0]; 334b8fb7609SAdeleke O. Bankole setup_context->ly = domain_size[1]; 335b8fb7609SAdeleke O. Bankole setup_context->lz = domain_size[2]; 336b8fb7609SAdeleke O. Bankole setup_context->time = 0; 337b8fb7609SAdeleke O. Bankole 338b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context)); 339b4c37c5cSJames Wright PetscCallCeed(ceed, 340b4c37c5cSJames Wright CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 341b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc)); 342b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, "evaluation time", offsetof(struct SetupContext_, time), 1, 343b4c37c5cSJames Wright "Time of evaluation")); 34415a3537eSJed Brown 345b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context)); 346b4c37c5cSJames Wright PetscCallCeed(ceed, 347b4c37c5cSJames Wright CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx)); 348b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, FreeContextPetsc)); 349b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1, 350b4c37c5cSJames Wright "Size of timestep, delta t")); 351b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift", 352b4c37c5cSJames Wright offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1, 353b4c37c5cSJames Wright "Shift for mass matrix in IJacobian")); 354b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1, 355b4c37c5cSJames Wright "Current solution time")); 356a8bca8acSJames Wright 35715a3537eSJed Brown problem->apply_vol_rhs.qfunction_context = newtonian_ig_context; 358b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ifunction.qfunction_context)); 359b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ijacobian.qfunction_context)); 360b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow.qfunction_context)); 361b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow_jacobian.qfunction_context)); 362f0b65372SJed Brown 3639ed3d70dSJames Wright if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 3649ed3d70dSJames Wright if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 3659ed3d70dSJames Wright if (bc->num_slip > 0) PetscCall(SlipBCSetup(problem, dm, ctx, newtonian_ig_context)); 3669ed3d70dSJames Wright 367f0b65372SJed Brown if (unit_tests) { 368f0b65372SJed Brown PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx)); 369f0b65372SJed Brown } 370d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3713a8779fbSJames Wright } 3723a8779fbSJames Wright 373991aef52SJames Wright PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx) { 3742d49c0afSJames Wright MPI_Comm comm = user->comm; 375b4c37c5cSJames Wright Ceed ceed = user->ceed; 37615a3537eSJed Brown NewtonianIdealGasContext newtonian_ctx; 37715a3537eSJed Brown 3783a8779fbSJames Wright PetscFunctionBeginUser; 379b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ctx)); 3802b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 38115a3537eSJed Brown " Problem:\n" 38215a3537eSJed Brown " Problem Name : %s\n" 38315a3537eSJed Brown " Stabilization : %s\n", 3842b916ea7SJeremy L Thompson app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization])); 385b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ctx)); 386d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3873a8779fbSJames Wright } 388