xref: /honee/problems/newtonian.c (revision a32db64d340db16914d4892be21e91c50f2a7cbd)
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", "ENTROPY", "StateVariable", "STATEVAR_", NULL};
20 
21 static PetscErrorCode CheckQWithTolerance(const CeedScalar Q_s[5], const CeedScalar Q_a[5], const CeedScalar Q_b[5], const char *name,
22                                           PetscReal rtol_0, PetscReal rtol_u, PetscReal rtol_4) {
23   CeedScalar relative_error[5];  // relative error
24   CeedScalar divisor_threshold = 10 * CEED_EPSILON;
25 
26   PetscFunctionBeginUser;
27   relative_error[0] = (Q_a[0] - Q_b[0]) / (fabs(Q_s[0]) > divisor_threshold ? Q_s[0] : 1);
28   relative_error[4] = (Q_a[4] - Q_b[4]) / (fabs(Q_s[4]) > divisor_threshold ? Q_s[4] : 1);
29 
30   CeedScalar u_magnitude = sqrt(Square(Q_s[1]) + Square(Q_s[2]) + Square(Q_s[3]));
31   CeedScalar u_divisor   = u_magnitude > divisor_threshold ? u_magnitude : 1;
32   for (int i = 1; i < 4; i++) {
33     relative_error[i] = (Q_a[i] - Q_b[i]) / u_divisor;
34   }
35 
36   if (fabs(relative_error[0]) >= rtol_0) {
37     printf("%s[0] error %g (expected %.10e, got %.10e)\n", name, relative_error[0], Q_s[0], Q_a[0]);
38   }
39   for (int i = 1; i < 4; i++) {
40     if (fabs(relative_error[i]) >= rtol_u) {
41       printf("%s[%d] error %g (expected %.10e, got %.10e)\n", name, i, relative_error[i], Q_s[i], Q_a[i]);
42     }
43   }
44   if (fabs(relative_error[4]) >= rtol_4) {
45     printf("%s[4] error %g (expected %.10e, got %.10e)\n", name, relative_error[4], Q_s[4], Q_a[4]);
46   }
47   PetscFunctionReturn(PETSC_SUCCESS);
48 }
49 
50 // @brief Verify `StateFromQ` by converting A0 -> B0 -> A0_test, where A0 should equal A0_test
51 static PetscErrorCode TestState(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5],
52                                 CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
53   CeedScalar        B0[5], A0_test[5];
54   char              buf[128];
55   const char *const StateVariables_Initial[] = {"U", "Y", "V"};
56 
57   PetscFunctionBeginUser;
58   const char *A_initial = StateVariables_Initial[state_var_A];
59   const char *B_initial = StateVariables_Initial[state_var_B];
60 
61   State state_A0 = StateFromQ(gas, A0, state_var_A);
62   StateToQ(gas, state_A0, B0, state_var_B);
63   State state_B0 = StateFromQ(gas, B0, state_var_B);
64   StateToQ(gas, state_B0, A0_test, state_var_A);
65 
66   snprintf(buf, sizeof buf, "%s->%s->%s: %s", A_initial, B_initial, A_initial, A_initial);
67   PetscCall(CheckQWithTolerance(A0, A0_test, A0, buf, rtol_0, rtol_u, rtol_4));
68   PetscFunctionReturn(PETSC_SUCCESS);
69 }
70 
71 // @brief Verify `StateFromQ_fwd` via a finite difference approximation
72 static PetscErrorCode TestState_fwd(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5],
73                                     CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
74   CeedScalar        eps = 4e-7;  // Finite difference step
75   char              buf[128];
76   const char *const StateVariables_Initial[] = {"U", "Y", "V"};
77 
78   PetscFunctionBeginUser;
79   const char *A_initial = StateVariables_Initial[state_var_A];
80   const char *B_initial = StateVariables_Initial[state_var_B];
81   State       state_0   = StateFromQ(gas, A0, state_var_A);
82 
83   for (int i = 0; i < 5; i++) {
84     CeedScalar dB[5] = {0.}, dB_fd[5] = {0.};
85     {  // Calculate dB using State functions
86       CeedScalar dA[5] = {0};
87 
88       dA[i]          = A0[i];
89       State dstate_0 = StateFromQ_fwd(gas, state_0, dA, state_var_A);
90       StateToQ_fwd(gas, state_0, dstate_0, dB, state_var_B);
91     }
92 
93     {  // Calculate dB_fd via finite difference approximation
94       CeedScalar A1[5], B0[5], B1[5];
95 
96       for (int j = 0; j < 5; j++) A1[j] = (1 + eps * (i == j)) * A0[j];
97       State state_1 = StateFromQ(gas, A1, state_var_A);
98       StateToQ(gas, state_0, B0, state_var_B);
99       StateToQ(gas, state_1, B1, state_var_B);
100       for (int j = 0; j < 5; j++) dB_fd[j] = (B1[j] - B0[j]) / eps;
101     }
102 
103     snprintf(buf, sizeof buf, "d%s->d%s: StateFrom%s_fwd i=%d: d%s", A_initial, B_initial, A_initial, i, B_initial);
104     PetscCall(CheckQWithTolerance(dB_fd, dB, dB_fd, buf, rtol_0, rtol_u, rtol_4));
105   }
106   PetscFunctionReturn(PETSC_SUCCESS);
107 }
108 
109 // @brief Test the Newtonian State transformation functions, `StateFrom*`
110 static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) {
111   Units            units = user->units;
112   const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, K = units->Kelvin;
113 
114   PetscFunctionBeginUser;
115   const CeedScalar T          = 200 * K;
116   const CeedScalar rho        = 1.2 * kg / Cube(m);
117   const CeedScalar P          = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T;
118   const CeedScalar u_base     = 40 * m / sec;
119   const CeedScalar u[3]       = {u_base, u_base * 1.1, u_base * 1.2};
120   const CeedScalar e_kinetic  = 0.5 * Dot3(u, u);
121   const CeedScalar e_internal = gas->cv * T;
122   const CeedScalar e_total    = e_kinetic + e_internal;
123   const CeedScalar gamma      = HeatCapacityRatio(gas);
124   const CeedScalar entropy    = log(P) - gamma * log(rho);
125   const CeedScalar rho_div_p  = rho / P;
126   const CeedScalar Y0[5]      = {P, u[0], u[1], u[2], T};
127   const CeedScalar U0[5]      = {rho, rho * u[0], rho * u[1], rho * u[2], rho * e_total};
128   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],
129                                  -rho_div_p};
130 
131   {
132     CeedScalar rtol = 20 * CEED_EPSILON;
133 
134     PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol));
135     PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol));
136     PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol));
137     PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, rtol, rtol, rtol));
138     PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, rtol, rtol, rtol));
139     PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, rtol, rtol, rtol));
140   }
141 
142   {
143     CeedScalar rtol = 5e-6;
144 
145     PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol));
146     PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol));
147     PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol));
148     PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, 10 * rtol, rtol, rtol));
149     PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, 5 * rtol, rtol, rtol));
150     PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, 5 * rtol, 5 * rtol, 5 * rtol));
151   }
152   PetscFunctionReturn(PETSC_SUCCESS);
153 }
154 
155 // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping
156 //
157 // Only used for SUPG stabilization
158 PetscErrorCode CreateKSPMassOperator_NewtonianStabilized(User user, CeedOperator *op_mass) {
159   Ceed                 ceed = user->ceed;
160   CeedInt              num_comp_q, q_data_size;
161   CeedQFunction        qf_mass;
162   CeedElemRestriction  elem_restr_q, elem_restr_qd_i;
163   CeedBasis            basis_q;
164   CeedVector           q_data;
165   CeedQFunctionContext qf_ctx = NULL;
166   PetscInt             dim    = 3;
167 
168   PetscFunctionBeginUser;
169   {  // Get restriction and basis from the RHS function
170     CeedOperator     *sub_ops;
171     CeedOperatorField field;
172     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
173 
174     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops));
175     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
176     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
177     PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q));
178 
179     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field));
180     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i));
181     PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data));
182 
183     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx));
184   }
185 
186   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
187   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size));
188 
189   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Newtonian_Conserv, MassFunction_Newtonian_Conserv_loc, &qf_mass));
190 
191   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx));
192   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
193   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
194   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
195   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
196   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
197   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
198 
199   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
200   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
201   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed));
202   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data));
203   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
204   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
205 
206   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
207   PetscFunctionReturn(PETSC_SUCCESS);
208 }
209 
210 PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
211   SetupContext             setup_context;
212   User                     user   = *(User *)ctx;
213   CeedInt                  degree = user->app_ctx->degree;
214   StabilizationType        stab;
215   StateVariable            state_var;
216   MPI_Comm                 comm = user->comm;
217   Ceed                     ceed = user->ceed;
218   PetscBool                implicit;
219   PetscBool                unit_tests;
220   NewtonianIdealGasContext newtonian_ig_ctx;
221   CeedQFunctionContext     newtonian_ig_context;
222 
223   PetscFunctionBeginUser;
224   PetscCall(PetscCalloc1(1, &setup_context));
225   PetscCall(PetscCalloc1(1, &newtonian_ig_ctx));
226 
227   // ------------------------------------------------------
228   //           Setup Generic Newtonian IG Problem
229   // ------------------------------------------------------
230   problem->dim                          = 3;
231   problem->jac_data_size_sur            = 11;
232   problem->compute_exact_solution_error = PETSC_FALSE;
233   problem->print_info                   = PRINT_NEWTONIAN;
234   problem->uses_newtonian               = PETSC_TRUE;
235 
236   // ------------------------------------------------------
237   //             Create the libCEED context
238   // ------------------------------------------------------
239   CeedScalar cv         = 717.;          // J/(kg K)
240   CeedScalar cp         = 1004.;         // J/(kg K)
241   CeedScalar g[3]       = {0, 0, 0};     // m/s^2
242   CeedScalar lambda     = -2. / 3.;      // -
243   CeedScalar mu         = 1.8e-5;        // Pa s, dynamic viscosity
244   CeedScalar k          = 0.02638;       // W/(m K)
245   CeedScalar c_tau      = 0.5 / degree;  // -
246   CeedScalar Ctau_t     = 1.0;           // -
247   CeedScalar Cv_func[3] = {36, 60, 128};
248   CeedScalar Ctau_v     = Cv_func[(CeedInt)Min(3, degree) - 1];
249   CeedScalar Ctau_C     = 0.25 / degree;
250   CeedScalar Ctau_M     = 0.25 / degree;
251   CeedScalar Ctau_E     = 0.125;
252   PetscReal  domain_min[3], domain_max[3], domain_size[3];
253   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
254   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
255 
256   StatePrimitive reference      = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15};
257   CeedScalar     idl_decay_time = -1, idl_start = 0, idl_length = 0, idl_pressure = reference.pressure;
258   PetscBool      idl_enable = PETSC_FALSE;
259 
260   // ------------------------------------------------------
261   //             Create the PETSc context
262   // ------------------------------------------------------
263   PetscScalar meter    = 1;  // 1 meter in scaled length units
264   PetscScalar kilogram = 1;  // 1 kilogram in scaled mass units
265   PetscScalar second   = 1;  // 1 second in scaled time units
266   PetscScalar Kelvin   = 1;  // 1 Kelvin in scaled temperature units
267   PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s;
268 
269   // ------------------------------------------------------
270   //              Command line Options
271   // ------------------------------------------------------
272   PetscBool given_option = PETSC_FALSE;
273   PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL);
274   // -- Conservative vs Primitive variables
275   PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE),
276                              (PetscEnum *)&state_var, NULL));
277 
278   switch (state_var) {
279     case STATEVAR_CONSERVATIVE:
280       problem->ics.qfunction                       = ICsNewtonianIG_Conserv;
281       problem->ics.qfunction_loc                   = ICsNewtonianIG_Conserv_loc;
282       problem->apply_vol_rhs.qfunction             = RHSFunction_Newtonian;
283       problem->apply_vol_rhs.qfunction_loc         = RHSFunction_Newtonian_loc;
284       problem->apply_vol_ifunction.qfunction       = IFunction_Newtonian_Conserv;
285       problem->apply_vol_ifunction.qfunction_loc   = IFunction_Newtonian_Conserv_loc;
286       problem->apply_vol_ijacobian.qfunction       = IJacobian_Newtonian_Conserv;
287       problem->apply_vol_ijacobian.qfunction_loc   = IJacobian_Newtonian_Conserv_loc;
288       problem->apply_inflow.qfunction              = BoundaryIntegral_Conserv;
289       problem->apply_inflow.qfunction_loc          = BoundaryIntegral_Conserv_loc;
290       problem->apply_inflow_jacobian.qfunction     = BoundaryIntegral_Jacobian_Conserv;
291       problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc;
292       break;
293 
294     case STATEVAR_PRIMITIVE:
295       problem->ics.qfunction                       = ICsNewtonianIG_Prim;
296       problem->ics.qfunction_loc                   = ICsNewtonianIG_Prim_loc;
297       problem->apply_vol_ifunction.qfunction       = IFunction_Newtonian_Prim;
298       problem->apply_vol_ifunction.qfunction_loc   = IFunction_Newtonian_Prim_loc;
299       problem->apply_vol_ijacobian.qfunction       = IJacobian_Newtonian_Prim;
300       problem->apply_vol_ijacobian.qfunction_loc   = IJacobian_Newtonian_Prim_loc;
301       problem->apply_inflow.qfunction              = BoundaryIntegral_Prim;
302       problem->apply_inflow.qfunction_loc          = BoundaryIntegral_Prim_loc;
303       problem->apply_inflow_jacobian.qfunction     = BoundaryIntegral_Jacobian_Prim;
304       problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc;
305       break;
306     case STATEVAR_ENTROPY:
307       problem->ics.qfunction                       = ICsNewtonianIG_Entropy;
308       problem->ics.qfunction_loc                   = ICsNewtonianIG_Entropy_loc;
309       problem->apply_vol_ifunction.qfunction       = IFunction_Newtonian_Entropy;
310       problem->apply_vol_ifunction.qfunction_loc   = IFunction_Newtonian_Entropy_loc;
311       problem->apply_vol_ijacobian.qfunction       = IJacobian_Newtonian_Entropy;
312       problem->apply_vol_ijacobian.qfunction_loc   = IJacobian_Newtonian_Entropy_loc;
313       problem->apply_inflow.qfunction              = BoundaryIntegral_Entropy;
314       problem->apply_inflow.qfunction_loc          = BoundaryIntegral_Entropy_loc;
315       problem->apply_inflow_jacobian.qfunction     = BoundaryIntegral_Jacobian_Entropy;
316       problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Entropy_loc;
317       break;
318   }
319 
320   // -- Physics
321   PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL));
322   PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL));
323   PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL));
324   PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL));
325   PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL));
326 
327   PetscInt dim = problem->dim;
328   PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL));
329   PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option));
330   if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim);
331 
332   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
333   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
334   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
335   PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL));
336   PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL));
337   PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL));
338   PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL));
339   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
340   PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL));
341 
342   dim = 3;
343   PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL));
344   PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL));
345   PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL));
346 
347   // -- Units
348   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
349   meter = fabs(meter);
350   PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL));
351   kilogram = fabs(kilogram);
352   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
353   second = fabs(second);
354   PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL));
355   Kelvin = fabs(Kelvin);
356 
357   // -- Warnings
358   PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP,
359              "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n");
360 
361   PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point",
362                                NULL, idl_decay_time, &idl_decay_time, &idl_enable));
363   PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero.");
364   if (idl_decay_time < 0) idl_enable = PETSC_FALSE;
365   PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL));
366   PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL));
367   idl_pressure = reference.pressure;
368   PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-reference_pressure`)", NULL, idl_pressure,
369                                &idl_pressure, NULL));
370   PetscOptionsEnd();
371 
372   if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized;
373 
374   // ------------------------------------------------------
375   //           Set up the PETSc context
376   // ------------------------------------------------------
377   // -- Define derived units
378   Pascal          = kilogram / (meter * PetscSqr(second));
379   J_per_kg_K      = PetscSqr(meter) / (PetscSqr(second) * Kelvin);
380   m_per_squared_s = meter / PetscSqr(second);
381   W_per_m_K       = kilogram * meter / (pow(second, 3) * Kelvin);
382 
383   user->units->meter           = meter;
384   user->units->kilogram        = kilogram;
385   user->units->second          = second;
386   user->units->Kelvin          = Kelvin;
387   user->units->Pascal          = Pascal;
388   user->units->J_per_kg_K      = J_per_kg_K;
389   user->units->m_per_squared_s = m_per_squared_s;
390   user->units->W_per_m_K       = W_per_m_K;
391 
392   // ------------------------------------------------------
393   //           Set up the libCEED context
394   // ------------------------------------------------------
395   // -- Scale variables to desired units
396   cv *= J_per_kg_K;
397   cp *= J_per_kg_K;
398   mu *= Pascal * second;
399   k *= W_per_m_K;
400   for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter;
401   for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s;
402   reference.pressure *= Pascal;
403   for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second;
404   reference.temperature *= Kelvin;
405   problem->dm_scale = meter;
406 
407   // -- Solver Settings
408   user->phys->implicit  = implicit;
409   user->phys->state_var = state_var;
410 
411   // -- QFunction Context
412   newtonian_ig_ctx->lambda        = lambda;
413   newtonian_ig_ctx->mu            = mu;
414   newtonian_ig_ctx->k             = k;
415   newtonian_ig_ctx->cv            = cv;
416   newtonian_ig_ctx->cp            = cp;
417   newtonian_ig_ctx->c_tau         = c_tau;
418   newtonian_ig_ctx->Ctau_t        = Ctau_t;
419   newtonian_ig_ctx->Ctau_v        = Ctau_v;
420   newtonian_ig_ctx->Ctau_C        = Ctau_C;
421   newtonian_ig_ctx->Ctau_M        = Ctau_M;
422   newtonian_ig_ctx->Ctau_E        = Ctau_E;
423   newtonian_ig_ctx->stabilization = stab;
424   newtonian_ig_ctx->is_implicit   = implicit;
425   newtonian_ig_ctx->state_var     = state_var;
426   newtonian_ig_ctx->idl_enable    = idl_enable;
427   newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * second);
428   newtonian_ig_ctx->idl_start     = idl_start * meter;
429   newtonian_ig_ctx->idl_length    = idl_length * meter;
430   newtonian_ig_ctx->idl_pressure  = idl_pressure;
431   PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3));
432 
433   // -- Setup Context
434   setup_context->reference = reference;
435   setup_context->gas       = *newtonian_ig_ctx;
436   setup_context->lx        = domain_size[0];
437   setup_context->ly        = domain_size[1];
438   setup_context->lz        = domain_size[2];
439   setup_context->time      = 0;
440 
441   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context));
442   PetscCallCeed(ceed,
443                 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
444   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc));
445   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, "evaluation time", offsetof(struct SetupContext_, time), 1,
446                                                          "Time of evaluation"));
447 
448   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context));
449   PetscCallCeed(ceed,
450                 CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx));
451   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, FreeContextPetsc));
452   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1,
453                                                          "Size of timestep, delta t"));
454   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift",
455                                                          offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1,
456                                                          "Shift for mass matrix in IJacobian"));
457   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1,
458                                                          "Current solution time"));
459 
460   problem->apply_vol_rhs.qfunction_context = newtonian_ig_context;
461   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ifunction.qfunction_context));
462   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ijacobian.qfunction_context));
463   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow.qfunction_context));
464   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow_jacobian.qfunction_context));
465 
466   if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
467   if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
468   if (bc->num_slip > 0) PetscCall(SlipBCSetup(problem, dm, ctx, newtonian_ig_context));
469 
470   if (unit_tests) {
471     PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx));
472   }
473   PetscFunctionReturn(PETSC_SUCCESS);
474 }
475 
476 PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx) {
477   MPI_Comm                 comm = user->comm;
478   Ceed                     ceed = user->ceed;
479   NewtonianIdealGasContext newtonian_ctx;
480 
481   PetscFunctionBeginUser;
482   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ctx));
483   PetscCall(PetscPrintf(comm,
484                         "  Problem:\n"
485                         "    Problem Name                       : %s\n"
486                         "    Stabilization                      : %s\n",
487                         app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization]));
488   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ctx));
489   PetscFunctionReturn(PETSC_SUCCESS);
490 }
491