xref: /libCEED/examples/fluids/problems/newtonian.c (revision d4cc18453651bd0f94c1a2e078b2646a92dafdcc)
1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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 
1349aac155SJeremy L Thompson #include <ceed.h>
1449aac155SJeremy L Thompson #include <petscdm.h>
1549aac155SJeremy L Thompson 
162b730f8bSJeremy L Thompson #include "../navierstokes.h"
176ef2784eSLeila Ghaffari 
18d64246eaSJed Brown // For use with PetscOptionsEnum
19fc39c77eSJames Wright static const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "ENTROPY", "StateVariable", "STATEVAR_", NULL};
20d64246eaSJed Brown 
CheckQWithTolerance(const CeedScalar Q_s[5],const CeedScalar Q_a[5],const CeedScalar Q_b[5],const char * name,PetscReal rtol_0,PetscReal rtol_u,PetscReal rtol_4)2102b29df7SJames Wright static PetscErrorCode CheckQWithTolerance(const CeedScalar Q_s[5], const CeedScalar Q_a[5], const CeedScalar Q_b[5], const char *name,
2202b29df7SJames Wright                                           PetscReal rtol_0, PetscReal rtol_u, PetscReal rtol_4) {
2302b29df7SJames Wright   CeedScalar relative_error[5];  // relative error
2402b29df7SJames Wright   CeedScalar divisor_threshold = 10 * CEED_EPSILON;
25f17d818dSJames Wright 
26f17d818dSJames Wright   PetscFunctionBeginUser;
2702b29df7SJames Wright   relative_error[0] = (Q_a[0] - Q_b[0]) / (fabs(Q_s[0]) > divisor_threshold ? Q_s[0] : 1);
2802b29df7SJames Wright   relative_error[4] = (Q_a[4] - Q_b[4]) / (fabs(Q_s[4]) > divisor_threshold ? Q_s[4] : 1);
2902b29df7SJames Wright 
3002b29df7SJames Wright   CeedScalar u_magnitude = sqrt(Square(Q_s[1]) + Square(Q_s[2]) + Square(Q_s[3]));
3102b29df7SJames Wright   CeedScalar u_divisor   = u_magnitude > divisor_threshold ? u_magnitude : 1;
3202b29df7SJames Wright   for (int i = 1; i < 4; i++) {
3302b29df7SJames Wright     relative_error[i] = (Q_a[i] - Q_b[i]) / u_divisor;
342b730f8bSJeremy L Thompson   }
3502b29df7SJames Wright 
3602b29df7SJames Wright   if (fabs(relative_error[0]) >= rtol_0) {
3702b29df7SJames Wright     printf("%s[0] error %g (expected %.10e, got %.10e)\n", name, relative_error[0], Q_s[0], Q_a[0]);
3802b29df7SJames Wright   }
3902b29df7SJames Wright   for (int i = 1; i < 4; i++) {
4002b29df7SJames Wright     if (fabs(relative_error[i]) >= rtol_u) {
4102b29df7SJames Wright       printf("%s[%d] error %g (expected %.10e, got %.10e)\n", name, i, relative_error[i], Q_s[i], Q_a[i]);
4202b29df7SJames Wright     }
4302b29df7SJames Wright   }
4402b29df7SJames Wright   if (fabs(relative_error[4]) >= rtol_4) {
4502b29df7SJames Wright     printf("%s[4] error %g (expected %.10e, got %.10e)\n", name, relative_error[4], Q_s[4], Q_a[4]);
4602b29df7SJames Wright   }
47ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
48e334ad8fSJed Brown }
49e334ad8fSJed Brown 
5002b29df7SJames Wright // @brief Verify `StateFromQ` by converting A0 -> B0 -> A0_test, where A0 should equal A0_test
TestState(StateVariable state_var_A,StateVariable state_var_B,NewtonianIdealGasContext gas,const CeedScalar A0[5],CeedScalar rtol_0,CeedScalar rtol_u,CeedScalar rtol_4)5102b29df7SJames Wright static PetscErrorCode TestState(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5],
5202b29df7SJames Wright                                 CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
5302b29df7SJames Wright   CeedScalar        B0[5], A0_test[5];
5402b29df7SJames Wright   char              buf[128];
5502b29df7SJames Wright   const char *const StateVariables_Initial[] = {"U", "Y", "V"};
5602b29df7SJames Wright 
5702b29df7SJames Wright   PetscFunctionBeginUser;
5802b29df7SJames Wright   const char *A_initial = StateVariables_Initial[state_var_A];
5902b29df7SJames Wright   const char *B_initial = StateVariables_Initial[state_var_B];
6002b29df7SJames Wright 
6102b29df7SJames Wright   State state_A0 = StateFromQ(gas, A0, state_var_A);
6202b29df7SJames Wright   StateToQ(gas, state_A0, B0, state_var_B);
6302b29df7SJames Wright   State state_B0 = StateFromQ(gas, B0, state_var_B);
6402b29df7SJames Wright   StateToQ(gas, state_B0, A0_test, state_var_A);
6502b29df7SJames Wright 
6602b29df7SJames Wright   snprintf(buf, sizeof buf, "%s->%s->%s: %s", A_initial, B_initial, A_initial, A_initial);
6702b29df7SJames Wright   PetscCall(CheckQWithTolerance(A0, A0_test, A0, buf, rtol_0, rtol_u, rtol_4));
6802b29df7SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
6902b29df7SJames Wright }
7002b29df7SJames Wright 
7102b29df7SJames Wright // @brief Verify `StateFromQ_fwd` via a finite difference approximation
TestState_fwd(StateVariable state_var_A,StateVariable state_var_B,NewtonianIdealGasContext gas,const CeedScalar A0[5],CeedScalar rtol_0,CeedScalar rtol_u,CeedScalar rtol_4)7202b29df7SJames Wright static PetscErrorCode TestState_fwd(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5],
7302b29df7SJames Wright                                     CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
7402b29df7SJames Wright   CeedScalar        eps = 4e-7;  // Finite difference step
7502b29df7SJames Wright   char              buf[128];
7602b29df7SJames Wright   const char *const StateVariables_Initial[] = {"U", "Y", "V"};
7702b29df7SJames Wright 
7802b29df7SJames Wright   PetscFunctionBeginUser;
7902b29df7SJames Wright   const char *A_initial = StateVariables_Initial[state_var_A];
8002b29df7SJames Wright   const char *B_initial = StateVariables_Initial[state_var_B];
8102b29df7SJames Wright   State       state_0   = StateFromQ(gas, A0, state_var_A);
8202b29df7SJames Wright 
8302b29df7SJames Wright   for (int i = 0; i < 5; i++) {
8402b29df7SJames Wright     CeedScalar dB[5] = {0.}, dB_fd[5] = {0.};
8502b29df7SJames Wright     {  // Calculate dB using State functions
8602b29df7SJames Wright       CeedScalar dA[5] = {0};
8702b29df7SJames Wright 
8802b29df7SJames Wright       dA[i]          = A0[i];
89eba67db8SJames Wright       State dstate_0 = StateFromQ_fwd(gas, state_0, dA, state_var_A);
90eba67db8SJames Wright       StateToQ_fwd(gas, state_0, dstate_0, dB, state_var_B);
9102b29df7SJames Wright     }
9202b29df7SJames Wright 
9302b29df7SJames Wright     {  // Calculate dB_fd via finite difference approximation
9402b29df7SJames Wright       CeedScalar A1[5], B0[5], B1[5];
9502b29df7SJames Wright 
9602b29df7SJames Wright       for (int j = 0; j < 5; j++) A1[j] = (1 + eps * (i == j)) * A0[j];
9702b29df7SJames Wright       State state_1 = StateFromQ(gas, A1, state_var_A);
9802b29df7SJames Wright       StateToQ(gas, state_0, B0, state_var_B);
9902b29df7SJames Wright       StateToQ(gas, state_1, B1, state_var_B);
10002b29df7SJames Wright       for (int j = 0; j < 5; j++) dB_fd[j] = (B1[j] - B0[j]) / eps;
10102b29df7SJames Wright     }
10202b29df7SJames Wright 
10302b29df7SJames Wright     snprintf(buf, sizeof buf, "d%s->d%s: StateFrom%s_fwd i=%d: d%s", A_initial, B_initial, A_initial, i, B_initial);
10402b29df7SJames Wright     PetscCall(CheckQWithTolerance(dB_fd, dB, dB_fd, buf, rtol_0, rtol_u, rtol_4));
10502b29df7SJames Wright   }
10602b29df7SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
10702b29df7SJames Wright }
10802b29df7SJames Wright 
10902b29df7SJames Wright // @brief Test the Newtonian State transformation functions, `StateFrom*`
UnitTests_Newtonian(User user,NewtonianIdealGasContext gas)1102b730f8bSJeremy L Thompson static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) {
111e334ad8fSJed Brown   Units            units = user->units;
11202b29df7SJames Wright   const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, K = units->Kelvin;
11302b29df7SJames Wright 
114e334ad8fSJed Brown   PetscFunctionBeginUser;
11502b29df7SJames Wright   const CeedScalar T          = 200 * K;
11602b29df7SJames Wright   const CeedScalar rho        = 1.2 * kg / Cube(m);
11702b29df7SJames Wright   const CeedScalar P          = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T;
11802b29df7SJames Wright   const CeedScalar u_base     = 40 * m / sec;
11902b29df7SJames Wright   const CeedScalar u[3]       = {u_base, u_base * 1.1, u_base * 1.2};
12002b29df7SJames Wright   const CeedScalar e_kinetic  = 0.5 * Dot3(u, u);
12102b29df7SJames Wright   const CeedScalar e_internal = gas->cv * T;
12202b29df7SJames Wright   const CeedScalar e_total    = e_kinetic + e_internal;
12302b29df7SJames Wright   const CeedScalar gamma      = HeatCapacityRatio(gas);
12402b29df7SJames Wright   const CeedScalar entropy    = log(P) - gamma * log(rho);
12502b29df7SJames Wright   const CeedScalar rho_div_p  = rho / P;
12602b29df7SJames Wright   const CeedScalar Y0[5]      = {P, u[0], u[1], u[2], T};
12702b29df7SJames Wright   const CeedScalar U0[5]      = {rho, rho * u[0], rho * u[1], rho * u[2], rho * e_total};
12802b29df7SJames Wright   const CeedScalar V0[5]      = {(gamma - entropy) / (gamma - 1) - rho_div_p * (e_kinetic), rho_div_p * u[0], rho_div_p * u[1], rho_div_p * u[2],
12902b29df7SJames Wright                                  -rho_div_p};
13002b29df7SJames Wright 
13102b29df7SJames Wright   {
13270952158SJames Wright     CeedScalar rtol = 40 * CEED_EPSILON;
13302b29df7SJames Wright 
13402b29df7SJames Wright     PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol));
13502b29df7SJames Wright     PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol));
13602b29df7SJames Wright     PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol));
13702b29df7SJames Wright     PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, rtol, rtol, rtol));
13802b29df7SJames Wright     PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, rtol, rtol, rtol));
13902b29df7SJames Wright     PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, rtol, rtol, rtol));
14002b29df7SJames Wright   }
14102b29df7SJames Wright 
14202b29df7SJames Wright   {
14302b29df7SJames Wright     CeedScalar rtol = 5e-6;
14402b29df7SJames Wright 
14502b29df7SJames Wright     PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol));
14602b29df7SJames Wright     PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol));
14702b29df7SJames Wright     PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol));
14802b29df7SJames Wright     PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, 10 * rtol, rtol, rtol));
14902b29df7SJames Wright     PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, 5 * rtol, rtol, rtol));
15002b29df7SJames Wright     PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, 5 * rtol, 5 * rtol, 5 * rtol));
151e334ad8fSJed Brown   }
152ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
153e334ad8fSJed Brown }
154e334ad8fSJed Brown 
1550fcbc436SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping
1560fcbc436SJames Wright //
1570fcbc436SJames Wright // Only used for SUPG stabilization
CreateKSPMassOperator_NewtonianStabilized(User user,CeedOperator * op_mass)1580fcbc436SJames Wright PetscErrorCode CreateKSPMassOperator_NewtonianStabilized(User user, CeedOperator *op_mass) {
1590fcbc436SJames Wright   Ceed                 ceed = user->ceed;
1600fcbc436SJames Wright   CeedInt              num_comp_q, q_data_size;
1610fcbc436SJames Wright   CeedQFunction        qf_mass;
1620fcbc436SJames Wright   CeedElemRestriction  elem_restr_q, elem_restr_qd_i;
1630fcbc436SJames Wright   CeedBasis            basis_q;
1640fcbc436SJames Wright   CeedVector           q_data;
1650fcbc436SJames Wright   CeedQFunctionContext qf_ctx = NULL;
1660fcbc436SJames Wright   PetscInt             dim    = 3;
1670fcbc436SJames Wright 
1680fcbc436SJames Wright   PetscFunctionBeginUser;
1690fcbc436SJames Wright   {  // Get restriction and basis from the RHS function
1700fcbc436SJames Wright     CeedOperator     *sub_ops;
1710fcbc436SJames Wright     CeedOperatorField field;
1720fcbc436SJames Wright     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
1730fcbc436SJames Wright 
174ed094490SJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(user->op_rhs_ctx->op, &sub_ops));
1750fcbc436SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
176681d0ea7SJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorFieldGetData(field, NULL, &elem_restr_q, &basis_q, NULL));
1770fcbc436SJames Wright 
1780fcbc436SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field));
179681d0ea7SJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorFieldGetData(field, NULL, &elem_restr_qd_i, NULL, &q_data));
1800fcbc436SJames Wright 
1810fcbc436SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx));
1820fcbc436SJames Wright   }
1830fcbc436SJames Wright 
1840fcbc436SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
1850fcbc436SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size));
1860fcbc436SJames Wright 
1870fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Newtonian_Conserv, MassFunction_Newtonian_Conserv_loc, &qf_mass));
1880fcbc436SJames Wright 
1890fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx));
1900fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
1910fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
1920fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
1930fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
1940fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
1950fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
1960fcbc436SJames Wright 
1970fcbc436SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
1980fcbc436SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
1990fcbc436SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed));
2000fcbc436SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data));
2010fcbc436SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
2020fcbc436SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
2030fcbc436SJames Wright 
204681d0ea7SJeremy L Thompson   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
205681d0ea7SJeremy L Thompson   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
206681d0ea7SJeremy L Thompson   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i));
207681d0ea7SJeremy L Thompson   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
20838690fecSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qf_ctx));
2090fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
2100fcbc436SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2110fcbc436SJames Wright }
212013a5551SJames Wright 
NS_NEWTONIAN_IG(ProblemData problem,DM dm,void * ctx,SimpleBC bc)213731c13d7SJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
214a0add3c9SJed Brown   SetupContext             setup_context;
21588b783a1SJames Wright   User                     user   = *(User *)ctx;
216ebd2ea64SLeila Ghaffari   CeedInt                  degree = user->app_ctx->degree;
21788b783a1SJames Wright   StabilizationType        stab;
21897baf651SJames Wright   StateVariable            state_var;
219a424bcd0SJames Wright   MPI_Comm                 comm = user->comm;
220a424bcd0SJames Wright   Ceed                     ceed = user->ceed;
22188b783a1SJames Wright   PetscBool                implicit;
222a2f6637eSJames Wright   PetscBool                unit_tests;
223841e4c73SJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
224841e4c73SJed Brown   CeedQFunctionContext     newtonian_ig_context;
22588b783a1SJames Wright 
226841e4c73SJed Brown   PetscFunctionBeginUser;
2272b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &setup_context));
2282b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &newtonian_ig_ctx));
22988b783a1SJames Wright 
23088b783a1SJames Wright   // ------------------------------------------------------
23188b783a1SJames Wright   //           Setup Generic Newtonian IG Problem
23288b783a1SJames Wright   // ------------------------------------------------------
23388b783a1SJames Wright   problem->dim                          = 3;
2340ec2498eSJames Wright   problem->jac_data_size_sur            = 11;
235de327db4SJames Wright   problem->compute_exact_solution_error = PETSC_FALSE;
236dc805cc4SLeila Ghaffari   problem->print_info                   = PRINT_NEWTONIAN;
237752a08a7SJames Wright   problem->uses_newtonian               = PETSC_TRUE;
23888b783a1SJames Wright 
23988b783a1SJames Wright   // ------------------------------------------------------
24088b783a1SJames Wright   //             Create the libCEED context
24188b783a1SJames Wright   // ------------------------------------------------------
24288b783a1SJames Wright   CeedScalar cv         = 717.;          // J/(kg K)
24388b783a1SJames Wright   CeedScalar cp         = 1004.;         // J/(kg K)
244a2726bdbSJames Wright   CeedScalar g[3]       = {0, 0, 0};     // m/s^2
24588b783a1SJames Wright   CeedScalar lambda     = -2. / 3.;      // -
24688626eedSJames Wright   CeedScalar mu         = 1.8e-5;        // Pa s, dynamic viscosity
24788b783a1SJames Wright   CeedScalar k          = 0.02638;       // W/(m K)
248ebd2ea64SLeila Ghaffari   CeedScalar c_tau      = 0.5 / degree;  // -
24988626eedSJames Wright   CeedScalar Ctau_t     = 1.0;           // -
25094c01735SLeila Ghaffari   CeedScalar Cv_func[3] = {36, 60, 128};
251c5dde687SLeila Ghaffari   CeedScalar Ctau_v     = Cv_func[(CeedInt)Min(3, degree) - 1];
25275538696SLeila Ghaffari   CeedScalar Ctau_C     = 0.25 / degree;
25375538696SLeila Ghaffari   CeedScalar Ctau_M     = 0.25 / degree;
25475538696SLeila Ghaffari   CeedScalar Ctau_E     = 0.125;
25588b783a1SJames Wright   PetscReal  domain_min[3], domain_max[3], domain_size[3];
2562b730f8bSJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
257ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
25888b783a1SJames Wright 
259d310b3d3SAdeleke O. Bankole   StatePrimitive reference      = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15};
2602249ac91SJames Wright   CeedScalar     idl_decay_time = -1, idl_start = 0, idl_length = 0, idl_pressure = reference.pressure;
261530ad8c4SKenneth E. Jansen   PetscBool      idl_enable = PETSC_FALSE;
262d310b3d3SAdeleke O. Bankole 
26388b783a1SJames Wright   // ------------------------------------------------------
26488b783a1SJames Wright   //             Create the PETSc context
26588b783a1SJames Wright   // ------------------------------------------------------
26688626eedSJames Wright   PetscScalar meter    = 1;  // 1 meter in scaled length units
26788626eedSJames Wright   PetscScalar kilogram = 1;  // 1 kilogram in scaled mass units
26888626eedSJames Wright   PetscScalar second   = 1;  // 1 second in scaled time units
26988b783a1SJames Wright   PetscScalar Kelvin   = 1;  // 1 Kelvin in scaled temperature units
27088b783a1SJames Wright   PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s;
27188b783a1SJames Wright 
27288b783a1SJames Wright   // ------------------------------------------------------
27388b783a1SJames Wright   //              Command line Options
27488b783a1SJames Wright   // ------------------------------------------------------
275a2726bdbSJames Wright   PetscBool given_option = PETSC_FALSE;
2762b730f8bSJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL);
277dc805cc4SLeila Ghaffari   // -- Conservative vs Primitive variables
2782b730f8bSJeremy L Thompson   PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE),
2792b730f8bSJeremy L Thompson                              (PetscEnum *)&state_var, NULL));
28097baf651SJames Wright 
28197baf651SJames Wright   switch (state_var) {
28297baf651SJames Wright     case STATEVAR_CONSERVATIVE:
283d310b3d3SAdeleke O. Bankole       problem->ics.qfunction                       = ICsNewtonianIG_Conserv;
284d310b3d3SAdeleke O. Bankole       problem->ics.qfunction_loc                   = ICsNewtonianIG_Conserv_loc;
285dc805cc4SLeila Ghaffari       problem->apply_vol_rhs.qfunction             = RHSFunction_Newtonian;
286dc805cc4SLeila Ghaffari       problem->apply_vol_rhs.qfunction_loc         = RHSFunction_Newtonian_loc;
2873d02368aSJames Wright       problem->apply_vol_ifunction.qfunction       = IFunction_Newtonian_Conserv;
2883d02368aSJames Wright       problem->apply_vol_ifunction.qfunction_loc   = IFunction_Newtonian_Conserv_loc;
2893d02368aSJames Wright       problem->apply_vol_ijacobian.qfunction       = IJacobian_Newtonian_Conserv;
2903d02368aSJames Wright       problem->apply_vol_ijacobian.qfunction_loc   = IJacobian_Newtonian_Conserv_loc;
29120840d50SJames Wright       problem->apply_inflow.qfunction              = BoundaryIntegral_Conserv;
29220840d50SJames Wright       problem->apply_inflow.qfunction_loc          = BoundaryIntegral_Conserv_loc;
29320840d50SJames Wright       problem->apply_inflow_jacobian.qfunction     = BoundaryIntegral_Jacobian_Conserv;
29420840d50SJames Wright       problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc;
29597baf651SJames Wright       break;
29697baf651SJames Wright     case STATEVAR_PRIMITIVE:
29797baf651SJames Wright       problem->ics.qfunction                       = ICsNewtonianIG_Prim;
29897baf651SJames Wright       problem->ics.qfunction_loc                   = ICsNewtonianIG_Prim_loc;
29997baf651SJames Wright       problem->apply_vol_ifunction.qfunction       = IFunction_Newtonian_Prim;
30097baf651SJames Wright       problem->apply_vol_ifunction.qfunction_loc   = IFunction_Newtonian_Prim_loc;
30197baf651SJames Wright       problem->apply_vol_ijacobian.qfunction       = IJacobian_Newtonian_Prim;
30297baf651SJames Wright       problem->apply_vol_ijacobian.qfunction_loc   = IJacobian_Newtonian_Prim_loc;
30397baf651SJames Wright       problem->apply_inflow.qfunction              = BoundaryIntegral_Prim;
30497baf651SJames Wright       problem->apply_inflow.qfunction_loc          = BoundaryIntegral_Prim_loc;
30597baf651SJames Wright       problem->apply_inflow_jacobian.qfunction     = BoundaryIntegral_Jacobian_Prim;
30697baf651SJames Wright       problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc;
30797baf651SJames Wright       break;
30802b29df7SJames Wright     case STATEVAR_ENTROPY:
309a2d72b6fSJames Wright       problem->ics.qfunction                       = ICsNewtonianIG_Entropy;
310a2d72b6fSJames Wright       problem->ics.qfunction_loc                   = ICsNewtonianIG_Entropy_loc;
311a2d72b6fSJames Wright       problem->apply_vol_ifunction.qfunction       = IFunction_Newtonian_Entropy;
312a2d72b6fSJames Wright       problem->apply_vol_ifunction.qfunction_loc   = IFunction_Newtonian_Entropy_loc;
313a2d72b6fSJames Wright       problem->apply_vol_ijacobian.qfunction       = IJacobian_Newtonian_Entropy;
314a2d72b6fSJames Wright       problem->apply_vol_ijacobian.qfunction_loc   = IJacobian_Newtonian_Entropy_loc;
315a2d72b6fSJames Wright       problem->apply_inflow.qfunction              = BoundaryIntegral_Entropy;
316a2d72b6fSJames Wright       problem->apply_inflow.qfunction_loc          = BoundaryIntegral_Entropy_loc;
317a2d72b6fSJames Wright       problem->apply_inflow_jacobian.qfunction     = BoundaryIntegral_Jacobian_Entropy;
318a2d72b6fSJames Wright       problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Entropy_loc;
31902b29df7SJames Wright       break;
320dc805cc4SLeila Ghaffari   }
32167490bc6SJeremy L Thompson 
32288b783a1SJames Wright   // -- Physics
3232b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL));
3242b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL));
3252b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL));
3262b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL));
3272b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL));
32888b783a1SJames Wright 
32988626eedSJames Wright   PetscInt dim = problem->dim;
330a2726bdbSJames Wright   PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL));
331a2726bdbSJames Wright   PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option));
332a2726bdbSJames Wright   if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim);
333a2726bdbSJames Wright 
3342b730f8bSJeremy L Thompson   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
3352b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
3362b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
3372b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL));
3382b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL));
3392b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL));
3402b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL));
3412b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
3422b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL));
34388b783a1SJames Wright 
344e30587cbSJames Wright   dim = 3;
345d310b3d3SAdeleke O. Bankole   PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL));
346d310b3d3SAdeleke O. Bankole   PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL));
347d310b3d3SAdeleke O. Bankole   PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL));
348d310b3d3SAdeleke O. Bankole 
34988b783a1SJames Wright   // -- Units
3502b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
35188b783a1SJames Wright   meter = fabs(meter);
3522b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL));
35388b783a1SJames Wright   kilogram = fabs(kilogram);
3542b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
35588b783a1SJames Wright   second = fabs(second);
3562b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL));
35788b783a1SJames Wright   Kelvin = fabs(Kelvin);
35888b783a1SJames Wright 
35988b783a1SJames Wright   // -- Warnings
3600e654f56SJames Wright   PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP,
3610e654f56SJames Wright              "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n");
362530ad8c4SKenneth E. Jansen 
363530ad8c4SKenneth E. Jansen   PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point",
364530ad8c4SKenneth E. Jansen                                NULL, idl_decay_time, &idl_decay_time, &idl_enable));
365fe1e732eSJames Wright   PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero.");
366fe1e732eSJames Wright   if (idl_decay_time < 0) idl_enable = PETSC_FALSE;
367530ad8c4SKenneth E. Jansen   PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL));
368530ad8c4SKenneth E. Jansen   PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL));
3692249ac91SJames Wright   idl_pressure = reference.pressure;
3702249ac91SJames Wright   PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-reference_pressure`)", NULL, idl_pressure,
3712249ac91SJames Wright                                &idl_pressure, NULL));
37267490bc6SJeremy L Thompson   PetscOptionsEnd();
37388b783a1SJames Wright 
3740fcbc436SJames Wright   if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized;
3750fcbc436SJames Wright 
37688b783a1SJames Wright   // ------------------------------------------------------
37788b783a1SJames Wright   //           Set up the PETSc context
37888b783a1SJames Wright   // ------------------------------------------------------
37988b783a1SJames Wright   // -- Define derived units
38088b783a1SJames Wright   Pascal          = kilogram / (meter * PetscSqr(second));
38188b783a1SJames Wright   J_per_kg_K      = PetscSqr(meter) / (PetscSqr(second) * Kelvin);
38288b783a1SJames Wright   m_per_squared_s = meter / PetscSqr(second);
38388b783a1SJames Wright   W_per_m_K       = kilogram * meter / (pow(second, 3) * Kelvin);
38488b783a1SJames Wright 
38588b783a1SJames Wright   user->units->meter           = meter;
38688b783a1SJames Wright   user->units->kilogram        = kilogram;
38788b783a1SJames Wright   user->units->second          = second;
38888b783a1SJames Wright   user->units->Kelvin          = Kelvin;
38988b783a1SJames Wright   user->units->Pascal          = Pascal;
39088b783a1SJames Wright   user->units->J_per_kg_K      = J_per_kg_K;
39188b783a1SJames Wright   user->units->m_per_squared_s = m_per_squared_s;
39288b783a1SJames Wright   user->units->W_per_m_K       = W_per_m_K;
39388b783a1SJames Wright 
39488b783a1SJames Wright   // ------------------------------------------------------
39588b783a1SJames Wright   //           Set up the libCEED context
39688b783a1SJames Wright   // ------------------------------------------------------
39788b783a1SJames Wright   // -- Scale variables to desired units
39888b783a1SJames Wright   cv *= J_per_kg_K;
39988b783a1SJames Wright   cp *= J_per_kg_K;
40088b783a1SJames Wright   mu *= Pascal * second;
40188b783a1SJames Wright   k *= W_per_m_K;
402ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter;
403ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s;
404d310b3d3SAdeleke O. Bankole   reference.pressure *= Pascal;
405d310b3d3SAdeleke O. Bankole   for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second;
406d310b3d3SAdeleke O. Bankole   reference.temperature *= Kelvin;
40788b783a1SJames Wright   problem->dm_scale = meter;
40888b783a1SJames Wright 
40988b783a1SJames Wright   // -- Solver Settings
41088b783a1SJames Wright   user->phys->implicit  = implicit;
41197baf651SJames Wright   user->phys->state_var = state_var;
41288b783a1SJames Wright 
41388b783a1SJames Wright   // -- QFunction Context
414841e4c73SJed Brown   newtonian_ig_ctx->lambda        = lambda;
415841e4c73SJed Brown   newtonian_ig_ctx->mu            = mu;
416841e4c73SJed Brown   newtonian_ig_ctx->k             = k;
417841e4c73SJed Brown   newtonian_ig_ctx->cv            = cv;
418841e4c73SJed Brown   newtonian_ig_ctx->cp            = cp;
419841e4c73SJed Brown   newtonian_ig_ctx->c_tau         = c_tau;
420841e4c73SJed Brown   newtonian_ig_ctx->Ctau_t        = Ctau_t;
421841e4c73SJed Brown   newtonian_ig_ctx->Ctau_v        = Ctau_v;
422841e4c73SJed Brown   newtonian_ig_ctx->Ctau_C        = Ctau_C;
423841e4c73SJed Brown   newtonian_ig_ctx->Ctau_M        = Ctau_M;
424841e4c73SJed Brown   newtonian_ig_ctx->Ctau_E        = Ctau_E;
425841e4c73SJed Brown   newtonian_ig_ctx->stabilization = stab;
42665dd5cafSJames Wright   newtonian_ig_ctx->is_implicit   = implicit;
42797baf651SJames Wright   newtonian_ig_ctx->state_var     = state_var;
428530ad8c4SKenneth E. Jansen   newtonian_ig_ctx->idl_enable    = idl_enable;
429530ad8c4SKenneth E. Jansen   newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * second);
430530ad8c4SKenneth E. Jansen   newtonian_ig_ctx->idl_start     = idl_start * meter;
431530ad8c4SKenneth E. Jansen   newtonian_ig_ctx->idl_length    = idl_length * meter;
4322249ac91SJames Wright   newtonian_ig_ctx->idl_pressure  = idl_pressure;
4332b730f8bSJeremy L Thompson   PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3));
43488b783a1SJames Wright 
435d310b3d3SAdeleke O. Bankole   // -- Setup Context
436d310b3d3SAdeleke O. Bankole   setup_context->reference = reference;
437d310b3d3SAdeleke O. Bankole   setup_context->gas       = *newtonian_ig_ctx;
438d310b3d3SAdeleke O. Bankole   setup_context->lx        = domain_size[0];
439d310b3d3SAdeleke O. Bankole   setup_context->ly        = domain_size[1];
440d310b3d3SAdeleke O. Bankole   setup_context->lz        = domain_size[2];
441d310b3d3SAdeleke O. Bankole   setup_context->time      = 0;
442d310b3d3SAdeleke O. Bankole 
443a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context));
444a424bcd0SJames Wright   PetscCallCeed(ceed,
445a424bcd0SJames Wright                 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
446a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc));
447a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, "evaluation time", offsetof(struct SetupContext_, time), 1,
448a424bcd0SJames Wright                                                          "Time of evaluation"));
449841e4c73SJed Brown 
450a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context));
451a424bcd0SJames Wright   PetscCallCeed(ceed,
452a424bcd0SJames Wright                 CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx));
453a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, FreeContextPetsc));
454a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1,
455a424bcd0SJames Wright                                                          "Size of timestep, delta t"));
456a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift",
457a424bcd0SJames Wright                                                          offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1,
458a424bcd0SJames Wright                                                          "Shift for mass matrix in IJacobian"));
459a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1,
460a424bcd0SJames Wright                                                          "Current solution time"));
461f179a332SJames Wright 
462841e4c73SJed Brown   problem->apply_vol_rhs.qfunction_context = newtonian_ig_context;
463a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ifunction.qfunction_context));
464a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ijacobian.qfunction_context));
465a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow.qfunction_context));
466a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow_jacobian.qfunction_context));
467e334ad8fSJed Brown 
4689f844368SJames Wright   if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
4699f844368SJames Wright   if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
4709f844368SJames Wright   if (bc->num_slip > 0) PetscCall(SlipBCSetup(problem, dm, ctx, newtonian_ig_context));
4719f844368SJames Wright 
472e334ad8fSJed Brown   if (unit_tests) {
473e334ad8fSJed Brown     PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx));
474e334ad8fSJed Brown   }
475ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
47688b783a1SJames Wright }
47788b783a1SJames Wright 
PRINT_NEWTONIAN(User user,ProblemData problem,AppCtx app_ctx)478731c13d7SJames Wright PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx) {
479367c849eSJames Wright   MPI_Comm                 comm = user->comm;
480a424bcd0SJames Wright   Ceed                     ceed = user->ceed;
481841e4c73SJed Brown   NewtonianIdealGasContext newtonian_ctx;
482841e4c73SJed Brown 
48388b783a1SJames Wright   PetscFunctionBeginUser;
484a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ctx));
4852b730f8bSJeremy L Thompson   PetscCall(PetscPrintf(comm,
486841e4c73SJed Brown                         "  Problem:\n"
487841e4c73SJed Brown                         "    Problem Name                       : %s\n"
488841e4c73SJed Brown                         "    Stabilization                      : %s\n",
4892b730f8bSJeremy L Thompson                         app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization]));
490a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ctx));
491ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
49288b783a1SJames Wright }
493