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