xref: /libCEED/examples/fluids/problems/newtonian.c (revision de327db462e8afadb04b2fd99bdb39cfb29a1d98)
1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Utility functions for setting up problems using the Newtonian Qfunction
10 
11 #include "../qfunctions/newtonian.h"
12 
13 #include <ceed.h>
14 #include <petscdm.h>
15 
16 #include "../navierstokes.h"
17 
18 // For use with PetscOptionsEnum
19 static const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "StateVariable", "STATEVAR_", NULL};
20 
21 // Compute relative error |a - b|/|s|
22 static PetscErrorCode CheckPrimitiveWithTolerance(StatePrimitive sY, StatePrimitive aY, StatePrimitive bY, const char *name, PetscReal rtol_pressure,
23                                                   PetscReal rtol_velocity, PetscReal rtol_temperature) {
24   StatePrimitive eY;  // relative error
25 
26   PetscFunctionBeginUser;
27   eY.pressure   = (aY.pressure - bY.pressure) / sY.pressure;
28   PetscScalar u = sqrt(Square(sY.velocity[0]) + Square(sY.velocity[1]) + Square(sY.velocity[2]));
29   for (int j = 0; j < 3; j++) eY.velocity[j] = (aY.velocity[j] - bY.velocity[j]) / u;
30   eY.temperature = (aY.temperature - bY.temperature) / sY.temperature;
31   if (fabs(eY.pressure) > rtol_pressure) printf("%s: pressure error %g\n", name, eY.pressure);
32   for (int j = 0; j < 3; j++) {
33     if (fabs(eY.velocity[j]) > rtol_velocity) printf("%s: velocity[%d] error %g\n", name, j, eY.velocity[j]);
34   }
35   if (fabs(eY.temperature) > rtol_temperature) printf("%s: temperature error %g\n", name, eY.temperature);
36   PetscFunctionReturn(PETSC_SUCCESS);
37 }
38 
39 static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) {
40   Units            units = user->units;
41   const CeedScalar eps   = 1e-6;
42   const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, Pascal = units->Pascal;
43   PetscFunctionBeginUser;
44   const CeedScalar rho = 1.2 * kg / (m * m * m), u = 40 * m / sec;
45   CeedScalar       U[5] = {rho, rho * u, rho * u * 1.1, rho * u * 1.2, 250e3 * Pascal + .5 * rho * u * u};
46   State            s    = StateFromU(gas, U);
47   for (int i = 0; i < 8; i++) {
48     CeedScalar dU[5] = {0};
49     if (i < 5) dU[i] = U[i];
50     State ds = StateFromU_fwd(gas, s, dU);
51     for (int j = 0; j < 5; j++) dU[j] = (1 + eps * (i == j)) * U[j];
52     State          t = StateFromU(gas, dU);
53     StatePrimitive dY;
54     dY.pressure = (t.Y.pressure - s.Y.pressure) / eps;
55     for (int j = 0; j < 3; j++) dY.velocity[j] = (t.Y.velocity[j] - s.Y.velocity[j]) / eps;
56     dY.temperature = (t.Y.temperature - s.Y.temperature) / eps;
57     char buf[128];
58     snprintf(buf, sizeof buf, "StateFromU_fwd i=%d", i);
59     PetscCall(CheckPrimitiveWithTolerance(dY, ds.Y, dY, buf, 5e-6, 1e-6, 1e-6));
60   }
61   PetscFunctionReturn(PETSC_SUCCESS);
62 }
63 
64 // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping
65 //
66 // Only used for SUPG stabilization
67 PetscErrorCode CreateKSPMassOperator_NewtonianStabilized(User user, CeedOperator *op_mass) {
68   Ceed                 ceed = user->ceed;
69   CeedInt              num_comp_q, q_data_size;
70   CeedQFunction        qf_mass;
71   CeedElemRestriction  elem_restr_q, elem_restr_qd_i;
72   CeedBasis            basis_q;
73   CeedVector           q_data;
74   CeedQFunctionContext qf_ctx = NULL;
75   PetscInt             dim    = 3;
76 
77   PetscFunctionBeginUser;
78   {  // Get restriction and basis from the RHS function
79     CeedOperator     *sub_ops;
80     CeedOperatorField field;
81     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
82 
83     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops));
84     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
85     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
86     PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q));
87 
88     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field));
89     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i));
90     PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data));
91 
92     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx));
93   }
94 
95   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
96   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size));
97 
98   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Newtonian_Conserv, MassFunction_Newtonian_Conserv_loc, &qf_mass));
99 
100   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx));
101   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
102   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
103   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
104   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
105   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
106   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
107 
108   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
109   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
110   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed));
111   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data));
112   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
113   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
114 
115   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
116   PetscFunctionReturn(PETSC_SUCCESS);
117 }
118 
119 PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
120   SetupContext             setup_context;
121   User                     user   = *(User *)ctx;
122   CeedInt                  degree = user->app_ctx->degree;
123   StabilizationType        stab;
124   StateVariable            state_var;
125   MPI_Comm                 comm = user->comm;
126   Ceed                     ceed = user->ceed;
127   PetscBool                implicit;
128   PetscBool                unit_tests;
129   NewtonianIdealGasContext newtonian_ig_ctx;
130   CeedQFunctionContext     newtonian_ig_context;
131 
132   PetscFunctionBeginUser;
133   PetscCall(PetscCalloc1(1, &setup_context));
134   PetscCall(PetscCalloc1(1, &newtonian_ig_ctx));
135 
136   // ------------------------------------------------------
137   //           Setup Generic Newtonian IG Problem
138   // ------------------------------------------------------
139   problem->dim                          = 3;
140   problem->jac_data_size_sur            = 11;
141   problem->compute_exact_solution_error = PETSC_FALSE;
142   problem->print_info                   = PRINT_NEWTONIAN;
143   problem->uses_newtonian               = PETSC_TRUE;
144 
145   // ------------------------------------------------------
146   //             Create the libCEED context
147   // ------------------------------------------------------
148   CeedScalar cv         = 717.;          // J/(kg K)
149   CeedScalar cp         = 1004.;         // J/(kg K)
150   CeedScalar g[3]       = {0, 0, 0};     // m/s^2
151   CeedScalar lambda     = -2. / 3.;      // -
152   CeedScalar mu         = 1.8e-5;        // Pa s, dynamic viscosity
153   CeedScalar k          = 0.02638;       // W/(m K)
154   CeedScalar c_tau      = 0.5 / degree;  // -
155   CeedScalar Ctau_t     = 1.0;           // -
156   CeedScalar Cv_func[3] = {36, 60, 128};
157   CeedScalar Ctau_v     = Cv_func[(CeedInt)Min(3, degree) - 1];
158   CeedScalar Ctau_C     = 0.25 / degree;
159   CeedScalar Ctau_M     = 0.25 / degree;
160   CeedScalar Ctau_E     = 0.125;
161   PetscReal  domain_min[3], domain_max[3], domain_size[3];
162   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
163   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
164 
165   StatePrimitive reference      = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15};
166   CeedScalar     idl_decay_time = -1, idl_start = 0, idl_length = 0, idl_pressure = reference.pressure;
167   PetscBool      idl_enable = PETSC_FALSE;
168 
169   // ------------------------------------------------------
170   //             Create the PETSc context
171   // ------------------------------------------------------
172   PetscScalar meter    = 1;  // 1 meter in scaled length units
173   PetscScalar kilogram = 1;  // 1 kilogram in scaled mass units
174   PetscScalar second   = 1;  // 1 second in scaled time units
175   PetscScalar Kelvin   = 1;  // 1 Kelvin in scaled temperature units
176   PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s;
177 
178   // ------------------------------------------------------
179   //              Command line Options
180   // ------------------------------------------------------
181   PetscBool given_option = PETSC_FALSE;
182   PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL);
183   // -- Conservative vs Primitive variables
184   PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE),
185                              (PetscEnum *)&state_var, NULL));
186 
187   switch (state_var) {
188     case STATEVAR_CONSERVATIVE:
189       problem->ics.qfunction                       = ICsNewtonianIG_Conserv;
190       problem->ics.qfunction_loc                   = ICsNewtonianIG_Conserv_loc;
191       problem->apply_vol_rhs.qfunction             = RHSFunction_Newtonian;
192       problem->apply_vol_rhs.qfunction_loc         = RHSFunction_Newtonian_loc;
193       problem->apply_vol_ifunction.qfunction       = IFunction_Newtonian_Conserv;
194       problem->apply_vol_ifunction.qfunction_loc   = IFunction_Newtonian_Conserv_loc;
195       problem->apply_vol_ijacobian.qfunction       = IJacobian_Newtonian_Conserv;
196       problem->apply_vol_ijacobian.qfunction_loc   = IJacobian_Newtonian_Conserv_loc;
197       problem->apply_inflow.qfunction              = BoundaryIntegral_Conserv;
198       problem->apply_inflow.qfunction_loc          = BoundaryIntegral_Conserv_loc;
199       problem->apply_inflow_jacobian.qfunction     = BoundaryIntegral_Jacobian_Conserv;
200       problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc;
201       break;
202 
203     case STATEVAR_PRIMITIVE:
204       problem->ics.qfunction                       = ICsNewtonianIG_Prim;
205       problem->ics.qfunction_loc                   = ICsNewtonianIG_Prim_loc;
206       problem->apply_vol_ifunction.qfunction       = IFunction_Newtonian_Prim;
207       problem->apply_vol_ifunction.qfunction_loc   = IFunction_Newtonian_Prim_loc;
208       problem->apply_vol_ijacobian.qfunction       = IJacobian_Newtonian_Prim;
209       problem->apply_vol_ijacobian.qfunction_loc   = IJacobian_Newtonian_Prim_loc;
210       problem->apply_inflow.qfunction              = BoundaryIntegral_Prim;
211       problem->apply_inflow.qfunction_loc          = BoundaryIntegral_Prim_loc;
212       problem->apply_inflow_jacobian.qfunction     = BoundaryIntegral_Jacobian_Prim;
213       problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc;
214       break;
215   }
216 
217   // -- Physics
218   PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL));
219   PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL));
220   PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL));
221   PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL));
222   PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL));
223 
224   PetscInt dim = problem->dim;
225   PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL));
226   PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option));
227   if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim);
228 
229   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
230   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
231   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
232   PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL));
233   PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL));
234   PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL));
235   PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL));
236   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
237   PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL));
238 
239   dim = 3;
240   PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL));
241   PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL));
242   PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL));
243 
244   // -- Units
245   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
246   meter = fabs(meter);
247   PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL));
248   kilogram = fabs(kilogram);
249   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
250   second = fabs(second);
251   PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL));
252   Kelvin = fabs(Kelvin);
253 
254   // -- Warnings
255   PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP,
256              "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n");
257 
258   PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point",
259                                NULL, idl_decay_time, &idl_decay_time, &idl_enable));
260   PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero.");
261   if (idl_decay_time < 0) idl_enable = PETSC_FALSE;
262   PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL));
263   PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL));
264   idl_pressure = reference.pressure;
265   PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-reference_pressure`)", NULL, idl_pressure,
266                                &idl_pressure, NULL));
267   PetscOptionsEnd();
268 
269   if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized;
270 
271   // ------------------------------------------------------
272   //           Set up the PETSc context
273   // ------------------------------------------------------
274   // -- Define derived units
275   Pascal          = kilogram / (meter * PetscSqr(second));
276   J_per_kg_K      = PetscSqr(meter) / (PetscSqr(second) * Kelvin);
277   m_per_squared_s = meter / PetscSqr(second);
278   W_per_m_K       = kilogram * meter / (pow(second, 3) * Kelvin);
279 
280   user->units->meter           = meter;
281   user->units->kilogram        = kilogram;
282   user->units->second          = second;
283   user->units->Kelvin          = Kelvin;
284   user->units->Pascal          = Pascal;
285   user->units->J_per_kg_K      = J_per_kg_K;
286   user->units->m_per_squared_s = m_per_squared_s;
287   user->units->W_per_m_K       = W_per_m_K;
288 
289   // ------------------------------------------------------
290   //           Set up the libCEED context
291   // ------------------------------------------------------
292   // -- Scale variables to desired units
293   cv *= J_per_kg_K;
294   cp *= J_per_kg_K;
295   mu *= Pascal * second;
296   k *= W_per_m_K;
297   for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter;
298   for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s;
299   reference.pressure *= Pascal;
300   for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second;
301   reference.temperature *= Kelvin;
302   problem->dm_scale = meter;
303 
304   // -- Solver Settings
305   user->phys->implicit  = implicit;
306   user->phys->state_var = state_var;
307 
308   // -- QFunction Context
309   newtonian_ig_ctx->lambda        = lambda;
310   newtonian_ig_ctx->mu            = mu;
311   newtonian_ig_ctx->k             = k;
312   newtonian_ig_ctx->cv            = cv;
313   newtonian_ig_ctx->cp            = cp;
314   newtonian_ig_ctx->c_tau         = c_tau;
315   newtonian_ig_ctx->Ctau_t        = Ctau_t;
316   newtonian_ig_ctx->Ctau_v        = Ctau_v;
317   newtonian_ig_ctx->Ctau_C        = Ctau_C;
318   newtonian_ig_ctx->Ctau_M        = Ctau_M;
319   newtonian_ig_ctx->Ctau_E        = Ctau_E;
320   newtonian_ig_ctx->stabilization = stab;
321   newtonian_ig_ctx->is_implicit   = implicit;
322   newtonian_ig_ctx->state_var     = state_var;
323   newtonian_ig_ctx->idl_enable    = idl_enable;
324   newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * second);
325   newtonian_ig_ctx->idl_start     = idl_start * meter;
326   newtonian_ig_ctx->idl_length    = idl_length * meter;
327   newtonian_ig_ctx->idl_pressure  = idl_pressure;
328   PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3));
329 
330   // -- Setup Context
331   setup_context->reference = reference;
332   setup_context->gas       = *newtonian_ig_ctx;
333   setup_context->lx        = domain_size[0];
334   setup_context->ly        = domain_size[1];
335   setup_context->lz        = domain_size[2];
336   setup_context->time      = 0;
337 
338   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context));
339   PetscCallCeed(ceed,
340                 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
341   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc));
342   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, "evaluation time", offsetof(struct SetupContext_, time), 1,
343                                                          "Time of evaluation"));
344 
345   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context));
346   PetscCallCeed(ceed,
347                 CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx));
348   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, FreeContextPetsc));
349   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1,
350                                                          "Size of timestep, delta t"));
351   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift",
352                                                          offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1,
353                                                          "Shift for mass matrix in IJacobian"));
354   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1,
355                                                          "Current solution time"));
356 
357   problem->apply_vol_rhs.qfunction_context = newtonian_ig_context;
358   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ifunction.qfunction_context));
359   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ijacobian.qfunction_context));
360   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow.qfunction_context));
361   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow_jacobian.qfunction_context));
362 
363   if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
364   if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
365   if (bc->num_slip > 0) PetscCall(SlipBCSetup(problem, dm, ctx, newtonian_ig_context));
366 
367   if (unit_tests) {
368     PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx));
369   }
370   PetscFunctionReturn(PETSC_SUCCESS);
371 }
372 
373 PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx) {
374   MPI_Comm                 comm = user->comm;
375   Ceed                     ceed = user->ceed;
376   NewtonianIdealGasContext newtonian_ctx;
377 
378   PetscFunctionBeginUser;
379   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ctx));
380   PetscCall(PetscPrintf(comm,
381                         "  Problem:\n"
382                         "    Problem Name                       : %s\n"
383                         "    Stabilization                      : %s\n",
384                         app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization]));
385   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ctx));
386   PetscFunctionReturn(PETSC_SUCCESS);
387 }
388