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 16f7534ca8SJed Brown // For use with PetscOptionsEnum 17f7534ca8SJed Brown static const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "StateVariable", "STATEVAR_", NULL}; 18f7534ca8SJed Brown 192b916ea7SJeremy L Thompson // Compute relative error |a - b|/|s| 202b916ea7SJeremy L Thompson static PetscErrorCode CheckPrimitiveWithTolerance(StatePrimitive sY, StatePrimitive aY, StatePrimitive bY, const char *name, PetscReal rtol_pressure, 212b916ea7SJeremy L Thompson PetscReal rtol_velocity, PetscReal rtol_temperature) { 22f0b65372SJed Brown PetscFunctionBeginUser; 23f0b65372SJed Brown StatePrimitive eY; // relative error 24f0b65372SJed Brown eY.pressure = (aY.pressure - bY.pressure) / sY.pressure; 252b916ea7SJeremy L Thompson PetscScalar u = sqrt(Square(sY.velocity[0]) + Square(sY.velocity[1]) + Square(sY.velocity[2])); 26f0b65372SJed Brown for (int j = 0; j < 3; j++) eY.velocity[j] = (aY.velocity[j] - bY.velocity[j]) / u; 27f0b65372SJed Brown eY.temperature = (aY.temperature - bY.temperature) / sY.temperature; 282b916ea7SJeremy L Thompson if (fabs(eY.pressure) > rtol_pressure) printf("%s: pressure error %g\n", name, eY.pressure); 292b916ea7SJeremy L Thompson for (int j = 0; j < 3; j++) { 302b916ea7SJeremy L Thompson if (fabs(eY.velocity[j]) > rtol_velocity) printf("%s: velocity[%d] error %g\n", name, j, eY.velocity[j]); 312b916ea7SJeremy L Thompson } 322b916ea7SJeremy L Thompson if (fabs(eY.temperature) > rtol_temperature) printf("%s: temperature error %g\n", name, eY.temperature); 33f0b65372SJed Brown PetscFunctionReturn(0); 34f0b65372SJed Brown } 35f0b65372SJed Brown 362b916ea7SJeremy L Thompson static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) { 37f0b65372SJed Brown Units units = user->units; 38f0b65372SJed Brown const CeedScalar eps = 1e-6; 392b916ea7SJeremy L Thompson const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, Pascal = units->Pascal; 40f0b65372SJed Brown PetscFunctionBeginUser; 412b916ea7SJeremy L Thompson const CeedScalar rho = 1.2 * kg / (m * m * m), u = 40 * m / sec; 42f0b65372SJed Brown CeedScalar U[5] = {rho, rho * u, rho * u * 1.1, rho * u * 1.2, 250e3 * Pascal + .5 * rho * u * u}; 43f0b65372SJed Brown const CeedScalar x[3] = {.1, .2, .3}; 44f0b65372SJed Brown State s = StateFromU(gas, U, x); 45f0b65372SJed Brown for (int i = 0; i < 8; i++) { 46f0b65372SJed Brown CeedScalar dU[5] = {0}, dx[3] = {0}; 47f0b65372SJed Brown if (i < 5) dU[i] = U[i]; 48f0b65372SJed Brown else dx[i - 5] = x[i - 5]; 49f0b65372SJed Brown State ds = StateFromU_fwd(gas, s, dU, x, dx); 50f0b65372SJed Brown for (int j = 0; j < 5; j++) dU[j] = (1 + eps * (i == j)) * U[j]; 51f0b65372SJed Brown for (int j = 0; j < 3; j++) dx[j] = (1 + eps * (i == 5 + j)) * x[j]; 52f0b65372SJed Brown State t = StateFromU(gas, dU, dx); 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 } 61f0b65372SJed Brown PetscFunctionReturn(0); 62f0b65372SJed Brown } 63f0b65372SJed Brown 64d1c51a42SJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) { 65b7f03f12SJed Brown SetupContext setup_context; 663a8779fbSJames Wright User user = *(User *)ctx; 674351d27fSLeila Ghaffari CeedInt degree = user->app_ctx->degree; 683a8779fbSJames Wright StabilizationType stab; 693636f6a4SJames Wright StateVariable state_var; 703a8779fbSJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 713a8779fbSJames Wright PetscBool implicit; 723636f6a4SJames Wright PetscBool has_curr_time = PETSC_FALSE, unit_tests; 7315a3537eSJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 7415a3537eSJed Brown CeedQFunctionContext newtonian_ig_context; 753a8779fbSJames Wright 7615a3537eSJed Brown PetscFunctionBeginUser; 772b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &setup_context)); 782b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &newtonian_ig_ctx)); 793a8779fbSJames Wright 803a8779fbSJames Wright // ------------------------------------------------------ 813a8779fbSJames Wright // Setup Generic Newtonian IG Problem 823a8779fbSJames Wright // ------------------------------------------------------ 833a8779fbSJames Wright problem->dim = 3; 843a8779fbSJames Wright problem->q_data_size_vol = 10; 85493642f1SJames Wright problem->q_data_size_sur = 10; 86b01ba163SJames Wright problem->jac_data_size_sur = 11; 879785fe93SJed Brown problem->setup_vol.qfunction = Setup; 889785fe93SJed Brown problem->setup_vol.qfunction_loc = Setup_loc; 899785fe93SJed Brown problem->setup_sur.qfunction = SetupBoundary; 909785fe93SJed Brown problem->setup_sur.qfunction_loc = SetupBoundary_loc; 91b7f03f12SJed Brown problem->bc = NULL; 92b7f03f12SJed Brown problem->bc_ctx = setup_context; 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) 101bb8a0c61SJames Wright CeedScalar g[3] = {0, 0, -9.81}; // 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}; 108*c176988bSLeila 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 1163a8779fbSJames Wright // ------------------------------------------------------ 1173a8779fbSJames Wright // Create the PETSc context 1183a8779fbSJames Wright // ------------------------------------------------------ 119bb8a0c61SJames Wright PetscScalar meter = 1; // 1 meter in scaled length units 120bb8a0c61SJames Wright PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 121bb8a0c61SJames Wright PetscScalar second = 1; // 1 second in scaled time units 1223a8779fbSJames Wright PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 1233a8779fbSJames Wright PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 1243a8779fbSJames Wright 1253a8779fbSJames Wright // ------------------------------------------------------ 1263a8779fbSJames Wright // Command line Options 1273a8779fbSJames Wright // ------------------------------------------------------ 1282b916ea7SJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL); 129cbe60e31SLeila Ghaffari // -- Conservative vs Primitive variables 1302b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE), 1312b916ea7SJeremy L Thompson (PetscEnum *)&state_var, NULL)); 1323636f6a4SJames Wright 1333636f6a4SJames Wright switch (state_var) { 1343636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 135cbe60e31SLeila Ghaffari problem->ics.qfunction = ICsNewtonianIG; 136cbe60e31SLeila Ghaffari problem->ics.qfunction_loc = ICsNewtonianIG_loc; 137cbe60e31SLeila Ghaffari problem->apply_vol_rhs.qfunction = RHSFunction_Newtonian; 138cbe60e31SLeila Ghaffari problem->apply_vol_rhs.qfunction_loc = RHSFunction_Newtonian_loc; 13976555becSJames Wright problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Conserv; 14076555becSJames Wright problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Conserv_loc; 14176555becSJames Wright problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Conserv; 14276555becSJames Wright problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Conserv_loc; 143d4559bbeSJames Wright problem->apply_inflow.qfunction = BoundaryIntegral_Conserv; 144d4559bbeSJames Wright problem->apply_inflow.qfunction_loc = BoundaryIntegral_Conserv_loc; 145d4559bbeSJames Wright problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Conserv; 146d4559bbeSJames Wright problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc; 147d4559bbeSJames Wright problem->apply_outflow.qfunction = PressureOutflow_Conserv; 148d4559bbeSJames Wright problem->apply_outflow.qfunction_loc = PressureOutflow_Conserv_loc; 149d4559bbeSJames Wright problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian_Conserv; 150d4559bbeSJames Wright problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Conserv_loc; 1513636f6a4SJames Wright break; 1523636f6a4SJames Wright 1533636f6a4SJames Wright case STATEVAR_PRIMITIVE: 1543636f6a4SJames Wright problem->ics.qfunction = ICsNewtonianIG_Prim; 1553636f6a4SJames Wright problem->ics.qfunction_loc = ICsNewtonianIG_Prim_loc; 1563636f6a4SJames Wright problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Prim; 1573636f6a4SJames Wright problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Prim_loc; 1583636f6a4SJames Wright problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Prim; 1593636f6a4SJames Wright problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Prim_loc; 1603636f6a4SJames Wright problem->apply_inflow.qfunction = BoundaryIntegral_Prim; 1613636f6a4SJames Wright problem->apply_inflow.qfunction_loc = BoundaryIntegral_Prim_loc; 1623636f6a4SJames Wright problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Prim; 1633636f6a4SJames Wright problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc; 1643636f6a4SJames Wright problem->apply_outflow.qfunction = PressureOutflow_Prim; 1653636f6a4SJames Wright problem->apply_outflow.qfunction_loc = PressureOutflow_Prim_loc; 1663636f6a4SJames Wright problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian_Prim; 1673636f6a4SJames Wright problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Prim_loc; 1683636f6a4SJames Wright break; 169cbe60e31SLeila Ghaffari } 1701485969bSJeremy L Thompson 1713a8779fbSJames Wright // -- Physics 1722b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL)); 1732b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL)); 1742b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL)); 1752b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL)); 1762b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL)); 1773a8779fbSJames Wright 178bb8a0c61SJames Wright PetscInt dim = problem->dim; 1792b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-g", "Gravitational acceleration", NULL, g, &dim, NULL)); 1802b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 1812b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 1822b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 1832b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL)); 1842b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL)); 1852b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL)); 1862b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL)); 1872b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 1882b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL)); 1893a8779fbSJames Wright 1903a8779fbSJames Wright // -- Units 1912b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 1923a8779fbSJames Wright meter = fabs(meter); 1932b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 1943a8779fbSJames Wright kilogram = fabs(kilogram); 1952b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 1963a8779fbSJames Wright second = fabs(second); 1972b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL)); 1983a8779fbSJames Wright Kelvin = fabs(Kelvin); 1993a8779fbSJames Wright 2003a8779fbSJames Wright // -- Warnings 2013a8779fbSJames Wright if (stab == STAB_SUPG && !implicit) { 2022b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n")); 2033a8779fbSJames Wright } 2043636f6a4SJames Wright if (state_var == STATEVAR_PRIMITIVE && !implicit) { 2052b916ea7SJeremy L Thompson SETERRQ(comm, PETSC_ERR_ARG_NULL, "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n"); 206cbe60e31SLeila Ghaffari } 2071485969bSJeremy L Thompson PetscOptionsEnd(); 2083a8779fbSJames Wright 2093a8779fbSJames Wright // ------------------------------------------------------ 2103a8779fbSJames Wright // Set up the PETSc context 2113a8779fbSJames Wright // ------------------------------------------------------ 2123a8779fbSJames Wright // -- Define derived units 2133a8779fbSJames Wright Pascal = kilogram / (meter * PetscSqr(second)); 2143a8779fbSJames Wright J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 2153a8779fbSJames Wright m_per_squared_s = meter / PetscSqr(second); 2163a8779fbSJames Wright W_per_m_K = kilogram * meter / (pow(second, 3) * Kelvin); 2173a8779fbSJames Wright 2183a8779fbSJames Wright user->units->meter = meter; 2193a8779fbSJames Wright user->units->kilogram = kilogram; 2203a8779fbSJames Wright user->units->second = second; 2213a8779fbSJames Wright user->units->Kelvin = Kelvin; 2223a8779fbSJames Wright user->units->Pascal = Pascal; 2233a8779fbSJames Wright user->units->J_per_kg_K = J_per_kg_K; 2243a8779fbSJames Wright user->units->m_per_squared_s = m_per_squared_s; 2253a8779fbSJames Wright user->units->W_per_m_K = W_per_m_K; 2263a8779fbSJames Wright 2273a8779fbSJames Wright // ------------------------------------------------------ 2283a8779fbSJames Wright // Set up the libCEED context 2293a8779fbSJames Wright // ------------------------------------------------------ 2303a8779fbSJames Wright // -- Scale variables to desired units 2313a8779fbSJames Wright cv *= J_per_kg_K; 2323a8779fbSJames Wright cp *= J_per_kg_K; 2333a8779fbSJames Wright mu *= Pascal * second; 2343a8779fbSJames Wright k *= W_per_m_K; 235493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter; 236493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s; 2373a8779fbSJames Wright problem->dm_scale = meter; 2383a8779fbSJames Wright 2393a8779fbSJames Wright // -- Setup Context 2403a8779fbSJames Wright setup_context->cv = cv; 2413a8779fbSJames Wright setup_context->cp = cp; 2423a8779fbSJames Wright setup_context->lx = domain_size[0]; 2433a8779fbSJames Wright setup_context->ly = domain_size[1]; 2443a8779fbSJames Wright setup_context->lz = domain_size[2]; 2453a8779fbSJames Wright setup_context->time = 0; 2462b916ea7SJeremy L Thompson PetscCall(PetscArraycpy(setup_context->g, g, 3)); 2473a8779fbSJames Wright 2483a8779fbSJames Wright // -- Solver Settings 2493a8779fbSJames Wright user->phys->stab = stab; 2503a8779fbSJames Wright user->phys->implicit = implicit; 2513636f6a4SJames Wright user->phys->state_var = state_var; 2523a8779fbSJames Wright user->phys->has_curr_time = has_curr_time; 2533a8779fbSJames Wright 2543a8779fbSJames Wright // -- QFunction Context 25515a3537eSJed Brown newtonian_ig_ctx->lambda = lambda; 25615a3537eSJed Brown newtonian_ig_ctx->mu = mu; 25715a3537eSJed Brown newtonian_ig_ctx->k = k; 25815a3537eSJed Brown newtonian_ig_ctx->cv = cv; 25915a3537eSJed Brown newtonian_ig_ctx->cp = cp; 26015a3537eSJed Brown newtonian_ig_ctx->c_tau = c_tau; 26115a3537eSJed Brown newtonian_ig_ctx->Ctau_t = Ctau_t; 26215a3537eSJed Brown newtonian_ig_ctx->Ctau_v = Ctau_v; 26315a3537eSJed Brown newtonian_ig_ctx->Ctau_C = Ctau_C; 26415a3537eSJed Brown newtonian_ig_ctx->Ctau_M = Ctau_M; 26515a3537eSJed Brown newtonian_ig_ctx->Ctau_E = Ctau_E; 26615a3537eSJed Brown newtonian_ig_ctx->stabilization = stab; 2678085925cSJames Wright newtonian_ig_ctx->is_implicit = implicit; 2683636f6a4SJames Wright newtonian_ig_ctx->state_var = state_var; 2692b916ea7SJeremy L Thompson PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3)); 2703a8779fbSJames Wright 2713e425e90SJames Wright if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx)); 2720bc3fea4SJames Wright 27315a3537eSJed Brown CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context); 2742b916ea7SJeremy L Thompson CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context); 2752b916ea7SJeremy L Thompson CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc); 2762b916ea7SJeremy L Thompson CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, "evaluation time", (char *)&setup_context->time - (char *)setup_context, 1, 2772b916ea7SJeremy L Thompson "Time of evaluation"); 27815a3537eSJed Brown 27915a3537eSJed Brown CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context); 2802b916ea7SJeremy L Thompson CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx); 2812b916ea7SJeremy L Thompson CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, FreeContextPetsc); 2822b916ea7SJeremy L Thompson CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1, 2832b916ea7SJeremy L Thompson "Size of timestep, delta t"); 2842b916ea7SJeremy L Thompson CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift", offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 2852b916ea7SJeremy L Thompson 1, "Shift for mass matrix in IJacobian"); 28615a3537eSJed Brown problem->apply_vol_rhs.qfunction_context = newtonian_ig_context; 2872b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ifunction.qfunction_context); 2882b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ijacobian.qfunction_context); 2892b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow.qfunction_context); 2902b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow_jacobian.qfunction_context); 2912b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_outflow.qfunction_context); 2922b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_outflow_jacobian.qfunction_context); 293f0b65372SJed Brown 294f0b65372SJed Brown if (unit_tests) { 295f0b65372SJed Brown PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx)); 296f0b65372SJed Brown } 2973a8779fbSJames Wright PetscFunctionReturn(0); 2983a8779fbSJames Wright } 2993a8779fbSJames Wright 300cbe60e31SLeila Ghaffari PetscErrorCode PRINT_NEWTONIAN(ProblemData *problem, AppCtx app_ctx) { 30115a3537eSJed Brown MPI_Comm comm = PETSC_COMM_WORLD; 30215a3537eSJed Brown NewtonianIdealGasContext newtonian_ctx; 30315a3537eSJed Brown 3043a8779fbSJames Wright PetscFunctionBeginUser; 3052b916ea7SJeremy L Thompson CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ctx); 3062b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 30715a3537eSJed Brown " Problem:\n" 30815a3537eSJed Brown " Problem Name : %s\n" 30915a3537eSJed Brown " Stabilization : %s\n", 3102b916ea7SJeremy L Thompson app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization])); 3112b916ea7SJeremy L Thompson CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ctx); 3123a8779fbSJames Wright PetscFunctionReturn(0); 3133a8779fbSJames Wright } 314