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