13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 388b783a1SJames Wright // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 588b783a1SJames Wright // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 788b783a1SJames Wright 888b783a1SJames Wright /// @file 988b783a1SJames Wright /// Utility functions for setting up problems using the Newtonian Qfunction 1088b783a1SJames Wright 1188b783a1SJames Wright #include "../qfunctions/newtonian.h" 1288b783a1SJames Wright 132b730f8bSJeremy L Thompson #include "../navierstokes.h" 142b730f8bSJeremy L Thompson #include "../qfunctions/setupgeo.h" 156ef2784eSLeila Ghaffari 16d64246eaSJed Brown // For use with PetscOptionsEnum 17d64246eaSJed Brown static const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "StateVariable", "STATEVAR_", NULL}; 18d64246eaSJed Brown 192b730f8bSJeremy L Thompson // Compute relative error |a - b|/|s| 202b730f8bSJeremy L Thompson static PetscErrorCode CheckPrimitiveWithTolerance(StatePrimitive sY, StatePrimitive aY, StatePrimitive bY, const char *name, PetscReal rtol_pressure, 212b730f8bSJeremy L Thompson PetscReal rtol_velocity, PetscReal rtol_temperature) { 22e334ad8fSJed Brown PetscFunctionBeginUser; 23e334ad8fSJed Brown StatePrimitive eY; // relative error 24e334ad8fSJed Brown eY.pressure = (aY.pressure - bY.pressure) / sY.pressure; 252b730f8bSJeremy L Thompson PetscScalar u = sqrt(Square(sY.velocity[0]) + Square(sY.velocity[1]) + Square(sY.velocity[2])); 26e334ad8fSJed Brown for (int j = 0; j < 3; j++) eY.velocity[j] = (aY.velocity[j] - bY.velocity[j]) / u; 27e334ad8fSJed Brown eY.temperature = (aY.temperature - bY.temperature) / sY.temperature; 282b730f8bSJeremy L Thompson if (fabs(eY.pressure) > rtol_pressure) printf("%s: pressure error %g\n", name, eY.pressure); 292b730f8bSJeremy L Thompson for (int j = 0; j < 3; j++) { 302b730f8bSJeremy L Thompson if (fabs(eY.velocity[j]) > rtol_velocity) printf("%s: velocity[%d] error %g\n", name, j, eY.velocity[j]); 312b730f8bSJeremy L Thompson } 322b730f8bSJeremy L Thompson if (fabs(eY.temperature) > rtol_temperature) printf("%s: temperature error %g\n", name, eY.temperature); 33e334ad8fSJed Brown PetscFunctionReturn(0); 34e334ad8fSJed Brown } 35e334ad8fSJed Brown 362b730f8bSJeremy L Thompson static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) { 37e334ad8fSJed Brown Units units = user->units; 38e334ad8fSJed Brown const CeedScalar eps = 1e-6; 392b730f8bSJeremy L Thompson const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, Pascal = units->Pascal; 40e334ad8fSJed Brown PetscFunctionBeginUser; 412b730f8bSJeremy L Thompson const CeedScalar rho = 1.2 * kg / (m * m * m), u = 40 * m / sec; 42e334ad8fSJed Brown CeedScalar U[5] = {rho, rho * u, rho * u * 1.1, rho * u * 1.2, 250e3 * Pascal + .5 * rho * u * u}; 43e334ad8fSJed Brown const CeedScalar x[3] = {.1, .2, .3}; 44e334ad8fSJed Brown State s = StateFromU(gas, U, x); 45e334ad8fSJed Brown for (int i = 0; i < 8; i++) { 46e334ad8fSJed Brown CeedScalar dU[5] = {0}, dx[3] = {0}; 47e334ad8fSJed Brown if (i < 5) dU[i] = U[i]; 48e334ad8fSJed Brown else dx[i - 5] = x[i - 5]; 49e334ad8fSJed Brown State ds = StateFromU_fwd(gas, s, dU, x, dx); 50e334ad8fSJed Brown for (int j = 0; j < 5; j++) dU[j] = (1 + eps * (i == j)) * U[j]; 51e334ad8fSJed Brown for (int j = 0; j < 3; j++) dx[j] = (1 + eps * (i == 5 + j)) * x[j]; 52e334ad8fSJed Brown State t = StateFromU(gas, dU, dx); 53e334ad8fSJed Brown StatePrimitive dY; 54e334ad8fSJed Brown dY.pressure = (t.Y.pressure - s.Y.pressure) / eps; 552b730f8bSJeremy L Thompson for (int j = 0; j < 3; j++) dY.velocity[j] = (t.Y.velocity[j] - s.Y.velocity[j]) / eps; 56e334ad8fSJed Brown dY.temperature = (t.Y.temperature - s.Y.temperature) / eps; 57e334ad8fSJed Brown char buf[128]; 58e334ad8fSJed Brown snprintf(buf, sizeof buf, "StateFromU_fwd i=%d", i); 59e334ad8fSJed Brown PetscCall(CheckPrimitiveWithTolerance(dY, ds.Y, dY, buf, 5e-6, 1e-6, 1e-6)); 60e334ad8fSJed Brown } 61e334ad8fSJed Brown PetscFunctionReturn(0); 62e334ad8fSJed Brown } 63e334ad8fSJed Brown 6446de7363SJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) { 65a0add3c9SJed Brown SetupContext setup_context; 6688b783a1SJames Wright User user = *(User *)ctx; 67*ebd2ea64SLeila Ghaffari CeedInt degree = user->app_ctx->degree; 6888b783a1SJames Wright StabilizationType stab; 6997baf651SJames Wright StateVariable state_var; 7088b783a1SJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 7188b783a1SJames Wright PetscBool implicit; 7297baf651SJames Wright PetscBool has_curr_time = PETSC_FALSE, unit_tests; 73841e4c73SJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 74841e4c73SJed Brown CeedQFunctionContext newtonian_ig_context; 7588b783a1SJames Wright 76841e4c73SJed Brown PetscFunctionBeginUser; 772b730f8bSJeremy L Thompson PetscCall(PetscCalloc1(1, &setup_context)); 782b730f8bSJeremy L Thompson PetscCall(PetscCalloc1(1, &newtonian_ig_ctx)); 7988b783a1SJames Wright 8088b783a1SJames Wright // ------------------------------------------------------ 8188b783a1SJames Wright // Setup Generic Newtonian IG Problem 8288b783a1SJames Wright // ------------------------------------------------------ 8388b783a1SJames Wright problem->dim = 3; 8488b783a1SJames Wright problem->q_data_size_vol = 10; 85ba6664aeSJames Wright problem->q_data_size_sur = 10; 860ec2498eSJames Wright problem->jac_data_size_sur = 11; 8791e5af17SJed Brown problem->setup_vol.qfunction = Setup; 8891e5af17SJed Brown problem->setup_vol.qfunction_loc = Setup_loc; 8991e5af17SJed Brown problem->setup_sur.qfunction = SetupBoundary; 9091e5af17SJed Brown problem->setup_sur.qfunction_loc = SetupBoundary_loc; 91a0add3c9SJed Brown problem->bc = NULL; 92a0add3c9SJed Brown problem->bc_ctx = setup_context; 9388b783a1SJames Wright problem->non_zero_time = PETSC_FALSE; 94dc805cc4SLeila Ghaffari problem->print_info = PRINT_NEWTONIAN; 9588b783a1SJames Wright 9688b783a1SJames Wright // ------------------------------------------------------ 9788b783a1SJames Wright // Create the libCEED context 9888b783a1SJames Wright // ------------------------------------------------------ 9988b783a1SJames Wright CeedScalar cv = 717.; // J/(kg K) 10088b783a1SJames Wright CeedScalar cp = 1004.; // J/(kg K) 10188626eedSJames Wright CeedScalar g[3] = {0, 0, -9.81}; // m/s^2 10288b783a1SJames Wright CeedScalar lambda = -2. / 3.; // - 10388626eedSJames Wright CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 10488b783a1SJames Wright CeedScalar k = 0.02638; // W/(m K) 105*ebd2ea64SLeila Ghaffari CeedScalar c_tau = 0.5 / degree; // - 10688626eedSJames Wright CeedScalar Ctau_t = 1.0; // - 10788626eedSJames Wright CeedScalar Ctau_v = 36.0; // TODO make function of degree 10888626eedSJames Wright CeedScalar Ctau_C = 1.0; // TODO make function of degree 10988626eedSJames Wright CeedScalar Ctau_M = 1.0; // TODO make function of degree 11088626eedSJames Wright CeedScalar Ctau_E = 1.0; // TODO make function of degree 11188b783a1SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 1122b730f8bSJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 113ba6664aeSJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 11488b783a1SJames Wright 11588b783a1SJames Wright // ------------------------------------------------------ 11688b783a1SJames Wright // Create the PETSc context 11788b783a1SJames Wright // ------------------------------------------------------ 11888626eedSJames Wright PetscScalar meter = 1; // 1 meter in scaled length units 11988626eedSJames Wright PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 12088626eedSJames Wright PetscScalar second = 1; // 1 second in scaled time units 12188b783a1SJames Wright PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 12288b783a1SJames Wright PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 12388b783a1SJames Wright 12488b783a1SJames Wright // ------------------------------------------------------ 12588b783a1SJames Wright // Command line Options 12688b783a1SJames Wright // ------------------------------------------------------ 1272b730f8bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL); 128dc805cc4SLeila Ghaffari // -- Conservative vs Primitive variables 1292b730f8bSJeremy L Thompson PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE), 1302b730f8bSJeremy L Thompson (PetscEnum *)&state_var, NULL)); 13197baf651SJames Wright 13297baf651SJames Wright switch (state_var) { 13397baf651SJames Wright case STATEVAR_CONSERVATIVE: 134dc805cc4SLeila Ghaffari problem->ics.qfunction = ICsNewtonianIG; 135dc805cc4SLeila Ghaffari problem->ics.qfunction_loc = ICsNewtonianIG_loc; 136dc805cc4SLeila Ghaffari problem->apply_vol_rhs.qfunction = RHSFunction_Newtonian; 137dc805cc4SLeila Ghaffari problem->apply_vol_rhs.qfunction_loc = RHSFunction_Newtonian_loc; 1383d02368aSJames Wright problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Conserv; 1393d02368aSJames Wright problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Conserv_loc; 1403d02368aSJames Wright problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Conserv; 1413d02368aSJames Wright problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Conserv_loc; 14220840d50SJames Wright problem->apply_inflow.qfunction = BoundaryIntegral_Conserv; 14320840d50SJames Wright problem->apply_inflow.qfunction_loc = BoundaryIntegral_Conserv_loc; 14420840d50SJames Wright problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Conserv; 14520840d50SJames Wright problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc; 14620840d50SJames Wright problem->apply_outflow.qfunction = PressureOutflow_Conserv; 14720840d50SJames Wright problem->apply_outflow.qfunction_loc = PressureOutflow_Conserv_loc; 14820840d50SJames Wright problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian_Conserv; 14920840d50SJames Wright problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Conserv_loc; 15097baf651SJames Wright break; 15197baf651SJames Wright 15297baf651SJames Wright case STATEVAR_PRIMITIVE: 15397baf651SJames Wright problem->ics.qfunction = ICsNewtonianIG_Prim; 15497baf651SJames Wright problem->ics.qfunction_loc = ICsNewtonianIG_Prim_loc; 15597baf651SJames Wright problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Prim; 15697baf651SJames Wright problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Prim_loc; 15797baf651SJames Wright problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Prim; 15897baf651SJames Wright problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Prim_loc; 15997baf651SJames Wright problem->apply_inflow.qfunction = BoundaryIntegral_Prim; 16097baf651SJames Wright problem->apply_inflow.qfunction_loc = BoundaryIntegral_Prim_loc; 16197baf651SJames Wright problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Prim; 16297baf651SJames Wright problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc; 16397baf651SJames Wright problem->apply_outflow.qfunction = PressureOutflow_Prim; 16497baf651SJames Wright problem->apply_outflow.qfunction_loc = PressureOutflow_Prim_loc; 16597baf651SJames Wright problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian_Prim; 16697baf651SJames Wright problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Prim_loc; 16797baf651SJames Wright break; 168dc805cc4SLeila Ghaffari } 16967490bc6SJeremy L Thompson 17088b783a1SJames Wright // -- Physics 1712b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL)); 1722b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL)); 1732b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL)); 1742b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL)); 1752b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL)); 17688b783a1SJames Wright 17788626eedSJames Wright PetscInt dim = problem->dim; 1782b730f8bSJeremy L Thompson PetscCall(PetscOptionsRealArray("-g", "Gravitational acceleration", NULL, g, &dim, NULL)); 1792b730f8bSJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 1802b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 1812b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 1822b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL)); 1832b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL)); 1842b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL)); 1852b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL)); 1862b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 1872b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL)); 18888b783a1SJames Wright 18988b783a1SJames Wright // -- Units 1902b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 19188b783a1SJames Wright meter = fabs(meter); 1922b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 19388b783a1SJames Wright kilogram = fabs(kilogram); 1942b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 19588b783a1SJames Wright second = fabs(second); 1962b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL)); 19788b783a1SJames Wright Kelvin = fabs(Kelvin); 19888b783a1SJames Wright 19988b783a1SJames Wright // -- Warnings 20088b783a1SJames Wright if (stab == STAB_SUPG && !implicit) { 2012b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n")); 20288b783a1SJames Wright } 20397baf651SJames Wright if (state_var == STATEVAR_PRIMITIVE && !implicit) { 2042b730f8bSJeremy L Thompson SETERRQ(comm, PETSC_ERR_ARG_NULL, "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n"); 205dc805cc4SLeila Ghaffari } 20667490bc6SJeremy L Thompson PetscOptionsEnd(); 20788b783a1SJames Wright 20888b783a1SJames Wright // ------------------------------------------------------ 20988b783a1SJames Wright // Set up the PETSc context 21088b783a1SJames Wright // ------------------------------------------------------ 21188b783a1SJames Wright // -- Define derived units 21288b783a1SJames Wright Pascal = kilogram / (meter * PetscSqr(second)); 21388b783a1SJames Wright J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 21488b783a1SJames Wright m_per_squared_s = meter / PetscSqr(second); 21588b783a1SJames Wright W_per_m_K = kilogram * meter / (pow(second, 3) * Kelvin); 21688b783a1SJames Wright 21788b783a1SJames Wright user->units->meter = meter; 21888b783a1SJames Wright user->units->kilogram = kilogram; 21988b783a1SJames Wright user->units->second = second; 22088b783a1SJames Wright user->units->Kelvin = Kelvin; 22188b783a1SJames Wright user->units->Pascal = Pascal; 22288b783a1SJames Wright user->units->J_per_kg_K = J_per_kg_K; 22388b783a1SJames Wright user->units->m_per_squared_s = m_per_squared_s; 22488b783a1SJames Wright user->units->W_per_m_K = W_per_m_K; 22588b783a1SJames Wright 22688b783a1SJames Wright // ------------------------------------------------------ 22788b783a1SJames Wright // Set up the libCEED context 22888b783a1SJames Wright // ------------------------------------------------------ 22988b783a1SJames Wright // -- Scale variables to desired units 23088b783a1SJames Wright cv *= J_per_kg_K; 23188b783a1SJames Wright cp *= J_per_kg_K; 23288b783a1SJames Wright mu *= Pascal * second; 23388b783a1SJames Wright k *= W_per_m_K; 234ba6664aeSJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter; 235ba6664aeSJames Wright for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s; 23688b783a1SJames Wright problem->dm_scale = meter; 23788b783a1SJames Wright 23888b783a1SJames Wright // -- Setup Context 23988b783a1SJames Wright setup_context->cv = cv; 24088b783a1SJames Wright setup_context->cp = cp; 24188b783a1SJames Wright setup_context->lx = domain_size[0]; 24288b783a1SJames Wright setup_context->ly = domain_size[1]; 24388b783a1SJames Wright setup_context->lz = domain_size[2]; 24488b783a1SJames Wright setup_context->time = 0; 2452b730f8bSJeremy L Thompson PetscCall(PetscArraycpy(setup_context->g, g, 3)); 24688b783a1SJames Wright 24788b783a1SJames Wright // -- Solver Settings 24888b783a1SJames Wright user->phys->stab = stab; 24988b783a1SJames Wright user->phys->implicit = implicit; 25097baf651SJames Wright user->phys->state_var = state_var; 25188b783a1SJames Wright user->phys->has_curr_time = has_curr_time; 25288b783a1SJames Wright 25388b783a1SJames Wright // -- QFunction Context 254841e4c73SJed Brown newtonian_ig_ctx->lambda = lambda; 255841e4c73SJed Brown newtonian_ig_ctx->mu = mu; 256841e4c73SJed Brown newtonian_ig_ctx->k = k; 257841e4c73SJed Brown newtonian_ig_ctx->cv = cv; 258841e4c73SJed Brown newtonian_ig_ctx->cp = cp; 259841e4c73SJed Brown newtonian_ig_ctx->c_tau = c_tau; 260841e4c73SJed Brown newtonian_ig_ctx->Ctau_t = Ctau_t; 261841e4c73SJed Brown newtonian_ig_ctx->Ctau_v = Ctau_v; 262841e4c73SJed Brown newtonian_ig_ctx->Ctau_C = Ctau_C; 263841e4c73SJed Brown newtonian_ig_ctx->Ctau_M = Ctau_M; 264841e4c73SJed Brown newtonian_ig_ctx->Ctau_E = Ctau_E; 265841e4c73SJed Brown newtonian_ig_ctx->stabilization = stab; 26665dd5cafSJames Wright newtonian_ig_ctx->is_implicit = implicit; 26797baf651SJames Wright newtonian_ig_ctx->state_var = state_var; 2682b730f8bSJeremy L Thompson PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3)); 26988b783a1SJames Wright 2701070106fSJames Wright if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx)); 271c5c39209SJames Wright 272841e4c73SJed Brown CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context); 2732b730f8bSJeremy L Thompson CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context); 2742b730f8bSJeremy L Thompson CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc); 2752b730f8bSJeremy L Thompson CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, "evaluation time", (char *)&setup_context->time - (char *)setup_context, 1, 2762b730f8bSJeremy L Thompson "Time of evaluation"); 277841e4c73SJed Brown 278841e4c73SJed Brown CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context); 2792b730f8bSJeremy L Thompson CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx); 2802b730f8bSJeremy L Thompson CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, FreeContextPetsc); 2812b730f8bSJeremy L Thompson CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1, 2822b730f8bSJeremy L Thompson "Size of timestep, delta t"); 2832b730f8bSJeremy L Thompson CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift", offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 2842b730f8bSJeremy L Thompson 1, "Shift for mass matrix in IJacobian"); 285841e4c73SJed Brown problem->apply_vol_rhs.qfunction_context = newtonian_ig_context; 2862b730f8bSJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ifunction.qfunction_context); 2872b730f8bSJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ijacobian.qfunction_context); 2882b730f8bSJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow.qfunction_context); 2892b730f8bSJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow_jacobian.qfunction_context); 2902b730f8bSJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_outflow.qfunction_context); 2912b730f8bSJeremy L Thompson CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_outflow_jacobian.qfunction_context); 292e334ad8fSJed Brown 293e334ad8fSJed Brown if (unit_tests) { 294e334ad8fSJed Brown PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx)); 295e334ad8fSJed Brown } 29688b783a1SJames Wright PetscFunctionReturn(0); 29788b783a1SJames Wright } 29888b783a1SJames Wright 299dc805cc4SLeila Ghaffari PetscErrorCode PRINT_NEWTONIAN(ProblemData *problem, AppCtx app_ctx) { 300841e4c73SJed Brown MPI_Comm comm = PETSC_COMM_WORLD; 301841e4c73SJed Brown NewtonianIdealGasContext newtonian_ctx; 302841e4c73SJed Brown 30388b783a1SJames Wright PetscFunctionBeginUser; 3042b730f8bSJeremy L Thompson CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ctx); 3052b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 306841e4c73SJed Brown " Problem:\n" 307841e4c73SJed Brown " Problem Name : %s\n" 308841e4c73SJed Brown " Stabilization : %s\n", 3092b730f8bSJeremy L Thompson app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization])); 3102b730f8bSJeremy L Thompson CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ctx); 31188b783a1SJames Wright PetscFunctionReturn(0); 31288b783a1SJames Wright } 313