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 132b916ea7SJeremy L Thompson #include "../navierstokes.h" 142b916ea7SJeremy L Thompson #include "../qfunctions/setupgeo.h" 156dd99bedSLeila Ghaffari 162b916ea7SJeremy L Thompson // Compute relative error |a - b|/|s| 172b916ea7SJeremy L Thompson static PetscErrorCode CheckPrimitiveWithTolerance(StatePrimitive sY, StatePrimitive aY, StatePrimitive bY, const char *name, PetscReal rtol_pressure, 182b916ea7SJeremy L Thompson PetscReal rtol_velocity, PetscReal rtol_temperature) { 19f0b65372SJed Brown PetscFunctionBeginUser; 20f0b65372SJed Brown StatePrimitive eY; // relative error 21f0b65372SJed Brown eY.pressure = (aY.pressure - bY.pressure) / sY.pressure; 222b916ea7SJeremy L Thompson PetscScalar u = sqrt(Square(sY.velocity[0]) + Square(sY.velocity[1]) + Square(sY.velocity[2])); 23f0b65372SJed Brown for (int j = 0; j < 3; j++) eY.velocity[j] = (aY.velocity[j] - bY.velocity[j]) / u; 24f0b65372SJed Brown eY.temperature = (aY.temperature - bY.temperature) / sY.temperature; 252b916ea7SJeremy L Thompson if (fabs(eY.pressure) > rtol_pressure) printf("%s: pressure error %g\n", name, eY.pressure); 262b916ea7SJeremy L Thompson for (int j = 0; j < 3; j++) { 272b916ea7SJeremy L Thompson if (fabs(eY.velocity[j]) > rtol_velocity) printf("%s: velocity[%d] error %g\n", name, j, eY.velocity[j]); 282b916ea7SJeremy L Thompson } 292b916ea7SJeremy L Thompson if (fabs(eY.temperature) > rtol_temperature) printf("%s: temperature error %g\n", name, eY.temperature); 30f0b65372SJed Brown PetscFunctionReturn(0); 31f0b65372SJed Brown } 32f0b65372SJed Brown 332b916ea7SJeremy L Thompson static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) { 34f0b65372SJed Brown Units units = user->units; 35f0b65372SJed Brown const CeedScalar eps = 1e-6; 362b916ea7SJeremy L Thompson const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, Pascal = units->Pascal; 37f0b65372SJed Brown PetscFunctionBeginUser; 382b916ea7SJeremy L Thompson const CeedScalar rho = 1.2 * kg / (m * m * m), u = 40 * m / sec; 39f0b65372SJed Brown CeedScalar U[5] = {rho, rho * u, rho * u * 1.1, rho * u * 1.2, 250e3 * Pascal + .5 * rho * u * u}; 40f0b65372SJed Brown const CeedScalar x[3] = {.1, .2, .3}; 41f0b65372SJed Brown State s = StateFromU(gas, U, x); 42f0b65372SJed Brown for (int i = 0; i < 8; i++) { 43f0b65372SJed Brown CeedScalar dU[5] = {0}, dx[3] = {0}; 44f0b65372SJed Brown if (i < 5) dU[i] = U[i]; 45f0b65372SJed Brown else dx[i - 5] = x[i - 5]; 46f0b65372SJed Brown State ds = StateFromU_fwd(gas, s, dU, x, dx); 47f0b65372SJed Brown for (int j = 0; j < 5; j++) dU[j] = (1 + eps * (i == j)) * U[j]; 48f0b65372SJed Brown for (int j = 0; j < 3; j++) dx[j] = (1 + eps * (i == 5 + j)) * x[j]; 49f0b65372SJed Brown State t = StateFromU(gas, dU, dx); 50f0b65372SJed Brown StatePrimitive dY; 51f0b65372SJed Brown dY.pressure = (t.Y.pressure - s.Y.pressure) / eps; 522b916ea7SJeremy L Thompson for (int j = 0; j < 3; j++) dY.velocity[j] = (t.Y.velocity[j] - s.Y.velocity[j]) / eps; 53f0b65372SJed Brown dY.temperature = (t.Y.temperature - s.Y.temperature) / eps; 54f0b65372SJed Brown char buf[128]; 55f0b65372SJed Brown snprintf(buf, sizeof buf, "StateFromU_fwd i=%d", i); 56f0b65372SJed Brown PetscCall(CheckPrimitiveWithTolerance(dY, ds.Y, dY, buf, 5e-6, 1e-6, 1e-6)); 57f0b65372SJed Brown } 58f0b65372SJed Brown PetscFunctionReturn(0); 59f0b65372SJed Brown } 60f0b65372SJed Brown 61b7f03f12SJed Brown PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx) { 62b7f03f12SJed Brown SetupContext setup_context; 633a8779fbSJames Wright User user = *(User *)ctx; 643a8779fbSJames Wright StabilizationType stab; 653636f6a4SJames Wright StateVariable state_var; 663a8779fbSJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 673a8779fbSJames Wright PetscBool implicit; 683636f6a4SJames Wright PetscBool has_curr_time = PETSC_FALSE, unit_tests; 6915a3537eSJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 7015a3537eSJed Brown CeedQFunctionContext newtonian_ig_context; 713a8779fbSJames Wright 7215a3537eSJed Brown PetscFunctionBeginUser; 732b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &setup_context)); 742b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &newtonian_ig_ctx)); 753a8779fbSJames Wright 763a8779fbSJames Wright // ------------------------------------------------------ 773a8779fbSJames Wright // Setup Generic Newtonian IG Problem 783a8779fbSJames Wright // ------------------------------------------------------ 793a8779fbSJames Wright problem->dim = 3; 803a8779fbSJames Wright problem->q_data_size_vol = 10; 81493642f1SJames Wright problem->q_data_size_sur = 10; 82b01ba163SJames Wright problem->jac_data_size_sur = 11; 839785fe93SJed Brown problem->setup_vol.qfunction = Setup; 849785fe93SJed Brown problem->setup_vol.qfunction_loc = Setup_loc; 859785fe93SJed Brown problem->setup_sur.qfunction = SetupBoundary; 869785fe93SJed Brown problem->setup_sur.qfunction_loc = SetupBoundary_loc; 87b7f03f12SJed Brown problem->bc = NULL; 88b7f03f12SJed Brown problem->bc_ctx = setup_context; 893a8779fbSJames Wright problem->non_zero_time = PETSC_FALSE; 90cbe60e31SLeila Ghaffari problem->print_info = PRINT_NEWTONIAN; 913a8779fbSJames Wright 923a8779fbSJames Wright // ------------------------------------------------------ 933a8779fbSJames Wright // Create the libCEED context 943a8779fbSJames Wright // ------------------------------------------------------ 953a8779fbSJames Wright CeedScalar cv = 717.; // J/(kg K) 963a8779fbSJames Wright CeedScalar cp = 1004.; // J/(kg K) 97bb8a0c61SJames Wright CeedScalar g[3] = {0, 0, -9.81}; // m/s^2 983a8779fbSJames Wright CeedScalar lambda = -2. / 3.; // - 99bb8a0c61SJames Wright CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 1003a8779fbSJames Wright CeedScalar k = 0.02638; // W/(m K) 1013a8779fbSJames Wright CeedScalar c_tau = 0.5; // - 102bb8a0c61SJames Wright CeedScalar Ctau_t = 1.0; // - 103bb8a0c61SJames Wright CeedScalar Ctau_v = 36.0; // TODO make function of degree 104bb8a0c61SJames Wright CeedScalar Ctau_C = 1.0; // TODO make function of degree 105bb8a0c61SJames Wright CeedScalar Ctau_M = 1.0; // TODO make function of degree 106bb8a0c61SJames Wright CeedScalar Ctau_E = 1.0; // TODO make function of degree 1073a8779fbSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 1082b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 109493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 1103a8779fbSJames Wright 1113a8779fbSJames Wright // ------------------------------------------------------ 1123a8779fbSJames Wright // Create the PETSc context 1133a8779fbSJames Wright // ------------------------------------------------------ 114bb8a0c61SJames Wright PetscScalar meter = 1; // 1 meter in scaled length units 115bb8a0c61SJames Wright PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 116bb8a0c61SJames Wright PetscScalar second = 1; // 1 second in scaled time units 1173a8779fbSJames Wright PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 1183a8779fbSJames Wright PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 1193a8779fbSJames Wright 1203a8779fbSJames Wright // ------------------------------------------------------ 1213a8779fbSJames Wright // Command line Options 1223a8779fbSJames Wright // ------------------------------------------------------ 1232b916ea7SJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL); 124cbe60e31SLeila Ghaffari // -- Conservative vs Primitive variables 1252b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE), 1262b916ea7SJeremy L Thompson (PetscEnum *)&state_var, NULL)); 1273636f6a4SJames Wright 128d4559bbeSJames Wright // *INDENT-OFF* 1293636f6a4SJames Wright switch (state_var) { 1303636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 131cbe60e31SLeila Ghaffari problem->ics.qfunction = ICsNewtonianIG; 132cbe60e31SLeila Ghaffari problem->ics.qfunction_loc = ICsNewtonianIG_loc; 133cbe60e31SLeila Ghaffari problem->apply_vol_rhs.qfunction = RHSFunction_Newtonian; 134cbe60e31SLeila Ghaffari problem->apply_vol_rhs.qfunction_loc = RHSFunction_Newtonian_loc; 13576555becSJames Wright problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Conserv; 13676555becSJames Wright problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Conserv_loc; 13776555becSJames Wright problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Conserv; 13876555becSJames Wright problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Conserv_loc; 139d4559bbeSJames Wright problem->apply_inflow.qfunction = BoundaryIntegral_Conserv; 140d4559bbeSJames Wright problem->apply_inflow.qfunction_loc = BoundaryIntegral_Conserv_loc; 141d4559bbeSJames Wright problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Conserv; 142d4559bbeSJames Wright problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc; 143d4559bbeSJames Wright problem->apply_outflow.qfunction = PressureOutflow_Conserv; 144d4559bbeSJames Wright problem->apply_outflow.qfunction_loc = PressureOutflow_Conserv_loc; 145d4559bbeSJames Wright problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian_Conserv; 146d4559bbeSJames Wright problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Conserv_loc; 1473636f6a4SJames Wright break; 1483636f6a4SJames Wright 1493636f6a4SJames Wright case STATEVAR_PRIMITIVE: 1503636f6a4SJames Wright problem->ics.qfunction = ICsNewtonianIG_Prim; 1513636f6a4SJames Wright problem->ics.qfunction_loc = ICsNewtonianIG_Prim_loc; 1523636f6a4SJames Wright problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Prim; 1533636f6a4SJames Wright problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Prim_loc; 1543636f6a4SJames Wright problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Prim; 1553636f6a4SJames Wright problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Prim_loc; 1563636f6a4SJames Wright problem->apply_inflow.qfunction = BoundaryIntegral_Prim; 1573636f6a4SJames Wright problem->apply_inflow.qfunction_loc = BoundaryIntegral_Prim_loc; 1583636f6a4SJames Wright problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Prim; 1593636f6a4SJames Wright problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc; 1603636f6a4SJames Wright problem->apply_outflow.qfunction = PressureOutflow_Prim; 1613636f6a4SJames Wright problem->apply_outflow.qfunction_loc = PressureOutflow_Prim_loc; 1623636f6a4SJames Wright problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian_Prim; 1633636f6a4SJames Wright problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Prim_loc; 1643636f6a4SJames Wright break; 165cbe60e31SLeila Ghaffari } 166d4559bbeSJames Wright // *INDENT-ON* 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; 1762b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-g", "Gravitational acceleration", NULL, g, &dim, NULL)); 1772b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 1782b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 1792b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 1802b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL)); 1812b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL)); 1822b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL)); 1832b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL)); 1842b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 1852b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL)); 1863a8779fbSJames Wright 1873a8779fbSJames Wright // -- Units 1882b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 1893a8779fbSJames Wright meter = fabs(meter); 1902b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 1913a8779fbSJames Wright kilogram = fabs(kilogram); 1922b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 1933a8779fbSJames Wright second = fabs(second); 1942b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL)); 1953a8779fbSJames Wright Kelvin = fabs(Kelvin); 1963a8779fbSJames Wright 1973a8779fbSJames Wright // -- Warnings 1983a8779fbSJames Wright if (stab == STAB_SUPG && !implicit) { 1992b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n")); 2003a8779fbSJames Wright } 2013636f6a4SJames Wright if (state_var == STATEVAR_PRIMITIVE && !implicit) { 2022b916ea7SJeremy L Thompson SETERRQ(comm, PETSC_ERR_ARG_NULL, "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n"); 203cbe60e31SLeila Ghaffari } 2041485969bSJeremy L Thompson PetscOptionsEnd(); 2053a8779fbSJames Wright 2063a8779fbSJames Wright // ------------------------------------------------------ 2073a8779fbSJames Wright // Set up the PETSc context 2083a8779fbSJames Wright // ------------------------------------------------------ 2093a8779fbSJames Wright // -- Define derived units 2103a8779fbSJames Wright Pascal = kilogram / (meter * PetscSqr(second)); 2113a8779fbSJames Wright J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 2123a8779fbSJames Wright m_per_squared_s = meter / PetscSqr(second); 2133a8779fbSJames Wright W_per_m_K = kilogram * meter / (pow(second, 3) * Kelvin); 2143a8779fbSJames Wright 2153a8779fbSJames Wright user->units->meter = meter; 2163a8779fbSJames Wright user->units->kilogram = kilogram; 2173a8779fbSJames Wright user->units->second = second; 2183a8779fbSJames Wright user->units->Kelvin = Kelvin; 2193a8779fbSJames Wright user->units->Pascal = Pascal; 2203a8779fbSJames Wright user->units->J_per_kg_K = J_per_kg_K; 2213a8779fbSJames Wright user->units->m_per_squared_s = m_per_squared_s; 2223a8779fbSJames Wright user->units->W_per_m_K = W_per_m_K; 2233a8779fbSJames Wright 2243a8779fbSJames Wright // ------------------------------------------------------ 2253a8779fbSJames Wright // Set up the libCEED context 2263a8779fbSJames Wright // ------------------------------------------------------ 2273a8779fbSJames Wright // -- Scale variables to desired units 2283a8779fbSJames Wright cv *= J_per_kg_K; 2293a8779fbSJames Wright cp *= J_per_kg_K; 2303a8779fbSJames Wright mu *= Pascal * second; 2313a8779fbSJames Wright k *= W_per_m_K; 232493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter; 233493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s; 2343a8779fbSJames Wright problem->dm_scale = meter; 2353a8779fbSJames Wright 2363a8779fbSJames Wright // -- Setup Context 2373a8779fbSJames Wright setup_context->cv = cv; 2383a8779fbSJames Wright setup_context->cp = cp; 2393a8779fbSJames Wright setup_context->lx = domain_size[0]; 2403a8779fbSJames Wright setup_context->ly = domain_size[1]; 2413a8779fbSJames Wright setup_context->lz = domain_size[2]; 2423a8779fbSJames Wright setup_context->time = 0; 2432b916ea7SJeremy L Thompson PetscCall(PetscArraycpy(setup_context->g, g, 3)); 2443a8779fbSJames Wright 2453a8779fbSJames Wright // -- Solver Settings 2463a8779fbSJames Wright user->phys->stab = stab; 2473a8779fbSJames Wright user->phys->implicit = implicit; 2483636f6a4SJames Wright user->phys->state_var = state_var; 2493a8779fbSJames Wright user->phys->has_curr_time = has_curr_time; 2503a8779fbSJames Wright 2513a8779fbSJames Wright // -- QFunction Context 25215a3537eSJed Brown newtonian_ig_ctx->lambda = lambda; 25315a3537eSJed Brown newtonian_ig_ctx->mu = mu; 25415a3537eSJed Brown newtonian_ig_ctx->k = k; 25515a3537eSJed Brown newtonian_ig_ctx->cv = cv; 25615a3537eSJed Brown newtonian_ig_ctx->cp = cp; 25715a3537eSJed Brown newtonian_ig_ctx->c_tau = c_tau; 25815a3537eSJed Brown newtonian_ig_ctx->Ctau_t = Ctau_t; 25915a3537eSJed Brown newtonian_ig_ctx->Ctau_v = Ctau_v; 26015a3537eSJed Brown newtonian_ig_ctx->Ctau_C = Ctau_C; 26115a3537eSJed Brown newtonian_ig_ctx->Ctau_M = Ctau_M; 26215a3537eSJed Brown newtonian_ig_ctx->Ctau_E = Ctau_E; 26315a3537eSJed Brown newtonian_ig_ctx->stabilization = stab; 2648085925cSJames Wright newtonian_ig_ctx->is_implicit = implicit; 2653636f6a4SJames Wright newtonian_ig_ctx->state_var = state_var; 2662b916ea7SJeremy L Thompson PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3)); 2673a8779fbSJames Wright 268*0bc3fea4SJames Wright PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx)); 269*0bc3fea4SJames Wright 27015a3537eSJed Brown CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context); 2712b916ea7SJeremy L Thompson CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context); 2722b916ea7SJeremy L Thompson CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc); 2732b916ea7SJeremy L Thompson CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, "evaluation time", (char *)&setup_context->time - (char *)setup_context, 1, 2742b916ea7SJeremy L Thompson "Time of evaluation"); 27515a3537eSJed Brown 27615a3537eSJed Brown CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context); 2772b916ea7SJeremy L Thompson CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx); 2782b916ea7SJeremy L Thompson CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, FreeContextPetsc); 2792b916ea7SJeremy L Thompson CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1, 2802b916ea7SJeremy L Thompson "Size of timestep, delta t"); 2812b916ea7SJeremy L Thompson CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift", offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 2822b916ea7SJeremy L Thompson 1, "Shift for mass matrix in IJacobian"); 28315a3537eSJed Brown problem->apply_vol_rhs.qfunction_context = newtonian_ig_context; 2842b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ifunction.qfunction_context); 2852b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ijacobian.qfunction_context); 2862b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow.qfunction_context); 2872b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow_jacobian.qfunction_context); 2882b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_outflow.qfunction_context); 2892b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_outflow_jacobian.qfunction_context); 290f0b65372SJed Brown 291f0b65372SJed Brown if (unit_tests) { 292f0b65372SJed Brown PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx)); 293f0b65372SJed Brown } 2943a8779fbSJames Wright PetscFunctionReturn(0); 2953a8779fbSJames Wright } 2963a8779fbSJames Wright 297cbe60e31SLeila Ghaffari PetscErrorCode PRINT_NEWTONIAN(ProblemData *problem, AppCtx app_ctx) { 29815a3537eSJed Brown MPI_Comm comm = PETSC_COMM_WORLD; 29915a3537eSJed Brown NewtonianIdealGasContext newtonian_ctx; 30015a3537eSJed Brown 3013a8779fbSJames Wright PetscFunctionBeginUser; 3022b916ea7SJeremy L Thompson CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ctx); 3032b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 30415a3537eSJed Brown " Problem:\n" 30515a3537eSJed Brown " Problem Name : %s\n" 30615a3537eSJed Brown " Stabilization : %s\n", 3072b916ea7SJeremy L Thompson app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization])); 3082b916ea7SJeremy L Thompson CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ctx); 3093a8779fbSJames Wright PetscFunctionReturn(0); 3103a8779fbSJames Wright } 311