1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, 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" 172b916ea7SJeremy L Thompson #include "../qfunctions/setupgeo.h" 186dd99bedSLeila Ghaffari 19f7534ca8SJed Brown // For use with PetscOptionsEnum 20f7534ca8SJed Brown static const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "StateVariable", "STATEVAR_", NULL}; 21f7534ca8SJed Brown 222b916ea7SJeremy L Thompson // Compute relative error |a - b|/|s| 232b916ea7SJeremy L Thompson static PetscErrorCode CheckPrimitiveWithTolerance(StatePrimitive sY, StatePrimitive aY, StatePrimitive bY, const char *name, PetscReal rtol_pressure, 242b916ea7SJeremy L Thompson PetscReal rtol_velocity, PetscReal rtol_temperature) { 25f0b65372SJed Brown StatePrimitive eY; // relative error 2606f41313SJames Wright 2706f41313SJames Wright PetscFunctionBeginUser; 28f0b65372SJed Brown eY.pressure = (aY.pressure - bY.pressure) / sY.pressure; 292b916ea7SJeremy L Thompson PetscScalar u = sqrt(Square(sY.velocity[0]) + Square(sY.velocity[1]) + Square(sY.velocity[2])); 30f0b65372SJed Brown for (int j = 0; j < 3; j++) eY.velocity[j] = (aY.velocity[j] - bY.velocity[j]) / u; 31f0b65372SJed Brown eY.temperature = (aY.temperature - bY.temperature) / sY.temperature; 322b916ea7SJeremy L Thompson if (fabs(eY.pressure) > rtol_pressure) printf("%s: pressure error %g\n", name, eY.pressure); 332b916ea7SJeremy L Thompson for (int j = 0; j < 3; j++) { 342b916ea7SJeremy L Thompson if (fabs(eY.velocity[j]) > rtol_velocity) printf("%s: velocity[%d] error %g\n", name, j, eY.velocity[j]); 352b916ea7SJeremy L Thompson } 362b916ea7SJeremy L Thompson if (fabs(eY.temperature) > rtol_temperature) printf("%s: temperature error %g\n", name, eY.temperature); 37d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 38f0b65372SJed Brown } 39f0b65372SJed Brown 402b916ea7SJeremy L Thompson static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) { 41f0b65372SJed Brown Units units = user->units; 42f0b65372SJed Brown const CeedScalar eps = 1e-6; 432b916ea7SJeremy L Thompson const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, Pascal = units->Pascal; 44f0b65372SJed Brown PetscFunctionBeginUser; 452b916ea7SJeremy L Thompson const CeedScalar rho = 1.2 * kg / (m * m * m), u = 40 * m / sec; 46f0b65372SJed Brown CeedScalar U[5] = {rho, rho * u, rho * u * 1.1, rho * u * 1.2, 250e3 * Pascal + .5 * rho * u * u}; 47*edcfef1bSKenneth E. Jansen State s = StateFromU(gas, U); 48f0b65372SJed Brown for (int i = 0; i < 8; i++) { 49f0b65372SJed Brown CeedScalar dU[5] = {0}, dx[3] = {0}; 50f0b65372SJed Brown if (i < 5) dU[i] = U[i]; 51*edcfef1bSKenneth E. Jansen State ds = StateFromU_fwd(gas, s, dU); 52f0b65372SJed Brown for (int j = 0; j < 5; j++) dU[j] = (1 + eps * (i == j)) * U[j]; 53*edcfef1bSKenneth E. Jansen State t = StateFromU(gas, dU); 54f0b65372SJed Brown StatePrimitive dY; 55f0b65372SJed Brown dY.pressure = (t.Y.pressure - s.Y.pressure) / eps; 562b916ea7SJeremy L Thompson for (int j = 0; j < 3; j++) dY.velocity[j] = (t.Y.velocity[j] - s.Y.velocity[j]) / eps; 57f0b65372SJed Brown dY.temperature = (t.Y.temperature - s.Y.temperature) / eps; 58f0b65372SJed Brown char buf[128]; 59f0b65372SJed Brown snprintf(buf, sizeof buf, "StateFromU_fwd i=%d", i); 60f0b65372SJed Brown PetscCall(CheckPrimitiveWithTolerance(dY, ds.Y, dY, buf, 5e-6, 1e-6, 1e-6)); 61f0b65372SJed Brown } 62d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 63f0b65372SJed Brown } 64f0b65372SJed Brown 65d1c51a42SJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) { 66b7f03f12SJed Brown SetupContext setup_context; 673a8779fbSJames Wright User user = *(User *)ctx; 684351d27fSLeila Ghaffari CeedInt degree = user->app_ctx->degree; 693a8779fbSJames Wright StabilizationType stab; 703636f6a4SJames Wright StateVariable state_var; 71b4c37c5cSJames Wright MPI_Comm comm = user->comm; 72b4c37c5cSJames Wright Ceed ceed = user->ceed; 733a8779fbSJames Wright PetscBool implicit; 743636f6a4SJames Wright PetscBool has_curr_time = PETSC_FALSE, unit_tests; 7515a3537eSJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 7615a3537eSJed Brown CeedQFunctionContext newtonian_ig_context; 773a8779fbSJames Wright 7815a3537eSJed Brown PetscFunctionBeginUser; 792b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &setup_context)); 802b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &newtonian_ig_ctx)); 813a8779fbSJames Wright 823a8779fbSJames Wright // ------------------------------------------------------ 833a8779fbSJames Wright // Setup Generic Newtonian IG Problem 843a8779fbSJames Wright // ------------------------------------------------------ 853a8779fbSJames Wright problem->dim = 3; 863a8779fbSJames Wright problem->q_data_size_vol = 10; 87493642f1SJames Wright problem->q_data_size_sur = 10; 88b01ba163SJames Wright problem->jac_data_size_sur = 11; 899785fe93SJed Brown problem->setup_vol.qfunction = Setup; 909785fe93SJed Brown problem->setup_vol.qfunction_loc = Setup_loc; 919785fe93SJed Brown problem->setup_sur.qfunction = SetupBoundary; 929785fe93SJed Brown problem->setup_sur.qfunction_loc = SetupBoundary_loc; 933a8779fbSJames Wright problem->non_zero_time = PETSC_FALSE; 94cbe60e31SLeila Ghaffari problem->print_info = PRINT_NEWTONIAN; 953a8779fbSJames Wright 963a8779fbSJames Wright // ------------------------------------------------------ 973a8779fbSJames Wright // Create the libCEED context 983a8779fbSJames Wright // ------------------------------------------------------ 993a8779fbSJames Wright CeedScalar cv = 717.; // J/(kg K) 1003a8779fbSJames Wright CeedScalar cp = 1004.; // J/(kg K) 101d9bb1cdbSJames Wright CeedScalar g[3] = {0, 0, 0}; // m/s^2 1023a8779fbSJames Wright CeedScalar lambda = -2. / 3.; // - 103bb8a0c61SJames Wright CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 1043a8779fbSJames Wright CeedScalar k = 0.02638; // W/(m K) 1054351d27fSLeila Ghaffari CeedScalar c_tau = 0.5 / degree; // - 106bb8a0c61SJames Wright CeedScalar Ctau_t = 1.0; // - 107b5786772SLeila Ghaffari CeedScalar Cv_func[3] = {36, 60, 128}; 108c176988bSLeila Ghaffari CeedScalar Ctau_v = Cv_func[(CeedInt)Min(3, degree) - 1]; 1092a0eee6eSLeila Ghaffari CeedScalar Ctau_C = 0.25 / degree; 1102a0eee6eSLeila Ghaffari CeedScalar Ctau_M = 0.25 / degree; 1112a0eee6eSLeila Ghaffari CeedScalar Ctau_E = 0.125; 1123a8779fbSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 1132b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 114493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 1153a8779fbSJames Wright 116b8fb7609SAdeleke O. Bankole StatePrimitive reference = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15}; 117e7754af5SKenneth E. Jansen CeedScalar idl_decay_time = -1, idl_start = 0, idl_length = 0; 118e7754af5SKenneth E. Jansen PetscBool idl_enable = PETSC_FALSE; 119b8fb7609SAdeleke O. Bankole 1203a8779fbSJames Wright // ------------------------------------------------------ 1213a8779fbSJames Wright // Create the PETSc context 1223a8779fbSJames Wright // ------------------------------------------------------ 123bb8a0c61SJames Wright PetscScalar meter = 1; // 1 meter in scaled length units 124bb8a0c61SJames Wright PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 125bb8a0c61SJames Wright PetscScalar second = 1; // 1 second in scaled time units 1263a8779fbSJames Wright PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 1273a8779fbSJames Wright PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 1283a8779fbSJames Wright 1293a8779fbSJames Wright // ------------------------------------------------------ 1303a8779fbSJames Wright // Command line Options 1313a8779fbSJames Wright // ------------------------------------------------------ 132d9bb1cdbSJames Wright PetscBool given_option = PETSC_FALSE; 1332b916ea7SJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL); 134cbe60e31SLeila Ghaffari // -- Conservative vs Primitive variables 1352b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE), 1362b916ea7SJeremy L Thompson (PetscEnum *)&state_var, NULL)); 1373636f6a4SJames Wright 1383636f6a4SJames Wright switch (state_var) { 1393636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 140b8fb7609SAdeleke O. Bankole problem->ics.qfunction = ICsNewtonianIG_Conserv; 141b8fb7609SAdeleke O. Bankole problem->ics.qfunction_loc = ICsNewtonianIG_Conserv_loc; 142cbe60e31SLeila Ghaffari problem->apply_vol_rhs.qfunction = RHSFunction_Newtonian; 143cbe60e31SLeila Ghaffari problem->apply_vol_rhs.qfunction_loc = RHSFunction_Newtonian_loc; 14476555becSJames Wright problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Conserv; 14576555becSJames Wright problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Conserv_loc; 14676555becSJames Wright problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Conserv; 14776555becSJames Wright problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Conserv_loc; 148d4559bbeSJames Wright problem->apply_inflow.qfunction = BoundaryIntegral_Conserv; 149d4559bbeSJames Wright problem->apply_inflow.qfunction_loc = BoundaryIntegral_Conserv_loc; 150d4559bbeSJames Wright problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Conserv; 151d4559bbeSJames Wright problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc; 1523636f6a4SJames Wright break; 1533636f6a4SJames Wright 1543636f6a4SJames Wright case STATEVAR_PRIMITIVE: 1553636f6a4SJames Wright problem->ics.qfunction = ICsNewtonianIG_Prim; 1563636f6a4SJames Wright problem->ics.qfunction_loc = ICsNewtonianIG_Prim_loc; 1573636f6a4SJames Wright problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Prim; 1583636f6a4SJames Wright problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Prim_loc; 1593636f6a4SJames Wright problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Prim; 1603636f6a4SJames Wright problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Prim_loc; 1613636f6a4SJames Wright problem->apply_inflow.qfunction = BoundaryIntegral_Prim; 1623636f6a4SJames Wright problem->apply_inflow.qfunction_loc = BoundaryIntegral_Prim_loc; 1633636f6a4SJames Wright problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Prim; 1643636f6a4SJames Wright problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc; 1653636f6a4SJames Wright break; 166cbe60e31SLeila Ghaffari } 1671485969bSJeremy L Thompson 1683a8779fbSJames Wright // -- Physics 1692b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL)); 1702b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL)); 1712b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL)); 1722b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL)); 1732b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL)); 1743a8779fbSJames Wright 175bb8a0c61SJames Wright PetscInt dim = problem->dim; 176d9bb1cdbSJames Wright PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL)); 177d9bb1cdbSJames Wright PetscCall(PetscOptionsRealArray("-g", "Gravitational acceleration vector", NULL, g, &dim, &given_option)); 178d9bb1cdbSJames Wright dim = problem->dim; 179d9bb1cdbSJames Wright PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option)); 180d9bb1cdbSJames Wright if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim); 181d9bb1cdbSJames Wright 1822b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 1832b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 1842b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 1852b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL)); 1862b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL)); 1872b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL)); 1882b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL)); 1892b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 1902b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL)); 1913a8779fbSJames Wright 192db95602eSJames Wright dim = 3; 193b8fb7609SAdeleke O. Bankole PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL)); 194b8fb7609SAdeleke O. Bankole PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL)); 195b8fb7609SAdeleke O. Bankole PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL)); 196b8fb7609SAdeleke O. Bankole 1973a8779fbSJames Wright // -- Units 1982b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 1993a8779fbSJames Wright meter = fabs(meter); 2002b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 2013a8779fbSJames Wright kilogram = fabs(kilogram); 2022b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 2033a8779fbSJames Wright second = fabs(second); 2042b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL)); 2053a8779fbSJames Wright Kelvin = fabs(Kelvin); 2063a8779fbSJames Wright 2073a8779fbSJames Wright // -- Warnings 2083a8779fbSJames Wright if (stab == STAB_SUPG && !implicit) { 2092b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n")); 2103a8779fbSJames Wright } 2115d28dccaSJames Wright PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP, 2125d28dccaSJames Wright "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n"); 213e7754af5SKenneth E. Jansen 214e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point", 215e7754af5SKenneth E. Jansen NULL, idl_decay_time, &idl_decay_time, &idl_enable)); 216e7754af5SKenneth E. Jansen if (idl_enable && idl_decay_time == 0) SETERRQ(comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero."); 217e7754af5SKenneth E. Jansen else if (idl_decay_time < 0) idl_enable = PETSC_FALSE; 218e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL)); 219e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL)); 2201485969bSJeremy L Thompson PetscOptionsEnd(); 2213a8779fbSJames Wright 2223a8779fbSJames Wright // ------------------------------------------------------ 2233a8779fbSJames Wright // Set up the PETSc context 2243a8779fbSJames Wright // ------------------------------------------------------ 2253a8779fbSJames Wright // -- Define derived units 2263a8779fbSJames Wright Pascal = kilogram / (meter * PetscSqr(second)); 2273a8779fbSJames Wright J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 2283a8779fbSJames Wright m_per_squared_s = meter / PetscSqr(second); 2293a8779fbSJames Wright W_per_m_K = kilogram * meter / (pow(second, 3) * Kelvin); 2303a8779fbSJames Wright 2313a8779fbSJames Wright user->units->meter = meter; 2323a8779fbSJames Wright user->units->kilogram = kilogram; 2333a8779fbSJames Wright user->units->second = second; 2343a8779fbSJames Wright user->units->Kelvin = Kelvin; 2353a8779fbSJames Wright user->units->Pascal = Pascal; 2363a8779fbSJames Wright user->units->J_per_kg_K = J_per_kg_K; 2373a8779fbSJames Wright user->units->m_per_squared_s = m_per_squared_s; 2383a8779fbSJames Wright user->units->W_per_m_K = W_per_m_K; 2393a8779fbSJames Wright 2403a8779fbSJames Wright // ------------------------------------------------------ 2413a8779fbSJames Wright // Set up the libCEED context 2423a8779fbSJames Wright // ------------------------------------------------------ 2433a8779fbSJames Wright // -- Scale variables to desired units 2443a8779fbSJames Wright cv *= J_per_kg_K; 2453a8779fbSJames Wright cp *= J_per_kg_K; 2463a8779fbSJames Wright mu *= Pascal * second; 2473a8779fbSJames Wright k *= W_per_m_K; 248493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter; 249493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s; 250b8fb7609SAdeleke O. Bankole reference.pressure *= Pascal; 251b8fb7609SAdeleke O. Bankole for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second; 252b8fb7609SAdeleke O. Bankole reference.temperature *= Kelvin; 2533a8779fbSJames Wright problem->dm_scale = meter; 2543a8779fbSJames Wright 2553a8779fbSJames Wright // -- Solver Settings 2563a8779fbSJames Wright user->phys->stab = stab; 2573a8779fbSJames Wright user->phys->implicit = implicit; 2583636f6a4SJames Wright user->phys->state_var = state_var; 2593a8779fbSJames Wright user->phys->has_curr_time = has_curr_time; 2603a8779fbSJames Wright 2613a8779fbSJames Wright // -- QFunction Context 26215a3537eSJed Brown newtonian_ig_ctx->lambda = lambda; 26315a3537eSJed Brown newtonian_ig_ctx->mu = mu; 26415a3537eSJed Brown newtonian_ig_ctx->k = k; 26515a3537eSJed Brown newtonian_ig_ctx->cv = cv; 26615a3537eSJed Brown newtonian_ig_ctx->cp = cp; 26715a3537eSJed Brown newtonian_ig_ctx->c_tau = c_tau; 26815a3537eSJed Brown newtonian_ig_ctx->Ctau_t = Ctau_t; 26915a3537eSJed Brown newtonian_ig_ctx->Ctau_v = Ctau_v; 27015a3537eSJed Brown newtonian_ig_ctx->Ctau_C = Ctau_C; 27115a3537eSJed Brown newtonian_ig_ctx->Ctau_M = Ctau_M; 27215a3537eSJed Brown newtonian_ig_ctx->Ctau_E = Ctau_E; 273e7754af5SKenneth E. Jansen newtonian_ig_ctx->P0 = reference.pressure; 27415a3537eSJed Brown newtonian_ig_ctx->stabilization = stab; 275c8c30d87SJed Brown newtonian_ig_ctx->P0 = reference.pressure; 2768085925cSJames Wright newtonian_ig_ctx->is_implicit = implicit; 2773636f6a4SJames Wright newtonian_ig_ctx->state_var = state_var; 278e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_enable = idl_enable; 279e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * second); 280e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_start = idl_start * meter; 281e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_length = idl_length * meter; 2822b916ea7SJeremy L Thompson PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3)); 2833a8779fbSJames Wright 284b8fb7609SAdeleke O. Bankole // -- Setup Context 285b8fb7609SAdeleke O. Bankole setup_context->reference = reference; 286b8fb7609SAdeleke O. Bankole setup_context->gas = *newtonian_ig_ctx; 287b8fb7609SAdeleke O. Bankole setup_context->lx = domain_size[0]; 288b8fb7609SAdeleke O. Bankole setup_context->ly = domain_size[1]; 289b8fb7609SAdeleke O. Bankole setup_context->lz = domain_size[2]; 290b8fb7609SAdeleke O. Bankole setup_context->time = 0; 291b8fb7609SAdeleke O. Bankole 292b8fb7609SAdeleke O. Bankole if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 293c8c30d87SJed Brown if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 2940bc3fea4SJames Wright 295b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context)); 296b4c37c5cSJames Wright PetscCallCeed(ceed, 297b4c37c5cSJames Wright CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 298b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc)); 299b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, "evaluation time", offsetof(struct SetupContext_, time), 1, 300b4c37c5cSJames Wright "Time of evaluation")); 30115a3537eSJed Brown 302b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context)); 303b4c37c5cSJames Wright PetscCallCeed(ceed, 304b4c37c5cSJames Wright CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx)); 305b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, FreeContextPetsc)); 306b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1, 307b4c37c5cSJames Wright "Size of timestep, delta t")); 308b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift", 309b4c37c5cSJames Wright offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1, 310b4c37c5cSJames Wright "Shift for mass matrix in IJacobian")); 311b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1, 312b4c37c5cSJames Wright "Current solution time")); 313a8bca8acSJames Wright 31415a3537eSJed Brown problem->apply_vol_rhs.qfunction_context = newtonian_ig_context; 315b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ifunction.qfunction_context)); 316b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ijacobian.qfunction_context)); 317b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow.qfunction_context)); 318b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow_jacobian.qfunction_context)); 319f0b65372SJed Brown 320f0b65372SJed Brown if (unit_tests) { 321f0b65372SJed Brown PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx)); 322f0b65372SJed Brown } 323d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3243a8779fbSJames Wright } 3253a8779fbSJames Wright 3262d49c0afSJames Wright PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData *problem, AppCtx app_ctx) { 3272d49c0afSJames Wright MPI_Comm comm = user->comm; 328b4c37c5cSJames Wright Ceed ceed = user->ceed; 32915a3537eSJed Brown NewtonianIdealGasContext newtonian_ctx; 33015a3537eSJed Brown 3313a8779fbSJames Wright PetscFunctionBeginUser; 332b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ctx)); 3332b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 33415a3537eSJed Brown " Problem:\n" 33515a3537eSJed Brown " Problem Name : %s\n" 33615a3537eSJed Brown " Stabilization : %s\n", 3372b916ea7SJeremy L Thompson app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization])); 338b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ctx)); 339d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3403a8779fbSJames Wright } 341