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 17c62f7daeSJames Wright static PetscErrorCode CheckQWithTolerance(const CeedScalar Q_s[5], const CeedScalar Q_a[5], const CeedScalar Q_b[5], const char *name, 18c62f7daeSJames Wright PetscReal rtol_0, PetscReal rtol_u, PetscReal rtol_4) { 19c62f7daeSJames Wright CeedScalar relative_error[5]; // relative error 20c62f7daeSJames Wright CeedScalar divisor_threshold = 10 * CEED_EPSILON; 2106f41313SJames Wright 2206f41313SJames Wright PetscFunctionBeginUser; 23c62f7daeSJames Wright relative_error[0] = (Q_a[0] - Q_b[0]) / (fabs(Q_s[0]) > divisor_threshold ? Q_s[0] : 1); 24c62f7daeSJames Wright relative_error[4] = (Q_a[4] - Q_b[4]) / (fabs(Q_s[4]) > divisor_threshold ? Q_s[4] : 1); 25c62f7daeSJames Wright 26c62f7daeSJames Wright CeedScalar u_magnitude = sqrt(Square(Q_s[1]) + Square(Q_s[2]) + Square(Q_s[3])); 27c62f7daeSJames Wright CeedScalar u_divisor = u_magnitude > divisor_threshold ? u_magnitude : 1; 28c62f7daeSJames Wright for (int i = 1; i < 4; i++) { 29c62f7daeSJames Wright relative_error[i] = (Q_a[i] - Q_b[i]) / u_divisor; 302b916ea7SJeremy L Thompson } 31c62f7daeSJames Wright 32c62f7daeSJames Wright if (fabs(relative_error[0]) >= rtol_0) { 33c62f7daeSJames Wright printf("%s[0] error %g (expected %.10e, got %.10e)\n", name, relative_error[0], Q_s[0], Q_a[0]); 34c62f7daeSJames Wright } 35c62f7daeSJames Wright for (int i = 1; i < 4; i++) { 36c62f7daeSJames Wright if (fabs(relative_error[i]) >= rtol_u) { 37c62f7daeSJames Wright printf("%s[%d] error %g (expected %.10e, got %.10e)\n", name, i, relative_error[i], Q_s[i], Q_a[i]); 38c62f7daeSJames Wright } 39c62f7daeSJames Wright } 40c62f7daeSJames Wright if (fabs(relative_error[4]) >= rtol_4) { 41c62f7daeSJames Wright printf("%s[4] error %g (expected %.10e, got %.10e)\n", name, relative_error[4], Q_s[4], Q_a[4]); 42c62f7daeSJames Wright } 43d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 44f0b65372SJed Brown } 45f0b65372SJed Brown 46c62f7daeSJames Wright // @brief Verify `StateFromQ` by converting A0 -> B0 -> A0_test, where A0 should equal A0_test 47c62f7daeSJames Wright static PetscErrorCode TestState(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5], 48c62f7daeSJames Wright CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) { 49c62f7daeSJames Wright CeedScalar B0[5], A0_test[5]; 50c62f7daeSJames Wright char buf[128]; 51c62f7daeSJames Wright const char *const StateVariables_Initial[] = {"U", "Y", "V"}; 52c62f7daeSJames Wright 53c62f7daeSJames Wright PetscFunctionBeginUser; 54c62f7daeSJames Wright const char *A_initial = StateVariables_Initial[state_var_A]; 55c62f7daeSJames Wright const char *B_initial = StateVariables_Initial[state_var_B]; 56c62f7daeSJames Wright 57c62f7daeSJames Wright State state_A0 = StateFromQ(gas, A0, state_var_A); 58c62f7daeSJames Wright StateToQ(gas, state_A0, B0, state_var_B); 59c62f7daeSJames Wright State state_B0 = StateFromQ(gas, B0, state_var_B); 60c62f7daeSJames Wright StateToQ(gas, state_B0, A0_test, state_var_A); 61c62f7daeSJames Wright 62c62f7daeSJames Wright snprintf(buf, sizeof buf, "%s->%s->%s: %s", A_initial, B_initial, A_initial, A_initial); 63c62f7daeSJames Wright PetscCall(CheckQWithTolerance(A0, A0_test, A0, buf, rtol_0, rtol_u, rtol_4)); 64c62f7daeSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 65c62f7daeSJames Wright } 66c62f7daeSJames Wright 67c62f7daeSJames Wright // @brief Verify `StateFromQ_fwd` via a finite difference approximation 68c62f7daeSJames Wright static PetscErrorCode TestState_fwd(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5], 69c62f7daeSJames Wright CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) { 70c62f7daeSJames Wright CeedScalar eps = 4e-7; // Finite difference step 71c62f7daeSJames Wright char buf[128]; 72c62f7daeSJames Wright const char *const StateVariables_Initial[] = {"U", "Y", "V"}; 73c62f7daeSJames Wright 74c62f7daeSJames Wright PetscFunctionBeginUser; 75c62f7daeSJames Wright const char *A_initial = StateVariables_Initial[state_var_A]; 76c62f7daeSJames Wright const char *B_initial = StateVariables_Initial[state_var_B]; 77c62f7daeSJames Wright State state_0 = StateFromQ(gas, A0, state_var_A); 78c62f7daeSJames Wright 79c62f7daeSJames Wright for (int i = 0; i < 5; i++) { 80c62f7daeSJames Wright CeedScalar dB[5] = {0.}, dB_fd[5] = {0.}; 81c62f7daeSJames Wright { // Calculate dB using State functions 82c62f7daeSJames Wright CeedScalar dA[5] = {0}; 83c62f7daeSJames Wright 84c62f7daeSJames Wright dA[i] = A0[i]; 85a6654c3eSJames Wright State dstate_0 = StateFromQ_fwd(gas, state_0, dA, state_var_A); 86a6654c3eSJames Wright StateToQ_fwd(gas, state_0, dstate_0, dB, state_var_B); 87c62f7daeSJames Wright } 88c62f7daeSJames Wright 89c62f7daeSJames Wright { // Calculate dB_fd via finite difference approximation 90c62f7daeSJames Wright CeedScalar A1[5], B0[5], B1[5]; 91c62f7daeSJames Wright 92c62f7daeSJames Wright for (int j = 0; j < 5; j++) A1[j] = (1 + eps * (i == j)) * A0[j]; 93c62f7daeSJames Wright State state_1 = StateFromQ(gas, A1, state_var_A); 94c62f7daeSJames Wright StateToQ(gas, state_0, B0, state_var_B); 95c62f7daeSJames Wright StateToQ(gas, state_1, B1, state_var_B); 96c62f7daeSJames Wright for (int j = 0; j < 5; j++) dB_fd[j] = (B1[j] - B0[j]) / eps; 97c62f7daeSJames Wright } 98c62f7daeSJames Wright 99c62f7daeSJames 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); 100c62f7daeSJames Wright PetscCall(CheckQWithTolerance(dB_fd, dB, dB_fd, buf, rtol_0, rtol_u, rtol_4)); 101c62f7daeSJames Wright } 102c62f7daeSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 103c62f7daeSJames Wright } 104c62f7daeSJames Wright 105c62f7daeSJames Wright // @brief Test the Newtonian State transformation functions, `StateFrom*` 1060c373b74SJames Wright static PetscErrorCode UnitTests_Newtonian(Honee honee, NewtonianIdealGasContext gas) { 1070c373b74SJames Wright Units units = honee->units; 108c62f7daeSJames Wright const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, K = units->Kelvin; 109c62f7daeSJames Wright 110f0b65372SJed Brown PetscFunctionBeginUser; 111c62f7daeSJames Wright const CeedScalar T = 200 * K; 112c62f7daeSJames Wright const CeedScalar rho = 1.2 * kg / Cube(m); 113c62f7daeSJames Wright const CeedScalar P = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T; 114c62f7daeSJames Wright const CeedScalar u_base = 40 * m / sec; 115c62f7daeSJames Wright const CeedScalar u[3] = {u_base, u_base * 1.1, u_base * 1.2}; 116c62f7daeSJames Wright const CeedScalar e_kinetic = 0.5 * Dot3(u, u); 117c62f7daeSJames Wright const CeedScalar e_internal = gas->cv * T; 118c62f7daeSJames Wright const CeedScalar e_total = e_kinetic + e_internal; 119c62f7daeSJames Wright const CeedScalar gamma = HeatCapacityRatio(gas); 120c62f7daeSJames Wright const CeedScalar entropy = log(P) - gamma * log(rho); 121c62f7daeSJames Wright const CeedScalar rho_div_p = rho / P; 122c62f7daeSJames Wright const CeedScalar Y0[5] = {P, u[0], u[1], u[2], T}; 123c62f7daeSJames Wright const CeedScalar U0[5] = {rho, rho * u[0], rho * u[1], rho * u[2], rho * e_total}; 124c62f7daeSJames 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], 125c62f7daeSJames Wright -rho_div_p}; 126c62f7daeSJames Wright 127c62f7daeSJames Wright { 128c62f7daeSJames Wright CeedScalar rtol = 20 * CEED_EPSILON; 129c62f7daeSJames Wright 130c62f7daeSJames Wright PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol)); 131c62f7daeSJames Wright PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol)); 132c62f7daeSJames Wright PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol)); 133c62f7daeSJames Wright PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, rtol, rtol, rtol)); 134c62f7daeSJames Wright PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, rtol, rtol, rtol)); 135c62f7daeSJames Wright PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, rtol, rtol, rtol)); 136c62f7daeSJames Wright } 137c62f7daeSJames Wright 138c62f7daeSJames Wright { 139c62f7daeSJames Wright CeedScalar rtol = 5e-6; 140c62f7daeSJames Wright 141c62f7daeSJames Wright PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol)); 142c62f7daeSJames Wright PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol)); 143c62f7daeSJames Wright PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol)); 144c62f7daeSJames Wright PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, 10 * rtol, rtol, rtol)); 145c62f7daeSJames Wright PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, 5 * rtol, rtol, rtol)); 146c62f7daeSJames Wright PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, 5 * rtol, 5 * rtol, 5 * rtol)); 147f0b65372SJed Brown } 148d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 149f0b65372SJed Brown } 150f0b65372SJed Brown 15165dee3d2SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping 15265dee3d2SJames Wright // 15365dee3d2SJames Wright // Only used for SUPG stabilization 1540c373b74SJames Wright PetscErrorCode CreateKSPMassOperator_NewtonianStabilized(Honee honee, CeedOperator *op_mass) { 1550c373b74SJames Wright Ceed ceed = honee->ceed; 15665dee3d2SJames Wright CeedInt num_comp_q, q_data_size; 15765dee3d2SJames Wright CeedQFunction qf_mass; 158*01e19bfaSJames Wright CeedElemRestriction elem_restr_q, elem_restr_qd; 15965dee3d2SJames Wright CeedBasis basis_q; 16065dee3d2SJames Wright CeedVector q_data; 161236a8ba6SJames Wright CeedQFunctionContext qfctx = NULL; 16265dee3d2SJames Wright PetscInt dim = 3; 16365dee3d2SJames Wright 16465dee3d2SJames Wright PetscFunctionBeginUser; 16565dee3d2SJames Wright { // Get restriction and basis from the RHS function 16665dee3d2SJames Wright CeedOperator *sub_ops; 167*01e19bfaSJames Wright CeedOperatorField op_field; 16865dee3d2SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 16965dee3d2SJames Wright 1700c373b74SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 171*01e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field)); 172*01e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL)); 173*01e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &op_field)); 174*01e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_qd, NULL, &q_data)); 17565dee3d2SJames Wright 176236a8ba6SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qfctx)); 17765dee3d2SJames Wright } 17865dee3d2SJames Wright 17965dee3d2SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 180*01e19bfaSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd, &q_data_size)); 18165dee3d2SJames Wright 18265dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Newtonian_Conserv, MassFunction_Newtonian_Conserv_loc, &qf_mass)); 18365dee3d2SJames Wright 184236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qfctx)); 18565dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0)); 18665dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP)); 18765dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP)); 18865dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE)); 18965dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP)); 19065dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD)); 19165dee3d2SJames Wright 19265dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 19365dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 1940c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, honee->q_ceed)); 195*01e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 19665dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 19765dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 19865dee3d2SJames Wright 199fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 200*01e19bfaSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 201fff85bd3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 202fff85bd3SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 203236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qfctx)); 20465dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 20565dee3d2SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 20665dee3d2SJames Wright } 2074c07ec22SJames Wright 20836038bbcSJames Wright /** 20936038bbcSJames Wright @brief Create RHS CeedOperator for direct projection of divergence of diffusive flux 21036038bbcSJames Wright 2110c373b74SJames Wright @param[in] honee `Honee` context 21236038bbcSJames Wright @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 21336038bbcSJames Wright @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 21436038bbcSJames Wright **/ 215e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_NS(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) { 2160c373b74SJames Wright Ceed ceed = honee->ceed; 21736038bbcSJames Wright NodalProjectionData projection = diff_flux_proj->projection; 21836038bbcSJames Wright CeedInt num_comp_q; 21936038bbcSJames Wright PetscInt dim, label_value = 0; 22036038bbcSJames Wright DMLabel domain_label = NULL; 221236a8ba6SJames Wright CeedQFunctionContext newtonian_qfctx = NULL; 22236038bbcSJames Wright 22336038bbcSJames Wright PetscFunctionBeginUser; 22436038bbcSJames Wright // -- Get Pre-requisite things 22536038bbcSJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 226e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q)); 22736038bbcSJames Wright 2280880fbb6SJames Wright { // Get newtonian QF context 22936038bbcSJames Wright CeedOperator *sub_ops; 23036038bbcSJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 23136038bbcSJames Wright 2320c373b74SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_ifunction, &sub_ops)); 233236a8ba6SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newtonian_qfctx)); 23436038bbcSJames Wright } 23536038bbcSJames Wright PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_rhs)); 23636038bbcSJames Wright { // Add the volume integral CeedOperator 23736038bbcSJames Wright CeedQFunction qf_rhs_volume; 23836038bbcSJames Wright CeedOperator op_rhs_volume; 23936038bbcSJames Wright CeedVector q_data; 2400880fbb6SJames Wright CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL; 2410880fbb6SJames Wright CeedBasis basis_diff_flux = NULL; 24236038bbcSJames Wright CeedInt q_data_size; 24336038bbcSJames Wright 2440880fbb6SJames Wright PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL)); 245e3663b90SJames 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, 246e3663b90SJames Wright &q_data_size)); 2470c373b74SJames Wright switch (honee->phys->state_var) { 24836038bbcSJames Wright case STATEVAR_PRIMITIVE: 24936038bbcSJames Wright PetscCallCeed(ceed, 25036038bbcSJames Wright CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Prim, DivDiffusiveFluxVolumeRHS_NS_Prim_loc, &qf_rhs_volume)); 25136038bbcSJames Wright break; 25236038bbcSJames Wright case STATEVAR_CONSERVATIVE: 25336038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Conserv, DivDiffusiveFluxVolumeRHS_NS_Conserv_loc, 25436038bbcSJames Wright &qf_rhs_volume)); 25536038bbcSJames Wright break; 25636038bbcSJames Wright case STATEVAR_ENTROPY: 25736038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_NS_Entropy, DivDiffusiveFluxVolumeRHS_NS_Entropy_loc, 25836038bbcSJames Wright &qf_rhs_volume)); 25936038bbcSJames Wright break; 26036038bbcSJames Wright } 26136038bbcSJames Wright 262236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, newtonian_qfctx)); 26336038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "q", num_comp_q, CEED_EVAL_INTERP)); 26436038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 26536038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE)); 26636038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD)); 26736038bbcSJames Wright 26836038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume)); 269e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 270e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 27136038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 27236038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE)); 27336038bbcSJames Wright 27436038bbcSJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_volume)); 27536038bbcSJames Wright 27636038bbcSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 27736038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 2780880fbb6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume)); 2790880fbb6SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 28036038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume)); 28136038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume)); 28236038bbcSJames Wright } 28336038bbcSJames Wright 28436038bbcSJames Wright { // Add the boundary integral CeedOperator 28536038bbcSJames Wright CeedQFunction qf_rhs_boundary; 28636038bbcSJames Wright DMLabel face_sets_label; 28736038bbcSJames Wright PetscInt num_face_set_values, *face_set_values; 28836038bbcSJames Wright CeedInt q_data_size; 28936038bbcSJames Wright 29036038bbcSJames Wright // -- Build RHS operator 2910c373b74SJames Wright switch (honee->phys->state_var) { 29236038bbcSJames Wright case STATEVAR_PRIMITIVE: 29336038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Prim, DivDiffusiveFluxBoundaryRHS_NS_Prim_loc, 29436038bbcSJames Wright &qf_rhs_boundary)); 29536038bbcSJames Wright break; 29636038bbcSJames Wright case STATEVAR_CONSERVATIVE: 29736038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Conserv, DivDiffusiveFluxBoundaryRHS_NS_Conserv_loc, 29836038bbcSJames Wright &qf_rhs_boundary)); 29936038bbcSJames Wright break; 30036038bbcSJames Wright case STATEVAR_ENTROPY: 30136038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_NS_Entropy, DivDiffusiveFluxBoundaryRHS_NS_Entropy_loc, 30236038bbcSJames Wright &qf_rhs_boundary)); 30336038bbcSJames Wright break; 30436038bbcSJames Wright } 30536038bbcSJames Wright 3060c373b74SJames Wright PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size)); 307236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, newtonian_qfctx)); 30836038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "q", num_comp_q, CEED_EVAL_INTERP)); 30936038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 31036038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE)); 31136038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP)); 31236038bbcSJames Wright 31336038bbcSJames Wright PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label)); 31436038bbcSJames Wright PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values)); 31536038bbcSJames Wright for (PetscInt f = 0; f < num_face_set_values; f++) { 31636038bbcSJames Wright DMLabel face_orientation_label; 31736038bbcSJames Wright PetscInt num_orientations_values, *orientation_values; 31836038bbcSJames Wright 31936038bbcSJames Wright { 32036038bbcSJames Wright char *face_orientation_label_name; 32136038bbcSJames Wright 32236038bbcSJames Wright PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name)); 32336038bbcSJames Wright PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label)); 32436038bbcSJames Wright PetscCall(PetscFree(face_orientation_label_name)); 32536038bbcSJames Wright } 32636038bbcSJames Wright PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values)); 32736038bbcSJames Wright for (PetscInt o = 0; o < num_orientations_values; o++) { 32836038bbcSJames Wright CeedOperator op_rhs_boundary; 32936038bbcSJames Wright CeedBasis basis_q, basis_diff_flux_boundary; 33036038bbcSJames Wright CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary; 33136038bbcSJames Wright CeedVector q_data; 33236038bbcSJames Wright CeedInt q_data_size; 33336038bbcSJames Wright PetscInt orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1; 33436038bbcSJames Wright 3350c373b74SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q)); 3360c373b74SJames Wright PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q)); 33736038bbcSJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, 33836038bbcSJames Wright &elem_restr_diff_flux_boundary)); 3390880fbb6SJames Wright PetscCall(CreateBasisFromPlex(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary)); 340e3663b90SJames Wright PetscCall( 341e3663b90SJames Wright QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, honee->x_coord, &elem_restr_qdata, &q_data, &q_data_size)); 34236038bbcSJames Wright 34336038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary)); 34436038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 34536038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 34636038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data)); 34736038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary, 34836038bbcSJames Wright CEED_VECTOR_ACTIVE)); 34936038bbcSJames Wright 35036038bbcSJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_boundary)); 35136038bbcSJames Wright 35236038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary)); 35336038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata)); 35436038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 35536038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary)); 35636038bbcSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 35736038bbcSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary)); 35836038bbcSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 35936038bbcSJames Wright } 36036038bbcSJames Wright PetscCall(PetscFree(orientation_values)); 36136038bbcSJames Wright } 36236038bbcSJames Wright PetscCall(PetscFree(face_set_values)); 36336038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary)); 36436038bbcSJames Wright } 36536038bbcSJames Wright 366236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newtonian_qfctx)); 36736038bbcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 36836038bbcSJames Wright } 36936038bbcSJames Wright 37036038bbcSJames Wright /** 37136038bbcSJames Wright @brief Create RHS CeedOperator for indirect projection of divergence of diffusive flux 37236038bbcSJames Wright 3730c373b74SJames Wright @param[in] honee `Honee` context 37436038bbcSJames Wright @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 37536038bbcSJames Wright @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 37636038bbcSJames Wright **/ 377e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Indirect_NS(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) { 3780c373b74SJames Wright Ceed ceed = honee->ceed; 37936038bbcSJames Wright NodalProjectionData projection = diff_flux_proj->projection; 38036038bbcSJames Wright CeedBasis basis_diff_flux; 38136038bbcSJames Wright CeedElemRestriction elem_restr_diff_flux, elem_restr_qd; 38236038bbcSJames Wright CeedVector q_data; 38336038bbcSJames Wright CeedInt num_comp_q, q_data_size; 38436038bbcSJames Wright PetscInt dim; 38536038bbcSJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 38636038bbcSJames Wright DMLabel domain_label = NULL; 38736038bbcSJames Wright CeedQFunction qf_rhs; 388236a8ba6SJames Wright CeedQFunctionContext newtonian_qfctx = NULL; 38936038bbcSJames Wright 39036038bbcSJames Wright PetscFunctionBeginUser; 39136038bbcSJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 392e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q)); 39336038bbcSJames Wright 39440b78511SJames Wright { // Get newtonian QF context 39536038bbcSJames Wright CeedOperator *sub_ops; 39636038bbcSJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 39736038bbcSJames Wright 3980c373b74SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_ifunction, &sub_ops)); 399236a8ba6SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &newtonian_qfctx)); 40036038bbcSJames Wright } 40136038bbcSJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux)); 40236038bbcSJames Wright PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux)); 403e3663b90SJames 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, 404e3663b90SJames Wright &q_data_size)); 40536038bbcSJames Wright 4060c373b74SJames Wright switch (honee->phys->state_var) { 40736038bbcSJames Wright case STATEVAR_PRIMITIVE: 40836038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Prim, DiffusiveFluxRHS_NS_Prim_loc, &qf_rhs)); 40936038bbcSJames Wright break; 41036038bbcSJames Wright case STATEVAR_CONSERVATIVE: 41136038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Conserv, DiffusiveFluxRHS_NS_Conserv_loc, &qf_rhs)); 41236038bbcSJames Wright break; 41336038bbcSJames Wright case STATEVAR_ENTROPY: 41436038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_NS_Entropy, DiffusiveFluxRHS_NS_Entropy_loc, &qf_rhs)); 41536038bbcSJames Wright break; 41636038bbcSJames Wright } 41736038bbcSJames Wright 418236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, newtonian_qfctx)); 41936038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "q", num_comp_q, CEED_EVAL_INTERP)); 42036038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 42136038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE)); 42236038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP)); 42336038bbcSJames Wright 42436038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, op_rhs)); 425e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 426e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 42736038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 42836038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 42936038bbcSJames Wright 43036038bbcSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs)); 431236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&newtonian_qfctx)); 43236038bbcSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 43336038bbcSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 43436038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 43536038bbcSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 43636038bbcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 43736038bbcSJames Wright } 43836038bbcSJames Wright 439991aef52SJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 440b7f03f12SJed Brown SetupContext setup_context; 4410c373b74SJames Wright Honee honee = *(Honee *)ctx; 4420c373b74SJames Wright CeedInt degree = honee->app_ctx->degree; 4433a8779fbSJames Wright StabilizationType stab; 4443636f6a4SJames Wright StateVariable state_var; 4450c373b74SJames Wright MPI_Comm comm = honee->comm; 4460c373b74SJames Wright Ceed ceed = honee->ceed; 4473a8779fbSJames Wright PetscBool implicit; 4480d8cd818SJames Wright PetscBool unit_tests; 44915a3537eSJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 450e07531f7SJames Wright CeedQFunctionContext newtonian_ig_qfctx; 4513a8779fbSJames Wright 45215a3537eSJed Brown PetscFunctionBeginUser; 4532b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &setup_context)); 4542b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &newtonian_ig_ctx)); 4553a8779fbSJames Wright 4563a8779fbSJames Wright // ------------------------------------------------------ 4573a8779fbSJames Wright // Setup Generic Newtonian IG Problem 4583a8779fbSJames Wright // ------------------------------------------------------ 45928160fc2SJames Wright problem->jac_data_size_vol = 14; 460b01ba163SJames Wright problem->jac_data_size_sur = 11; 46158ce1233SJames Wright problem->compute_exact_solution_error = PETSC_FALSE; 462cbe60e31SLeila Ghaffari problem->print_info = PRINT_NEWTONIAN; 4633a8779fbSJames Wright 4640c373b74SJames Wright PetscCall(DivDiffFluxProjectionCreate(honee, 4, &honee->diff_flux_proj)); 4650c373b74SJames Wright if (honee->diff_flux_proj) { 4660c373b74SJames Wright DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj; 46736038bbcSJames Wright NodalProjectionData projection = diff_flux_proj->projection; 46836038bbcSJames Wright 46936038bbcSJames Wright diff_flux_proj->CreateRHSOperator_Direct = DivDiffFluxProjectionCreateRHS_Direct_NS; 47036038bbcSJames Wright diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_NS; 47136038bbcSJames Wright 4720c373b74SJames Wright switch (honee->diff_flux_proj->method) { 47336038bbcSJames Wright case DIV_DIFF_FLUX_PROJ_DIRECT: { 47436038bbcSJames Wright PetscSection section; 47536038bbcSJames Wright 47636038bbcSJames Wright PetscCall(DMGetLocalSection(projection->dm, §ion)); 47736038bbcSJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 47836038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_MomentumX")); 47936038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "DivDiffusiveFlux_MomentumY")); 48036038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "DivDiffusiveFlux_MomentumZ")); 48136038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "DivDiffusiveFlux_Energy")); 48236038bbcSJames Wright } break; 48336038bbcSJames Wright case DIV_DIFF_FLUX_PROJ_INDIRECT: { 48436038bbcSJames Wright PetscSection section; 48536038bbcSJames Wright 48636038bbcSJames Wright PetscCall(DMGetLocalSection(projection->dm, §ion)); 48736038bbcSJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 48836038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_MomentumXX")); 48936038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_MomentumXY")); 49036038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_MomentumXZ")); 49136038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "DiffusiveFlux_MomentumYX")); 49236038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "DiffusiveFlux_MomentumYY")); 49336038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 5, "DiffusiveFlux_MomentumYZ")); 49436038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 6, "DiffusiveFlux_MomentumZX")); 49536038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 7, "DiffusiveFlux_MomentumZY")); 49636038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 8, "DiffusiveFlux_MomentumZZ")); 49736038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 9, "DiffusiveFlux_EnergyX")); 49836038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 10, "DiffusiveFlux_EnergyY")); 49936038bbcSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 11, "DiffusiveFlux_EnergyZ")); 50036038bbcSJames Wright } break; 50136038bbcSJames Wright case DIV_DIFF_FLUX_PROJ_NONE: 5020c373b74SJames Wright SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 5030c373b74SJames Wright DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]); 50436038bbcSJames Wright break; 50536038bbcSJames Wright } 50636038bbcSJames Wright } 50736038bbcSJames Wright 5083a8779fbSJames Wright // ------------------------------------------------------ 5093a8779fbSJames Wright // Create the libCEED context 5103a8779fbSJames Wright // ------------------------------------------------------ 5113a8779fbSJames Wright CeedScalar cv = 717.; // J/(kg K) 5123a8779fbSJames Wright CeedScalar cp = 1004.; // J/(kg K) 513d9bb1cdbSJames Wright CeedScalar g[3] = {0, 0, 0}; // m/s^2 5143a8779fbSJames Wright CeedScalar lambda = -2. / 3.; // - 515bb8a0c61SJames Wright CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 5163a8779fbSJames Wright CeedScalar k = 0.02638; // W/(m K) 5174351d27fSLeila Ghaffari CeedScalar c_tau = 0.5 / degree; // - 518bb8a0c61SJames Wright CeedScalar Ctau_t = 1.0; // - 519b5786772SLeila Ghaffari CeedScalar Cv_func[3] = {36, 60, 128}; 520c176988bSLeila Ghaffari CeedScalar Ctau_v = Cv_func[(CeedInt)Min(3, degree) - 1]; 5212a0eee6eSLeila Ghaffari CeedScalar Ctau_C = 0.25 / degree; 5222a0eee6eSLeila Ghaffari CeedScalar Ctau_M = 0.25 / degree; 5232a0eee6eSLeila Ghaffari CeedScalar Ctau_E = 0.125; 5243a8779fbSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 5252b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 526493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 5273a8779fbSJames Wright 528b8fb7609SAdeleke O. Bankole StatePrimitive reference = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15}; 529a213b8aaSJames Wright CeedScalar idl_decay_time = -1, idl_start = 0, idl_length = 0, idl_pressure = reference.pressure; 530e7754af5SKenneth E. Jansen PetscBool idl_enable = PETSC_FALSE; 531b8fb7609SAdeleke O. Bankole 5323a8779fbSJames Wright // ------------------------------------------------------ 5333a8779fbSJames Wright // Create the PETSc context 5343a8779fbSJames Wright // ------------------------------------------------------ 535bb8a0c61SJames Wright PetscScalar meter = 1; // 1 meter in scaled length units 536bb8a0c61SJames Wright PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 537bb8a0c61SJames Wright PetscScalar second = 1; // 1 second in scaled time units 5383a8779fbSJames Wright PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 5393a8779fbSJames Wright PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 5403a8779fbSJames Wright 5413a8779fbSJames Wright // ------------------------------------------------------ 5423a8779fbSJames Wright // Command line Options 5433a8779fbSJames Wright // ------------------------------------------------------ 544d9bb1cdbSJames Wright PetscBool given_option = PETSC_FALSE; 5452b916ea7SJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL); 546cbe60e31SLeila Ghaffari // -- Conservative vs Primitive variables 5472b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE), 5482b916ea7SJeremy L Thompson (PetscEnum *)&state_var, NULL)); 5493636f6a4SJames Wright 5503636f6a4SJames Wright switch (state_var) { 5513636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 552e07531f7SJames Wright problem->ics.qf_func_ptr = ICsNewtonianIG_Conserv; 553e07531f7SJames Wright problem->ics.qf_loc = ICsNewtonianIG_Conserv_loc; 554e07531f7SJames Wright problem->apply_vol_rhs.qf_func_ptr = RHSFunction_Newtonian; 555e07531f7SJames Wright problem->apply_vol_rhs.qf_loc = RHSFunction_Newtonian_loc; 556e07531f7SJames Wright problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Conserv; 557e07531f7SJames Wright problem->apply_vol_ifunction.qf_loc = IFunction_Newtonian_Conserv_loc; 558e07531f7SJames Wright problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Conserv; 559e07531f7SJames Wright problem->apply_vol_ijacobian.qf_loc = IJacobian_Newtonian_Conserv_loc; 560e07531f7SJames Wright problem->apply_inflow.qf_func_ptr = BoundaryIntegral_Conserv; 561e07531f7SJames Wright problem->apply_inflow.qf_loc = BoundaryIntegral_Conserv_loc; 562e07531f7SJames Wright problem->apply_inflow_jacobian.qf_func_ptr = BoundaryIntegral_Jacobian_Conserv; 563e07531f7SJames Wright problem->apply_inflow_jacobian.qf_loc = BoundaryIntegral_Jacobian_Conserv_loc; 5643636f6a4SJames Wright break; 5653636f6a4SJames Wright case STATEVAR_PRIMITIVE: 566e07531f7SJames Wright problem->ics.qf_func_ptr = ICsNewtonianIG_Prim; 567e07531f7SJames Wright problem->ics.qf_loc = ICsNewtonianIG_Prim_loc; 568e07531f7SJames Wright problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Prim; 569e07531f7SJames Wright problem->apply_vol_ifunction.qf_loc = IFunction_Newtonian_Prim_loc; 570e07531f7SJames Wright problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Prim; 571e07531f7SJames Wright problem->apply_vol_ijacobian.qf_loc = IJacobian_Newtonian_Prim_loc; 572e07531f7SJames Wright problem->apply_inflow.qf_func_ptr = BoundaryIntegral_Prim; 573e07531f7SJames Wright problem->apply_inflow.qf_loc = BoundaryIntegral_Prim_loc; 574e07531f7SJames Wright problem->apply_inflow_jacobian.qf_func_ptr = BoundaryIntegral_Jacobian_Prim; 575e07531f7SJames Wright problem->apply_inflow_jacobian.qf_loc = BoundaryIntegral_Jacobian_Prim_loc; 5763636f6a4SJames Wright break; 577c62f7daeSJames Wright case STATEVAR_ENTROPY: 578e07531f7SJames Wright problem->ics.qf_func_ptr = ICsNewtonianIG_Entropy; 579e07531f7SJames Wright problem->ics.qf_loc = ICsNewtonianIG_Entropy_loc; 580e07531f7SJames Wright problem->apply_vol_ifunction.qf_func_ptr = IFunction_Newtonian_Entropy; 581e07531f7SJames Wright problem->apply_vol_ifunction.qf_loc = IFunction_Newtonian_Entropy_loc; 582e07531f7SJames Wright problem->apply_vol_ijacobian.qf_func_ptr = IJacobian_Newtonian_Entropy; 583e07531f7SJames Wright problem->apply_vol_ijacobian.qf_loc = IJacobian_Newtonian_Entropy_loc; 584e07531f7SJames Wright problem->apply_inflow.qf_func_ptr = BoundaryIntegral_Entropy; 585e07531f7SJames Wright problem->apply_inflow.qf_loc = BoundaryIntegral_Entropy_loc; 586e07531f7SJames Wright problem->apply_inflow_jacobian.qf_func_ptr = BoundaryIntegral_Jacobian_Entropy; 587e07531f7SJames Wright problem->apply_inflow_jacobian.qf_loc = BoundaryIntegral_Jacobian_Entropy_loc; 588c62f7daeSJames Wright break; 589cbe60e31SLeila Ghaffari } 5901485969bSJeremy L Thompson 5913a8779fbSJames Wright // -- Physics 5922b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL)); 5932b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL)); 5942b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL)); 5952b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL)); 5962b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL)); 5973a8779fbSJames Wright 59866d54740SJames Wright PetscInt dim = 3; 599d9bb1cdbSJames Wright PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL)); 600d9bb1cdbSJames Wright PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option)); 601d9bb1cdbSJames Wright if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim); 602d9bb1cdbSJames Wright 6032b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 6042b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 6052b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 6062b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL)); 6072b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL)); 6082b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL)); 6092b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL)); 6102b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 6112b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL)); 6123a8779fbSJames Wright 613db95602eSJames Wright dim = 3; 614b8fb7609SAdeleke O. Bankole PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL)); 615b8fb7609SAdeleke O. Bankole PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL)); 616b8fb7609SAdeleke O. Bankole PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL)); 617b8fb7609SAdeleke O. Bankole 6183a8779fbSJames Wright // -- Units 6192b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 6203a8779fbSJames Wright meter = fabs(meter); 6212b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 6223a8779fbSJames Wright kilogram = fabs(kilogram); 6232b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 6243a8779fbSJames Wright second = fabs(second); 6252b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL)); 6263a8779fbSJames Wright Kelvin = fabs(Kelvin); 6273a8779fbSJames Wright 6283a8779fbSJames Wright // -- Warnings 6295d28dccaSJames Wright PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP, 6305d28dccaSJames Wright "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n"); 6310c373b74SJames Wright PetscCheck(!(honee->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE && !implicit), comm, PETSC_ERR_SUP, 6325f952e8dSJames Wright "Projection of divergence of diffusive flux is not implemented for Navier-Stokes with explicit timestepping"); 633e7754af5SKenneth E. Jansen 634e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point", 635e7754af5SKenneth E. Jansen NULL, idl_decay_time, &idl_decay_time, &idl_enable)); 6369cbdf780SJames Wright PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero."); 6379cbdf780SJames Wright if (idl_decay_time < 0) idl_enable = PETSC_FALSE; 63828160fc2SJames Wright if (idl_enable) problem->jac_data_size_vol++; 639e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL)); 640e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL)); 641a213b8aaSJames Wright idl_pressure = reference.pressure; 642a213b8aaSJames Wright PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-reference_pressure`)", NULL, idl_pressure, 643a213b8aaSJames Wright &idl_pressure, NULL)); 6441485969bSJeremy L Thompson PetscOptionsEnd(); 6453a8779fbSJames Wright 64665dee3d2SJames Wright if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized; 64765dee3d2SJames Wright 6483a8779fbSJames Wright // ------------------------------------------------------ 6493a8779fbSJames Wright // Set up the PETSc context 6503a8779fbSJames Wright // ------------------------------------------------------ 6513a8779fbSJames Wright // -- Define derived units 6523a8779fbSJames Wright Pascal = kilogram / (meter * PetscSqr(second)); 6533a8779fbSJames Wright J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 6543a8779fbSJames Wright m_per_squared_s = meter / PetscSqr(second); 6553a8779fbSJames Wright W_per_m_K = kilogram * meter / (pow(second, 3) * Kelvin); 6563a8779fbSJames Wright 6570c373b74SJames Wright honee->units->meter = meter; 6580c373b74SJames Wright honee->units->kilogram = kilogram; 6590c373b74SJames Wright honee->units->second = second; 6600c373b74SJames Wright honee->units->Kelvin = Kelvin; 6610c373b74SJames Wright honee->units->Pascal = Pascal; 6620c373b74SJames Wright honee->units->J_per_kg_K = J_per_kg_K; 6630c373b74SJames Wright honee->units->m_per_squared_s = m_per_squared_s; 6640c373b74SJames Wright honee->units->W_per_m_K = W_per_m_K; 6653a8779fbSJames Wright 6663a8779fbSJames Wright // ------------------------------------------------------ 6673a8779fbSJames Wright // Set up the libCEED context 6683a8779fbSJames Wright // ------------------------------------------------------ 6693a8779fbSJames Wright // -- Scale variables to desired units 6703a8779fbSJames Wright cv *= J_per_kg_K; 6713a8779fbSJames Wright cp *= J_per_kg_K; 6723a8779fbSJames Wright mu *= Pascal * second; 6733a8779fbSJames Wright k *= W_per_m_K; 674493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter; 675493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s; 676b8fb7609SAdeleke O. Bankole reference.pressure *= Pascal; 677b8fb7609SAdeleke O. Bankole for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second; 678b8fb7609SAdeleke O. Bankole reference.temperature *= Kelvin; 6793a8779fbSJames Wright 6803a8779fbSJames Wright // -- Solver Settings 6810c373b74SJames Wright honee->phys->implicit = implicit; 6820c373b74SJames Wright honee->phys->state_var = state_var; 6833a8779fbSJames Wright 6843a8779fbSJames Wright // -- QFunction Context 68515a3537eSJed Brown newtonian_ig_ctx->lambda = lambda; 68615a3537eSJed Brown newtonian_ig_ctx->mu = mu; 68715a3537eSJed Brown newtonian_ig_ctx->k = k; 68815a3537eSJed Brown newtonian_ig_ctx->cv = cv; 68915a3537eSJed Brown newtonian_ig_ctx->cp = cp; 69015a3537eSJed Brown newtonian_ig_ctx->c_tau = c_tau; 69115a3537eSJed Brown newtonian_ig_ctx->Ctau_t = Ctau_t; 69215a3537eSJed Brown newtonian_ig_ctx->Ctau_v = Ctau_v; 69315a3537eSJed Brown newtonian_ig_ctx->Ctau_C = Ctau_C; 69415a3537eSJed Brown newtonian_ig_ctx->Ctau_M = Ctau_M; 69515a3537eSJed Brown newtonian_ig_ctx->Ctau_E = Ctau_E; 69615a3537eSJed Brown newtonian_ig_ctx->stabilization = stab; 6978085925cSJames Wright newtonian_ig_ctx->is_implicit = implicit; 6983636f6a4SJames Wright newtonian_ig_ctx->state_var = state_var; 699e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_enable = idl_enable; 700e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * second); 701e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_start = idl_start * meter; 702e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_length = idl_length * meter; 703a213b8aaSJames Wright newtonian_ig_ctx->idl_pressure = idl_pressure; 7040c373b74SJames Wright newtonian_ig_ctx->divFdiff_method = honee->app_ctx->divFdiffproj_method; 7052b916ea7SJeremy L Thompson PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3)); 7063a8779fbSJames Wright 707b8fb7609SAdeleke O. Bankole // -- Setup Context 708b8fb7609SAdeleke O. Bankole setup_context->reference = reference; 709b8fb7609SAdeleke O. Bankole setup_context->gas = *newtonian_ig_ctx; 710b8fb7609SAdeleke O. Bankole setup_context->lx = domain_size[0]; 711b8fb7609SAdeleke O. Bankole setup_context->ly = domain_size[1]; 712b8fb7609SAdeleke O. Bankole setup_context->lz = domain_size[2]; 713b8fb7609SAdeleke O. Bankole setup_context->time = 0; 714b8fb7609SAdeleke O. Bankole 7150c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &problem->ics.qfctx)); 716e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 717e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc)); 718e07531f7SJames Wright PetscCallCeed( 719e07531f7SJames Wright ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfctx, "evaluation time", offsetof(struct SetupContext_, time), 1, "Time of evaluation")); 72015a3537eSJed Brown 7210c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &newtonian_ig_qfctx)); 722e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(newtonian_ig_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx)); 723e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 724e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1, 725b4c37c5cSJames Wright "Size of timestep, delta t")); 726e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "ijacobian time shift", 727b4c37c5cSJames Wright offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1, 728b4c37c5cSJames Wright "Shift for mass matrix in IJacobian")); 729e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_qfctx, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1, 730b4c37c5cSJames Wright "Current solution time")); 731a8bca8acSJames Wright 732e07531f7SJames Wright problem->apply_vol_rhs.qfctx = newtonian_ig_qfctx; 733e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ifunction.qfctx)); 734e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_vol_ijacobian.qfctx)); 735e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_inflow.qfctx)); 736e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_qfctx, &problem->apply_inflow_jacobian.qfctx)); 737f0b65372SJed Brown 7389ed3d70dSJames Wright if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 7399ed3d70dSJames Wright if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 740e07531f7SJames Wright if (bc->num_slip > 0) PetscCall(SlipBCSetup(problem, dm, ctx, newtonian_ig_qfctx)); 7419ed3d70dSJames Wright 742f0b65372SJed Brown if (unit_tests) { 7430c373b74SJames Wright PetscCall(UnitTests_Newtonian(honee, newtonian_ig_ctx)); 744f0b65372SJed Brown } 745d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 7463a8779fbSJames Wright } 7473a8779fbSJames Wright 7480c373b74SJames Wright PetscErrorCode PRINT_NEWTONIAN(Honee honee, ProblemData problem, AppCtx app_ctx) { 7490c373b74SJames Wright MPI_Comm comm = honee->comm; 7500c373b74SJames Wright Ceed ceed = honee->ceed; 75115a3537eSJed Brown NewtonianIdealGasContext newtonian_ctx; 75215a3537eSJed Brown 7533a8779fbSJames Wright PetscFunctionBeginUser; 754e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ctx)); 7552b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 75615a3537eSJed Brown " Problem:\n" 75715a3537eSJed Brown " Problem Name : %s\n" 75815a3537eSJed Brown " Stabilization : %s\n", 7592b916ea7SJeremy L Thompson app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization])); 760e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ctx)); 761d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 7623a8779fbSJames Wright } 763