xref: /honee/problems/newtonian.c (revision 01e19bfa2194a1f80d34d4c13eafa992c6e37e28)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 
4 /// @file
5 /// Utility functions for setting up problems using the Newtonian Qfunction
6 
7 #include "../qfunctions/newtonian.h"
8 
9 #include <ceed.h>
10 #include <petscdm.h>
11 
12 #include <navierstokes.h>
13 
14 // For use with PetscOptionsEnum
15 static const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "ENTROPY", "StateVariable", "STATEVAR_", NULL};
16 
17 static PetscErrorCode CheckQWithTolerance(const CeedScalar Q_s[5], const CeedScalar Q_a[5], const CeedScalar Q_b[5], const char *name,
18                                           PetscReal rtol_0, PetscReal rtol_u, PetscReal rtol_4) {
19   CeedScalar relative_error[5];  // relative error
20   CeedScalar divisor_threshold = 10 * CEED_EPSILON;
21 
22   PetscFunctionBeginUser;
23   relative_error[0] = (Q_a[0] - Q_b[0]) / (fabs(Q_s[0]) > divisor_threshold ? Q_s[0] : 1);
24   relative_error[4] = (Q_a[4] - Q_b[4]) / (fabs(Q_s[4]) > divisor_threshold ? Q_s[4] : 1);
25 
26   CeedScalar u_magnitude = sqrt(Square(Q_s[1]) + Square(Q_s[2]) + Square(Q_s[3]));
27   CeedScalar u_divisor   = u_magnitude > divisor_threshold ? u_magnitude : 1;
28   for (int i = 1; i < 4; i++) {
29     relative_error[i] = (Q_a[i] - Q_b[i]) / u_divisor;
30   }
31 
32   if (fabs(relative_error[0]) >= rtol_0) {
33     printf("%s[0] error %g (expected %.10e, got %.10e)\n", name, relative_error[0], Q_s[0], Q_a[0]);
34   }
35   for (int i = 1; i < 4; i++) {
36     if (fabs(relative_error[i]) >= rtol_u) {
37       printf("%s[%d] error %g (expected %.10e, got %.10e)\n", name, i, relative_error[i], Q_s[i], Q_a[i]);
38     }
39   }
40   if (fabs(relative_error[4]) >= rtol_4) {
41     printf("%s[4] error %g (expected %.10e, got %.10e)\n", name, relative_error[4], Q_s[4], Q_a[4]);
42   }
43   PetscFunctionReturn(PETSC_SUCCESS);
44 }
45 
46 // @brief Verify `StateFromQ` by converting A0 -> B0 -> A0_test, where A0 should equal A0_test
47 static PetscErrorCode TestState(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5],
48                                 CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
49   CeedScalar        B0[5], A0_test[5];
50   char              buf[128];
51   const char *const StateVariables_Initial[] = {"U", "Y", "V"};
52 
53   PetscFunctionBeginUser;
54   const char *A_initial = StateVariables_Initial[state_var_A];
55   const char *B_initial = StateVariables_Initial[state_var_B];
56 
57   State state_A0 = StateFromQ(gas, A0, state_var_A);
58   StateToQ(gas, state_A0, B0, state_var_B);
59   State state_B0 = StateFromQ(gas, B0, state_var_B);
60   StateToQ(gas, state_B0, A0_test, state_var_A);
61 
62   snprintf(buf, sizeof buf, "%s->%s->%s: %s", A_initial, B_initial, A_initial, A_initial);
63   PetscCall(CheckQWithTolerance(A0, A0_test, A0, buf, rtol_0, rtol_u, rtol_4));
64   PetscFunctionReturn(PETSC_SUCCESS);
65 }
66 
67 // @brief Verify `StateFromQ_fwd` via a finite difference approximation
68 static PetscErrorCode TestState_fwd(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5],
69                                     CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
70   CeedScalar        eps = 4e-7;  // Finite difference step
71   char              buf[128];
72   const char *const StateVariables_Initial[] = {"U", "Y", "V"};
73 
74   PetscFunctionBeginUser;
75   const char *A_initial = StateVariables_Initial[state_var_A];
76   const char *B_initial = StateVariables_Initial[state_var_B];
77   State       state_0   = StateFromQ(gas, A0, state_var_A);
78 
79   for (int i = 0; i < 5; i++) {
80     CeedScalar dB[5] = {0.}, dB_fd[5] = {0.};
81     {  // Calculate dB using State functions
82       CeedScalar dA[5] = {0};
83 
84       dA[i]          = A0[i];
85       State dstate_0 = StateFromQ_fwd(gas, state_0, dA, state_var_A);
86       StateToQ_fwd(gas, state_0, dstate_0, dB, state_var_B);
87     }
88 
89     {  // Calculate dB_fd via finite difference approximation
90       CeedScalar A1[5], B0[5], B1[5];
91 
92       for (int j = 0; j < 5; j++) A1[j] = (1 + eps * (i == j)) * A0[j];
93       State state_1 = StateFromQ(gas, A1, state_var_A);
94       StateToQ(gas, state_0, B0, state_var_B);
95       StateToQ(gas, state_1, B1, state_var_B);
96       for (int j = 0; j < 5; j++) dB_fd[j] = (B1[j] - B0[j]) / eps;
97     }
98 
99     snprintf(buf, sizeof buf, "d%s->d%s: StateFrom%s_fwd i=%d: d%s", A_initial, B_initial, A_initial, i, B_initial);
100     PetscCall(CheckQWithTolerance(dB_fd, dB, dB_fd, buf, rtol_0, rtol_u, rtol_4));
101   }
102   PetscFunctionReturn(PETSC_SUCCESS);
103 }
104 
105 // @brief Test the Newtonian State transformation functions, `StateFrom*`
106 static PetscErrorCode UnitTests_Newtonian(Honee honee, NewtonianIdealGasContext gas) {
107   Units            units = honee->units;
108   const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, K = units->Kelvin;
109 
110   PetscFunctionBeginUser;
111   const CeedScalar T          = 200 * K;
112   const CeedScalar rho        = 1.2 * kg / Cube(m);
113   const CeedScalar P          = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T;
114   const CeedScalar u_base     = 40 * m / sec;
115   const CeedScalar u[3]       = {u_base, u_base * 1.1, u_base * 1.2};
116   const CeedScalar e_kinetic  = 0.5 * Dot3(u, u);
117   const CeedScalar e_internal = gas->cv * T;
118   const CeedScalar e_total    = e_kinetic + e_internal;
119   const CeedScalar gamma      = HeatCapacityRatio(gas);
120   const CeedScalar entropy    = log(P) - gamma * log(rho);
121   const CeedScalar rho_div_p  = rho / P;
122   const CeedScalar Y0[5]      = {P, u[0], u[1], u[2], T};
123   const CeedScalar U0[5]      = {rho, rho * u[0], rho * u[1], rho * u[2], rho * e_total};
124   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],
125                                  -rho_div_p};
126 
127   {
128     CeedScalar rtol = 20 * CEED_EPSILON;
129 
130     PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol));
131     PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol));
132     PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol));
133     PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, rtol, rtol, rtol));
134     PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, rtol, rtol, rtol));
135     PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, rtol, rtol, rtol));
136   }
137 
138   {
139     CeedScalar rtol = 5e-6;
140 
141     PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol));
142     PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol));
143     PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol));
144     PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, 10 * rtol, rtol, rtol));
145     PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, 5 * rtol, rtol, rtol));
146     PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, 5 * rtol, 5 * rtol, 5 * rtol));
147   }
148   PetscFunctionReturn(PETSC_SUCCESS);
149 }
150 
151 // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping
152 //
153 // Only used for SUPG stabilization
154 PetscErrorCode CreateKSPMassOperator_NewtonianStabilized(Honee honee, CeedOperator *op_mass) {
155   Ceed                 ceed = honee->ceed;
156   CeedInt              num_comp_q, q_data_size;
157   CeedQFunction        qf_mass;
158   CeedElemRestriction  elem_restr_q, elem_restr_qd;
159   CeedBasis            basis_q;
160   CeedVector           q_data;
161   CeedQFunctionContext qfctx = NULL;
162   PetscInt             dim   = 3;
163 
164   PetscFunctionBeginUser;
165   {  // Get restriction and basis from the RHS function
166     CeedOperator     *sub_ops;
167     CeedOperatorField op_field;
168     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
169 
170     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops));
171     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field));
172     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL));
173     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &op_field));
174     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_qd, NULL, &q_data));
175 
176     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qfctx));
177   }
178 
179   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
180   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd, &q_data_size));
181 
182   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Newtonian_Conserv, MassFunction_Newtonian_Conserv_loc, &qf_mass));
183 
184   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qfctx));
185   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
186   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
187   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
188   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
189   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
190   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
191 
192   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
193   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
194   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, honee->q_ceed));
195   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
196   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
197   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
198 
199   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
200   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
201   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
202   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
203   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qfctx));
204   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
205   PetscFunctionReturn(PETSC_SUCCESS);
206 }
207 
208 /**
209   @brief Create RHS CeedOperator for direct projection of divergence of diffusive flux
210 
211   @param[in]  honee          `Honee` context
212   @param[in]  diff_flux_proj `DivDiffFluxProjectionData` object
213   @param[out] op_rhs         Operator to calculate the RHS of the L^2 projection
214 **/
215 static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_NS(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) {
216   Ceed                 ceed       = honee->ceed;
217   NodalProjectionData  projection = diff_flux_proj->projection;
218   CeedInt              num_comp_q;
219   PetscInt             dim, label_value = 0;
220   DMLabel              domain_label    = NULL;
221   CeedQFunctionContext newtonian_qfctx = NULL;
222 
223   PetscFunctionBeginUser;
224   // -- Get Pre-requisite things
225   PetscCall(DMGetDimension(projection->dm, &dim));
226   PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q));
227 
228   {  // Get newtonian QF context
229     CeedOperator *sub_ops;
230     PetscInt      sub_op_index = 0;  // will be 0 for the volume op
231 
232     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_ifunction, &sub_ops));
233     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newtonian_qfctx));
234   }
235   PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_rhs));
236   {  // Add the volume integral CeedOperator
237     CeedQFunction       qf_rhs_volume;
238     CeedOperator        op_rhs_volume;
239     CeedVector          q_data;
240     CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL;
241     CeedBasis           basis_diff_flux = NULL;
242     CeedInt             q_data_size;
243 
244     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL));
245     PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data,
246                        &q_data_size));
247     switch (honee->phys->state_var) {
248       case STATEVAR_PRIMITIVE:
249         PetscCallCeed(ceed,
250                       CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Prim, DivDiffusiveFluxVolumeRHS_NS_Prim_loc, &qf_rhs_volume));
251         break;
252       case STATEVAR_CONSERVATIVE:
253         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Conserv, DivDiffusiveFluxVolumeRHS_NS_Conserv_loc,
254                                                         &qf_rhs_volume));
255         break;
256       case STATEVAR_ENTROPY:
257         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Entropy, DivDiffusiveFluxVolumeRHS_NS_Entropy_loc,
258                                                         &qf_rhs_volume));
259         break;
260     }
261 
262     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, newtonian_qfctx));
263     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "q", num_comp_q, CEED_EVAL_INTERP));
264     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
265     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE));
266     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD));
267 
268     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume));
269     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
270     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
271     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
272     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE));
273 
274     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_volume));
275 
276     PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
277     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
278     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume));
279     PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
280     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume));
281     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume));
282   }
283 
284   {  // Add the boundary integral CeedOperator
285     CeedQFunction qf_rhs_boundary;
286     DMLabel       face_sets_label;
287     PetscInt      num_face_set_values, *face_set_values;
288     CeedInt       q_data_size;
289 
290     // -- Build RHS operator
291     switch (honee->phys->state_var) {
292       case STATEVAR_PRIMITIVE:
293         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Prim, DivDiffusiveFluxBoundaryRHS_NS_Prim_loc,
294                                                         &qf_rhs_boundary));
295         break;
296       case STATEVAR_CONSERVATIVE:
297         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Conserv, DivDiffusiveFluxBoundaryRHS_NS_Conserv_loc,
298                                                         &qf_rhs_boundary));
299         break;
300       case STATEVAR_ENTROPY:
301         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Entropy, DivDiffusiveFluxBoundaryRHS_NS_Entropy_loc,
302                                                         &qf_rhs_boundary));
303         break;
304     }
305 
306     PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size));
307     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, newtonian_qfctx));
308     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "q", num_comp_q, CEED_EVAL_INTERP));
309     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
310     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE));
311     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP));
312 
313     PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label));
314     PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values));
315     for (PetscInt f = 0; f < num_face_set_values; f++) {
316       DMLabel  face_orientation_label;
317       PetscInt num_orientations_values, *orientation_values;
318 
319       {
320         char *face_orientation_label_name;
321 
322         PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name));
323         PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label));
324         PetscCall(PetscFree(face_orientation_label_name));
325       }
326       PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values));
327       for (PetscInt o = 0; o < num_orientations_values; o++) {
328         CeedOperator        op_rhs_boundary;
329         CeedBasis           basis_q, basis_diff_flux_boundary;
330         CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary;
331         CeedVector          q_data;
332         CeedInt             q_data_size;
333         PetscInt            orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1;
334 
335         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q));
336         PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q));
337         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0,
338                                                   &elem_restr_diff_flux_boundary));
339         PetscCall(CreateBasisFromPlex(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary));
340         PetscCall(
341             QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, honee->x_coord, &elem_restr_qdata, &q_data, &q_data_size));
342 
343         PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary));
344         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
345         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
346         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data));
347         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary,
348                                                  CEED_VECTOR_ACTIVE));
349 
350         PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_boundary));
351 
352         PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary));
353         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata));
354         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
355         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary));
356         PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
357         PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary));
358         PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
359       }
360       PetscCall(PetscFree(orientation_values));
361     }
362     PetscCall(PetscFree(face_set_values));
363     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary));
364   }
365 
366   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newtonian_qfctx));
367   PetscFunctionReturn(PETSC_SUCCESS);
368 }
369 
370 /**
371   @brief Create RHS CeedOperator for indirect projection of divergence of diffusive flux
372 
373   @param[in]  honee          `Honee` context
374   @param[in]  diff_flux_proj `DivDiffFluxProjectionData` object
375   @param[out] op_rhs         Operator to calculate the RHS of the L^2 projection
376 **/
377 static PetscErrorCode DivDiffFluxProjectionCreateRHS_Indirect_NS(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) {
378   Ceed                 ceed       = honee->ceed;
379   NodalProjectionData  projection = diff_flux_proj->projection;
380   CeedBasis            basis_diff_flux;
381   CeedElemRestriction  elem_restr_diff_flux, elem_restr_qd;
382   CeedVector           q_data;
383   CeedInt              num_comp_q, q_data_size;
384   PetscInt             dim;
385   PetscInt             label_value = 0, height = 0, dm_field = 0;
386   DMLabel              domain_label = NULL;
387   CeedQFunction        qf_rhs;
388   CeedQFunctionContext newtonian_qfctx = NULL;
389 
390   PetscFunctionBeginUser;
391   PetscCall(DMGetDimension(projection->dm, &dim));
392   PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q));
393 
394   {  // Get newtonian QF context
395     CeedOperator *sub_ops;
396     PetscInt      sub_op_index = 0;  // will be 0 for the volume op
397 
398     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_ifunction, &sub_ops));
399     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newtonian_qfctx));
400   }
401   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux));
402   PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux));
403   PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data,
404                      &q_data_size));
405 
406   switch (honee->phys->state_var) {
407     case STATEVAR_PRIMITIVE:
408       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Prim, DiffusiveFluxRHS_NS_Prim_loc, &qf_rhs));
409       break;
410     case STATEVAR_CONSERVATIVE:
411       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Conserv, DiffusiveFluxRHS_NS_Conserv_loc, &qf_rhs));
412       break;
413     case STATEVAR_ENTROPY:
414       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Entropy, DiffusiveFluxRHS_NS_Entropy_loc, &qf_rhs));
415       break;
416   }
417 
418   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, newtonian_qfctx));
419   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "q", num_comp_q, CEED_EVAL_INTERP));
420   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
421   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE));
422   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP));
423 
424   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, op_rhs));
425   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
426   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
427   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
428   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
429 
430   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs));
431   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newtonian_qfctx));
432   PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
433   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
434   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
435   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux));
436   PetscFunctionReturn(PETSC_SUCCESS);
437 }
438 
439 PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
440   SetupContext             setup_context;
441   Honee                    honee  = *(Honee *)ctx;
442   CeedInt                  degree = honee->app_ctx->degree;
443   StabilizationType        stab;
444   StateVariable            state_var;
445   MPI_Comm                 comm = honee->comm;
446   Ceed                     ceed = honee->ceed;
447   PetscBool                implicit;
448   PetscBool                unit_tests;
449   NewtonianIdealGasContext newtonian_ig_ctx;
450   CeedQFunctionContext     newtonian_ig_qfctx;
451 
452   PetscFunctionBeginUser;
453   PetscCall(PetscCalloc1(1, &setup_context));
454   PetscCall(PetscCalloc1(1, &newtonian_ig_ctx));
455 
456   // ------------------------------------------------------
457   //           Setup Generic Newtonian IG Problem
458   // ------------------------------------------------------
459   problem->jac_data_size_vol            = 14;
460   problem->jac_data_size_sur            = 11;
461   problem->compute_exact_solution_error = PETSC_FALSE;
462   problem->print_info                   = PRINT_NEWTONIAN;
463 
464   PetscCall(DivDiffFluxProjectionCreate(honee, 4, &honee->diff_flux_proj));
465   if (honee->diff_flux_proj) {
466     DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj;
467     NodalProjectionData       projection     = diff_flux_proj->projection;
468 
469     diff_flux_proj->CreateRHSOperator_Direct   = DivDiffFluxProjectionCreateRHS_Direct_NS;
470     diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_NS;
471 
472     switch (honee->diff_flux_proj->method) {
473       case DIV_DIFF_FLUX_PROJ_DIRECT: {
474         PetscSection section;
475 
476         PetscCall(DMGetLocalSection(projection->dm, &section));
477         PetscCall(PetscSectionSetFieldName(section, 0, ""));
478         PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_MomentumX"));
479         PetscCall(PetscSectionSetComponentName(section, 0, 1, "DivDiffusiveFlux_MomentumY"));
480         PetscCall(PetscSectionSetComponentName(section, 0, 2, "DivDiffusiveFlux_MomentumZ"));
481         PetscCall(PetscSectionSetComponentName(section, 0, 3, "DivDiffusiveFlux_Energy"));
482       } break;
483       case DIV_DIFF_FLUX_PROJ_INDIRECT: {
484         PetscSection section;
485 
486         PetscCall(DMGetLocalSection(projection->dm, &section));
487         PetscCall(PetscSectionSetFieldName(section, 0, ""));
488         PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_MomentumXX"));
489         PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_MomentumXY"));
490         PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_MomentumXZ"));
491         PetscCall(PetscSectionSetComponentName(section, 0, 3, "DiffusiveFlux_MomentumYX"));
492         PetscCall(PetscSectionSetComponentName(section, 0, 4, "DiffusiveFlux_MomentumYY"));
493         PetscCall(PetscSectionSetComponentName(section, 0, 5, "DiffusiveFlux_MomentumYZ"));
494         PetscCall(PetscSectionSetComponentName(section, 0, 6, "DiffusiveFlux_MomentumZX"));
495         PetscCall(PetscSectionSetComponentName(section, 0, 7, "DiffusiveFlux_MomentumZY"));
496         PetscCall(PetscSectionSetComponentName(section, 0, 8, "DiffusiveFlux_MomentumZZ"));
497         PetscCall(PetscSectionSetComponentName(section, 0, 9, "DiffusiveFlux_EnergyX"));
498         PetscCall(PetscSectionSetComponentName(section, 0, 10, "DiffusiveFlux_EnergyY"));
499         PetscCall(PetscSectionSetComponentName(section, 0, 11, "DiffusiveFlux_EnergyZ"));
500       } break;
501       case DIV_DIFF_FLUX_PROJ_NONE:
502         SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
503                 DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]);
504         break;
505     }
506   }
507 
508   // ------------------------------------------------------
509   //             Create the libCEED context
510   // ------------------------------------------------------
511   CeedScalar cv         = 717.;          // J/(kg K)
512   CeedScalar cp         = 1004.;         // J/(kg K)
513   CeedScalar g[3]       = {0, 0, 0};     // m/s^2
514   CeedScalar lambda     = -2. / 3.;      // -
515   CeedScalar mu         = 1.8e-5;        // Pa s, dynamic viscosity
516   CeedScalar k          = 0.02638;       // W/(m K)
517   CeedScalar c_tau      = 0.5 / degree;  // -
518   CeedScalar Ctau_t     = 1.0;           // -
519   CeedScalar Cv_func[3] = {36, 60, 128};
520   CeedScalar Ctau_v     = Cv_func[(CeedInt)Min(3, degree) - 1];
521   CeedScalar Ctau_C     = 0.25 / degree;
522   CeedScalar Ctau_M     = 0.25 / degree;
523   CeedScalar Ctau_E     = 0.125;
524   PetscReal  domain_min[3], domain_max[3], domain_size[3];
525   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
526   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
527 
528   StatePrimitive reference      = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15};
529   CeedScalar     idl_decay_time = -1, idl_start = 0, idl_length = 0, idl_pressure = reference.pressure;
530   PetscBool      idl_enable = PETSC_FALSE;
531 
532   // ------------------------------------------------------
533   //             Create the PETSc context
534   // ------------------------------------------------------
535   PetscScalar meter    = 1;  // 1 meter in scaled length units
536   PetscScalar kilogram = 1;  // 1 kilogram in scaled mass units
537   PetscScalar second   = 1;  // 1 second in scaled time units
538   PetscScalar Kelvin   = 1;  // 1 Kelvin in scaled temperature units
539   PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s;
540 
541   // ------------------------------------------------------
542   //              Command line Options
543   // ------------------------------------------------------
544   PetscBool given_option = PETSC_FALSE;
545   PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL);
546   // -- Conservative vs Primitive variables
547   PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE),
548                              (PetscEnum *)&state_var, NULL));
549 
550   switch (state_var) {
551     case STATEVAR_CONSERVATIVE:
552       problem->ics.qf_func_ptr                   = ICsNewtonianIG_Conserv;
553       problem->ics.qf_loc                        = ICsNewtonianIG_Conserv_loc;
554       problem->apply_vol_rhs.qf_func_ptr         = RHSFunction_Newtonian;
555       problem->apply_vol_rhs.qf_loc              = RHSFunction_Newtonian_loc;
556       problem->apply_vol_ifunction.qf_func_ptr   = IFunction_Newtonian_Conserv;
557       problem->apply_vol_ifunction.qf_loc        = IFunction_Newtonian_Conserv_loc;
558       problem->apply_vol_ijacobian.qf_func_ptr   = IJacobian_Newtonian_Conserv;
559       problem->apply_vol_ijacobian.qf_loc        = IJacobian_Newtonian_Conserv_loc;
560       problem->apply_inflow.qf_func_ptr          = BoundaryIntegral_Conserv;
561       problem->apply_inflow.qf_loc               = BoundaryIntegral_Conserv_loc;
562       problem->apply_inflow_jacobian.qf_func_ptr = BoundaryIntegral_Jacobian_Conserv;
563       problem->apply_inflow_jacobian.qf_loc      = BoundaryIntegral_Jacobian_Conserv_loc;
564       break;
565     case STATEVAR_PRIMITIVE:
566       problem->ics.qf_func_ptr                   = ICsNewtonianIG_Prim;
567       problem->ics.qf_loc                        = ICsNewtonianIG_Prim_loc;
568       problem->apply_vol_ifunction.qf_func_ptr   = IFunction_Newtonian_Prim;
569       problem->apply_vol_ifunction.qf_loc        = IFunction_Newtonian_Prim_loc;
570       problem->apply_vol_ijacobian.qf_func_ptr   = IJacobian_Newtonian_Prim;
571       problem->apply_vol_ijacobian.qf_loc        = IJacobian_Newtonian_Prim_loc;
572       problem->apply_inflow.qf_func_ptr          = BoundaryIntegral_Prim;
573       problem->apply_inflow.qf_loc               = BoundaryIntegral_Prim_loc;
574       problem->apply_inflow_jacobian.qf_func_ptr = BoundaryIntegral_Jacobian_Prim;
575       problem->apply_inflow_jacobian.qf_loc      = BoundaryIntegral_Jacobian_Prim_loc;
576       break;
577     case STATEVAR_ENTROPY:
578       problem->ics.qf_func_ptr                   = ICsNewtonianIG_Entropy;
579       problem->ics.qf_loc                        = ICsNewtonianIG_Entropy_loc;
580       problem->apply_vol_ifunction.qf_func_ptr   = IFunction_Newtonian_Entropy;
581       problem->apply_vol_ifunction.qf_loc        = IFunction_Newtonian_Entropy_loc;
582       problem->apply_vol_ijacobian.qf_func_ptr   = IJacobian_Newtonian_Entropy;
583       problem->apply_vol_ijacobian.qf_loc        = IJacobian_Newtonian_Entropy_loc;
584       problem->apply_inflow.qf_func_ptr          = BoundaryIntegral_Entropy;
585       problem->apply_inflow.qf_loc               = BoundaryIntegral_Entropy_loc;
586       problem->apply_inflow_jacobian.qf_func_ptr = BoundaryIntegral_Jacobian_Entropy;
587       problem->apply_inflow_jacobian.qf_loc      = BoundaryIntegral_Jacobian_Entropy_loc;
588       break;
589   }
590 
591   // -- Physics
592   PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL));
593   PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL));
594   PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL));
595   PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL));
596   PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL));
597 
598   PetscInt dim = 3;
599   PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL));
600   PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option));
601   if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim);
602 
603   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
604   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
605   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
606   PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL));
607   PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL));
608   PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL));
609   PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL));
610   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
611   PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL));
612 
613   dim = 3;
614   PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL));
615   PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL));
616   PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL));
617 
618   // -- Units
619   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
620   meter = fabs(meter);
621   PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL));
622   kilogram = fabs(kilogram);
623   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
624   second = fabs(second);
625   PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL));
626   Kelvin = fabs(Kelvin);
627 
628   // -- Warnings
629   PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP,
630              "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n");
631   PetscCheck(!(honee->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE && !implicit), comm, PETSC_ERR_SUP,
632              "Projection of divergence of diffusive flux is not implemented for Navier-Stokes with explicit timestepping");
633 
634   PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point",
635                                NULL, idl_decay_time, &idl_decay_time, &idl_enable));
636   PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero.");
637   if (idl_decay_time < 0) idl_enable = PETSC_FALSE;
638   if (idl_enable) problem->jac_data_size_vol++;
639   PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL));
640   PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL));
641   idl_pressure = reference.pressure;
642   PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-reference_pressure`)", NULL, idl_pressure,
643                                &idl_pressure, NULL));
644   PetscOptionsEnd();
645 
646   if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized;
647 
648   // ------------------------------------------------------
649   //           Set up the PETSc context
650   // ------------------------------------------------------
651   // -- Define derived units
652   Pascal          = kilogram / (meter * PetscSqr(second));
653   J_per_kg_K      = PetscSqr(meter) / (PetscSqr(second) * Kelvin);
654   m_per_squared_s = meter / PetscSqr(second);
655   W_per_m_K       = kilogram * meter / (pow(second, 3) * Kelvin);
656 
657   honee->units->meter           = meter;
658   honee->units->kilogram        = kilogram;
659   honee->units->second          = second;
660   honee->units->Kelvin          = Kelvin;
661   honee->units->Pascal          = Pascal;
662   honee->units->J_per_kg_K      = J_per_kg_K;
663   honee->units->m_per_squared_s = m_per_squared_s;
664   honee->units->W_per_m_K       = W_per_m_K;
665 
666   // ------------------------------------------------------
667   //           Set up the libCEED context
668   // ------------------------------------------------------
669   // -- Scale variables to desired units
670   cv *= J_per_kg_K;
671   cp *= J_per_kg_K;
672   mu *= Pascal * second;
673   k *= W_per_m_K;
674   for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter;
675   for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s;
676   reference.pressure *= Pascal;
677   for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second;
678   reference.temperature *= Kelvin;
679 
680   // -- Solver Settings
681   honee->phys->implicit  = implicit;
682   honee->phys->state_var = state_var;
683 
684   // -- QFunction Context
685   newtonian_ig_ctx->lambda          = lambda;
686   newtonian_ig_ctx->mu              = mu;
687   newtonian_ig_ctx->k               = k;
688   newtonian_ig_ctx->cv              = cv;
689   newtonian_ig_ctx->cp              = cp;
690   newtonian_ig_ctx->c_tau           = c_tau;
691   newtonian_ig_ctx->Ctau_t          = Ctau_t;
692   newtonian_ig_ctx->Ctau_v          = Ctau_v;
693   newtonian_ig_ctx->Ctau_C          = Ctau_C;
694   newtonian_ig_ctx->Ctau_M          = Ctau_M;
695   newtonian_ig_ctx->Ctau_E          = Ctau_E;
696   newtonian_ig_ctx->stabilization   = stab;
697   newtonian_ig_ctx->is_implicit     = implicit;
698   newtonian_ig_ctx->state_var       = state_var;
699   newtonian_ig_ctx->idl_enable      = idl_enable;
700   newtonian_ig_ctx->idl_amplitude   = 1 / (idl_decay_time * second);
701   newtonian_ig_ctx->idl_start       = idl_start * meter;
702   newtonian_ig_ctx->idl_length      = idl_length * meter;
703   newtonian_ig_ctx->idl_pressure    = idl_pressure;
704   newtonian_ig_ctx->divFdiff_method = honee->app_ctx->divFdiffproj_method;
705   PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3));
706 
707   // -- Setup Context
708   setup_context->reference = reference;
709   setup_context->gas       = *newtonian_ig_ctx;
710   setup_context->lx        = domain_size[0];
711   setup_context->ly        = domain_size[1];
712   setup_context->lz        = domain_size[2];
713   setup_context->time      = 0;
714 
715   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &problem->ics.qfctx));
716   PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
717   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc));
718   PetscCallCeed(
719       ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfctx, "evaluation time", offsetof(struct SetupContext_, time), 1, "Time of evaluation"));
720 
721   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &newtonian_ig_qfctx));
722   PetscCallCeed(ceed, CeedQFunctionContextSetData(newtonian_ig_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx));
723   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_qfctx, CEED_MEM_HOST, FreeContextPetsc));
724   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1,
725                                                          "Size of timestep, delta t"));
726   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "ijacobian time shift",
727                                                          offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1,
728                                                          "Shift for mass matrix in IJacobian"));
729   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1,
730                                                          "Current solution time"));
731 
732   problem->apply_vol_rhs.qfctx = newtonian_ig_qfctx;
733   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ifunction.qfctx));
734   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ijacobian.qfctx));
735   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_inflow.qfctx));
736   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_inflow_jacobian.qfctx));
737 
738   if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
739   if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
740   if (bc->num_slip > 0) PetscCall(SlipBCSetup(problem, dm, ctx, newtonian_ig_qfctx));
741 
742   if (unit_tests) {
743     PetscCall(UnitTests_Newtonian(honee, newtonian_ig_ctx));
744   }
745   PetscFunctionReturn(PETSC_SUCCESS);
746 }
747 
748 PetscErrorCode PRINT_NEWTONIAN(Honee honee, ProblemData problem, AppCtx app_ctx) {
749   MPI_Comm                 comm = honee->comm;
750   Ceed                     ceed = honee->ceed;
751   NewtonianIdealGasContext newtonian_ctx;
752 
753   PetscFunctionBeginUser;
754   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ctx));
755   PetscCall(PetscPrintf(comm,
756                         "  Problem:\n"
757                         "    Problem Name                       : %s\n"
758                         "    Stabilization                      : %s\n",
759                         app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization]));
760   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ctx));
761   PetscFunctionReturn(PETSC_SUCCESS);
762 }
763