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