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