xref: /honee/problems/newtonian.c (revision 68ae065aca9a501369b2d5dac99bb38ff65b6b2d)
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 "../navierstokes.h"
123a8779fbSJames Wright #include "../qfunctions/setupgeo.h"
133a8779fbSJames Wright #include "../qfunctions/newtonian.h"
143a8779fbSJames Wright 
15f0b65372SJed Brown // Compute relative error |a - b|/|s|
16f0b65372SJed Brown static PetscErrorCode CheckPrimitiveWithTolerance(StatePrimitive sY,
17f0b65372SJed Brown     StatePrimitive aY, StatePrimitive bY, const char *name, PetscReal rtol_pressure,
18f0b65372SJed Brown     PetscReal rtol_velocity, PetscReal rtol_temperature) {
196dd99bedSLeila Ghaffari 
20f0b65372SJed Brown   PetscFunctionBeginUser;
21f0b65372SJed Brown   StatePrimitive eY; // relative error
22f0b65372SJed Brown   eY.pressure = (aY.pressure - bY.pressure) / sY.pressure;
23f0b65372SJed Brown   PetscScalar u = sqrt(Square(sY.velocity[0]) + Square(sY.velocity[1]) + Square(
24f0b65372SJed Brown                          sY.velocity[2]));
25f0b65372SJed Brown   for (int j=0; j<3; j++) eY.velocity[j] = (aY.velocity[j] - bY.velocity[j]) / u;
26f0b65372SJed Brown   eY.temperature = (aY.temperature - bY.temperature) / sY.temperature;
27edd152dcSJed Brown   if (fabs(eY.pressure) > rtol_pressure)
28edd152dcSJed Brown     printf("%s: pressure error %g\n", name, eY.pressure);
29edd152dcSJed Brown   for (int j=0; j<3; j++)
30edd152dcSJed Brown     if (fabs(eY.velocity[j]) > rtol_velocity)
31f0b65372SJed Brown       printf("%s: velocity[%d] error %g\n", name, j, eY.velocity[j]);
32f0b65372SJed Brown   if (fabs(eY.temperature) > rtol_temperature)
33f0b65372SJed Brown     printf("%s: temperature error %g\n", name, eY.temperature);
34f0b65372SJed Brown   PetscFunctionReturn(0);
35f0b65372SJed Brown }
36f0b65372SJed Brown 
37f0b65372SJed Brown static PetscErrorCode UnitTests_Newtonian(User user,
38f0b65372SJed Brown     NewtonianIdealGasContext gas) {
396dd99bedSLeila Ghaffari 
40f0b65372SJed Brown   Units units = user->units;
41f0b65372SJed Brown   const CeedScalar eps = 1e-6;
42f0b65372SJed Brown   const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second,
43f0b65372SJed Brown                    Pascal = units->Pascal;
44f0b65372SJed Brown 
45f0b65372SJed Brown   PetscFunctionBeginUser;
46f0b65372SJed Brown   const CeedScalar rho = 1.2 * kg / (m*m*m), u = 40 * m/sec;
47f0b65372SJed Brown   CeedScalar U[5] = {rho, rho*u, rho *u*1.1, rho *u*1.2, 250e3*Pascal + .5*rho *u*u};
48f0b65372SJed Brown   const CeedScalar x[3] = {.1, .2, .3};
49f0b65372SJed Brown   State s = StateFromU(gas, U, x);
50f0b65372SJed Brown   for (int i=0; i<8; i++) {
51f0b65372SJed Brown     CeedScalar dU[5] = {0}, dx[3] = {0};
52f0b65372SJed Brown     if (i < 5) dU[i] = U[i];
53f0b65372SJed Brown     else dx[i-5] = x[i-5];
54f0b65372SJed Brown     State ds = StateFromU_fwd(gas, s, dU, x, dx);
55f0b65372SJed Brown     for (int j=0; j<5; j++) dU[j] = (1 + eps * (i == j)) * U[j];
56f0b65372SJed Brown     for (int j=0; j<3; j++) dx[j] = (1 + eps * (i == 5 + j)) * x[j];
57f0b65372SJed Brown     State t = StateFromU(gas, dU, dx);
58f0b65372SJed Brown     StatePrimitive dY;
59f0b65372SJed Brown     dY.pressure = (t.Y.pressure - s.Y.pressure) / eps;
60f0b65372SJed Brown     for (int j=0; j<3; j++)
61f0b65372SJed Brown       dY.velocity[j] = (t.Y.velocity[j] - s.Y.velocity[j]) / eps;
62f0b65372SJed Brown     dY.temperature = (t.Y.temperature - s.Y.temperature) / eps;
63f0b65372SJed Brown     char buf[128];
64f0b65372SJed Brown     snprintf(buf, sizeof buf, "StateFromU_fwd i=%d", i);
65f0b65372SJed Brown     PetscCall(CheckPrimitiveWithTolerance(dY, ds.Y, dY, buf, 5e-6, 1e-6, 1e-6));
66f0b65372SJed Brown   }
67f0b65372SJed Brown   PetscFunctionReturn(0);
68f0b65372SJed Brown }
69f0b65372SJed Brown 
70b7f03f12SJed Brown PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx) {
716dd99bedSLeila Ghaffari 
72b7f03f12SJed Brown   SetupContext      setup_context;
733a8779fbSJames Wright   User              user = *(User *)ctx;
743a8779fbSJames Wright   StabilizationType stab;
753a8779fbSJames Wright   MPI_Comm          comm = PETSC_COMM_WORLD;
763a8779fbSJames Wright   PetscBool         implicit;
77f0b65372SJed Brown   PetscBool         has_curr_time = PETSC_FALSE, unit_tests;
783a8779fbSJames Wright   PetscInt          ierr;
7915a3537eSJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
8015a3537eSJed Brown   CeedQFunctionContext newtonian_ig_context;
813a8779fbSJames Wright 
8215a3537eSJed Brown   PetscFunctionBeginUser;
83b7f03f12SJed Brown   ierr = PetscCalloc1(1, &setup_context); CHKERRQ(ierr);
8415a3537eSJed Brown   ierr = PetscCalloc1(1, &newtonian_ig_ctx); CHKERRQ(ierr);
853a8779fbSJames Wright 
863a8779fbSJames Wright   // ------------------------------------------------------
873a8779fbSJames Wright   //           Setup Generic Newtonian IG Problem
883a8779fbSJames Wright   // ------------------------------------------------------
893a8779fbSJames Wright   problem->dim                                  = 3;
903a8779fbSJames Wright   problem->q_data_size_vol                      = 10;
91493642f1SJames Wright   problem->q_data_size_sur                      = 10;
92b01ba163SJames Wright   problem->jac_data_size_sur                    = 11;
939785fe93SJed Brown   problem->setup_vol.qfunction                  = Setup;
949785fe93SJed Brown   problem->setup_vol.qfunction_loc              = Setup_loc;
959785fe93SJed Brown   problem->ics.qfunction                        = ICsNewtonianIG;
969785fe93SJed Brown   problem->ics.qfunction_loc                    = ICsNewtonianIG_loc;
979785fe93SJed Brown   problem->setup_sur.qfunction                  = SetupBoundary;
989785fe93SJed Brown   problem->setup_sur.qfunction_loc              = SetupBoundary_loc;
99c1a52365SJed Brown   problem->apply_vol_rhs.qfunction              = RHSFunction_Newtonian;
100c1a52365SJed Brown   problem->apply_vol_rhs.qfunction_loc          = RHSFunction_Newtonian_loc;
1019785fe93SJed Brown   problem->apply_vol_ifunction.qfunction        = IFunction_Newtonian;
1029785fe93SJed Brown   problem->apply_vol_ifunction.qfunction_loc    = IFunction_Newtonian_loc;
103f0b65372SJed Brown   problem->apply_vol_ijacobian.qfunction        = IJacobian_Newtonian;
104f0b65372SJed Brown   problem->apply_vol_ijacobian.qfunction_loc    = IJacobian_Newtonian_loc;
1058085925cSJames Wright   problem->apply_inflow.qfunction               = BoundaryIntegral;
1068085925cSJames Wright   problem->apply_inflow.qfunction_loc           = BoundaryIntegral_loc;
107*68ae065aSJames Wright   problem->apply_inflow_jacobian.qfunction      = BoundaryIntegral_Jacobian;
108*68ae065aSJames Wright   problem->apply_inflow_jacobian.qfunction_loc  = BoundaryIntegral_Jacobian_loc;
10904b9037bSJames Wright   problem->apply_outflow.qfunction              = PressureOutflow;
11004b9037bSJames Wright   problem->apply_outflow.qfunction_loc          = PressureOutflow_loc;
11104b9037bSJames Wright   problem->apply_outflow_jacobian.qfunction     = PressureOutflow_Jacobian;
11204b9037bSJames Wright   problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_loc;
113b7f03f12SJed Brown   problem->bc                                   = NULL;
114b7f03f12SJed Brown   problem->bc_ctx                               = setup_context;
1153a8779fbSJames Wright   problem->non_zero_time                        = PETSC_FALSE;
1163a8779fbSJames Wright   problem->print_info                           = PRINT_DENSITY_CURRENT;
1173a8779fbSJames Wright 
1183a8779fbSJames Wright   // ------------------------------------------------------
1193a8779fbSJames Wright   //             Create the libCEED context
1203a8779fbSJames Wright   // ------------------------------------------------------
1213a8779fbSJames Wright   CeedScalar cv     = 717.;          // J/(kg K)
1223a8779fbSJames Wright   CeedScalar cp     = 1004.;         // J/(kg K)
123bb8a0c61SJames Wright   CeedScalar g[3]   = {0, 0, -9.81}; // m/s^2
1243a8779fbSJames Wright   CeedScalar lambda = -2./3.;        // -
125bb8a0c61SJames Wright   CeedScalar mu     = 1.8e-5;        // Pa s, dynamic viscosity
1263a8779fbSJames Wright   CeedScalar k      = 0.02638;       // W/(m K)
1273a8779fbSJames Wright   CeedScalar c_tau  = 0.5;           // -
128bb8a0c61SJames Wright   CeedScalar Ctau_t  = 1.0;          // -
129bb8a0c61SJames Wright   CeedScalar Ctau_v  = 36.0;         // TODO make function of degree
130bb8a0c61SJames Wright   CeedScalar Ctau_C  = 1.0;          // TODO make function of degree
131bb8a0c61SJames Wright   CeedScalar Ctau_M  = 1.0;          // TODO make function of degree
132bb8a0c61SJames Wright   CeedScalar Ctau_E  = 1.0;          // TODO make function of degree
1333a8779fbSJames Wright   PetscReal domain_min[3], domain_max[3], domain_size[3];
1343a8779fbSJames Wright   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
135493642f1SJames Wright   for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
1363a8779fbSJames Wright 
1373a8779fbSJames Wright   // ------------------------------------------------------
1383a8779fbSJames Wright   //             Create the PETSc context
1393a8779fbSJames Wright   // ------------------------------------------------------
140bb8a0c61SJames Wright   PetscScalar meter    = 1;  // 1 meter in scaled length units
141bb8a0c61SJames Wright   PetscScalar kilogram = 1;  // 1 kilogram in scaled mass units
142bb8a0c61SJames Wright   PetscScalar second   = 1;  // 1 second in scaled time units
1433a8779fbSJames Wright   PetscScalar Kelvin   = 1;     // 1 Kelvin in scaled temperature units
1443a8779fbSJames Wright   PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s;
1453a8779fbSJames Wright 
1463a8779fbSJames Wright   // ------------------------------------------------------
1473a8779fbSJames Wright   //              Command line Options
1483a8779fbSJames Wright   // ------------------------------------------------------
1491485969bSJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem",
1501485969bSJeremy L Thompson                     NULL);
1511485969bSJeremy L Thompson 
1523a8779fbSJames Wright   // -- Physics
1533a8779fbSJames Wright   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
1543a8779fbSJames Wright                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
1553a8779fbSJames Wright   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
1563a8779fbSJames Wright                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
1573a8779fbSJames Wright   ierr = PetscOptionsScalar("-lambda",
1583a8779fbSJames Wright                             "Stokes hypothesis second viscosity coefficient",
1593a8779fbSJames Wright                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
1603a8779fbSJames Wright   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
1613a8779fbSJames Wright                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
1623a8779fbSJames Wright   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
1633a8779fbSJames Wright                             NULL, k, &k, NULL); CHKERRQ(ierr);
1643a8779fbSJames Wright 
165bb8a0c61SJames Wright   PetscInt dim = problem->dim;
166bb8a0c61SJames Wright   ierr = PetscOptionsRealArray("-g", "Gravitational acceleration",
167bb8a0c61SJames Wright                                NULL, g, &dim, NULL); CHKERRQ(ierr);
1683a8779fbSJames Wright   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
1693a8779fbSJames Wright                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
1703a8779fbSJames Wright                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
1713a8779fbSJames Wright   ierr = PetscOptionsScalar("-c_tau", "Stabilization constant",
1723a8779fbSJames Wright                             NULL, c_tau, &c_tau, NULL); CHKERRQ(ierr);
173bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-Ctau_t", "Stabilization time constant",
174bb8a0c61SJames Wright                             NULL, Ctau_t, &Ctau_t, NULL); CHKERRQ(ierr);
175bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant",
176bb8a0c61SJames Wright                             NULL, Ctau_v, &Ctau_v, NULL); CHKERRQ(ierr);
177bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant",
178bb8a0c61SJames Wright                             NULL, Ctau_C, &Ctau_C, NULL); CHKERRQ(ierr);
179bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant",
180bb8a0c61SJames Wright                             NULL, Ctau_M, &Ctau_M, NULL); CHKERRQ(ierr);
181bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-Ctau_E", "Stabilization energy constant",
182bb8a0c61SJames Wright                             NULL, Ctau_E, &Ctau_E, NULL); CHKERRQ(ierr);
1833a8779fbSJames Wright   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
1843a8779fbSJames Wright                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
1853a8779fbSJames Wright   CHKERRQ(ierr);
186f0b65372SJed Brown   ierr = PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests",
187f0b65372SJed Brown                           NULL, unit_tests=PETSC_FALSE, &unit_tests, NULL);
188f0b65372SJed Brown   CHKERRQ(ierr);
1893a8779fbSJames Wright 
1903a8779fbSJames Wright   // -- Units
1913a8779fbSJames Wright   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
1923a8779fbSJames Wright                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
1933a8779fbSJames Wright   meter = fabs(meter);
1943a8779fbSJames Wright   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
1953a8779fbSJames Wright                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
1963a8779fbSJames Wright   kilogram = fabs(kilogram);
1973a8779fbSJames Wright   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
1983a8779fbSJames Wright                             NULL, second, &second, NULL); CHKERRQ(ierr);
1993a8779fbSJames Wright   second = fabs(second);
2003a8779fbSJames Wright   ierr = PetscOptionsScalar("-units_Kelvin",
2013a8779fbSJames Wright                             "1 Kelvin in scaled temperature units",
2023a8779fbSJames Wright                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
2033a8779fbSJames Wright   Kelvin = fabs(Kelvin);
2043a8779fbSJames Wright 
2053a8779fbSJames Wright   // -- Warnings
2063a8779fbSJames Wright   if (stab == STAB_SUPG && !implicit) {
2073a8779fbSJames Wright     ierr = PetscPrintf(comm,
2083a8779fbSJames Wright                        "Warning! Use -stab supg only with -implicit\n");
2093a8779fbSJames Wright     CHKERRQ(ierr);
2103a8779fbSJames Wright   }
2111485969bSJeremy L Thompson   PetscOptionsEnd();
2123a8779fbSJames Wright 
2133a8779fbSJames Wright   // ------------------------------------------------------
2143a8779fbSJames Wright   //           Set up the PETSc context
2153a8779fbSJames Wright   // ------------------------------------------------------
2163a8779fbSJames Wright   // -- Define derived units
2173a8779fbSJames Wright   Pascal          = kilogram / (meter * PetscSqr(second));
2183a8779fbSJames Wright   J_per_kg_K      =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
2193a8779fbSJames Wright   m_per_squared_s = meter / PetscSqr(second);
2203a8779fbSJames Wright   W_per_m_K       = kilogram * meter / (pow(second,3) * Kelvin);
2213a8779fbSJames Wright 
2223a8779fbSJames Wright   user->units->meter           = meter;
2233a8779fbSJames Wright   user->units->kilogram        = kilogram;
2243a8779fbSJames Wright   user->units->second          = second;
2253a8779fbSJames Wright   user->units->Kelvin          = Kelvin;
2263a8779fbSJames Wright   user->units->Pascal          = Pascal;
2273a8779fbSJames Wright   user->units->J_per_kg_K      = J_per_kg_K;
2283a8779fbSJames Wright   user->units->m_per_squared_s = m_per_squared_s;
2293a8779fbSJames Wright   user->units->W_per_m_K       = W_per_m_K;
2303a8779fbSJames Wright 
2313a8779fbSJames Wright   // ------------------------------------------------------
2323a8779fbSJames Wright   //           Set up the libCEED context
2333a8779fbSJames Wright   // ------------------------------------------------------
2343a8779fbSJames Wright   // -- Scale variables to desired units
2353a8779fbSJames Wright   cv     *= J_per_kg_K;
2363a8779fbSJames Wright   cp     *= J_per_kg_K;
2373a8779fbSJames Wright   mu     *= Pascal * second;
2383a8779fbSJames Wright   k      *= W_per_m_K;
239493642f1SJames Wright   for (PetscInt i=0; i<3; i++) domain_size[i] *= meter;
240493642f1SJames Wright   for (PetscInt i=0; i<3; i++) g[i]           *= m_per_squared_s;
2413a8779fbSJames Wright   problem->dm_scale = meter;
2423a8779fbSJames Wright 
2433a8779fbSJames Wright   // -- Setup Context
2443a8779fbSJames Wright   setup_context->cv         = cv;
2453a8779fbSJames Wright   setup_context->cp         = cp;
2463a8779fbSJames Wright   setup_context->lx         = domain_size[0];
2473a8779fbSJames Wright   setup_context->ly         = domain_size[1];
2483a8779fbSJames Wright   setup_context->lz         = domain_size[2];
2493a8779fbSJames Wright   setup_context->time       = 0;
250bb8a0c61SJames Wright   ierr = PetscArraycpy(setup_context->g, g, 3); CHKERRQ(ierr);
2513a8779fbSJames Wright 
2523a8779fbSJames Wright   // -- Solver Settings
2533a8779fbSJames Wright   user->phys->stab          = stab;
2543a8779fbSJames Wright   user->phys->implicit      = implicit;
2553a8779fbSJames Wright   user->phys->has_curr_time = has_curr_time;
2563a8779fbSJames Wright 
2573a8779fbSJames Wright   // -- QFunction Context
25815a3537eSJed Brown   newtonian_ig_ctx->lambda        = lambda;
25915a3537eSJed Brown   newtonian_ig_ctx->mu            = mu;
26015a3537eSJed Brown   newtonian_ig_ctx->k             = k;
26115a3537eSJed Brown   newtonian_ig_ctx->cv            = cv;
26215a3537eSJed Brown   newtonian_ig_ctx->cp            = cp;
26315a3537eSJed Brown   newtonian_ig_ctx->c_tau         = c_tau;
26415a3537eSJed Brown   newtonian_ig_ctx->Ctau_t        = Ctau_t;
26515a3537eSJed Brown   newtonian_ig_ctx->Ctau_v        = Ctau_v;
26615a3537eSJed Brown   newtonian_ig_ctx->Ctau_C        = Ctau_C;
26715a3537eSJed Brown   newtonian_ig_ctx->Ctau_M        = Ctau_M;
26815a3537eSJed Brown   newtonian_ig_ctx->Ctau_E        = Ctau_E;
26915a3537eSJed Brown   newtonian_ig_ctx->stabilization = stab;
2708085925cSJames Wright   newtonian_ig_ctx->is_implicit   = implicit;
27115a3537eSJed Brown   ierr = PetscArraycpy(newtonian_ig_ctx->g, g, 3); CHKERRQ(ierr);
2723a8779fbSJames Wright 
27315a3537eSJed Brown   CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context);
27415a3537eSJed Brown   CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST,
27515a3537eSJed Brown                               CEED_USE_POINTER, sizeof(*setup_context), setup_context);
276270bbb13SJeremy L Thompson   CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context,
277270bbb13SJeremy L Thompson                                      CEED_MEM_HOST,
278270bbb13SJeremy L Thompson                                      FreeContextPetsc);
27915a3537eSJed Brown   CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context,
28015a3537eSJed Brown                                      "evaluation time",
28115a3537eSJed Brown                                      (char *)&setup_context->time - (char *)setup_context, 1, "Time of evaluation");
28215a3537eSJed Brown 
28315a3537eSJed Brown   CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context);
28415a3537eSJed Brown   CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST,
28515a3537eSJed Brown                               CEED_USE_POINTER,
28615a3537eSJed Brown                               sizeof(*newtonian_ig_ctx), newtonian_ig_ctx);
28715a3537eSJed Brown   CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST,
28815a3537eSJed Brown                                      FreeContextPetsc);
28915a3537eSJed Brown   CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size",
29015a3537eSJed Brown                                      offsetof(struct NewtonianIdealGasContext_, dt), 1, "Size of timestep, delta t");
291f0b65372SJed Brown   CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift",
292f0b65372SJed Brown                                      offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1,
293f0b65372SJed Brown                                      "Shift for mass matrix in IJacobian");
29415a3537eSJed Brown   problem->apply_vol_rhs.qfunction_context = newtonian_ig_context;
29515a3537eSJed Brown   CeedQFunctionContextReferenceCopy(newtonian_ig_context,
29615a3537eSJed Brown                                     &problem->apply_vol_ifunction.qfunction_context);
297f0b65372SJed Brown   CeedQFunctionContextReferenceCopy(newtonian_ig_context,
298f0b65372SJed Brown                                     &problem->apply_vol_ijacobian.qfunction_context);
2998085925cSJames Wright   CeedQFunctionContextReferenceCopy(newtonian_ig_context,
3008085925cSJames Wright                                     &problem->apply_inflow.qfunction_context);
30104b9037bSJames Wright   CeedQFunctionContextReferenceCopy(newtonian_ig_context,
30204b9037bSJames Wright                                     &problem->apply_outflow.qfunction_context);
30304b9037bSJames Wright   CeedQFunctionContextReferenceCopy(newtonian_ig_context,
30404b9037bSJames Wright                                     &problem->apply_outflow_jacobian.qfunction_context);
305f0b65372SJed Brown 
306f0b65372SJed Brown   if (unit_tests) {
307f0b65372SJed Brown     PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx));
308f0b65372SJed Brown   }
3093a8779fbSJames Wright   PetscFunctionReturn(0);
3103a8779fbSJames Wright }
3113a8779fbSJames Wright 
31215a3537eSJed Brown PetscErrorCode PRINT_DENSITY_CURRENT(ProblemData *problem,
31315a3537eSJed Brown                                      AppCtx app_ctx) {
31415a3537eSJed Brown   MPI_Comm comm = PETSC_COMM_WORLD;
31515a3537eSJed Brown   PetscErrorCode ierr;
31615a3537eSJed Brown   NewtonianIdealGasContext newtonian_ctx;
31715a3537eSJed Brown 
3183a8779fbSJames Wright   PetscFunctionBeginUser;
31915a3537eSJed Brown   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context,
32015a3537eSJed Brown                               CEED_MEM_HOST, &newtonian_ctx);
32115a3537eSJed Brown   ierr = PetscPrintf(comm,
32215a3537eSJed Brown                      "  Problem:\n"
32315a3537eSJed Brown                      "    Problem Name                       : %s\n"
32415a3537eSJed Brown                      "    Stabilization                      : %s\n",
32515a3537eSJed Brown                      app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization]);
32615a3537eSJed Brown   CHKERRQ(ierr);
32715a3537eSJed Brown   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context,
32815a3537eSJed Brown                                   &newtonian_ctx);
3293a8779fbSJames Wright   PetscFunctionReturn(0);
3303a8779fbSJames Wright }
331