xref: /honee/problems/newtonian.c (revision bff204bcbef0cb2ac109c0e2b33a4eb8a241421f)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
33a8779fbSJames Wright 
43a8779fbSJames Wright /// @file
53a8779fbSJames Wright /// Utility functions for setting up problems using the Newtonian Qfunction
63a8779fbSJames Wright 
73a8779fbSJames Wright #include "../qfunctions/newtonian.h"
83a8779fbSJames Wright 
9e419654dSJeremy L Thompson #include <ceed.h>
10e419654dSJeremy L Thompson #include <petscdm.h>
11e419654dSJeremy L Thompson 
12149fb536SJames Wright #include <navierstokes.h>
136dd99bedSLeila Ghaffari 
14e5a8cae0SJames Wright const char *const StateVariables[]     = {"CONSERVATIVE", "PRIMITIVE", "ENTROPY", "StateVariable", "STATEVAR_", NULL};
15e5a8cae0SJames Wright const char *const StabilizationTypes[] = {"NONE", "SU", "SUPG", "StabilizationType", "STAB_", NULL};
16e5a8cae0SJames Wright 
1748417fb7SJames Wright static PetscErrorCode UnitTests_Newtonian(Honee honee, NewtonianIdealGasContext gas);
18f0b65372SJed Brown 
199b05e62eSJames Wright static PetscErrorCode PRINT_NEWTONIAN(Honee honee, ProblemData problem, AppCtx app_ctx) {
209b05e62eSJames Wright   MPI_Comm                 comm = honee->comm;
219b05e62eSJames Wright   Ceed                     ceed = honee->ceed;
229b05e62eSJames Wright   NewtonianIdealGasContext newtonian_ctx;
239b05e62eSJames Wright 
249b05e62eSJames Wright   PetscFunctionBeginUser;
259b05e62eSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ctx));
269b05e62eSJames Wright   PetscCall(PetscPrintf(comm,
279b05e62eSJames Wright                         "  Problem:\n"
289b05e62eSJames Wright                         "    Problem Name                       : %s\n"
299b05e62eSJames Wright                         "    Stabilization                      : %s\n",
309b05e62eSJames Wright                         app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization]));
319b05e62eSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ctx));
329b05e62eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
339b05e62eSJames Wright }
349b05e62eSJames Wright 
3565dee3d2SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping
3665dee3d2SJames Wright //
3765dee3d2SJames Wright // Only used for SUPG stabilization
380c373b74SJames Wright PetscErrorCode CreateKSPMassOperator_NewtonianStabilized(Honee honee, CeedOperator *op_mass) {
390c373b74SJames Wright   Ceed                 ceed = honee->ceed;
4065dee3d2SJames Wright   CeedInt              num_comp_q, q_data_size;
4165dee3d2SJames Wright   CeedQFunction        qf_mass;
4201e19bfaSJames Wright   CeedElemRestriction  elem_restr_q, elem_restr_qd;
4365dee3d2SJames Wright   CeedBasis            basis_q;
4465dee3d2SJames Wright   CeedVector           q_data;
45236a8ba6SJames Wright   CeedQFunctionContext qfctx = NULL;
4665dee3d2SJames Wright   PetscInt             dim   = 3;
4765dee3d2SJames Wright 
4865dee3d2SJames Wright   PetscFunctionBeginUser;
4965dee3d2SJames Wright   {  // Get restriction and basis from the RHS function
5065dee3d2SJames Wright     CeedOperator     *sub_ops;
5101e19bfaSJames Wright     CeedOperatorField op_field;
5265dee3d2SJames Wright     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
5365dee3d2SJames Wright 
540c373b74SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops));
5501e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field));
5601e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL));
5701e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &op_field));
5801e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_qd, NULL, &q_data));
5965dee3d2SJames Wright 
60236a8ba6SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qfctx));
6165dee3d2SJames Wright   }
6265dee3d2SJames Wright 
6365dee3d2SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
6401e19bfaSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd, &q_data_size));
6565dee3d2SJames Wright 
6665dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Newtonian_Conserv, MassFunction_Newtonian_Conserv_loc, &qf_mass));
6765dee3d2SJames Wright 
68236a8ba6SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qfctx));
6965dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
7065dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
7165dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
7265dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
7365dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
7465dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
7565dee3d2SJames Wright 
7665dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
7771f2ed29SJames Wright   PetscCallCeed(ceed, CeedOperatorSetName(*op_mass, "RHS Mass Operator, Newtonian Stabilized"));
7865dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
790c373b74SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, honee->q_ceed));
8001e19bfaSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
8165dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
8265dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
8365dee3d2SJames Wright 
84fff85bd3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
8501e19bfaSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
86fff85bd3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
87fff85bd3SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
88236a8ba6SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qfctx));
8965dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
9065dee3d2SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
9165dee3d2SJames Wright }
924c07ec22SJames Wright 
9336038bbcSJames Wright /**
9436038bbcSJames Wright   @brief Create RHS CeedOperator for direct projection of divergence of diffusive flux
9536038bbcSJames Wright 
960c373b74SJames Wright   @param[in]  honee          `Honee` context
9736038bbcSJames Wright   @param[in]  diff_flux_proj `DivDiffFluxProjectionData` object
9836038bbcSJames Wright   @param[out] op_rhs         Operator to calculate the RHS of the L^2 projection
9936038bbcSJames Wright **/
100e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_NS(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) {
1010c373b74SJames Wright   Ceed                 ceed       = honee->ceed;
10236038bbcSJames Wright   NodalProjectionData  projection = diff_flux_proj->projection;
10336038bbcSJames Wright   CeedInt              num_comp_q;
10436038bbcSJames Wright   PetscInt             dim, label_value = 0;
10536038bbcSJames Wright   DMLabel              domain_label    = NULL;
106236a8ba6SJames Wright   CeedQFunctionContext newtonian_qfctx = NULL;
10736038bbcSJames Wright 
10836038bbcSJames Wright   PetscFunctionBeginUser;
10936038bbcSJames Wright   // -- Get Pre-requisite things
11036038bbcSJames Wright   PetscCall(DMGetDimension(projection->dm, &dim));
111e3663b90SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q));
11236038bbcSJames Wright 
1130880fbb6SJames Wright   {  // Get newtonian QF context
114b3b24828SJames Wright     CeedOperator *sub_ops, main_op = honee->op_ifunction ? honee->op_ifunction : honee->op_rhs_ctx->op;
11536038bbcSJames Wright     PetscInt      sub_op_index = 0;  // will be 0 for the volume op
11636038bbcSJames Wright 
117b3b24828SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(main_op, &sub_ops));
118236a8ba6SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newtonian_qfctx));
11936038bbcSJames Wright   }
12036038bbcSJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_rhs));
12136038bbcSJames Wright   {  // Add the volume integral CeedOperator
12236038bbcSJames Wright     CeedQFunction       qf_rhs_volume;
12336038bbcSJames Wright     CeedOperator        op_rhs_volume;
12436038bbcSJames Wright     CeedVector          q_data;
1250880fbb6SJames Wright     CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL;
1260880fbb6SJames Wright     CeedBasis           basis_diff_flux = NULL;
12736038bbcSJames Wright     CeedInt             q_data_size;
12836038bbcSJames Wright 
1290880fbb6SJames Wright     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL));
130e3663b90SJames Wright     PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data,
131e3663b90SJames Wright                        &q_data_size));
1320c373b74SJames Wright     switch (honee->phys->state_var) {
13336038bbcSJames Wright       case STATEVAR_PRIMITIVE:
13436038bbcSJames Wright         PetscCallCeed(ceed,
13536038bbcSJames Wright                       CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Prim, DivDiffusiveFluxVolumeRHS_NS_Prim_loc, &qf_rhs_volume));
13636038bbcSJames Wright         break;
13736038bbcSJames Wright       case STATEVAR_CONSERVATIVE:
13836038bbcSJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Conserv, DivDiffusiveFluxVolumeRHS_NS_Conserv_loc,
13936038bbcSJames Wright                                                         &qf_rhs_volume));
14036038bbcSJames Wright         break;
14136038bbcSJames Wright       case STATEVAR_ENTROPY:
14236038bbcSJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Entropy, DivDiffusiveFluxVolumeRHS_NS_Entropy_loc,
14336038bbcSJames Wright                                                         &qf_rhs_volume));
14436038bbcSJames Wright         break;
14536038bbcSJames Wright     }
14636038bbcSJames Wright 
147236a8ba6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, newtonian_qfctx));
14836038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "q", num_comp_q, CEED_EVAL_INTERP));
14936038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
15036038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE));
15136038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD));
15236038bbcSJames Wright 
15336038bbcSJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume));
154e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
155e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
15636038bbcSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
15736038bbcSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE));
15836038bbcSJames Wright 
15936038bbcSJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_volume));
16036038bbcSJames Wright 
16136038bbcSJames Wright     PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
16236038bbcSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
1630880fbb6SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume));
1640880fbb6SJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
16536038bbcSJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume));
16636038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume));
16736038bbcSJames Wright   }
16836038bbcSJames Wright 
16936038bbcSJames Wright   {  // Add the boundary integral CeedOperator
17036038bbcSJames Wright     CeedQFunction qf_rhs_boundary;
17136038bbcSJames Wright     DMLabel       face_sets_label;
17236038bbcSJames Wright     PetscInt      num_face_set_values, *face_set_values;
17336038bbcSJames Wright     CeedInt       q_data_size;
17436038bbcSJames Wright 
17536038bbcSJames Wright     // -- Build RHS operator
1760c373b74SJames Wright     switch (honee->phys->state_var) {
17736038bbcSJames Wright       case STATEVAR_PRIMITIVE:
17836038bbcSJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Prim, DivDiffusiveFluxBoundaryRHS_NS_Prim_loc,
17936038bbcSJames Wright                                                         &qf_rhs_boundary));
18036038bbcSJames Wright         break;
18136038bbcSJames Wright       case STATEVAR_CONSERVATIVE:
18236038bbcSJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Conserv, DivDiffusiveFluxBoundaryRHS_NS_Conserv_loc,
18336038bbcSJames Wright                                                         &qf_rhs_boundary));
18436038bbcSJames Wright         break;
18536038bbcSJames Wright       case STATEVAR_ENTROPY:
18636038bbcSJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Entropy, DivDiffusiveFluxBoundaryRHS_NS_Entropy_loc,
18736038bbcSJames Wright                                                         &qf_rhs_boundary));
18836038bbcSJames Wright         break;
18936038bbcSJames Wright     }
19036038bbcSJames Wright 
1910c373b74SJames Wright     PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size));
192236a8ba6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, newtonian_qfctx));
19336038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "q", num_comp_q, CEED_EVAL_INTERP));
19436038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
19536038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE));
19636038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP));
19736038bbcSJames Wright 
19836038bbcSJames Wright     PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label));
19936038bbcSJames Wright     PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values));
20036038bbcSJames Wright     for (PetscInt f = 0; f < num_face_set_values; f++) {
20136038bbcSJames Wright       DMLabel  face_orientation_label;
20236038bbcSJames Wright       PetscInt num_orientations_values, *orientation_values;
20336038bbcSJames Wright 
20436038bbcSJames Wright       {
20536038bbcSJames Wright         char *face_orientation_label_name;
20636038bbcSJames Wright 
20736038bbcSJames Wright         PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name));
20836038bbcSJames Wright         PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label));
20936038bbcSJames Wright         PetscCall(PetscFree(face_orientation_label_name));
21036038bbcSJames Wright       }
21136038bbcSJames Wright       PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values));
21236038bbcSJames Wright       for (PetscInt o = 0; o < num_orientations_values; o++) {
21336038bbcSJames Wright         CeedOperator        op_rhs_boundary;
21436038bbcSJames Wright         CeedBasis           basis_q, basis_diff_flux_boundary;
21536038bbcSJames Wright         CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary;
21636038bbcSJames Wright         CeedVector          q_data;
21736038bbcSJames Wright         CeedInt             q_data_size;
21836038bbcSJames Wright         PetscInt            orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1;
21936038bbcSJames Wright 
2200c373b74SJames Wright         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q));
2210c373b74SJames Wright         PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q));
22236038bbcSJames Wright         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0,
22336038bbcSJames Wright                                                   &elem_restr_diff_flux_boundary));
2240880fbb6SJames Wright         PetscCall(CreateBasisFromPlex(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary));
225e3663b90SJames Wright         PetscCall(
226e3663b90SJames Wright             QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, honee->x_coord, &elem_restr_qdata, &q_data, &q_data_size));
22736038bbcSJames Wright 
22836038bbcSJames Wright         PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary));
22936038bbcSJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
23036038bbcSJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
23136038bbcSJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data));
23236038bbcSJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary,
23336038bbcSJames Wright                                                  CEED_VECTOR_ACTIVE));
23436038bbcSJames Wright 
23536038bbcSJames Wright         PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_boundary));
23636038bbcSJames Wright 
23736038bbcSJames Wright         PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary));
23836038bbcSJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata));
23936038bbcSJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
24036038bbcSJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary));
24136038bbcSJames Wright         PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
24236038bbcSJames Wright         PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary));
24336038bbcSJames Wright         PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
24436038bbcSJames Wright       }
24536038bbcSJames Wright       PetscCall(PetscFree(orientation_values));
24636038bbcSJames Wright     }
24736038bbcSJames Wright     PetscCall(PetscFree(face_set_values));
24836038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary));
24936038bbcSJames Wright   }
25036038bbcSJames Wright 
251236a8ba6SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newtonian_qfctx));
25236038bbcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
25336038bbcSJames Wright }
25436038bbcSJames Wright 
25536038bbcSJames Wright /**
25636038bbcSJames Wright   @brief Create RHS CeedOperator for indirect projection of divergence of diffusive flux
25736038bbcSJames Wright 
2580c373b74SJames Wright   @param[in]  honee          `Honee` context
25936038bbcSJames Wright   @param[in]  diff_flux_proj `DivDiffFluxProjectionData` object
26036038bbcSJames Wright   @param[out] op_rhs         Operator to calculate the RHS of the L^2 projection
26136038bbcSJames Wright **/
262e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Indirect_NS(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) {
2630c373b74SJames Wright   Ceed                 ceed       = honee->ceed;
26436038bbcSJames Wright   NodalProjectionData  projection = diff_flux_proj->projection;
26536038bbcSJames Wright   CeedBasis            basis_diff_flux;
26636038bbcSJames Wright   CeedElemRestriction  elem_restr_diff_flux, elem_restr_qd;
26736038bbcSJames Wright   CeedVector           q_data;
26836038bbcSJames Wright   CeedInt              num_comp_q, q_data_size;
26936038bbcSJames Wright   PetscInt             dim;
27036038bbcSJames Wright   PetscInt             label_value = 0, height = 0, dm_field = 0;
27136038bbcSJames Wright   DMLabel              domain_label = NULL;
27236038bbcSJames Wright   CeedQFunction        qf_rhs;
273236a8ba6SJames Wright   CeedQFunctionContext newtonian_qfctx = NULL;
27436038bbcSJames Wright 
27536038bbcSJames Wright   PetscFunctionBeginUser;
27636038bbcSJames Wright   PetscCall(DMGetDimension(projection->dm, &dim));
277e3663b90SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q));
27836038bbcSJames Wright 
27940b78511SJames Wright   {  // Get newtonian QF context
280b3b24828SJames Wright     CeedOperator *sub_ops, main_op = honee->op_ifunction ? honee->op_ifunction : honee->op_rhs_ctx->op;
28136038bbcSJames Wright     PetscInt      sub_op_index = 0;  // will be 0 for the volume op
28236038bbcSJames Wright 
283b3b24828SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(main_op, &sub_ops));
284236a8ba6SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newtonian_qfctx));
28536038bbcSJames Wright   }
28636038bbcSJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux));
28736038bbcSJames Wright   PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux));
288e3663b90SJames Wright   PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data,
289e3663b90SJames Wright                      &q_data_size));
29036038bbcSJames Wright 
2910c373b74SJames Wright   switch (honee->phys->state_var) {
29236038bbcSJames Wright     case STATEVAR_PRIMITIVE:
29336038bbcSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Prim, DiffusiveFluxRHS_NS_Prim_loc, &qf_rhs));
29436038bbcSJames Wright       break;
29536038bbcSJames Wright     case STATEVAR_CONSERVATIVE:
29636038bbcSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Conserv, DiffusiveFluxRHS_NS_Conserv_loc, &qf_rhs));
29736038bbcSJames Wright       break;
29836038bbcSJames Wright     case STATEVAR_ENTROPY:
29936038bbcSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Entropy, DiffusiveFluxRHS_NS_Entropy_loc, &qf_rhs));
30036038bbcSJames Wright       break;
30136038bbcSJames Wright   }
30236038bbcSJames Wright 
303236a8ba6SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, newtonian_qfctx));
30436038bbcSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "q", num_comp_q, CEED_EVAL_INTERP));
30536038bbcSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
30636038bbcSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE));
30736038bbcSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP));
30836038bbcSJames Wright 
30936038bbcSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, op_rhs));
310e3663b90SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
311e3663b90SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
31236038bbcSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
31336038bbcSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
31436038bbcSJames Wright 
31536038bbcSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs));
316236a8ba6SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newtonian_qfctx));
31736038bbcSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
31836038bbcSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
31936038bbcSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
32036038bbcSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux));
32136038bbcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
32236038bbcSJames Wright }
32336038bbcSJames Wright 
324d6cac220SJames Wright static PetscErrorCode BoundaryIntegralBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) {
325d6cac220SJames Wright   HoneeBCStruct honee_bc;
326d6cac220SJames Wright 
327d6cac220SJames Wright   PetscFunctionBeginUser;
328d6cac220SJames Wright   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
329d6cac220SJames Wright   Honee honee = honee_bc->honee;
330d6cac220SJames Wright 
331d6cac220SJames Wright   switch (honee->phys->state_var) {
332d6cac220SJames Wright     case STATEVAR_CONSERVATIVE:
333d6cac220SJames Wright       PetscCall(HoneeBCCreateIFunctionQF(bc_def, BoundaryIntegral_Conserv, BoundaryIntegral_Conserv_loc, honee_bc->qfctx, qf));
334d6cac220SJames Wright       break;
335d6cac220SJames Wright     case STATEVAR_PRIMITIVE:
336d6cac220SJames Wright       PetscCall(HoneeBCCreateIFunctionQF(bc_def, BoundaryIntegral_Prim, BoundaryIntegral_Prim_loc, honee_bc->qfctx, qf));
337d6cac220SJames Wright       break;
338d6cac220SJames Wright     case STATEVAR_ENTROPY:
339d6cac220SJames Wright       PetscCall(HoneeBCCreateIFunctionQF(bc_def, BoundaryIntegral_Entropy, BoundaryIntegral_Entropy_loc, honee_bc->qfctx, qf));
340d6cac220SJames Wright       break;
341d6cac220SJames Wright   }
342d6cac220SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
343d6cac220SJames Wright }
344d6cac220SJames Wright 
345d6cac220SJames Wright static PetscErrorCode BoundaryIntegralBCSetup_CreateIJacobianQF(BCDefinition bc_def, CeedQFunction *qf) {
346d6cac220SJames Wright   HoneeBCStruct honee_bc;
347d6cac220SJames Wright 
348d6cac220SJames Wright   PetscFunctionBeginUser;
349d6cac220SJames Wright   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
350d6cac220SJames Wright   Honee honee = honee_bc->honee;
351d6cac220SJames Wright 
352d6cac220SJames Wright   switch (honee->phys->state_var) {
353d6cac220SJames Wright     case STATEVAR_CONSERVATIVE:
354d6cac220SJames Wright       PetscCall(HoneeBCCreateIJacobianQF(bc_def, BoundaryIntegral_Jacobian_Conserv, BoundaryIntegral_Jacobian_Conserv_loc, honee_bc->qfctx, qf));
355d6cac220SJames Wright       break;
356d6cac220SJames Wright     case STATEVAR_PRIMITIVE:
357d6cac220SJames Wright       PetscCall(HoneeBCCreateIJacobianQF(bc_def, BoundaryIntegral_Jacobian_Prim, BoundaryIntegral_Jacobian_Prim_loc, honee_bc->qfctx, qf));
358d6cac220SJames Wright       break;
359d6cac220SJames Wright     case STATEVAR_ENTROPY:
360d6cac220SJames Wright       PetscCall(HoneeBCCreateIJacobianQF(bc_def, BoundaryIntegral_Jacobian_Entropy, BoundaryIntegral_Jacobian_Entropy_loc, honee_bc->qfctx, qf));
361d6cac220SJames Wright       break;
362d6cac220SJames Wright   }
363d6cac220SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
364d6cac220SJames Wright }
365d6cac220SJames Wright 
366d3c60affSJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx) {
367b7f03f12SJed Brown   SetupContext             setup_context;
3680c373b74SJames Wright   Honee                    honee  = *(Honee *)ctx;
3690c373b74SJames Wright   CeedInt                  degree = honee->app_ctx->degree;
3703a8779fbSJames Wright   StabilizationType        stab;
3713636f6a4SJames Wright   StateVariable            state_var;
3720c373b74SJames Wright   MPI_Comm                 comm = honee->comm;
3730c373b74SJames Wright   Ceed                     ceed = honee->ceed;
3743a8779fbSJames Wright   PetscBool                implicit;
3750d8cd818SJames Wright   PetscBool                unit_tests;
37615a3537eSJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
377e07531f7SJames Wright   CeedQFunctionContext     newtonian_ig_qfctx;
3783a8779fbSJames Wright 
37915a3537eSJed Brown   PetscFunctionBeginUser;
380*bff204bcSJames Wright   // Option Defaults
381*bff204bcSJames Wright   StatePrimitive reference      = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15};
382*bff204bcSJames Wright   CeedScalar     idl_decay_time = -1;
383*bff204bcSJames Wright   {
384*bff204bcSJames Wright     PetscCall(PetscNew(&newtonian_ig_ctx));
385*bff204bcSJames Wright     const CeedScalar                 Cv_func[3]                = {36, 60, 128};
386*bff204bcSJames Wright     struct NewtonianIdealGasContext_ newtonian_ig_ctx_defaults = {
387*bff204bcSJames Wright         .cv           = 717., // J/(kg K)
388*bff204bcSJames Wright         .cp           = 1004., // J/(kg K)
389*bff204bcSJames Wright         .g            = {0, 0, 0}, // m/s^2
390*bff204bcSJames Wright         .lambda       = -2. / 3., // -
391*bff204bcSJames Wright         .mu           = 1.8e-5, // Pa s, dynamic viscosity
392*bff204bcSJames Wright         .k            = 0.02638, // W/(m K)
393*bff204bcSJames Wright         .c_tau        = 0.5 / degree, // -
394*bff204bcSJames Wright         .Ctau_t       = 1.0, // -
395*bff204bcSJames Wright         .Ctau_v       = Cv_func[(CeedInt)Min(3, degree) - 1],
396*bff204bcSJames Wright         .Ctau_C       = 0.25 / degree,
397*bff204bcSJames Wright         .Ctau_M       = 0.25 / degree,
398*bff204bcSJames Wright         .Ctau_E       = 0.125,
399*bff204bcSJames Wright         .idl_start    = 0,
400*bff204bcSJames Wright         .idl_length   = 0,
401*bff204bcSJames Wright         .idl_pressure = reference.pressure,
402*bff204bcSJames Wright         .idl_enable   = PETSC_FALSE,
403*bff204bcSJames Wright     };
404*bff204bcSJames Wright     *newtonian_ig_ctx = newtonian_ig_ctx_defaults;
405*bff204bcSJames Wright   }
406*bff204bcSJames Wright 
4078d8b1f0cSJames Wright   PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL);
4088d8b1f0cSJames Wright   PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE),
4098d8b1f0cSJames Wright                              (PetscEnum *)&state_var, NULL));
4103a8779fbSJames Wright 
4118d8b1f0cSJames Wright   // Newtonian fluid properties
412*bff204bcSJames Wright   PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, newtonian_ig_ctx->cv, &newtonian_ig_ctx->cv, NULL));
413*bff204bcSJames Wright   PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, newtonian_ig_ctx->cp, &newtonian_ig_ctx->cp, NULL));
414*bff204bcSJames Wright   PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, newtonian_ig_ctx->lambda, &newtonian_ig_ctx->lambda,
415*bff204bcSJames Wright                                NULL));
416*bff204bcSJames Wright   PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, newtonian_ig_ctx->mu, &newtonian_ig_ctx->mu, NULL));
417*bff204bcSJames Wright   PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, newtonian_ig_ctx->k, &newtonian_ig_ctx->k, NULL));
4188d8b1f0cSJames Wright 
4198d8b1f0cSJames Wright   PetscInt  dim          = 3;
420*bff204bcSJames Wright   PetscBool given_option = PETSC_FALSE;
4218d8b1f0cSJames Wright   PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL));
422*bff204bcSJames Wright   PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, newtonian_ig_ctx->g, &dim, &given_option));
4238d8b1f0cSJames Wright   if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim);
4248d8b1f0cSJames Wright 
4258d8b1f0cSJames Wright   // Stabilization parameters
4268d8b1f0cSJames Wright   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
427*bff204bcSJames Wright   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, newtonian_ig_ctx->c_tau, &newtonian_ig_ctx->c_tau, NULL));
428*bff204bcSJames Wright   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, newtonian_ig_ctx->Ctau_t, &newtonian_ig_ctx->Ctau_t, NULL));
429*bff204bcSJames Wright   PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, newtonian_ig_ctx->Ctau_v, &newtonian_ig_ctx->Ctau_v, NULL));
430*bff204bcSJames Wright   PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, newtonian_ig_ctx->Ctau_C, &newtonian_ig_ctx->Ctau_C, NULL));
431*bff204bcSJames Wright   PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, newtonian_ig_ctx->Ctau_M, &newtonian_ig_ctx->Ctau_M, NULL));
432*bff204bcSJames Wright   PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, newtonian_ig_ctx->Ctau_E, &newtonian_ig_ctx->Ctau_E, NULL));
4338d8b1f0cSJames Wright 
4348d8b1f0cSJames Wright   dim = 3;
4358d8b1f0cSJames Wright   PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL));
4368d8b1f0cSJames Wright   PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL));
4378d8b1f0cSJames Wright   PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL));
4388d8b1f0cSJames Wright 
4398d8b1f0cSJames Wright   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
4408d8b1f0cSJames Wright   PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL));
4418d8b1f0cSJames Wright   PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP,
4428d8b1f0cSJames Wright              "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n");
4438d8b1f0cSJames Wright 
4448d8b1f0cSJames Wright   // IDL Settings
445*bff204bcSJames Wright   {
446*bff204bcSJames Wright     PetscBool idl_enable = (PetscBool)newtonian_ig_ctx->idl_enable;  // Need PetscBool variable to read in from PetscOptionsScalar()
4478d8b1f0cSJames Wright     PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point",
4488d8b1f0cSJames Wright                                  NULL, idl_decay_time, &idl_decay_time, &idl_enable));
4498d8b1f0cSJames Wright     PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero.");
4508d8b1f0cSJames Wright     if (idl_decay_time < 0) idl_enable = PETSC_FALSE;
451*bff204bcSJames Wright     newtonian_ig_ctx->idl_enable = idl_enable;
452*bff204bcSJames Wright 
453*bff204bcSJames Wright     PetscCall(
454*bff204bcSJames Wright         PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, newtonian_ig_ctx->idl_start, &newtonian_ig_ctx->idl_start, NULL));
455*bff204bcSJames Wright     PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, newtonian_ig_ctx->idl_length,
456*bff204bcSJames Wright                                  &newtonian_ig_ctx->idl_length, NULL));
457*bff204bcSJames Wright     newtonian_ig_ctx->idl_pressure = reference.pressure;
458*bff204bcSJames Wright     PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-reference_pressure`)", NULL,
459*bff204bcSJames Wright                                  newtonian_ig_ctx->idl_pressure, &newtonian_ig_ctx->idl_pressure, NULL));
460*bff204bcSJames Wright   }
4618d8b1f0cSJames Wright   PetscOptionsEnd();
4628d8b1f0cSJames Wright 
4638d8b1f0cSJames Wright   // Set problem information
4641abc2837SJames Wright   problem->num_comps_jac_data = 14;
465*bff204bcSJames Wright   if (newtonian_ig_ctx->idl_enable) problem->num_comps_jac_data += 1;
46658ce1233SJames Wright   problem->compute_exact_solution_error = PETSC_FALSE;
467cbe60e31SLeila Ghaffari   problem->print_info                   = PRINT_NEWTONIAN;
468f27c2204SJames Wright   problem->num_components               = 5;
469f27c2204SJames Wright   PetscCall(PetscMalloc1(problem->num_components, &problem->component_names));
4708d8b1f0cSJames Wright   static const char *const conservative_component_names[] = {"Density", "MomentumX", "MomentumY", "MomentumZ", "TotalEnergy"};
4718d8b1f0cSJames Wright   static const char *const primitive_component_names[]    = {"Pressure", "VelocityX", "VelocityY", "VelocityZ", "Temperature"};
4728d8b1f0cSJames Wright   static const char *const entropy_component_names[]      = {"EntropyDensity", "EntropyMomentumX", "EntropyMomentumY", "EntropyMomentumZ",
4738d8b1f0cSJames Wright                                                              "EntropyTotalEnergy"};
4748d8b1f0cSJames Wright 
4758d8b1f0cSJames Wright   switch (state_var) {
476f27c2204SJames Wright     case STATEVAR_CONSERVATIVE:
4778d8b1f0cSJames Wright       problem->ics.qf_func_ptr                 = ICsNewtonianIG_Conserv;
4788d8b1f0cSJames Wright       problem->ics.qf_loc                      = ICsNewtonianIG_Conserv_loc;
4798d8b1f0cSJames Wright       problem->apply_vol_rhs.qf_func_ptr       = RHSFunction_Newtonian;
4808d8b1f0cSJames Wright       problem->apply_vol_rhs.qf_loc            = RHSFunction_Newtonian_loc;
4818d8b1f0cSJames Wright       problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Conserv;
4828d8b1f0cSJames Wright       problem->apply_vol_ifunction.qf_loc      = IFunction_Newtonian_Conserv_loc;
4838d8b1f0cSJames Wright       problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Conserv;
4848d8b1f0cSJames Wright       problem->apply_vol_ijacobian.qf_loc      = IJacobian_Newtonian_Conserv_loc;
485f27c2204SJames Wright       for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(conservative_component_names[i], &problem->component_names[i]));
486f27c2204SJames Wright       break;
487f27c2204SJames Wright     case STATEVAR_PRIMITIVE:
4888d8b1f0cSJames Wright       problem->ics.qf_func_ptr                 = ICsNewtonianIG_Prim;
4898d8b1f0cSJames Wright       problem->ics.qf_loc                      = ICsNewtonianIG_Prim_loc;
4908d8b1f0cSJames Wright       problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Prim;
4918d8b1f0cSJames Wright       problem->apply_vol_ifunction.qf_loc      = IFunction_Newtonian_Prim_loc;
4928d8b1f0cSJames Wright       problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Prim;
4938d8b1f0cSJames Wright       problem->apply_vol_ijacobian.qf_loc      = IJacobian_Newtonian_Prim_loc;
494f27c2204SJames Wright       for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(primitive_component_names[i], &problem->component_names[i]));
495f27c2204SJames Wright       break;
496f27c2204SJames Wright     case STATEVAR_ENTROPY:
4978d8b1f0cSJames Wright       problem->ics.qf_func_ptr                 = ICsNewtonianIG_Entropy;
4988d8b1f0cSJames Wright       problem->ics.qf_loc                      = ICsNewtonianIG_Entropy_loc;
4998d8b1f0cSJames Wright       problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Entropy;
5008d8b1f0cSJames Wright       problem->apply_vol_ifunction.qf_loc      = IFunction_Newtonian_Entropy_loc;
5018d8b1f0cSJames Wright       problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Entropy;
5028d8b1f0cSJames Wright       problem->apply_vol_ijacobian.qf_loc      = IJacobian_Newtonian_Entropy_loc;
503f27c2204SJames Wright       for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(entropy_component_names[i], &problem->component_names[i]));
504f27c2204SJames Wright       break;
505f27c2204SJames Wright   }
506f27c2204SJames Wright 
5078d8b1f0cSJames Wright   if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized;
5088d8b1f0cSJames Wright 
5090c373b74SJames Wright   PetscCall(DivDiffFluxProjectionCreate(honee, 4, &honee->diff_flux_proj));
5100c373b74SJames Wright   if (honee->diff_flux_proj) {
5110c373b74SJames Wright     DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj;
51236038bbcSJames Wright     NodalProjectionData       projection     = diff_flux_proj->projection;
51336038bbcSJames Wright 
51436038bbcSJames Wright     diff_flux_proj->CreateRHSOperator_Direct   = DivDiffFluxProjectionCreateRHS_Direct_NS;
51536038bbcSJames Wright     diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_NS;
51636038bbcSJames Wright 
5170c373b74SJames Wright     switch (honee->diff_flux_proj->method) {
51836038bbcSJames Wright       case DIV_DIFF_FLUX_PROJ_DIRECT: {
51936038bbcSJames Wright         PetscSection section;
52036038bbcSJames Wright 
52136038bbcSJames Wright         PetscCall(DMGetLocalSection(projection->dm, &section));
52236038bbcSJames Wright         PetscCall(PetscSectionSetFieldName(section, 0, ""));
52336038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_MomentumX"));
52436038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 1, "DivDiffusiveFlux_MomentumY"));
52536038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 2, "DivDiffusiveFlux_MomentumZ"));
52636038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 3, "DivDiffusiveFlux_Energy"));
52736038bbcSJames Wright       } break;
52836038bbcSJames Wright       case DIV_DIFF_FLUX_PROJ_INDIRECT: {
52936038bbcSJames Wright         PetscSection section;
53036038bbcSJames Wright 
53136038bbcSJames Wright         PetscCall(DMGetLocalSection(projection->dm, &section));
53236038bbcSJames Wright         PetscCall(PetscSectionSetFieldName(section, 0, ""));
53336038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_MomentumXX"));
53436038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_MomentumXY"));
53536038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_MomentumXZ"));
53636038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 3, "DiffusiveFlux_MomentumYX"));
53736038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 4, "DiffusiveFlux_MomentumYY"));
53836038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 5, "DiffusiveFlux_MomentumYZ"));
53936038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 6, "DiffusiveFlux_MomentumZX"));
54036038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 7, "DiffusiveFlux_MomentumZY"));
54136038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 8, "DiffusiveFlux_MomentumZZ"));
54236038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 9, "DiffusiveFlux_EnergyX"));
54336038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 10, "DiffusiveFlux_EnergyY"));
54436038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 11, "DiffusiveFlux_EnergyZ"));
54536038bbcSJames Wright       } break;
54636038bbcSJames Wright       case DIV_DIFF_FLUX_PROJ_NONE:
5470c373b74SJames Wright         SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
5480c373b74SJames Wright                 DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]);
54936038bbcSJames Wright         break;
55036038bbcSJames Wright     }
55136038bbcSJames Wright   }
55236038bbcSJames Wright 
5533a8779fbSJames Wright   // ------------------------------------------------------
554ea615d4cSJames Wright   //           Set up the QFunction context
5553a8779fbSJames Wright   // ------------------------------------------------------
5563a8779fbSJames Wright   // -- Scale variables to desired units
557c9f37605SMohammed Amin   Units units = honee->units;
558*bff204bcSJames Wright   newtonian_ig_ctx->cv *= units->J_per_kg_K;
559*bff204bcSJames Wright   newtonian_ig_ctx->cp *= units->J_per_kg_K;
560*bff204bcSJames Wright   newtonian_ig_ctx->mu *= units->Pascal * units->second;
561*bff204bcSJames Wright   newtonian_ig_ctx->k *= units->W_per_m_K;
562*bff204bcSJames Wright   for (PetscInt i = 0; i < 3; i++) newtonian_ig_ctx->g[i] *= units->m_per_squared_s;
563c9f37605SMohammed Amin   reference.pressure *= units->Pascal;
564c9f37605SMohammed Amin   for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= units->meter / units->second;
565c9f37605SMohammed Amin   reference.temperature *= units->Kelvin;
5663a8779fbSJames Wright 
5678d8b1f0cSJames Wright   PetscReal domain_min[3], domain_max[3], domain_size[3];
5688d8b1f0cSJames Wright   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
5698d8b1f0cSJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = (domain_max[i] - domain_min[i]) * units->meter;
5708d8b1f0cSJames Wright 
5713a8779fbSJames Wright   // -- Solver Settings
5720c373b74SJames Wright   honee->phys->implicit  = implicit;
5730c373b74SJames Wright   honee->phys->state_var = state_var;
5743a8779fbSJames Wright 
5753a8779fbSJames Wright   // -- QFunction Context
57615a3537eSJed Brown   newtonian_ig_ctx->stabilization   = stab;
5778085925cSJames Wright   newtonian_ig_ctx->is_implicit     = implicit;
5783636f6a4SJames Wright   newtonian_ig_ctx->state_var       = state_var;
579c9f37605SMohammed Amin   newtonian_ig_ctx->idl_amplitude   = 1 / (idl_decay_time * units->second);
5800c373b74SJames Wright   newtonian_ig_ctx->divFdiff_method = honee->app_ctx->divFdiffproj_method;
5813a8779fbSJames Wright 
582b8fb7609SAdeleke O. Bankole   // -- Setup Context
583*bff204bcSJames Wright   PetscCall(PetscNew(&setup_context));
584b8fb7609SAdeleke O. Bankole   setup_context->reference = reference;
585b8fb7609SAdeleke O. Bankole   setup_context->gas       = *newtonian_ig_ctx;
586b8fb7609SAdeleke O. Bankole   setup_context->lx        = domain_size[0];
587b8fb7609SAdeleke O. Bankole   setup_context->ly        = domain_size[1];
588b8fb7609SAdeleke O. Bankole   setup_context->lz        = domain_size[2];
589b8fb7609SAdeleke O. Bankole   setup_context->time      = 0;
590b8fb7609SAdeleke O. Bankole 
5910c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &problem->ics.qfctx));
592e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
593e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc));
594e07531f7SJames Wright   PetscCallCeed(
595e07531f7SJames Wright       ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfctx, "evaluation time", offsetof(struct SetupContext_, time), 1, "Time of evaluation"));
59615a3537eSJed Brown 
5970c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &newtonian_ig_qfctx));
598e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(newtonian_ig_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx));
599e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_qfctx, CEED_MEM_HOST, FreeContextPetsc));
600e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1,
601b4c37c5cSJames Wright                                                          "Size of timestep, delta t"));
602e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "ijacobian time shift",
603b4c37c5cSJames Wright                                                          offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1,
604b4c37c5cSJames Wright                                                          "Shift for mass matrix in IJacobian"));
605e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1,
606b4c37c5cSJames Wright                                                          "Current solution time"));
607a8bca8acSJames Wright 
608e07531f7SJames Wright   problem->apply_vol_rhs.qfctx = newtonian_ig_qfctx;
609e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ifunction.qfctx));
610e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ijacobian.qfctx));
611f0b65372SJed Brown 
6125e79d562SJames Wright   for (PetscCount b = 0; b < problem->num_bc_defs; b++) {
6135e79d562SJames Wright     BCDefinition bc_def = problem->bc_defs[b];
6145e79d562SJames Wright     const char  *name;
6155e79d562SJames Wright 
6165e79d562SJames Wright     PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL));
6175e79d562SJames Wright     if (!strcmp(name, "slip")) {
6185e79d562SJames Wright       PetscCall(SlipBCSetup(bc_def, problem, dm, ctx, newtonian_ig_qfctx));
619d4e423e7SJames Wright     } else if (!strcmp(name, "freestream")) {
620d4e423e7SJames Wright       PetscCall(FreestreamBCSetup(bc_def, problem, dm, ctx, newtonian_ig_ctx, &reference));
621f978755dSJames Wright     } else if (!strcmp(name, "outflow")) {
622f978755dSJames Wright       PetscCall(OutflowBCSetup(bc_def, problem, dm, ctx, newtonian_ig_ctx, &reference));
623d6cac220SJames Wright     } else if (!strcmp(name, "inflow")) {
624d6cac220SJames Wright       HoneeBCStruct honee_bc;
625d6cac220SJames Wright 
626d6cac220SJames Wright       PetscCall(PetscNew(&honee_bc));
627d6cac220SJames Wright       PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &honee_bc->qfctx));
628d6cac220SJames Wright       honee_bc->honee              = honee;
6291abc2837SJames Wright       honee_bc->num_comps_jac_data = honee->phys->implicit ? 11 : 0;
630d6cac220SJames Wright       PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc));
631d6cac220SJames Wright 
632d6cac220SJames Wright       PetscCall(BCDefinitionSetIFunction(bc_def, BoundaryIntegralBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
633d6cac220SJames Wright       PetscCall(BCDefinitionSetIJacobian(bc_def, BoundaryIntegralBCSetup_CreateIJacobianQF, HoneeBCAddIJacobianOp));
6345e79d562SJames Wright     }
6355e79d562SJames Wright   }
6369ed3d70dSJames Wright 
6378d8b1f0cSJames Wright   if (unit_tests) PetscCall(UnitTests_Newtonian(honee, newtonian_ig_ctx));
638d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
6393a8779fbSJames Wright }
6403a8779fbSJames Wright 
6419b05e62eSJames Wright //------------------------------------
6429b05e62eSJames Wright // Unit test functions
6439b05e62eSJames Wright //------------------------------------
64448417fb7SJames Wright 
64548417fb7SJames Wright static PetscErrorCode CheckQWithTolerance(const CeedScalar Q_s[5], const CeedScalar Q_a[5], const CeedScalar Q_b[5], const char *name,
64648417fb7SJames Wright                                           PetscReal rtol_0, PetscReal rtol_u, PetscReal rtol_4) {
64748417fb7SJames Wright   CeedScalar relative_error[5];  // relative error
64848417fb7SJames Wright   CeedScalar divisor_threshold = 10 * CEED_EPSILON;
64948417fb7SJames Wright 
65048417fb7SJames Wright   PetscFunctionBeginUser;
65148417fb7SJames Wright   relative_error[0] = (Q_a[0] - Q_b[0]) / (fabs(Q_s[0]) > divisor_threshold ? Q_s[0] : 1);
65248417fb7SJames Wright   relative_error[4] = (Q_a[4] - Q_b[4]) / (fabs(Q_s[4]) > divisor_threshold ? Q_s[4] : 1);
65348417fb7SJames Wright 
65448417fb7SJames Wright   CeedScalar u_magnitude = sqrt(Square(Q_s[1]) + Square(Q_s[2]) + Square(Q_s[3]));
65548417fb7SJames Wright   CeedScalar u_divisor   = u_magnitude > divisor_threshold ? u_magnitude : 1;
65648417fb7SJames Wright   for (int i = 1; i < 4; i++) {
65748417fb7SJames Wright     relative_error[i] = (Q_a[i] - Q_b[i]) / u_divisor;
65848417fb7SJames Wright   }
65948417fb7SJames Wright 
66048417fb7SJames Wright   if (fabs(relative_error[0]) >= rtol_0) {
66148417fb7SJames Wright     printf("%s[0] error %g (expected %.10e, got %.10e)\n", name, relative_error[0], Q_s[0], Q_a[0]);
66248417fb7SJames Wright   }
66348417fb7SJames Wright   for (int i = 1; i < 4; i++) {
66448417fb7SJames Wright     if (fabs(relative_error[i]) >= rtol_u) {
66548417fb7SJames Wright       printf("%s[%d] error %g (expected %.10e, got %.10e)\n", name, i, relative_error[i], Q_s[i], Q_a[i]);
66648417fb7SJames Wright     }
66748417fb7SJames Wright   }
66848417fb7SJames Wright   if (fabs(relative_error[4]) >= rtol_4) {
66948417fb7SJames Wright     printf("%s[4] error %g (expected %.10e, got %.10e)\n", name, relative_error[4], Q_s[4], Q_a[4]);
67048417fb7SJames Wright   }
67148417fb7SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
67248417fb7SJames Wright }
67348417fb7SJames Wright 
67448417fb7SJames Wright // @brief Verify `StateFromQ` by converting A0 -> B0 -> A0_test, where A0 should equal A0_test
67548417fb7SJames Wright static PetscErrorCode TestState(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5],
67648417fb7SJames Wright                                 CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
67748417fb7SJames Wright   CeedScalar        B0[5], A0_test[5];
67848417fb7SJames Wright   char              buf[128];
67948417fb7SJames Wright   const char *const StateVariables_Initial[] = {"U", "Y", "V"};
68048417fb7SJames Wright 
68148417fb7SJames Wright   PetscFunctionBeginUser;
68248417fb7SJames Wright   const char *A_initial = StateVariables_Initial[state_var_A];
68348417fb7SJames Wright   const char *B_initial = StateVariables_Initial[state_var_B];
68448417fb7SJames Wright 
68548417fb7SJames Wright   State state_A0 = StateFromQ(gas, A0, state_var_A);
68648417fb7SJames Wright   StateToQ(gas, state_A0, B0, state_var_B);
68748417fb7SJames Wright   State state_B0 = StateFromQ(gas, B0, state_var_B);
68848417fb7SJames Wright   StateToQ(gas, state_B0, A0_test, state_var_A);
68948417fb7SJames Wright 
69048417fb7SJames Wright   snprintf(buf, sizeof buf, "%s->%s->%s: %s", A_initial, B_initial, A_initial, A_initial);
69148417fb7SJames Wright   PetscCall(CheckQWithTolerance(A0, A0_test, A0, buf, rtol_0, rtol_u, rtol_4));
69248417fb7SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
69348417fb7SJames Wright }
69448417fb7SJames Wright 
69548417fb7SJames Wright // @brief Verify `StateFromQ_fwd` via a finite difference approximation
69648417fb7SJames Wright static PetscErrorCode TestState_fwd(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5],
69748417fb7SJames Wright                                     CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
69848417fb7SJames Wright   CeedScalar        eps = 4e-7;  // Finite difference step
69948417fb7SJames Wright   char              buf[128];
70048417fb7SJames Wright   const char *const StateVariables_Initial[] = {"U", "Y", "V"};
70148417fb7SJames Wright 
70248417fb7SJames Wright   PetscFunctionBeginUser;
70348417fb7SJames Wright   const char *A_initial = StateVariables_Initial[state_var_A];
70448417fb7SJames Wright   const char *B_initial = StateVariables_Initial[state_var_B];
70548417fb7SJames Wright   State       state_0   = StateFromQ(gas, A0, state_var_A);
70648417fb7SJames Wright 
70748417fb7SJames Wright   for (int i = 0; i < 5; i++) {
70848417fb7SJames Wright     CeedScalar dB[5] = {0.}, dB_fd[5] = {0.};
70948417fb7SJames Wright     {  // Calculate dB using State functions
71048417fb7SJames Wright       CeedScalar dA[5] = {0};
71148417fb7SJames Wright 
71248417fb7SJames Wright       dA[i]          = A0[i];
71348417fb7SJames Wright       State dstate_0 = StateFromQ_fwd(gas, state_0, dA, state_var_A);
71448417fb7SJames Wright       StateToQ_fwd(gas, state_0, dstate_0, dB, state_var_B);
71548417fb7SJames Wright     }
71648417fb7SJames Wright 
71748417fb7SJames Wright     {  // Calculate dB_fd via finite difference approximation
71848417fb7SJames Wright       CeedScalar A1[5], B0[5], B1[5];
71948417fb7SJames Wright 
72048417fb7SJames Wright       for (int j = 0; j < 5; j++) A1[j] = (1 + eps * (i == j)) * A0[j];
72148417fb7SJames Wright       State state_1 = StateFromQ(gas, A1, state_var_A);
72248417fb7SJames Wright       StateToQ(gas, state_0, B0, state_var_B);
72348417fb7SJames Wright       StateToQ(gas, state_1, B1, state_var_B);
72448417fb7SJames Wright       for (int j = 0; j < 5; j++) dB_fd[j] = (B1[j] - B0[j]) / eps;
72548417fb7SJames Wright     }
72648417fb7SJames Wright 
72748417fb7SJames Wright     snprintf(buf, sizeof buf, "d%s->d%s: StateFrom%s_fwd i=%d: d%s", A_initial, B_initial, A_initial, i, B_initial);
72848417fb7SJames Wright     PetscCall(CheckQWithTolerance(dB_fd, dB, dB_fd, buf, rtol_0, rtol_u, rtol_4));
72948417fb7SJames Wright   }
73048417fb7SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
73148417fb7SJames Wright }
73248417fb7SJames Wright 
73348417fb7SJames Wright // @brief Test the Newtonian State transformation functions, `StateFrom*`
73448417fb7SJames Wright static PetscErrorCode UnitTests_Newtonian(Honee honee, NewtonianIdealGasContext gas) {
73548417fb7SJames Wright   Units            units = honee->units;
73648417fb7SJames Wright   const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, K = units->Kelvin;
73748417fb7SJames Wright   CeedScalar       rtol;
73848417fb7SJames Wright 
73948417fb7SJames Wright   PetscFunctionBeginUser;
74048417fb7SJames Wright   const CeedScalar T          = 200 * K;
74148417fb7SJames Wright   const CeedScalar rho        = 1.2 * kg / Cube(m);
74248417fb7SJames Wright   const CeedScalar P          = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T;
74348417fb7SJames Wright   const CeedScalar u_base     = 40 * m / sec;
74448417fb7SJames Wright   const CeedScalar u[3]       = {u_base, u_base * 1.1, u_base * 1.2};
74548417fb7SJames Wright   const CeedScalar e_kinetic  = 0.5 * Dot3(u, u);
74648417fb7SJames Wright   const CeedScalar e_internal = gas->cv * T;
74748417fb7SJames Wright   const CeedScalar e_total    = e_kinetic + e_internal;
74848417fb7SJames Wright   const CeedScalar gamma      = HeatCapacityRatio(gas);
74948417fb7SJames Wright   const CeedScalar entropy    = log(P) - gamma * log(rho);
75048417fb7SJames Wright   const CeedScalar rho_div_p  = rho / P;
75148417fb7SJames Wright   const CeedScalar Y0[5]      = {P, u[0], u[1], u[2], T};
75248417fb7SJames Wright   const CeedScalar U0[5]      = {rho, rho * u[0], rho * u[1], rho * u[2], rho * e_total};
75348417fb7SJames Wright   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],
75448417fb7SJames Wright                                  -rho_div_p};
75548417fb7SJames Wright 
75648417fb7SJames Wright   rtol = 20 * CEED_EPSILON;
75748417fb7SJames Wright   PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol));
75848417fb7SJames Wright   PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol));
75948417fb7SJames Wright   PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol));
76048417fb7SJames Wright   PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, rtol, rtol, rtol));
76148417fb7SJames Wright   PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, rtol, rtol, rtol));
76248417fb7SJames Wright   PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, rtol, rtol, rtol));
76348417fb7SJames Wright 
76448417fb7SJames Wright   rtol = 5e-6;
76548417fb7SJames Wright   PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol));
76648417fb7SJames Wright   PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol));
76748417fb7SJames Wright   PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol));
76848417fb7SJames Wright   PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, 10 * rtol, rtol, rtol));
76948417fb7SJames Wright   PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, 5 * rtol, rtol, rtol));
77048417fb7SJames Wright   PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, 5 * rtol, 5 * rtol, 5 * rtol));
77148417fb7SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
77248417fb7SJames Wright }
773