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 14f7534ca8SJed Brown // For use with PetscOptionsEnum 156dfcbb05SJames Wright static const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "ENTROPY", "StateVariable", "STATEVAR_", NULL}; 16f7534ca8SJed Brown 17*48417fb7SJames 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)); 6165dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 620c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, honee->q_ceed)); 6301e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 6465dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 6565dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 6665dee3d2SJames Wright 67fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 6801e19bfaSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 69fff85bd3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 70fff85bd3SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 71236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qfctx)); 7265dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 7365dee3d2SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 7465dee3d2SJames Wright } 754c07ec22SJames Wright 7636038bbcSJames Wright /** 7736038bbcSJames Wright @brief Create RHS CeedOperator for direct projection of divergence of diffusive flux 7836038bbcSJames Wright 790c373b74SJames Wright @param[in] honee `Honee` context 8036038bbcSJames Wright @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 8136038bbcSJames Wright @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 8236038bbcSJames Wright **/ 83e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_NS(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) { 840c373b74SJames Wright Ceed ceed = honee->ceed; 8536038bbcSJames Wright NodalProjectionData projection = diff_flux_proj->projection; 8636038bbcSJames Wright CeedInt num_comp_q; 8736038bbcSJames Wright PetscInt dim, label_value = 0; 8836038bbcSJames Wright DMLabel domain_label = NULL; 89236a8ba6SJames Wright CeedQFunctionContext newtonian_qfctx = NULL; 9036038bbcSJames Wright 9136038bbcSJames Wright PetscFunctionBeginUser; 9236038bbcSJames Wright // -- Get Pre-requisite things 9336038bbcSJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 94e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q)); 9536038bbcSJames Wright 960880fbb6SJames Wright { // Get newtonian QF context 9736038bbcSJames Wright CeedOperator *sub_ops; 9836038bbcSJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 9936038bbcSJames Wright 1000c373b74SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_ifunction, &sub_ops)); 101236a8ba6SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newtonian_qfctx)); 10236038bbcSJames Wright } 10336038bbcSJames Wright PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_rhs)); 10436038bbcSJames Wright { // Add the volume integral CeedOperator 10536038bbcSJames Wright CeedQFunction qf_rhs_volume; 10636038bbcSJames Wright CeedOperator op_rhs_volume; 10736038bbcSJames Wright CeedVector q_data; 1080880fbb6SJames Wright CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL; 1090880fbb6SJames Wright CeedBasis basis_diff_flux = NULL; 11036038bbcSJames Wright CeedInt q_data_size; 11136038bbcSJames Wright 1120880fbb6SJames Wright PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL)); 113e3663b90SJames 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, 114e3663b90SJames Wright &q_data_size)); 1150c373b74SJames Wright switch (honee->phys->state_var) { 11636038bbcSJames Wright case STATEVAR_PRIMITIVE: 11736038bbcSJames Wright PetscCallCeed(ceed, 11836038bbcSJames Wright CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Prim, DivDiffusiveFluxVolumeRHS_NS_Prim_loc, &qf_rhs_volume)); 11936038bbcSJames Wright break; 12036038bbcSJames Wright case STATEVAR_CONSERVATIVE: 12136038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Conserv, DivDiffusiveFluxVolumeRHS_NS_Conserv_loc, 12236038bbcSJames Wright &qf_rhs_volume)); 12336038bbcSJames Wright break; 12436038bbcSJames Wright case STATEVAR_ENTROPY: 12536038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Entropy, DivDiffusiveFluxVolumeRHS_NS_Entropy_loc, 12636038bbcSJames Wright &qf_rhs_volume)); 12736038bbcSJames Wright break; 12836038bbcSJames Wright } 12936038bbcSJames Wright 130236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, newtonian_qfctx)); 13136038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "q", num_comp_q, CEED_EVAL_INTERP)); 13236038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 13336038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE)); 13436038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD)); 13536038bbcSJames Wright 13636038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume)); 137e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 138e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 13936038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 14036038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE)); 14136038bbcSJames Wright 14236038bbcSJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_volume)); 14336038bbcSJames Wright 14436038bbcSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 14536038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 1460880fbb6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume)); 1470880fbb6SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 14836038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume)); 14936038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume)); 15036038bbcSJames Wright } 15136038bbcSJames Wright 15236038bbcSJames Wright { // Add the boundary integral CeedOperator 15336038bbcSJames Wright CeedQFunction qf_rhs_boundary; 15436038bbcSJames Wright DMLabel face_sets_label; 15536038bbcSJames Wright PetscInt num_face_set_values, *face_set_values; 15636038bbcSJames Wright CeedInt q_data_size; 15736038bbcSJames Wright 15836038bbcSJames Wright // -- Build RHS operator 1590c373b74SJames Wright switch (honee->phys->state_var) { 16036038bbcSJames Wright case STATEVAR_PRIMITIVE: 16136038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Prim, DivDiffusiveFluxBoundaryRHS_NS_Prim_loc, 16236038bbcSJames Wright &qf_rhs_boundary)); 16336038bbcSJames Wright break; 16436038bbcSJames Wright case STATEVAR_CONSERVATIVE: 16536038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Conserv, DivDiffusiveFluxBoundaryRHS_NS_Conserv_loc, 16636038bbcSJames Wright &qf_rhs_boundary)); 16736038bbcSJames Wright break; 16836038bbcSJames Wright case STATEVAR_ENTROPY: 16936038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Entropy, DivDiffusiveFluxBoundaryRHS_NS_Entropy_loc, 17036038bbcSJames Wright &qf_rhs_boundary)); 17136038bbcSJames Wright break; 17236038bbcSJames Wright } 17336038bbcSJames Wright 1740c373b74SJames Wright PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size)); 175236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, newtonian_qfctx)); 17636038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "q", num_comp_q, CEED_EVAL_INTERP)); 17736038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 17836038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE)); 17936038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP)); 18036038bbcSJames Wright 18136038bbcSJames Wright PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label)); 18236038bbcSJames Wright PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values)); 18336038bbcSJames Wright for (PetscInt f = 0; f < num_face_set_values; f++) { 18436038bbcSJames Wright DMLabel face_orientation_label; 18536038bbcSJames Wright PetscInt num_orientations_values, *orientation_values; 18636038bbcSJames Wright 18736038bbcSJames Wright { 18836038bbcSJames Wright char *face_orientation_label_name; 18936038bbcSJames Wright 19036038bbcSJames Wright PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name)); 19136038bbcSJames Wright PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label)); 19236038bbcSJames Wright PetscCall(PetscFree(face_orientation_label_name)); 19336038bbcSJames Wright } 19436038bbcSJames Wright PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values)); 19536038bbcSJames Wright for (PetscInt o = 0; o < num_orientations_values; o++) { 19636038bbcSJames Wright CeedOperator op_rhs_boundary; 19736038bbcSJames Wright CeedBasis basis_q, basis_diff_flux_boundary; 19836038bbcSJames Wright CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary; 19936038bbcSJames Wright CeedVector q_data; 20036038bbcSJames Wright CeedInt q_data_size; 20136038bbcSJames Wright PetscInt orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1; 20236038bbcSJames Wright 2030c373b74SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q)); 2040c373b74SJames Wright PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q)); 20536038bbcSJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, 20636038bbcSJames Wright &elem_restr_diff_flux_boundary)); 2070880fbb6SJames Wright PetscCall(CreateBasisFromPlex(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary)); 208e3663b90SJames Wright PetscCall( 209e3663b90SJames Wright QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, honee->x_coord, &elem_restr_qdata, &q_data, &q_data_size)); 21036038bbcSJames Wright 21136038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary)); 21236038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 21336038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 21436038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data)); 21536038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary, 21636038bbcSJames Wright CEED_VECTOR_ACTIVE)); 21736038bbcSJames Wright 21836038bbcSJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_boundary)); 21936038bbcSJames Wright 22036038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary)); 22136038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata)); 22236038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 22336038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary)); 22436038bbcSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 22536038bbcSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary)); 22636038bbcSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 22736038bbcSJames Wright } 22836038bbcSJames Wright PetscCall(PetscFree(orientation_values)); 22936038bbcSJames Wright } 23036038bbcSJames Wright PetscCall(PetscFree(face_set_values)); 23136038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary)); 23236038bbcSJames Wright } 23336038bbcSJames Wright 234236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newtonian_qfctx)); 23536038bbcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 23636038bbcSJames Wright } 23736038bbcSJames Wright 23836038bbcSJames Wright /** 23936038bbcSJames Wright @brief Create RHS CeedOperator for indirect projection of divergence of diffusive flux 24036038bbcSJames Wright 2410c373b74SJames Wright @param[in] honee `Honee` context 24236038bbcSJames Wright @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 24336038bbcSJames Wright @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 24436038bbcSJames Wright **/ 245e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Indirect_NS(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) { 2460c373b74SJames Wright Ceed ceed = honee->ceed; 24736038bbcSJames Wright NodalProjectionData projection = diff_flux_proj->projection; 24836038bbcSJames Wright CeedBasis basis_diff_flux; 24936038bbcSJames Wright CeedElemRestriction elem_restr_diff_flux, elem_restr_qd; 25036038bbcSJames Wright CeedVector q_data; 25136038bbcSJames Wright CeedInt num_comp_q, q_data_size; 25236038bbcSJames Wright PetscInt dim; 25336038bbcSJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 25436038bbcSJames Wright DMLabel domain_label = NULL; 25536038bbcSJames Wright CeedQFunction qf_rhs; 256236a8ba6SJames Wright CeedQFunctionContext newtonian_qfctx = NULL; 25736038bbcSJames Wright 25836038bbcSJames Wright PetscFunctionBeginUser; 25936038bbcSJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 260e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q)); 26136038bbcSJames Wright 26240b78511SJames Wright { // Get newtonian QF context 26336038bbcSJames Wright CeedOperator *sub_ops; 26436038bbcSJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 26536038bbcSJames Wright 2660c373b74SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_ifunction, &sub_ops)); 267236a8ba6SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newtonian_qfctx)); 26836038bbcSJames Wright } 26936038bbcSJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux)); 27036038bbcSJames Wright PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux)); 271e3663b90SJames 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, 272e3663b90SJames Wright &q_data_size)); 27336038bbcSJames Wright 2740c373b74SJames Wright switch (honee->phys->state_var) { 27536038bbcSJames Wright case STATEVAR_PRIMITIVE: 27636038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Prim, DiffusiveFluxRHS_NS_Prim_loc, &qf_rhs)); 27736038bbcSJames Wright break; 27836038bbcSJames Wright case STATEVAR_CONSERVATIVE: 27936038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Conserv, DiffusiveFluxRHS_NS_Conserv_loc, &qf_rhs)); 28036038bbcSJames Wright break; 28136038bbcSJames Wright case STATEVAR_ENTROPY: 28236038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Entropy, DiffusiveFluxRHS_NS_Entropy_loc, &qf_rhs)); 28336038bbcSJames Wright break; 28436038bbcSJames Wright } 28536038bbcSJames Wright 286236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, newtonian_qfctx)); 28736038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "q", num_comp_q, CEED_EVAL_INTERP)); 28836038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 28936038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE)); 29036038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP)); 29136038bbcSJames Wright 29236038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, op_rhs)); 293e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 294e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 29536038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 29636038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 29736038bbcSJames Wright 29836038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs)); 299236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newtonian_qfctx)); 30036038bbcSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 30136038bbcSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 30236038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 30336038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 30436038bbcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 30536038bbcSJames Wright } 30636038bbcSJames Wright 307991aef52SJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 308b7f03f12SJed Brown SetupContext setup_context; 3090c373b74SJames Wright Honee honee = *(Honee *)ctx; 3100c373b74SJames Wright CeedInt degree = honee->app_ctx->degree; 3113a8779fbSJames Wright StabilizationType stab; 3123636f6a4SJames Wright StateVariable state_var; 3130c373b74SJames Wright MPI_Comm comm = honee->comm; 3140c373b74SJames Wright Ceed ceed = honee->ceed; 3153a8779fbSJames Wright PetscBool implicit; 3160d8cd818SJames Wright PetscBool unit_tests; 31715a3537eSJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 318e07531f7SJames Wright CeedQFunctionContext newtonian_ig_qfctx; 3193a8779fbSJames Wright 32015a3537eSJed Brown PetscFunctionBeginUser; 3212b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &setup_context)); 3222b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &newtonian_ig_ctx)); 3233a8779fbSJames Wright 3243a8779fbSJames Wright // ------------------------------------------------------ 3253a8779fbSJames Wright // Setup Generic Newtonian IG Problem 3263a8779fbSJames Wright // ------------------------------------------------------ 32728160fc2SJames Wright problem->jac_data_size_vol = 14; 328b01ba163SJames Wright problem->jac_data_size_sur = 11; 32958ce1233SJames Wright problem->compute_exact_solution_error = PETSC_FALSE; 330cbe60e31SLeila Ghaffari problem->print_info = PRINT_NEWTONIAN; 3313a8779fbSJames Wright 3320c373b74SJames Wright PetscCall(DivDiffFluxProjectionCreate(honee, 4, &honee->diff_flux_proj)); 3330c373b74SJames Wright if (honee->diff_flux_proj) { 3340c373b74SJames Wright DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj; 33536038bbcSJames Wright NodalProjectionData projection = diff_flux_proj->projection; 33636038bbcSJames Wright 33736038bbcSJames Wright diff_flux_proj->CreateRHSOperator_Direct = DivDiffFluxProjectionCreateRHS_Direct_NS; 33836038bbcSJames Wright diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_NS; 33936038bbcSJames Wright 3400c373b74SJames Wright switch (honee->diff_flux_proj->method) { 34136038bbcSJames Wright case DIV_DIFF_FLUX_PROJ_DIRECT: { 34236038bbcSJames Wright PetscSection section; 34336038bbcSJames Wright 34436038bbcSJames Wright PetscCall(DMGetLocalSection(projection->dm, §ion)); 34536038bbcSJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 34636038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_MomentumX")); 34736038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "DivDiffusiveFlux_MomentumY")); 34836038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "DivDiffusiveFlux_MomentumZ")); 34936038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "DivDiffusiveFlux_Energy")); 35036038bbcSJames Wright } break; 35136038bbcSJames Wright case DIV_DIFF_FLUX_PROJ_INDIRECT: { 35236038bbcSJames Wright PetscSection section; 35336038bbcSJames Wright 35436038bbcSJames Wright PetscCall(DMGetLocalSection(projection->dm, §ion)); 35536038bbcSJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 35636038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_MomentumXX")); 35736038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_MomentumXY")); 35836038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_MomentumXZ")); 35936038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "DiffusiveFlux_MomentumYX")); 36036038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "DiffusiveFlux_MomentumYY")); 36136038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 5, "DiffusiveFlux_MomentumYZ")); 36236038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 6, "DiffusiveFlux_MomentumZX")); 36336038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 7, "DiffusiveFlux_MomentumZY")); 36436038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 8, "DiffusiveFlux_MomentumZZ")); 36536038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 9, "DiffusiveFlux_EnergyX")); 36636038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 10, "DiffusiveFlux_EnergyY")); 36736038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 11, "DiffusiveFlux_EnergyZ")); 36836038bbcSJames Wright } break; 36936038bbcSJames Wright case DIV_DIFF_FLUX_PROJ_NONE: 3700c373b74SJames Wright SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 3710c373b74SJames Wright DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]); 37236038bbcSJames Wright break; 37336038bbcSJames Wright } 37436038bbcSJames Wright } 37536038bbcSJames Wright 3763a8779fbSJames Wright // ------------------------------------------------------ 3773a8779fbSJames Wright // Create the libCEED context 3783a8779fbSJames Wright // ------------------------------------------------------ 3793a8779fbSJames Wright CeedScalar cv = 717.; // J/(kg K) 3803a8779fbSJames Wright CeedScalar cp = 1004.; // J/(kg K) 381d9bb1cdbSJames Wright CeedScalar g[3] = {0, 0, 0}; // m/s^2 3823a8779fbSJames Wright CeedScalar lambda = -2. / 3.; // - 383bb8a0c61SJames Wright CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 3843a8779fbSJames Wright CeedScalar k = 0.02638; // W/(m K) 3854351d27fSLeila Ghaffari CeedScalar c_tau = 0.5 / degree; // - 386bb8a0c61SJames Wright CeedScalar Ctau_t = 1.0; // - 387b5786772SLeila Ghaffari CeedScalar Cv_func[3] = {36, 60, 128}; 388c176988bSLeila Ghaffari CeedScalar Ctau_v = Cv_func[(CeedInt)Min(3, degree) - 1]; 3892a0eee6eSLeila Ghaffari CeedScalar Ctau_C = 0.25 / degree; 3902a0eee6eSLeila Ghaffari CeedScalar Ctau_M = 0.25 / degree; 3912a0eee6eSLeila Ghaffari CeedScalar Ctau_E = 0.125; 3923a8779fbSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 3932b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 394493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 3953a8779fbSJames Wright 396b8fb7609SAdeleke O. Bankole StatePrimitive reference = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15}; 397a213b8aaSJames Wright CeedScalar idl_decay_time = -1, idl_start = 0, idl_length = 0, idl_pressure = reference.pressure; 398e7754af5SKenneth E. Jansen PetscBool idl_enable = PETSC_FALSE; 399b8fb7609SAdeleke O. Bankole 4003a8779fbSJames Wright // ------------------------------------------------------ 4013a8779fbSJames Wright // Create the PETSc context 4023a8779fbSJames Wright // ------------------------------------------------------ 403bb8a0c61SJames Wright PetscScalar meter = 1; // 1 meter in scaled length units 404bb8a0c61SJames Wright PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 405bb8a0c61SJames Wright PetscScalar second = 1; // 1 second in scaled time units 4063a8779fbSJames Wright PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 4073a8779fbSJames Wright PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 4083a8779fbSJames Wright 4093a8779fbSJames Wright // ------------------------------------------------------ 4103a8779fbSJames Wright // Command line Options 4113a8779fbSJames Wright // ------------------------------------------------------ 412d9bb1cdbSJames Wright PetscBool given_option = PETSC_FALSE; 4132b916ea7SJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL); 414cbe60e31SLeila Ghaffari // -- Conservative vs Primitive variables 4152b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE), 4162b916ea7SJeremy L Thompson (PetscEnum *)&state_var, NULL)); 4173636f6a4SJames Wright 4183636f6a4SJames Wright switch (state_var) { 4193636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 420e07531f7SJames Wright problem->ics.qf_func_ptr = ICsNewtonianIG_Conserv; 421e07531f7SJames Wright problem->ics.qf_loc = ICsNewtonianIG_Conserv_loc; 422e07531f7SJames Wright problem->apply_vol_rhs.qf_func_ptr = RHSFunction_Newtonian; 423e07531f7SJames Wright problem->apply_vol_rhs.qf_loc = RHSFunction_Newtonian_loc; 424e07531f7SJames Wright problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Conserv; 425e07531f7SJames Wright problem->apply_vol_ifunction.qf_loc = IFunction_Newtonian_Conserv_loc; 426e07531f7SJames Wright problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Conserv; 427e07531f7SJames Wright problem->apply_vol_ijacobian.qf_loc = IJacobian_Newtonian_Conserv_loc; 428e07531f7SJames Wright problem->apply_inflow.qf_func_ptr = BoundaryIntegral_Conserv; 429e07531f7SJames Wright problem->apply_inflow.qf_loc = BoundaryIntegral_Conserv_loc; 430e07531f7SJames Wright problem->apply_inflow_jacobian.qf_func_ptr = BoundaryIntegral_Jacobian_Conserv; 431e07531f7SJames Wright problem->apply_inflow_jacobian.qf_loc = BoundaryIntegral_Jacobian_Conserv_loc; 4323636f6a4SJames Wright break; 4333636f6a4SJames Wright case STATEVAR_PRIMITIVE: 434e07531f7SJames Wright problem->ics.qf_func_ptr = ICsNewtonianIG_Prim; 435e07531f7SJames Wright problem->ics.qf_loc = ICsNewtonianIG_Prim_loc; 436e07531f7SJames Wright problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Prim; 437e07531f7SJames Wright problem->apply_vol_ifunction.qf_loc = IFunction_Newtonian_Prim_loc; 438e07531f7SJames Wright problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Prim; 439e07531f7SJames Wright problem->apply_vol_ijacobian.qf_loc = IJacobian_Newtonian_Prim_loc; 440e07531f7SJames Wright problem->apply_inflow.qf_func_ptr = BoundaryIntegral_Prim; 441e07531f7SJames Wright problem->apply_inflow.qf_loc = BoundaryIntegral_Prim_loc; 442e07531f7SJames Wright problem->apply_inflow_jacobian.qf_func_ptr = BoundaryIntegral_Jacobian_Prim; 443e07531f7SJames Wright problem->apply_inflow_jacobian.qf_loc = BoundaryIntegral_Jacobian_Prim_loc; 4443636f6a4SJames Wright break; 445c62f7daeSJames Wright case STATEVAR_ENTROPY: 446e07531f7SJames Wright problem->ics.qf_func_ptr = ICsNewtonianIG_Entropy; 447e07531f7SJames Wright problem->ics.qf_loc = ICsNewtonianIG_Entropy_loc; 448e07531f7SJames Wright problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Entropy; 449e07531f7SJames Wright problem->apply_vol_ifunction.qf_loc = IFunction_Newtonian_Entropy_loc; 450e07531f7SJames Wright problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Entropy; 451e07531f7SJames Wright problem->apply_vol_ijacobian.qf_loc = IJacobian_Newtonian_Entropy_loc; 452e07531f7SJames Wright problem->apply_inflow.qf_func_ptr = BoundaryIntegral_Entropy; 453e07531f7SJames Wright problem->apply_inflow.qf_loc = BoundaryIntegral_Entropy_loc; 454e07531f7SJames Wright problem->apply_inflow_jacobian.qf_func_ptr = BoundaryIntegral_Jacobian_Entropy; 455e07531f7SJames Wright problem->apply_inflow_jacobian.qf_loc = BoundaryIntegral_Jacobian_Entropy_loc; 456c62f7daeSJames Wright break; 457cbe60e31SLeila Ghaffari } 4581485969bSJeremy L Thompson 4593a8779fbSJames Wright // -- Physics 4602b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL)); 4612b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL)); 4622b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL)); 4632b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL)); 4642b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL)); 4653a8779fbSJames Wright 46666d54740SJames Wright PetscInt dim = 3; 467d9bb1cdbSJames Wright PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL)); 468d9bb1cdbSJames Wright PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option)); 469d9bb1cdbSJames Wright if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim); 470d9bb1cdbSJames Wright 4712b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 4722b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 4732b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 4742b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL)); 4752b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL)); 4762b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL)); 4772b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL)); 4782b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 4792b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL)); 4803a8779fbSJames Wright 481db95602eSJames Wright dim = 3; 482b8fb7609SAdeleke O. Bankole PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL)); 483b8fb7609SAdeleke O. Bankole PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL)); 484b8fb7609SAdeleke O. Bankole PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL)); 485b8fb7609SAdeleke O. Bankole 4863a8779fbSJames Wright // -- Units 4872b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 4883a8779fbSJames Wright meter = fabs(meter); 4892b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 4903a8779fbSJames Wright kilogram = fabs(kilogram); 4912b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 4923a8779fbSJames Wright second = fabs(second); 4932b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL)); 4943a8779fbSJames Wright Kelvin = fabs(Kelvin); 4953a8779fbSJames Wright 4963a8779fbSJames Wright // -- Warnings 4975d28dccaSJames Wright PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP, 4985d28dccaSJames Wright "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n"); 4990c373b74SJames Wright PetscCheck(!(honee->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE && !implicit), comm, PETSC_ERR_SUP, 5005f952e8dSJames Wright "Projection of divergence of diffusive flux is not implemented for Navier-Stokes with explicit timestepping"); 501e7754af5SKenneth E. Jansen 502e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point", 503e7754af5SKenneth E. Jansen NULL, idl_decay_time, &idl_decay_time, &idl_enable)); 5049cbdf780SJames Wright PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero."); 5059cbdf780SJames Wright if (idl_decay_time < 0) idl_enable = PETSC_FALSE; 50628160fc2SJames Wright if (idl_enable) problem->jac_data_size_vol++; 507e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL)); 508e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL)); 509a213b8aaSJames Wright idl_pressure = reference.pressure; 510a213b8aaSJames Wright PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-reference_pressure`)", NULL, idl_pressure, 511a213b8aaSJames Wright &idl_pressure, NULL)); 5121485969bSJeremy L Thompson PetscOptionsEnd(); 5133a8779fbSJames Wright 51465dee3d2SJames Wright if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized; 51565dee3d2SJames Wright 5163a8779fbSJames Wright // ------------------------------------------------------ 5173a8779fbSJames Wright // Set up the PETSc context 5183a8779fbSJames Wright // ------------------------------------------------------ 5193a8779fbSJames Wright // -- Define derived units 5203a8779fbSJames Wright Pascal = kilogram / (meter * PetscSqr(second)); 5213a8779fbSJames Wright J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 5223a8779fbSJames Wright m_per_squared_s = meter / PetscSqr(second); 5233a8779fbSJames Wright W_per_m_K = kilogram * meter / (pow(second, 3) * Kelvin); 5243a8779fbSJames Wright 5250c373b74SJames Wright honee->units->meter = meter; 5260c373b74SJames Wright honee->units->kilogram = kilogram; 5270c373b74SJames Wright honee->units->second = second; 5280c373b74SJames Wright honee->units->Kelvin = Kelvin; 5290c373b74SJames Wright honee->units->Pascal = Pascal; 5300c373b74SJames Wright honee->units->J_per_kg_K = J_per_kg_K; 5310c373b74SJames Wright honee->units->m_per_squared_s = m_per_squared_s; 5320c373b74SJames Wright honee->units->W_per_m_K = W_per_m_K; 5333a8779fbSJames Wright 5343a8779fbSJames Wright // ------------------------------------------------------ 5353a8779fbSJames Wright // Set up the libCEED context 5363a8779fbSJames Wright // ------------------------------------------------------ 5373a8779fbSJames Wright // -- Scale variables to desired units 5383a8779fbSJames Wright cv *= J_per_kg_K; 5393a8779fbSJames Wright cp *= J_per_kg_K; 5403a8779fbSJames Wright mu *= Pascal * second; 5413a8779fbSJames Wright k *= W_per_m_K; 542493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter; 543493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s; 544b8fb7609SAdeleke O. Bankole reference.pressure *= Pascal; 545b8fb7609SAdeleke O. Bankole for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second; 546b8fb7609SAdeleke O. Bankole reference.temperature *= Kelvin; 5473a8779fbSJames Wright 5483a8779fbSJames Wright // -- Solver Settings 5490c373b74SJames Wright honee->phys->implicit = implicit; 5500c373b74SJames Wright honee->phys->state_var = state_var; 5513a8779fbSJames Wright 5523a8779fbSJames Wright // -- QFunction Context 55315a3537eSJed Brown newtonian_ig_ctx->lambda = lambda; 55415a3537eSJed Brown newtonian_ig_ctx->mu = mu; 55515a3537eSJed Brown newtonian_ig_ctx->k = k; 55615a3537eSJed Brown newtonian_ig_ctx->cv = cv; 55715a3537eSJed Brown newtonian_ig_ctx->cp = cp; 55815a3537eSJed Brown newtonian_ig_ctx->c_tau = c_tau; 55915a3537eSJed Brown newtonian_ig_ctx->Ctau_t = Ctau_t; 56015a3537eSJed Brown newtonian_ig_ctx->Ctau_v = Ctau_v; 56115a3537eSJed Brown newtonian_ig_ctx->Ctau_C = Ctau_C; 56215a3537eSJed Brown newtonian_ig_ctx->Ctau_M = Ctau_M; 56315a3537eSJed Brown newtonian_ig_ctx->Ctau_E = Ctau_E; 56415a3537eSJed Brown newtonian_ig_ctx->stabilization = stab; 5658085925cSJames Wright newtonian_ig_ctx->is_implicit = implicit; 5663636f6a4SJames Wright newtonian_ig_ctx->state_var = state_var; 567e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_enable = idl_enable; 568e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * second); 569e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_start = idl_start * meter; 570e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_length = idl_length * meter; 571a213b8aaSJames Wright newtonian_ig_ctx->idl_pressure = idl_pressure; 5720c373b74SJames Wright newtonian_ig_ctx->divFdiff_method = honee->app_ctx->divFdiffproj_method; 5732b916ea7SJeremy L Thompson PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3)); 5743a8779fbSJames Wright 575b8fb7609SAdeleke O. Bankole // -- Setup Context 576b8fb7609SAdeleke O. Bankole setup_context->reference = reference; 577b8fb7609SAdeleke O. Bankole setup_context->gas = *newtonian_ig_ctx; 578b8fb7609SAdeleke O. Bankole setup_context->lx = domain_size[0]; 579b8fb7609SAdeleke O. Bankole setup_context->ly = domain_size[1]; 580b8fb7609SAdeleke O. Bankole setup_context->lz = domain_size[2]; 581b8fb7609SAdeleke O. Bankole setup_context->time = 0; 582b8fb7609SAdeleke O. Bankole 5830c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &problem->ics.qfctx)); 584e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 585e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc)); 586e07531f7SJames Wright PetscCallCeed( 587e07531f7SJames Wright ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfctx, "evaluation time", offsetof(struct SetupContext_, time), 1, "Time of evaluation")); 58815a3537eSJed Brown 5890c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &newtonian_ig_qfctx)); 590e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(newtonian_ig_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx)); 591e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 592e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1, 593b4c37c5cSJames Wright "Size of timestep, delta t")); 594e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "ijacobian time shift", 595b4c37c5cSJames Wright offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1, 596b4c37c5cSJames Wright "Shift for mass matrix in IJacobian")); 597e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1, 598b4c37c5cSJames Wright "Current solution time")); 599a8bca8acSJames Wright 600e07531f7SJames Wright problem->apply_vol_rhs.qfctx = newtonian_ig_qfctx; 601e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ifunction.qfctx)); 602e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ijacobian.qfctx)); 603e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_inflow.qfctx)); 604e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_inflow_jacobian.qfctx)); 605f0b65372SJed Brown 6069ed3d70dSJames Wright if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 6079ed3d70dSJames Wright if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 608e07531f7SJames Wright if (bc->num_slip > 0) PetscCall(SlipBCSetup(problem, dm, ctx, newtonian_ig_qfctx)); 6099ed3d70dSJames Wright 610f0b65372SJed Brown if (unit_tests) { 6110c373b74SJames Wright PetscCall(UnitTests_Newtonian(honee, newtonian_ig_ctx)); 612f0b65372SJed Brown } 613d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 6143a8779fbSJames Wright } 6153a8779fbSJames Wright 6160c373b74SJames Wright PetscErrorCode PRINT_NEWTONIAN(Honee honee, ProblemData problem, AppCtx app_ctx) { 6170c373b74SJames Wright MPI_Comm comm = honee->comm; 6180c373b74SJames Wright Ceed ceed = honee->ceed; 61915a3537eSJed Brown NewtonianIdealGasContext newtonian_ctx; 62015a3537eSJed Brown 6213a8779fbSJames Wright PetscFunctionBeginUser; 622e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ctx)); 6232b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 62415a3537eSJed Brown " Problem:\n" 62515a3537eSJed Brown " Problem Name : %s\n" 62615a3537eSJed Brown " Stabilization : %s\n", 6272b916ea7SJeremy L Thompson app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization])); 628e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ctx)); 629d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 6303a8779fbSJames Wright } 631*48417fb7SJames Wright 632*48417fb7SJames Wright static PetscErrorCode CheckQWithTolerance(const CeedScalar Q_s[5], const CeedScalar Q_a[5], const CeedScalar Q_b[5], const char *name, 633*48417fb7SJames Wright PetscReal rtol_0, PetscReal rtol_u, PetscReal rtol_4) { 634*48417fb7SJames Wright CeedScalar relative_error[5]; // relative error 635*48417fb7SJames Wright CeedScalar divisor_threshold = 10 * CEED_EPSILON; 636*48417fb7SJames Wright 637*48417fb7SJames Wright PetscFunctionBeginUser; 638*48417fb7SJames Wright relative_error[0] = (Q_a[0] - Q_b[0]) / (fabs(Q_s[0]) > divisor_threshold ? Q_s[0] : 1); 639*48417fb7SJames Wright relative_error[4] = (Q_a[4] - Q_b[4]) / (fabs(Q_s[4]) > divisor_threshold ? Q_s[4] : 1); 640*48417fb7SJames Wright 641*48417fb7SJames Wright CeedScalar u_magnitude = sqrt(Square(Q_s[1]) + Square(Q_s[2]) + Square(Q_s[3])); 642*48417fb7SJames Wright CeedScalar u_divisor = u_magnitude > divisor_threshold ? u_magnitude : 1; 643*48417fb7SJames Wright for (int i = 1; i < 4; i++) { 644*48417fb7SJames Wright relative_error[i] = (Q_a[i] - Q_b[i]) / u_divisor; 645*48417fb7SJames Wright } 646*48417fb7SJames Wright 647*48417fb7SJames Wright if (fabs(relative_error[0]) >= rtol_0) { 648*48417fb7SJames Wright printf("%s[0] error %g (expected %.10e, got %.10e)\n", name, relative_error[0], Q_s[0], Q_a[0]); 649*48417fb7SJames Wright } 650*48417fb7SJames Wright for (int i = 1; i < 4; i++) { 651*48417fb7SJames Wright if (fabs(relative_error[i]) >= rtol_u) { 652*48417fb7SJames Wright printf("%s[%d] error %g (expected %.10e, got %.10e)\n", name, i, relative_error[i], Q_s[i], Q_a[i]); 653*48417fb7SJames Wright } 654*48417fb7SJames Wright } 655*48417fb7SJames Wright if (fabs(relative_error[4]) >= rtol_4) { 656*48417fb7SJames Wright printf("%s[4] error %g (expected %.10e, got %.10e)\n", name, relative_error[4], Q_s[4], Q_a[4]); 657*48417fb7SJames Wright } 658*48417fb7SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 659*48417fb7SJames Wright } 660*48417fb7SJames Wright 661*48417fb7SJames Wright // @brief Verify `StateFromQ` by converting A0 -> B0 -> A0_test, where A0 should equal A0_test 662*48417fb7SJames Wright static PetscErrorCode TestState(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5], 663*48417fb7SJames Wright CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) { 664*48417fb7SJames Wright CeedScalar B0[5], A0_test[5]; 665*48417fb7SJames Wright char buf[128]; 666*48417fb7SJames Wright const char *const StateVariables_Initial[] = {"U", "Y", "V"}; 667*48417fb7SJames Wright 668*48417fb7SJames Wright PetscFunctionBeginUser; 669*48417fb7SJames Wright const char *A_initial = StateVariables_Initial[state_var_A]; 670*48417fb7SJames Wright const char *B_initial = StateVariables_Initial[state_var_B]; 671*48417fb7SJames Wright 672*48417fb7SJames Wright State state_A0 = StateFromQ(gas, A0, state_var_A); 673*48417fb7SJames Wright StateToQ(gas, state_A0, B0, state_var_B); 674*48417fb7SJames Wright State state_B0 = StateFromQ(gas, B0, state_var_B); 675*48417fb7SJames Wright StateToQ(gas, state_B0, A0_test, state_var_A); 676*48417fb7SJames Wright 677*48417fb7SJames Wright snprintf(buf, sizeof buf, "%s->%s->%s: %s", A_initial, B_initial, A_initial, A_initial); 678*48417fb7SJames Wright PetscCall(CheckQWithTolerance(A0, A0_test, A0, buf, rtol_0, rtol_u, rtol_4)); 679*48417fb7SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 680*48417fb7SJames Wright } 681*48417fb7SJames Wright 682*48417fb7SJames Wright // @brief Verify `StateFromQ_fwd` via a finite difference approximation 683*48417fb7SJames Wright static PetscErrorCode TestState_fwd(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5], 684*48417fb7SJames Wright CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) { 685*48417fb7SJames Wright CeedScalar eps = 4e-7; // Finite difference step 686*48417fb7SJames Wright char buf[128]; 687*48417fb7SJames Wright const char *const StateVariables_Initial[] = {"U", "Y", "V"}; 688*48417fb7SJames Wright 689*48417fb7SJames Wright PetscFunctionBeginUser; 690*48417fb7SJames Wright const char *A_initial = StateVariables_Initial[state_var_A]; 691*48417fb7SJames Wright const char *B_initial = StateVariables_Initial[state_var_B]; 692*48417fb7SJames Wright State state_0 = StateFromQ(gas, A0, state_var_A); 693*48417fb7SJames Wright 694*48417fb7SJames Wright for (int i = 0; i < 5; i++) { 695*48417fb7SJames Wright CeedScalar dB[5] = {0.}, dB_fd[5] = {0.}; 696*48417fb7SJames Wright { // Calculate dB using State functions 697*48417fb7SJames Wright CeedScalar dA[5] = {0}; 698*48417fb7SJames Wright 699*48417fb7SJames Wright dA[i] = A0[i]; 700*48417fb7SJames Wright State dstate_0 = StateFromQ_fwd(gas, state_0, dA, state_var_A); 701*48417fb7SJames Wright StateToQ_fwd(gas, state_0, dstate_0, dB, state_var_B); 702*48417fb7SJames Wright } 703*48417fb7SJames Wright 704*48417fb7SJames Wright { // Calculate dB_fd via finite difference approximation 705*48417fb7SJames Wright CeedScalar A1[5], B0[5], B1[5]; 706*48417fb7SJames Wright 707*48417fb7SJames Wright for (int j = 0; j < 5; j++) A1[j] = (1 + eps * (i == j)) * A0[j]; 708*48417fb7SJames Wright State state_1 = StateFromQ(gas, A1, state_var_A); 709*48417fb7SJames Wright StateToQ(gas, state_0, B0, state_var_B); 710*48417fb7SJames Wright StateToQ(gas, state_1, B1, state_var_B); 711*48417fb7SJames Wright for (int j = 0; j < 5; j++) dB_fd[j] = (B1[j] - B0[j]) / eps; 712*48417fb7SJames Wright } 713*48417fb7SJames Wright 714*48417fb7SJames 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); 715*48417fb7SJames Wright PetscCall(CheckQWithTolerance(dB_fd, dB, dB_fd, buf, rtol_0, rtol_u, rtol_4)); 716*48417fb7SJames Wright } 717*48417fb7SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 718*48417fb7SJames Wright } 719*48417fb7SJames Wright 720*48417fb7SJames Wright // @brief Test the Newtonian State transformation functions, `StateFrom*` 721*48417fb7SJames Wright static PetscErrorCode UnitTests_Newtonian(Honee honee, NewtonianIdealGasContext gas) { 722*48417fb7SJames Wright Units units = honee->units; 723*48417fb7SJames Wright const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, K = units->Kelvin; 724*48417fb7SJames Wright CeedScalar rtol; 725*48417fb7SJames Wright 726*48417fb7SJames Wright PetscFunctionBeginUser; 727*48417fb7SJames Wright const CeedScalar T = 200 * K; 728*48417fb7SJames Wright const CeedScalar rho = 1.2 * kg / Cube(m); 729*48417fb7SJames Wright const CeedScalar P = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T; 730*48417fb7SJames Wright const CeedScalar u_base = 40 * m / sec; 731*48417fb7SJames Wright const CeedScalar u[3] = {u_base, u_base * 1.1, u_base * 1.2}; 732*48417fb7SJames Wright const CeedScalar e_kinetic = 0.5 * Dot3(u, u); 733*48417fb7SJames Wright const CeedScalar e_internal = gas->cv * T; 734*48417fb7SJames Wright const CeedScalar e_total = e_kinetic + e_internal; 735*48417fb7SJames Wright const CeedScalar gamma = HeatCapacityRatio(gas); 736*48417fb7SJames Wright const CeedScalar entropy = log(P) - gamma * log(rho); 737*48417fb7SJames Wright const CeedScalar rho_div_p = rho / P; 738*48417fb7SJames Wright const CeedScalar Y0[5] = {P, u[0], u[1], u[2], T}; 739*48417fb7SJames Wright const CeedScalar U0[5] = {rho, rho * u[0], rho * u[1], rho * u[2], rho * e_total}; 740*48417fb7SJames 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], 741*48417fb7SJames Wright -rho_div_p}; 742*48417fb7SJames Wright 743*48417fb7SJames Wright rtol = 20 * CEED_EPSILON; 744*48417fb7SJames Wright PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol)); 745*48417fb7SJames Wright PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol)); 746*48417fb7SJames Wright PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol)); 747*48417fb7SJames Wright PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, rtol, rtol, rtol)); 748*48417fb7SJames Wright PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, rtol, rtol, rtol)); 749*48417fb7SJames Wright PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, rtol, rtol, rtol)); 750*48417fb7SJames Wright 751*48417fb7SJames Wright rtol = 5e-6; 752*48417fb7SJames Wright PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol)); 753*48417fb7SJames Wright PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol)); 754*48417fb7SJames Wright PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol)); 755*48417fb7SJames Wright PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, 10 * rtol, rtol, rtol)); 756*48417fb7SJames Wright PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, 5 * rtol, rtol, rtol)); 757*48417fb7SJames Wright PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, 5 * rtol, 5 * rtol, 5 * rtol)); 758*48417fb7SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 759*48417fb7SJames Wright } 760