xref: /honee/problems/newtonian.c (revision e5a8cae06a62eeac7fdfa27f57fc80d5b628dcec)
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 
14*e5a8cae0SJames Wright const char *const StateVariables[]     = {"CONSERVATIVE", "PRIMITIVE", "ENTROPY", "StateVariable", "STATEVAR_", NULL};
15*e5a8cae0SJames Wright const char *const StabilizationTypes[] = {"NONE", "SU", "SUPG", "StabilizationType", "STAB_", NULL};
16*e5a8cae0SJames Wright 
1748417fb7SJames Wright static PetscErrorCode UnitTests_Newtonian(Honee honee, NewtonianIdealGasContext gas);
18f0b65372SJed Brown 
1965dee3d2SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping
2065dee3d2SJames Wright //
2165dee3d2SJames Wright // Only used for SUPG stabilization
220c373b74SJames Wright PetscErrorCode CreateKSPMassOperator_NewtonianStabilized(Honee honee, CeedOperator *op_mass) {
230c373b74SJames Wright   Ceed                 ceed = honee->ceed;
2465dee3d2SJames Wright   CeedInt              num_comp_q, q_data_size;
2565dee3d2SJames Wright   CeedQFunction        qf_mass;
2601e19bfaSJames Wright   CeedElemRestriction  elem_restr_q, elem_restr_qd;
2765dee3d2SJames Wright   CeedBasis            basis_q;
2865dee3d2SJames Wright   CeedVector           q_data;
29236a8ba6SJames Wright   CeedQFunctionContext qfctx = NULL;
3065dee3d2SJames Wright   PetscInt             dim   = 3;
3165dee3d2SJames Wright 
3265dee3d2SJames Wright   PetscFunctionBeginUser;
3365dee3d2SJames Wright   {  // Get restriction and basis from the RHS function
3465dee3d2SJames Wright     CeedOperator     *sub_ops;
3501e19bfaSJames Wright     CeedOperatorField op_field;
3665dee3d2SJames Wright     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
3765dee3d2SJames Wright 
380c373b74SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops));
3901e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field));
4001e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL));
4101e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &op_field));
4201e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_qd, NULL, &q_data));
4365dee3d2SJames Wright 
44236a8ba6SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qfctx));
4565dee3d2SJames Wright   }
4665dee3d2SJames Wright 
4765dee3d2SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
4801e19bfaSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd, &q_data_size));
4965dee3d2SJames Wright 
5065dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Newtonian_Conserv, MassFunction_Newtonian_Conserv_loc, &qf_mass));
5165dee3d2SJames Wright 
52236a8ba6SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qfctx));
5365dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
5465dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
5565dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
5665dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
5765dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
5865dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
5965dee3d2SJames Wright 
6065dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
6171f2ed29SJames Wright   PetscCallCeed(ceed, CeedOperatorSetName(*op_mass, "RHS Mass Operator, Newtonian Stabilized"));
6265dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
630c373b74SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, honee->q_ceed));
6401e19bfaSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
6565dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
6665dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
6765dee3d2SJames Wright 
68fff85bd3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
6901e19bfaSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
70fff85bd3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
71fff85bd3SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
72236a8ba6SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qfctx));
7365dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
7465dee3d2SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
7565dee3d2SJames Wright }
764c07ec22SJames Wright 
7736038bbcSJames Wright /**
7836038bbcSJames Wright   @brief Create RHS CeedOperator for direct projection of divergence of diffusive flux
7936038bbcSJames Wright 
800c373b74SJames Wright   @param[in]  honee          `Honee` context
8136038bbcSJames Wright   @param[in]  diff_flux_proj `DivDiffFluxProjectionData` object
8236038bbcSJames Wright   @param[out] op_rhs         Operator to calculate the RHS of the L^2 projection
8336038bbcSJames Wright **/
84e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_NS(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) {
850c373b74SJames Wright   Ceed                 ceed       = honee->ceed;
8636038bbcSJames Wright   NodalProjectionData  projection = diff_flux_proj->projection;
8736038bbcSJames Wright   CeedInt              num_comp_q;
8836038bbcSJames Wright   PetscInt             dim, label_value = 0;
8936038bbcSJames Wright   DMLabel              domain_label    = NULL;
90236a8ba6SJames Wright   CeedQFunctionContext newtonian_qfctx = NULL;
9136038bbcSJames Wright 
9236038bbcSJames Wright   PetscFunctionBeginUser;
9336038bbcSJames Wright   // -- Get Pre-requisite things
9436038bbcSJames Wright   PetscCall(DMGetDimension(projection->dm, &dim));
95e3663b90SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q));
9636038bbcSJames Wright 
970880fbb6SJames Wright   {  // Get newtonian QF context
98b3b24828SJames Wright     CeedOperator *sub_ops, main_op = honee->op_ifunction ? honee->op_ifunction : honee->op_rhs_ctx->op;
9936038bbcSJames Wright     PetscInt      sub_op_index = 0;  // will be 0 for the volume op
10036038bbcSJames Wright 
101b3b24828SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(main_op, &sub_ops));
102236a8ba6SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newtonian_qfctx));
10336038bbcSJames Wright   }
10436038bbcSJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_rhs));
10536038bbcSJames Wright   {  // Add the volume integral CeedOperator
10636038bbcSJames Wright     CeedQFunction       qf_rhs_volume;
10736038bbcSJames Wright     CeedOperator        op_rhs_volume;
10836038bbcSJames Wright     CeedVector          q_data;
1090880fbb6SJames Wright     CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL;
1100880fbb6SJames Wright     CeedBasis           basis_diff_flux = NULL;
11136038bbcSJames Wright     CeedInt             q_data_size;
11236038bbcSJames Wright 
1130880fbb6SJames Wright     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL));
114e3663b90SJames 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,
115e3663b90SJames Wright                        &q_data_size));
1160c373b74SJames Wright     switch (honee->phys->state_var) {
11736038bbcSJames Wright       case STATEVAR_PRIMITIVE:
11836038bbcSJames Wright         PetscCallCeed(ceed,
11936038bbcSJames Wright                       CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Prim, DivDiffusiveFluxVolumeRHS_NS_Prim_loc, &qf_rhs_volume));
12036038bbcSJames Wright         break;
12136038bbcSJames Wright       case STATEVAR_CONSERVATIVE:
12236038bbcSJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Conserv, DivDiffusiveFluxVolumeRHS_NS_Conserv_loc,
12336038bbcSJames Wright                                                         &qf_rhs_volume));
12436038bbcSJames Wright         break;
12536038bbcSJames Wright       case STATEVAR_ENTROPY:
12636038bbcSJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Entropy, DivDiffusiveFluxVolumeRHS_NS_Entropy_loc,
12736038bbcSJames Wright                                                         &qf_rhs_volume));
12836038bbcSJames Wright         break;
12936038bbcSJames Wright     }
13036038bbcSJames Wright 
131236a8ba6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, newtonian_qfctx));
13236038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "q", num_comp_q, CEED_EVAL_INTERP));
13336038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
13436038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE));
13536038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD));
13636038bbcSJames Wright 
13736038bbcSJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume));
138e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
139e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
14036038bbcSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
14136038bbcSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE));
14236038bbcSJames Wright 
14336038bbcSJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_volume));
14436038bbcSJames Wright 
14536038bbcSJames Wright     PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
14636038bbcSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
1470880fbb6SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume));
1480880fbb6SJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
14936038bbcSJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume));
15036038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume));
15136038bbcSJames Wright   }
15236038bbcSJames Wright 
15336038bbcSJames Wright   {  // Add the boundary integral CeedOperator
15436038bbcSJames Wright     CeedQFunction qf_rhs_boundary;
15536038bbcSJames Wright     DMLabel       face_sets_label;
15636038bbcSJames Wright     PetscInt      num_face_set_values, *face_set_values;
15736038bbcSJames Wright     CeedInt       q_data_size;
15836038bbcSJames Wright 
15936038bbcSJames Wright     // -- Build RHS operator
1600c373b74SJames Wright     switch (honee->phys->state_var) {
16136038bbcSJames Wright       case STATEVAR_PRIMITIVE:
16236038bbcSJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Prim, DivDiffusiveFluxBoundaryRHS_NS_Prim_loc,
16336038bbcSJames Wright                                                         &qf_rhs_boundary));
16436038bbcSJames Wright         break;
16536038bbcSJames Wright       case STATEVAR_CONSERVATIVE:
16636038bbcSJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Conserv, DivDiffusiveFluxBoundaryRHS_NS_Conserv_loc,
16736038bbcSJames Wright                                                         &qf_rhs_boundary));
16836038bbcSJames Wright         break;
16936038bbcSJames Wright       case STATEVAR_ENTROPY:
17036038bbcSJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Entropy, DivDiffusiveFluxBoundaryRHS_NS_Entropy_loc,
17136038bbcSJames Wright                                                         &qf_rhs_boundary));
17236038bbcSJames Wright         break;
17336038bbcSJames Wright     }
17436038bbcSJames Wright 
1750c373b74SJames Wright     PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size));
176236a8ba6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, newtonian_qfctx));
17736038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "q", num_comp_q, CEED_EVAL_INTERP));
17836038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
17936038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE));
18036038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP));
18136038bbcSJames Wright 
18236038bbcSJames Wright     PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label));
18336038bbcSJames Wright     PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values));
18436038bbcSJames Wright     for (PetscInt f = 0; f < num_face_set_values; f++) {
18536038bbcSJames Wright       DMLabel  face_orientation_label;
18636038bbcSJames Wright       PetscInt num_orientations_values, *orientation_values;
18736038bbcSJames Wright 
18836038bbcSJames Wright       {
18936038bbcSJames Wright         char *face_orientation_label_name;
19036038bbcSJames Wright 
19136038bbcSJames Wright         PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name));
19236038bbcSJames Wright         PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label));
19336038bbcSJames Wright         PetscCall(PetscFree(face_orientation_label_name));
19436038bbcSJames Wright       }
19536038bbcSJames Wright       PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values));
19636038bbcSJames Wright       for (PetscInt o = 0; o < num_orientations_values; o++) {
19736038bbcSJames Wright         CeedOperator        op_rhs_boundary;
19836038bbcSJames Wright         CeedBasis           basis_q, basis_diff_flux_boundary;
19936038bbcSJames Wright         CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary;
20036038bbcSJames Wright         CeedVector          q_data;
20136038bbcSJames Wright         CeedInt             q_data_size;
20236038bbcSJames Wright         PetscInt            orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1;
20336038bbcSJames Wright 
2040c373b74SJames Wright         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q));
2050c373b74SJames Wright         PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q));
20636038bbcSJames Wright         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0,
20736038bbcSJames Wright                                                   &elem_restr_diff_flux_boundary));
2080880fbb6SJames Wright         PetscCall(CreateBasisFromPlex(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary));
209e3663b90SJames Wright         PetscCall(
210e3663b90SJames Wright             QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, honee->x_coord, &elem_restr_qdata, &q_data, &q_data_size));
21136038bbcSJames Wright 
21236038bbcSJames Wright         PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary));
21336038bbcSJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
21436038bbcSJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
21536038bbcSJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data));
21636038bbcSJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary,
21736038bbcSJames Wright                                                  CEED_VECTOR_ACTIVE));
21836038bbcSJames Wright 
21936038bbcSJames Wright         PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_boundary));
22036038bbcSJames Wright 
22136038bbcSJames Wright         PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary));
22236038bbcSJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata));
22336038bbcSJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
22436038bbcSJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary));
22536038bbcSJames Wright         PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
22636038bbcSJames Wright         PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary));
22736038bbcSJames Wright         PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
22836038bbcSJames Wright       }
22936038bbcSJames Wright       PetscCall(PetscFree(orientation_values));
23036038bbcSJames Wright     }
23136038bbcSJames Wright     PetscCall(PetscFree(face_set_values));
23236038bbcSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary));
23336038bbcSJames Wright   }
23436038bbcSJames Wright 
235236a8ba6SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newtonian_qfctx));
23636038bbcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
23736038bbcSJames Wright }
23836038bbcSJames Wright 
23936038bbcSJames Wright /**
24036038bbcSJames Wright   @brief Create RHS CeedOperator for indirect projection of divergence of diffusive flux
24136038bbcSJames Wright 
2420c373b74SJames Wright   @param[in]  honee          `Honee` context
24336038bbcSJames Wright   @param[in]  diff_flux_proj `DivDiffFluxProjectionData` object
24436038bbcSJames Wright   @param[out] op_rhs         Operator to calculate the RHS of the L^2 projection
24536038bbcSJames Wright **/
246e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Indirect_NS(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) {
2470c373b74SJames Wright   Ceed                 ceed       = honee->ceed;
24836038bbcSJames Wright   NodalProjectionData  projection = diff_flux_proj->projection;
24936038bbcSJames Wright   CeedBasis            basis_diff_flux;
25036038bbcSJames Wright   CeedElemRestriction  elem_restr_diff_flux, elem_restr_qd;
25136038bbcSJames Wright   CeedVector           q_data;
25236038bbcSJames Wright   CeedInt              num_comp_q, q_data_size;
25336038bbcSJames Wright   PetscInt             dim;
25436038bbcSJames Wright   PetscInt             label_value = 0, height = 0, dm_field = 0;
25536038bbcSJames Wright   DMLabel              domain_label = NULL;
25636038bbcSJames Wright   CeedQFunction        qf_rhs;
257236a8ba6SJames Wright   CeedQFunctionContext newtonian_qfctx = NULL;
25836038bbcSJames Wright 
25936038bbcSJames Wright   PetscFunctionBeginUser;
26036038bbcSJames Wright   PetscCall(DMGetDimension(projection->dm, &dim));
261e3663b90SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q));
26236038bbcSJames Wright 
26340b78511SJames Wright   {  // Get newtonian QF context
264b3b24828SJames Wright     CeedOperator *sub_ops, main_op = honee->op_ifunction ? honee->op_ifunction : honee->op_rhs_ctx->op;
26536038bbcSJames Wright     PetscInt      sub_op_index = 0;  // will be 0 for the volume op
26636038bbcSJames Wright 
267b3b24828SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(main_op, &sub_ops));
268236a8ba6SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newtonian_qfctx));
26936038bbcSJames Wright   }
27036038bbcSJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux));
27136038bbcSJames Wright   PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux));
272e3663b90SJames 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,
273e3663b90SJames Wright                      &q_data_size));
27436038bbcSJames Wright 
2750c373b74SJames Wright   switch (honee->phys->state_var) {
27636038bbcSJames Wright     case STATEVAR_PRIMITIVE:
27736038bbcSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Prim, DiffusiveFluxRHS_NS_Prim_loc, &qf_rhs));
27836038bbcSJames Wright       break;
27936038bbcSJames Wright     case STATEVAR_CONSERVATIVE:
28036038bbcSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Conserv, DiffusiveFluxRHS_NS_Conserv_loc, &qf_rhs));
28136038bbcSJames Wright       break;
28236038bbcSJames Wright     case STATEVAR_ENTROPY:
28336038bbcSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Entropy, DiffusiveFluxRHS_NS_Entropy_loc, &qf_rhs));
28436038bbcSJames Wright       break;
28536038bbcSJames Wright   }
28636038bbcSJames Wright 
287236a8ba6SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, newtonian_qfctx));
28836038bbcSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "q", num_comp_q, CEED_EVAL_INTERP));
28936038bbcSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
29036038bbcSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE));
29136038bbcSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP));
29236038bbcSJames Wright 
29336038bbcSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, op_rhs));
294e3663b90SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
295e3663b90SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
29636038bbcSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
29736038bbcSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
29836038bbcSJames Wright 
29936038bbcSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs));
300236a8ba6SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newtonian_qfctx));
30136038bbcSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
30236038bbcSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
30336038bbcSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
30436038bbcSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux));
30536038bbcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
30636038bbcSJames Wright }
30736038bbcSJames Wright 
308d6cac220SJames Wright static PetscErrorCode BoundaryIntegralBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) {
309d6cac220SJames Wright   HoneeBCStruct honee_bc;
310d6cac220SJames Wright 
311d6cac220SJames Wright   PetscFunctionBeginUser;
312d6cac220SJames Wright   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
313d6cac220SJames Wright   Honee honee = honee_bc->honee;
314d6cac220SJames Wright 
315d6cac220SJames Wright   switch (honee->phys->state_var) {
316d6cac220SJames Wright     case STATEVAR_CONSERVATIVE:
317d6cac220SJames Wright       PetscCall(HoneeBCCreateIFunctionQF(bc_def, BoundaryIntegral_Conserv, BoundaryIntegral_Conserv_loc, honee_bc->qfctx, qf));
318d6cac220SJames Wright       break;
319d6cac220SJames Wright     case STATEVAR_PRIMITIVE:
320d6cac220SJames Wright       PetscCall(HoneeBCCreateIFunctionQF(bc_def, BoundaryIntegral_Prim, BoundaryIntegral_Prim_loc, honee_bc->qfctx, qf));
321d6cac220SJames Wright       break;
322d6cac220SJames Wright     case STATEVAR_ENTROPY:
323d6cac220SJames Wright       PetscCall(HoneeBCCreateIFunctionQF(bc_def, BoundaryIntegral_Entropy, BoundaryIntegral_Entropy_loc, honee_bc->qfctx, qf));
324d6cac220SJames Wright       break;
325d6cac220SJames Wright   }
326d6cac220SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
327d6cac220SJames Wright }
328d6cac220SJames Wright 
329d6cac220SJames Wright static PetscErrorCode BoundaryIntegralBCSetup_CreateIJacobianQF(BCDefinition bc_def, CeedQFunction *qf) {
330d6cac220SJames Wright   HoneeBCStruct honee_bc;
331d6cac220SJames Wright 
332d6cac220SJames Wright   PetscFunctionBeginUser;
333d6cac220SJames Wright   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
334d6cac220SJames Wright   Honee honee = honee_bc->honee;
335d6cac220SJames Wright 
336d6cac220SJames Wright   switch (honee->phys->state_var) {
337d6cac220SJames Wright     case STATEVAR_CONSERVATIVE:
338d6cac220SJames Wright       PetscCall(HoneeBCCreateIJacobianQF(bc_def, BoundaryIntegral_Jacobian_Conserv, BoundaryIntegral_Jacobian_Conserv_loc, honee_bc->qfctx, qf));
339d6cac220SJames Wright       break;
340d6cac220SJames Wright     case STATEVAR_PRIMITIVE:
341d6cac220SJames Wright       PetscCall(HoneeBCCreateIJacobianQF(bc_def, BoundaryIntegral_Jacobian_Prim, BoundaryIntegral_Jacobian_Prim_loc, honee_bc->qfctx, qf));
342d6cac220SJames Wright       break;
343d6cac220SJames Wright     case STATEVAR_ENTROPY:
344d6cac220SJames Wright       PetscCall(HoneeBCCreateIJacobianQF(bc_def, BoundaryIntegral_Jacobian_Entropy, BoundaryIntegral_Jacobian_Entropy_loc, honee_bc->qfctx, qf));
345d6cac220SJames Wright       break;
346d6cac220SJames Wright   }
347d6cac220SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
348d6cac220SJames Wright }
349d6cac220SJames Wright 
350d3c60affSJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx) {
351b7f03f12SJed Brown   SetupContext             setup_context;
3520c373b74SJames Wright   Honee                    honee  = *(Honee *)ctx;
3530c373b74SJames Wright   CeedInt                  degree = honee->app_ctx->degree;
3543a8779fbSJames Wright   StabilizationType        stab;
3553636f6a4SJames Wright   StateVariable            state_var;
3560c373b74SJames Wright   MPI_Comm                 comm = honee->comm;
3570c373b74SJames Wright   Ceed                     ceed = honee->ceed;
3583a8779fbSJames Wright   PetscBool                implicit;
3590d8cd818SJames Wright   PetscBool                unit_tests;
36015a3537eSJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
361e07531f7SJames Wright   CeedQFunctionContext     newtonian_ig_qfctx;
3623a8779fbSJames Wright 
36315a3537eSJed Brown   PetscFunctionBeginUser;
3642b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &setup_context));
3652b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &newtonian_ig_ctx));
3663a8779fbSJames Wright 
3673a8779fbSJames Wright   // ------------------------------------------------------
3683a8779fbSJames Wright   //           Setup Generic Newtonian IG Problem
3693a8779fbSJames Wright   // ------------------------------------------------------
3701abc2837SJames Wright   problem->num_comps_jac_data           = 14;
37158ce1233SJames Wright   problem->compute_exact_solution_error = PETSC_FALSE;
372cbe60e31SLeila Ghaffari   problem->print_info                   = PRINT_NEWTONIAN;
3733a8779fbSJames Wright 
3740c373b74SJames Wright   PetscCall(DivDiffFluxProjectionCreate(honee, 4, &honee->diff_flux_proj));
3750c373b74SJames Wright   if (honee->diff_flux_proj) {
3760c373b74SJames Wright     DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj;
37736038bbcSJames Wright     NodalProjectionData       projection     = diff_flux_proj->projection;
37836038bbcSJames Wright 
37936038bbcSJames Wright     diff_flux_proj->CreateRHSOperator_Direct   = DivDiffFluxProjectionCreateRHS_Direct_NS;
38036038bbcSJames Wright     diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_NS;
38136038bbcSJames Wright 
3820c373b74SJames Wright     switch (honee->diff_flux_proj->method) {
38336038bbcSJames Wright       case DIV_DIFF_FLUX_PROJ_DIRECT: {
38436038bbcSJames Wright         PetscSection section;
38536038bbcSJames Wright 
38636038bbcSJames Wright         PetscCall(DMGetLocalSection(projection->dm, &section));
38736038bbcSJames Wright         PetscCall(PetscSectionSetFieldName(section, 0, ""));
38836038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_MomentumX"));
38936038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 1, "DivDiffusiveFlux_MomentumY"));
39036038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 2, "DivDiffusiveFlux_MomentumZ"));
39136038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 3, "DivDiffusiveFlux_Energy"));
39236038bbcSJames Wright       } break;
39336038bbcSJames Wright       case DIV_DIFF_FLUX_PROJ_INDIRECT: {
39436038bbcSJames Wright         PetscSection section;
39536038bbcSJames Wright 
39636038bbcSJames Wright         PetscCall(DMGetLocalSection(projection->dm, &section));
39736038bbcSJames Wright         PetscCall(PetscSectionSetFieldName(section, 0, ""));
39836038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_MomentumXX"));
39936038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_MomentumXY"));
40036038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_MomentumXZ"));
40136038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 3, "DiffusiveFlux_MomentumYX"));
40236038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 4, "DiffusiveFlux_MomentumYY"));
40336038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 5, "DiffusiveFlux_MomentumYZ"));
40436038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 6, "DiffusiveFlux_MomentumZX"));
40536038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 7, "DiffusiveFlux_MomentumZY"));
40636038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 8, "DiffusiveFlux_MomentumZZ"));
40736038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 9, "DiffusiveFlux_EnergyX"));
40836038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 10, "DiffusiveFlux_EnergyY"));
40936038bbcSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 11, "DiffusiveFlux_EnergyZ"));
41036038bbcSJames Wright       } break;
41136038bbcSJames Wright       case DIV_DIFF_FLUX_PROJ_NONE:
4120c373b74SJames Wright         SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
4130c373b74SJames Wright                 DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]);
41436038bbcSJames Wright         break;
41536038bbcSJames Wright     }
41636038bbcSJames Wright   }
41736038bbcSJames Wright 
4183a8779fbSJames Wright   // ------------------------------------------------------
419ea615d4cSJames Wright   //             Create the QFunction context
4203a8779fbSJames Wright   // ------------------------------------------------------
4213a8779fbSJames Wright   CeedScalar cv         = 717.;          // J/(kg K)
4223a8779fbSJames Wright   CeedScalar cp         = 1004.;         // J/(kg K)
423d9bb1cdbSJames Wright   CeedScalar g[3]       = {0, 0, 0};     // m/s^2
4243a8779fbSJames Wright   CeedScalar lambda     = -2. / 3.;      // -
425bb8a0c61SJames Wright   CeedScalar mu         = 1.8e-5;        // Pa s, dynamic viscosity
4263a8779fbSJames Wright   CeedScalar k          = 0.02638;       // W/(m K)
4274351d27fSLeila Ghaffari   CeedScalar c_tau      = 0.5 / degree;  // -
428bb8a0c61SJames Wright   CeedScalar Ctau_t     = 1.0;           // -
429b5786772SLeila Ghaffari   CeedScalar Cv_func[3] = {36, 60, 128};
430c176988bSLeila Ghaffari   CeedScalar Ctau_v     = Cv_func[(CeedInt)Min(3, degree) - 1];
4312a0eee6eSLeila Ghaffari   CeedScalar Ctau_C     = 0.25 / degree;
4322a0eee6eSLeila Ghaffari   CeedScalar Ctau_M     = 0.25 / degree;
4332a0eee6eSLeila Ghaffari   CeedScalar Ctau_E     = 0.125;
4343a8779fbSJames Wright   PetscReal  domain_min[3], domain_max[3], domain_size[3];
4352b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
436493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
4373a8779fbSJames Wright 
438b8fb7609SAdeleke O. Bankole   StatePrimitive reference      = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15};
439a213b8aaSJames Wright   CeedScalar     idl_decay_time = -1, idl_start = 0, idl_length = 0, idl_pressure = reference.pressure;
440e7754af5SKenneth E. Jansen   PetscBool      idl_enable = PETSC_FALSE;
441b8fb7609SAdeleke O. Bankole 
4423a8779fbSJames Wright   // ------------------------------------------------------
4433a8779fbSJames Wright   //             Create the PETSc context
4443a8779fbSJames Wright   // ------------------------------------------------------
445bb8a0c61SJames Wright   PetscScalar meter    = 1;  // 1 meter in scaled length units
446bb8a0c61SJames Wright   PetscScalar kilogram = 1;  // 1 kilogram in scaled mass units
447bb8a0c61SJames Wright   PetscScalar second   = 1;  // 1 second in scaled time units
4483a8779fbSJames Wright   PetscScalar Kelvin   = 1;  // 1 Kelvin in scaled temperature units
4493a8779fbSJames Wright   PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s;
4503a8779fbSJames Wright 
4513a8779fbSJames Wright   // ------------------------------------------------------
4523a8779fbSJames Wright   //              Command line Options
4533a8779fbSJames Wright   // ------------------------------------------------------
454d9bb1cdbSJames Wright   PetscBool given_option = PETSC_FALSE;
4552b916ea7SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL);
456cbe60e31SLeila Ghaffari   // -- Conservative vs Primitive variables
4572b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE),
4582b916ea7SJeremy L Thompson                              (PetscEnum *)&state_var, NULL));
4593636f6a4SJames Wright 
4603636f6a4SJames Wright   switch (state_var) {
4613636f6a4SJames Wright     case STATEVAR_CONSERVATIVE:
462e07531f7SJames Wright       problem->ics.qf_func_ptr                 = ICsNewtonianIG_Conserv;
463e07531f7SJames Wright       problem->ics.qf_loc                      = ICsNewtonianIG_Conserv_loc;
464e07531f7SJames Wright       problem->apply_vol_rhs.qf_func_ptr       = RHSFunction_Newtonian;
465e07531f7SJames Wright       problem->apply_vol_rhs.qf_loc            = RHSFunction_Newtonian_loc;
466e07531f7SJames Wright       problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Conserv;
467e07531f7SJames Wright       problem->apply_vol_ifunction.qf_loc      = IFunction_Newtonian_Conserv_loc;
468e07531f7SJames Wright       problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Conserv;
469e07531f7SJames Wright       problem->apply_vol_ijacobian.qf_loc      = IJacobian_Newtonian_Conserv_loc;
4703636f6a4SJames Wright       break;
4713636f6a4SJames Wright     case STATEVAR_PRIMITIVE:
472e07531f7SJames Wright       problem->ics.qf_func_ptr                 = ICsNewtonianIG_Prim;
473e07531f7SJames Wright       problem->ics.qf_loc                      = ICsNewtonianIG_Prim_loc;
474e07531f7SJames Wright       problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Prim;
475e07531f7SJames Wright       problem->apply_vol_ifunction.qf_loc      = IFunction_Newtonian_Prim_loc;
476e07531f7SJames Wright       problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Prim;
477e07531f7SJames Wright       problem->apply_vol_ijacobian.qf_loc      = IJacobian_Newtonian_Prim_loc;
4783636f6a4SJames Wright       break;
479c62f7daeSJames Wright     case STATEVAR_ENTROPY:
480e07531f7SJames Wright       problem->ics.qf_func_ptr                 = ICsNewtonianIG_Entropy;
481e07531f7SJames Wright       problem->ics.qf_loc                      = ICsNewtonianIG_Entropy_loc;
482e07531f7SJames Wright       problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Entropy;
483e07531f7SJames Wright       problem->apply_vol_ifunction.qf_loc      = IFunction_Newtonian_Entropy_loc;
484e07531f7SJames Wright       problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Entropy;
485e07531f7SJames Wright       problem->apply_vol_ijacobian.qf_loc      = IJacobian_Newtonian_Entropy_loc;
486c62f7daeSJames Wright       break;
487cbe60e31SLeila Ghaffari   }
4881485969bSJeremy L Thompson 
4893a8779fbSJames Wright   // -- Physics
4902b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL));
4912b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL));
4922b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL));
4932b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL));
4942b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL));
4953a8779fbSJames Wright 
49666d54740SJames Wright   PetscInt dim = 3;
497d9bb1cdbSJames Wright   PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL));
498d9bb1cdbSJames Wright   PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option));
499d9bb1cdbSJames Wright   if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim);
500d9bb1cdbSJames Wright 
5012b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
5022b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
5032b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
5042b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL));
5052b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL));
5062b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL));
5072b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL));
5082b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
5092b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL));
5103a8779fbSJames Wright 
511db95602eSJames Wright   dim = 3;
512b8fb7609SAdeleke O. Bankole   PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL));
513b8fb7609SAdeleke O. Bankole   PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL));
514b8fb7609SAdeleke O. Bankole   PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL));
515b8fb7609SAdeleke O. Bankole 
5163a8779fbSJames Wright   // -- Units
5172b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
5183a8779fbSJames Wright   meter = fabs(meter);
5192b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL));
5203a8779fbSJames Wright   kilogram = fabs(kilogram);
5212b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
5223a8779fbSJames Wright   second = fabs(second);
5232b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL));
5243a8779fbSJames Wright   Kelvin = fabs(Kelvin);
5253a8779fbSJames Wright 
5263a8779fbSJames Wright   // -- Warnings
5275d28dccaSJames Wright   PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP,
5285d28dccaSJames Wright              "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n");
529e7754af5SKenneth E. Jansen 
530e7754af5SKenneth E. Jansen   PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point",
531e7754af5SKenneth E. Jansen                                NULL, idl_decay_time, &idl_decay_time, &idl_enable));
5329cbdf780SJames Wright   PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero.");
5339cbdf780SJames Wright   if (idl_decay_time < 0) idl_enable = PETSC_FALSE;
5341abc2837SJames Wright   if (idl_enable) problem->num_comps_jac_data++;
535e7754af5SKenneth E. Jansen   PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL));
536e7754af5SKenneth E. Jansen   PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL));
537a213b8aaSJames Wright   idl_pressure = reference.pressure;
538a213b8aaSJames Wright   PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-reference_pressure`)", NULL, idl_pressure,
539a213b8aaSJames Wright                                &idl_pressure, NULL));
5401485969bSJeremy L Thompson   PetscOptionsEnd();
5413a8779fbSJames Wright 
54265dee3d2SJames Wright   if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized;
54365dee3d2SJames Wright 
5443a8779fbSJames Wright   // ------------------------------------------------------
5453a8779fbSJames Wright   //           Set up the PETSc context
5463a8779fbSJames Wright   // ------------------------------------------------------
5473a8779fbSJames Wright   // -- Define derived units
5483a8779fbSJames Wright   Pascal          = kilogram / (meter * PetscSqr(second));
5493a8779fbSJames Wright   J_per_kg_K      = PetscSqr(meter) / (PetscSqr(second) * Kelvin);
5503a8779fbSJames Wright   m_per_squared_s = meter / PetscSqr(second);
5513a8779fbSJames Wright   W_per_m_K       = kilogram * meter / (pow(second, 3) * Kelvin);
5523a8779fbSJames Wright 
5530c373b74SJames Wright   honee->units->meter           = meter;
5540c373b74SJames Wright   honee->units->kilogram        = kilogram;
5550c373b74SJames Wright   honee->units->second          = second;
5560c373b74SJames Wright   honee->units->Kelvin          = Kelvin;
5570c373b74SJames Wright   honee->units->Pascal          = Pascal;
5580c373b74SJames Wright   honee->units->J_per_kg_K      = J_per_kg_K;
5590c373b74SJames Wright   honee->units->m_per_squared_s = m_per_squared_s;
5600c373b74SJames Wright   honee->units->W_per_m_K       = W_per_m_K;
5613a8779fbSJames Wright 
5623a8779fbSJames Wright   // ------------------------------------------------------
563ea615d4cSJames Wright   //           Set up the QFunction context
5643a8779fbSJames Wright   // ------------------------------------------------------
5653a8779fbSJames Wright   // -- Scale variables to desired units
5663a8779fbSJames Wright   cv *= J_per_kg_K;
5673a8779fbSJames Wright   cp *= J_per_kg_K;
5683a8779fbSJames Wright   mu *= Pascal * second;
5693a8779fbSJames Wright   k *= W_per_m_K;
570493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter;
571493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s;
572b8fb7609SAdeleke O. Bankole   reference.pressure *= Pascal;
573b8fb7609SAdeleke O. Bankole   for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second;
574b8fb7609SAdeleke O. Bankole   reference.temperature *= Kelvin;
5753a8779fbSJames Wright 
5763a8779fbSJames Wright   // -- Solver Settings
5770c373b74SJames Wright   honee->phys->implicit  = implicit;
5780c373b74SJames Wright   honee->phys->state_var = state_var;
5793a8779fbSJames Wright 
5803a8779fbSJames Wright   // -- QFunction Context
58115a3537eSJed Brown   newtonian_ig_ctx->lambda          = lambda;
58215a3537eSJed Brown   newtonian_ig_ctx->mu              = mu;
58315a3537eSJed Brown   newtonian_ig_ctx->k               = k;
58415a3537eSJed Brown   newtonian_ig_ctx->cv              = cv;
58515a3537eSJed Brown   newtonian_ig_ctx->cp              = cp;
58615a3537eSJed Brown   newtonian_ig_ctx->c_tau           = c_tau;
58715a3537eSJed Brown   newtonian_ig_ctx->Ctau_t          = Ctau_t;
58815a3537eSJed Brown   newtonian_ig_ctx->Ctau_v          = Ctau_v;
58915a3537eSJed Brown   newtonian_ig_ctx->Ctau_C          = Ctau_C;
59015a3537eSJed Brown   newtonian_ig_ctx->Ctau_M          = Ctau_M;
59115a3537eSJed Brown   newtonian_ig_ctx->Ctau_E          = Ctau_E;
59215a3537eSJed Brown   newtonian_ig_ctx->stabilization   = stab;
5938085925cSJames Wright   newtonian_ig_ctx->is_implicit     = implicit;
5943636f6a4SJames Wright   newtonian_ig_ctx->state_var       = state_var;
595e7754af5SKenneth E. Jansen   newtonian_ig_ctx->idl_enable      = idl_enable;
596e7754af5SKenneth E. Jansen   newtonian_ig_ctx->idl_amplitude   = 1 / (idl_decay_time * second);
597e7754af5SKenneth E. Jansen   newtonian_ig_ctx->idl_start       = idl_start * meter;
598e7754af5SKenneth E. Jansen   newtonian_ig_ctx->idl_length      = idl_length * meter;
599a213b8aaSJames Wright   newtonian_ig_ctx->idl_pressure    = idl_pressure;
6000c373b74SJames Wright   newtonian_ig_ctx->divFdiff_method = honee->app_ctx->divFdiffproj_method;
6012b916ea7SJeremy L Thompson   PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3));
6023a8779fbSJames Wright 
603b8fb7609SAdeleke O. Bankole   // -- Setup Context
604b8fb7609SAdeleke O. Bankole   setup_context->reference = reference;
605b8fb7609SAdeleke O. Bankole   setup_context->gas       = *newtonian_ig_ctx;
606b8fb7609SAdeleke O. Bankole   setup_context->lx        = domain_size[0];
607b8fb7609SAdeleke O. Bankole   setup_context->ly        = domain_size[1];
608b8fb7609SAdeleke O. Bankole   setup_context->lz        = domain_size[2];
609b8fb7609SAdeleke O. Bankole   setup_context->time      = 0;
610b8fb7609SAdeleke O. Bankole 
6110c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &problem->ics.qfctx));
612e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
613e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc));
614e07531f7SJames Wright   PetscCallCeed(
615e07531f7SJames Wright       ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfctx, "evaluation time", offsetof(struct SetupContext_, time), 1, "Time of evaluation"));
61615a3537eSJed Brown 
6170c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &newtonian_ig_qfctx));
618e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(newtonian_ig_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx));
619e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_qfctx, CEED_MEM_HOST, FreeContextPetsc));
620e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1,
621b4c37c5cSJames Wright                                                          "Size of timestep, delta t"));
622e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "ijacobian time shift",
623b4c37c5cSJames Wright                                                          offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1,
624b4c37c5cSJames Wright                                                          "Shift for mass matrix in IJacobian"));
625e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1,
626b4c37c5cSJames Wright                                                          "Current solution time"));
627a8bca8acSJames Wright 
628e07531f7SJames Wright   problem->apply_vol_rhs.qfctx = newtonian_ig_qfctx;
629e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ifunction.qfctx));
630e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ijacobian.qfctx));
631f0b65372SJed Brown 
6325e79d562SJames Wright   for (PetscCount b = 0; b < problem->num_bc_defs; b++) {
6335e79d562SJames Wright     BCDefinition bc_def = problem->bc_defs[b];
6345e79d562SJames Wright     const char  *name;
6355e79d562SJames Wright 
6365e79d562SJames Wright     PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL));
6375e79d562SJames Wright     if (!strcmp(name, "slip")) {
6385e79d562SJames Wright       PetscCall(SlipBCSetup(bc_def, problem, dm, ctx, newtonian_ig_qfctx));
639d4e423e7SJames Wright     } else if (!strcmp(name, "freestream")) {
640d4e423e7SJames Wright       PetscCall(FreestreamBCSetup(bc_def, problem, dm, ctx, newtonian_ig_ctx, &reference));
641f978755dSJames Wright     } else if (!strcmp(name, "outflow")) {
642f978755dSJames Wright       PetscCall(OutflowBCSetup(bc_def, problem, dm, ctx, newtonian_ig_ctx, &reference));
643d6cac220SJames Wright     } else if (!strcmp(name, "inflow")) {
644d6cac220SJames Wright       HoneeBCStruct honee_bc;
645d6cac220SJames Wright 
646d6cac220SJames Wright       PetscCall(PetscNew(&honee_bc));
647d6cac220SJames Wright       PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &honee_bc->qfctx));
648d6cac220SJames Wright       honee_bc->honee              = honee;
6491abc2837SJames Wright       honee_bc->num_comps_jac_data = honee->phys->implicit ? 11 : 0;
650d6cac220SJames Wright       PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc));
651d6cac220SJames Wright 
652d6cac220SJames Wright       PetscCall(BCDefinitionSetIFunction(bc_def, BoundaryIntegralBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
653d6cac220SJames Wright       PetscCall(BCDefinitionSetIJacobian(bc_def, BoundaryIntegralBCSetup_CreateIJacobianQF, HoneeBCAddIJacobianOp));
6545e79d562SJames Wright     }
6555e79d562SJames Wright   }
6569ed3d70dSJames Wright 
657f0b65372SJed Brown   if (unit_tests) {
6580c373b74SJames Wright     PetscCall(UnitTests_Newtonian(honee, newtonian_ig_ctx));
659f0b65372SJed Brown   }
660d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
6613a8779fbSJames Wright }
6623a8779fbSJames Wright 
6630c373b74SJames Wright PetscErrorCode PRINT_NEWTONIAN(Honee honee, ProblemData problem, AppCtx app_ctx) {
6640c373b74SJames Wright   MPI_Comm                 comm = honee->comm;
6650c373b74SJames Wright   Ceed                     ceed = honee->ceed;
66615a3537eSJed Brown   NewtonianIdealGasContext newtonian_ctx;
66715a3537eSJed Brown 
6683a8779fbSJames Wright   PetscFunctionBeginUser;
669e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ctx));
6702b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(comm,
67115a3537eSJed Brown                         "  Problem:\n"
67215a3537eSJed Brown                         "    Problem Name                       : %s\n"
67315a3537eSJed Brown                         "    Stabilization                      : %s\n",
6742b916ea7SJeremy L Thompson                         app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization]));
675e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ctx));
676d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
6773a8779fbSJames Wright }
67848417fb7SJames Wright 
67948417fb7SJames Wright static PetscErrorCode CheckQWithTolerance(const CeedScalar Q_s[5], const CeedScalar Q_a[5], const CeedScalar Q_b[5], const char *name,
68048417fb7SJames Wright                                           PetscReal rtol_0, PetscReal rtol_u, PetscReal rtol_4) {
68148417fb7SJames Wright   CeedScalar relative_error[5];  // relative error
68248417fb7SJames Wright   CeedScalar divisor_threshold = 10 * CEED_EPSILON;
68348417fb7SJames Wright 
68448417fb7SJames Wright   PetscFunctionBeginUser;
68548417fb7SJames Wright   relative_error[0] = (Q_a[0] - Q_b[0]) / (fabs(Q_s[0]) > divisor_threshold ? Q_s[0] : 1);
68648417fb7SJames Wright   relative_error[4] = (Q_a[4] - Q_b[4]) / (fabs(Q_s[4]) > divisor_threshold ? Q_s[4] : 1);
68748417fb7SJames Wright 
68848417fb7SJames Wright   CeedScalar u_magnitude = sqrt(Square(Q_s[1]) + Square(Q_s[2]) + Square(Q_s[3]));
68948417fb7SJames Wright   CeedScalar u_divisor   = u_magnitude > divisor_threshold ? u_magnitude : 1;
69048417fb7SJames Wright   for (int i = 1; i < 4; i++) {
69148417fb7SJames Wright     relative_error[i] = (Q_a[i] - Q_b[i]) / u_divisor;
69248417fb7SJames Wright   }
69348417fb7SJames Wright 
69448417fb7SJames Wright   if (fabs(relative_error[0]) >= rtol_0) {
69548417fb7SJames Wright     printf("%s[0] error %g (expected %.10e, got %.10e)\n", name, relative_error[0], Q_s[0], Q_a[0]);
69648417fb7SJames Wright   }
69748417fb7SJames Wright   for (int i = 1; i < 4; i++) {
69848417fb7SJames Wright     if (fabs(relative_error[i]) >= rtol_u) {
69948417fb7SJames Wright       printf("%s[%d] error %g (expected %.10e, got %.10e)\n", name, i, relative_error[i], Q_s[i], Q_a[i]);
70048417fb7SJames Wright     }
70148417fb7SJames Wright   }
70248417fb7SJames Wright   if (fabs(relative_error[4]) >= rtol_4) {
70348417fb7SJames Wright     printf("%s[4] error %g (expected %.10e, got %.10e)\n", name, relative_error[4], Q_s[4], Q_a[4]);
70448417fb7SJames Wright   }
70548417fb7SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
70648417fb7SJames Wright }
70748417fb7SJames Wright 
70848417fb7SJames Wright // @brief Verify `StateFromQ` by converting A0 -> B0 -> A0_test, where A0 should equal A0_test
70948417fb7SJames Wright static PetscErrorCode TestState(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5],
71048417fb7SJames Wright                                 CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
71148417fb7SJames Wright   CeedScalar        B0[5], A0_test[5];
71248417fb7SJames Wright   char              buf[128];
71348417fb7SJames Wright   const char *const StateVariables_Initial[] = {"U", "Y", "V"};
71448417fb7SJames Wright 
71548417fb7SJames Wright   PetscFunctionBeginUser;
71648417fb7SJames Wright   const char *A_initial = StateVariables_Initial[state_var_A];
71748417fb7SJames Wright   const char *B_initial = StateVariables_Initial[state_var_B];
71848417fb7SJames Wright 
71948417fb7SJames Wright   State state_A0 = StateFromQ(gas, A0, state_var_A);
72048417fb7SJames Wright   StateToQ(gas, state_A0, B0, state_var_B);
72148417fb7SJames Wright   State state_B0 = StateFromQ(gas, B0, state_var_B);
72248417fb7SJames Wright   StateToQ(gas, state_B0, A0_test, state_var_A);
72348417fb7SJames Wright 
72448417fb7SJames Wright   snprintf(buf, sizeof buf, "%s->%s->%s: %s", A_initial, B_initial, A_initial, A_initial);
72548417fb7SJames Wright   PetscCall(CheckQWithTolerance(A0, A0_test, A0, buf, rtol_0, rtol_u, rtol_4));
72648417fb7SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
72748417fb7SJames Wright }
72848417fb7SJames Wright 
72948417fb7SJames Wright // @brief Verify `StateFromQ_fwd` via a finite difference approximation
73048417fb7SJames Wright static PetscErrorCode TestState_fwd(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5],
73148417fb7SJames Wright                                     CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
73248417fb7SJames Wright   CeedScalar        eps = 4e-7;  // Finite difference step
73348417fb7SJames Wright   char              buf[128];
73448417fb7SJames Wright   const char *const StateVariables_Initial[] = {"U", "Y", "V"};
73548417fb7SJames Wright 
73648417fb7SJames Wright   PetscFunctionBeginUser;
73748417fb7SJames Wright   const char *A_initial = StateVariables_Initial[state_var_A];
73848417fb7SJames Wright   const char *B_initial = StateVariables_Initial[state_var_B];
73948417fb7SJames Wright   State       state_0   = StateFromQ(gas, A0, state_var_A);
74048417fb7SJames Wright 
74148417fb7SJames Wright   for (int i = 0; i < 5; i++) {
74248417fb7SJames Wright     CeedScalar dB[5] = {0.}, dB_fd[5] = {0.};
74348417fb7SJames Wright     {  // Calculate dB using State functions
74448417fb7SJames Wright       CeedScalar dA[5] = {0};
74548417fb7SJames Wright 
74648417fb7SJames Wright       dA[i]          = A0[i];
74748417fb7SJames Wright       State dstate_0 = StateFromQ_fwd(gas, state_0, dA, state_var_A);
74848417fb7SJames Wright       StateToQ_fwd(gas, state_0, dstate_0, dB, state_var_B);
74948417fb7SJames Wright     }
75048417fb7SJames Wright 
75148417fb7SJames Wright     {  // Calculate dB_fd via finite difference approximation
75248417fb7SJames Wright       CeedScalar A1[5], B0[5], B1[5];
75348417fb7SJames Wright 
75448417fb7SJames Wright       for (int j = 0; j < 5; j++) A1[j] = (1 + eps * (i == j)) * A0[j];
75548417fb7SJames Wright       State state_1 = StateFromQ(gas, A1, state_var_A);
75648417fb7SJames Wright       StateToQ(gas, state_0, B0, state_var_B);
75748417fb7SJames Wright       StateToQ(gas, state_1, B1, state_var_B);
75848417fb7SJames Wright       for (int j = 0; j < 5; j++) dB_fd[j] = (B1[j] - B0[j]) / eps;
75948417fb7SJames Wright     }
76048417fb7SJames Wright 
76148417fb7SJames 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);
76248417fb7SJames Wright     PetscCall(CheckQWithTolerance(dB_fd, dB, dB_fd, buf, rtol_0, rtol_u, rtol_4));
76348417fb7SJames Wright   }
76448417fb7SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
76548417fb7SJames Wright }
76648417fb7SJames Wright 
76748417fb7SJames Wright // @brief Test the Newtonian State transformation functions, `StateFrom*`
76848417fb7SJames Wright static PetscErrorCode UnitTests_Newtonian(Honee honee, NewtonianIdealGasContext gas) {
76948417fb7SJames Wright   Units            units = honee->units;
77048417fb7SJames Wright   const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, K = units->Kelvin;
77148417fb7SJames Wright   CeedScalar       rtol;
77248417fb7SJames Wright 
77348417fb7SJames Wright   PetscFunctionBeginUser;
77448417fb7SJames Wright   const CeedScalar T          = 200 * K;
77548417fb7SJames Wright   const CeedScalar rho        = 1.2 * kg / Cube(m);
77648417fb7SJames Wright   const CeedScalar P          = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T;
77748417fb7SJames Wright   const CeedScalar u_base     = 40 * m / sec;
77848417fb7SJames Wright   const CeedScalar u[3]       = {u_base, u_base * 1.1, u_base * 1.2};
77948417fb7SJames Wright   const CeedScalar e_kinetic  = 0.5 * Dot3(u, u);
78048417fb7SJames Wright   const CeedScalar e_internal = gas->cv * T;
78148417fb7SJames Wright   const CeedScalar e_total    = e_kinetic + e_internal;
78248417fb7SJames Wright   const CeedScalar gamma      = HeatCapacityRatio(gas);
78348417fb7SJames Wright   const CeedScalar entropy    = log(P) - gamma * log(rho);
78448417fb7SJames Wright   const CeedScalar rho_div_p  = rho / P;
78548417fb7SJames Wright   const CeedScalar Y0[5]      = {P, u[0], u[1], u[2], T};
78648417fb7SJames Wright   const CeedScalar U0[5]      = {rho, rho * u[0], rho * u[1], rho * u[2], rho * e_total};
78748417fb7SJames 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],
78848417fb7SJames Wright                                  -rho_div_p};
78948417fb7SJames Wright 
79048417fb7SJames Wright   rtol = 20 * CEED_EPSILON;
79148417fb7SJames Wright   PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol));
79248417fb7SJames Wright   PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol));
79348417fb7SJames Wright   PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol));
79448417fb7SJames Wright   PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, rtol, rtol, rtol));
79548417fb7SJames Wright   PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, rtol, rtol, rtol));
79648417fb7SJames Wright   PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, rtol, rtol, rtol));
79748417fb7SJames Wright 
79848417fb7SJames Wright   rtol = 5e-6;
79948417fb7SJames Wright   PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol));
80048417fb7SJames Wright   PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol));
80148417fb7SJames Wright   PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol));
80248417fb7SJames Wright   PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, 10 * rtol, rtol, rtol));
80348417fb7SJames Wright   PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, 5 * rtol, rtol, rtol));
80448417fb7SJames Wright   PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, 5 * rtol, 5 * rtol, 5 * rtol));
80548417fb7SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
80648417fb7SJames Wright }
807