1*ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2*ae2b091fSJames 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*` 1062b916ea7SJeremy L Thompson static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) { 107f0b65372SJed Brown Units units = user->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 15465dee3d2SJames Wright PetscErrorCode CreateKSPMassOperator_NewtonianStabilized(User user, CeedOperator *op_mass) { 15565dee3d2SJames Wright Ceed ceed = user->ceed; 15665dee3d2SJames Wright CeedInt num_comp_q, q_data_size; 15765dee3d2SJames Wright CeedQFunction qf_mass; 15865dee3d2SJames Wright CeedElemRestriction elem_restr_q, elem_restr_qd_i; 15965dee3d2SJames Wright CeedBasis basis_q; 16065dee3d2SJames Wright CeedVector q_data; 16165dee3d2SJames Wright CeedQFunctionContext qf_ctx = 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; 16765dee3d2SJames Wright CeedOperatorField field; 16865dee3d2SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 16965dee3d2SJames Wright 17065dee3d2SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops)); 17165dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field)); 17265dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q)); 17365dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q)); 17465dee3d2SJames Wright 17565dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field)); 17665dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i)); 17765dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data)); 17865dee3d2SJames Wright 17965dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx)); 18065dee3d2SJames Wright } 18165dee3d2SJames Wright 18265dee3d2SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 18365dee3d2SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size)); 18465dee3d2SJames Wright 18565dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Newtonian_Conserv, MassFunction_Newtonian_Conserv_loc, &qf_mass)); 18665dee3d2SJames Wright 18765dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx)); 18865dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0)); 18965dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP)); 19065dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP)); 19165dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE)); 19265dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP)); 19365dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD)); 19465dee3d2SJames Wright 19565dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 19665dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 19765dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed)); 19865dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data)); 19965dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 20065dee3d2SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 20165dee3d2SJames Wright 202d438a78dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qf_ctx)); 20365dee3d2SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 20465dee3d2SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 20565dee3d2SJames Wright } 2064c07ec22SJames Wright 207991aef52SJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 208b7f03f12SJed Brown SetupContext setup_context; 2093a8779fbSJames Wright User user = *(User *)ctx; 2104351d27fSLeila Ghaffari CeedInt degree = user->app_ctx->degree; 2113a8779fbSJames Wright StabilizationType stab; 2123636f6a4SJames Wright StateVariable state_var; 213b4c37c5cSJames Wright MPI_Comm comm = user->comm; 214b4c37c5cSJames Wright Ceed ceed = user->ceed; 2153a8779fbSJames Wright PetscBool implicit; 2160d8cd818SJames Wright PetscBool unit_tests; 21715a3537eSJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 21815a3537eSJed Brown CeedQFunctionContext newtonian_ig_context; 2193a8779fbSJames Wright 22015a3537eSJed Brown PetscFunctionBeginUser; 2212b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &setup_context)); 2222b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &newtonian_ig_ctx)); 2233a8779fbSJames Wright 2243a8779fbSJames Wright // ------------------------------------------------------ 2253a8779fbSJames Wright // Setup Generic Newtonian IG Problem 2263a8779fbSJames Wright // ------------------------------------------------------ 2273a8779fbSJames Wright problem->dim = 3; 228b01ba163SJames Wright problem->jac_data_size_sur = 11; 22958ce1233SJames Wright problem->compute_exact_solution_error = PETSC_FALSE; 230cbe60e31SLeila Ghaffari problem->print_info = PRINT_NEWTONIAN; 2317d8a615bSJames Wright problem->uses_newtonian = PETSC_TRUE; 2323a8779fbSJames Wright 2333a8779fbSJames Wright // ------------------------------------------------------ 2343a8779fbSJames Wright // Create the libCEED context 2353a8779fbSJames Wright // ------------------------------------------------------ 2363a8779fbSJames Wright CeedScalar cv = 717.; // J/(kg K) 2373a8779fbSJames Wright CeedScalar cp = 1004.; // J/(kg K) 238d9bb1cdbSJames Wright CeedScalar g[3] = {0, 0, 0}; // m/s^2 2393a8779fbSJames Wright CeedScalar lambda = -2. / 3.; // - 240bb8a0c61SJames Wright CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 2413a8779fbSJames Wright CeedScalar k = 0.02638; // W/(m K) 2424351d27fSLeila Ghaffari CeedScalar c_tau = 0.5 / degree; // - 243bb8a0c61SJames Wright CeedScalar Ctau_t = 1.0; // - 244b5786772SLeila Ghaffari CeedScalar Cv_func[3] = {36, 60, 128}; 245c176988bSLeila Ghaffari CeedScalar Ctau_v = Cv_func[(CeedInt)Min(3, degree) - 1]; 2462a0eee6eSLeila Ghaffari CeedScalar Ctau_C = 0.25 / degree; 2472a0eee6eSLeila Ghaffari CeedScalar Ctau_M = 0.25 / degree; 2482a0eee6eSLeila Ghaffari CeedScalar Ctau_E = 0.125; 2493a8779fbSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 2502b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 251493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 2523a8779fbSJames Wright 253b8fb7609SAdeleke O. Bankole StatePrimitive reference = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15}; 254a213b8aaSJames Wright CeedScalar idl_decay_time = -1, idl_start = 0, idl_length = 0, idl_pressure = reference.pressure; 255e7754af5SKenneth E. Jansen PetscBool idl_enable = PETSC_FALSE; 256b8fb7609SAdeleke O. Bankole 2573a8779fbSJames Wright // ------------------------------------------------------ 2583a8779fbSJames Wright // Create the PETSc context 2593a8779fbSJames Wright // ------------------------------------------------------ 260bb8a0c61SJames Wright PetscScalar meter = 1; // 1 meter in scaled length units 261bb8a0c61SJames Wright PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 262bb8a0c61SJames Wright PetscScalar second = 1; // 1 second in scaled time units 2633a8779fbSJames Wright PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 2643a8779fbSJames Wright PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 2653a8779fbSJames Wright 2663a8779fbSJames Wright // ------------------------------------------------------ 2673a8779fbSJames Wright // Command line Options 2683a8779fbSJames Wright // ------------------------------------------------------ 269d9bb1cdbSJames Wright PetscBool given_option = PETSC_FALSE; 2702b916ea7SJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL); 271cbe60e31SLeila Ghaffari // -- Conservative vs Primitive variables 2722b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE), 2732b916ea7SJeremy L Thompson (PetscEnum *)&state_var, NULL)); 2743636f6a4SJames Wright 2753636f6a4SJames Wright switch (state_var) { 2763636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 277b8fb7609SAdeleke O. Bankole problem->ics.qfunction = ICsNewtonianIG_Conserv; 278b8fb7609SAdeleke O. Bankole problem->ics.qfunction_loc = ICsNewtonianIG_Conserv_loc; 279cbe60e31SLeila Ghaffari problem->apply_vol_rhs.qfunction = RHSFunction_Newtonian; 280cbe60e31SLeila Ghaffari problem->apply_vol_rhs.qfunction_loc = RHSFunction_Newtonian_loc; 28176555becSJames Wright problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Conserv; 28276555becSJames Wright problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Conserv_loc; 28376555becSJames Wright problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Conserv; 28476555becSJames Wright problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Conserv_loc; 285d4559bbeSJames Wright problem->apply_inflow.qfunction = BoundaryIntegral_Conserv; 286d4559bbeSJames Wright problem->apply_inflow.qfunction_loc = BoundaryIntegral_Conserv_loc; 287d4559bbeSJames Wright problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Conserv; 288d4559bbeSJames Wright problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc; 2893636f6a4SJames Wright break; 2903636f6a4SJames Wright case STATEVAR_PRIMITIVE: 2913636f6a4SJames Wright problem->ics.qfunction = ICsNewtonianIG_Prim; 2923636f6a4SJames Wright problem->ics.qfunction_loc = ICsNewtonianIG_Prim_loc; 2933636f6a4SJames Wright problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Prim; 2943636f6a4SJames Wright problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Prim_loc; 2953636f6a4SJames Wright problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Prim; 2963636f6a4SJames Wright problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Prim_loc; 2973636f6a4SJames Wright problem->apply_inflow.qfunction = BoundaryIntegral_Prim; 2983636f6a4SJames Wright problem->apply_inflow.qfunction_loc = BoundaryIntegral_Prim_loc; 2993636f6a4SJames Wright problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Prim; 3003636f6a4SJames Wright problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc; 3013636f6a4SJames Wright break; 302c62f7daeSJames Wright case STATEVAR_ENTROPY: 3039b103f75SJames Wright problem->ics.qfunction = ICsNewtonianIG_Entropy; 3049b103f75SJames Wright problem->ics.qfunction_loc = ICsNewtonianIG_Entropy_loc; 3059b103f75SJames Wright problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Entropy; 3069b103f75SJames Wright problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Entropy_loc; 3079b103f75SJames Wright problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Entropy; 3089b103f75SJames Wright problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Entropy_loc; 3099b103f75SJames Wright problem->apply_inflow.qfunction = BoundaryIntegral_Entropy; 3109b103f75SJames Wright problem->apply_inflow.qfunction_loc = BoundaryIntegral_Entropy_loc; 3119b103f75SJames Wright problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Entropy; 3129b103f75SJames Wright problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Entropy_loc; 313c62f7daeSJames Wright break; 314cbe60e31SLeila Ghaffari } 3151485969bSJeremy L Thompson 3163a8779fbSJames Wright // -- Physics 3172b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL)); 3182b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL)); 3192b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL)); 3202b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL)); 3212b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL)); 3223a8779fbSJames Wright 323bb8a0c61SJames Wright PetscInt dim = problem->dim; 324d9bb1cdbSJames Wright PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL)); 325d9bb1cdbSJames Wright PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option)); 326d9bb1cdbSJames Wright if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim); 327d9bb1cdbSJames Wright 3282b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 3292b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 3302b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 3312b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL)); 3322b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL)); 3332b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL)); 3342b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL)); 3352b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 3362b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL)); 3373a8779fbSJames Wright 338db95602eSJames Wright dim = 3; 339b8fb7609SAdeleke O. Bankole PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL)); 340b8fb7609SAdeleke O. Bankole PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL)); 341b8fb7609SAdeleke O. Bankole PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL)); 342b8fb7609SAdeleke O. Bankole 3433a8779fbSJames Wright // -- Units 3442b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 3453a8779fbSJames Wright meter = fabs(meter); 3462b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 3473a8779fbSJames Wright kilogram = fabs(kilogram); 3482b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 3493a8779fbSJames Wright second = fabs(second); 3502b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL)); 3513a8779fbSJames Wright Kelvin = fabs(Kelvin); 3523a8779fbSJames Wright 3533a8779fbSJames Wright // -- Warnings 3545d28dccaSJames Wright PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP, 3555d28dccaSJames Wright "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n"); 356e7754af5SKenneth E. Jansen 357e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point", 358e7754af5SKenneth E. Jansen NULL, idl_decay_time, &idl_decay_time, &idl_enable)); 3599cbdf780SJames Wright PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero."); 3609cbdf780SJames Wright if (idl_decay_time < 0) idl_enable = PETSC_FALSE; 361e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL)); 362e7754af5SKenneth E. Jansen PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL)); 363a213b8aaSJames Wright idl_pressure = reference.pressure; 364a213b8aaSJames Wright PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-reference_pressure`)", NULL, idl_pressure, 365a213b8aaSJames Wright &idl_pressure, NULL)); 3661485969bSJeremy L Thompson PetscOptionsEnd(); 3673a8779fbSJames Wright 36865dee3d2SJames Wright if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized; 36965dee3d2SJames Wright 3703a8779fbSJames Wright // ------------------------------------------------------ 3713a8779fbSJames Wright // Set up the PETSc context 3723a8779fbSJames Wright // ------------------------------------------------------ 3733a8779fbSJames Wright // -- Define derived units 3743a8779fbSJames Wright Pascal = kilogram / (meter * PetscSqr(second)); 3753a8779fbSJames Wright J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 3763a8779fbSJames Wright m_per_squared_s = meter / PetscSqr(second); 3773a8779fbSJames Wright W_per_m_K = kilogram * meter / (pow(second, 3) * Kelvin); 3783a8779fbSJames Wright 3793a8779fbSJames Wright user->units->meter = meter; 3803a8779fbSJames Wright user->units->kilogram = kilogram; 3813a8779fbSJames Wright user->units->second = second; 3823a8779fbSJames Wright user->units->Kelvin = Kelvin; 3833a8779fbSJames Wright user->units->Pascal = Pascal; 3843a8779fbSJames Wright user->units->J_per_kg_K = J_per_kg_K; 3853a8779fbSJames Wright user->units->m_per_squared_s = m_per_squared_s; 3863a8779fbSJames Wright user->units->W_per_m_K = W_per_m_K; 3873a8779fbSJames Wright 3883a8779fbSJames Wright // ------------------------------------------------------ 3893a8779fbSJames Wright // Set up the libCEED context 3903a8779fbSJames Wright // ------------------------------------------------------ 3913a8779fbSJames Wright // -- Scale variables to desired units 3923a8779fbSJames Wright cv *= J_per_kg_K; 3933a8779fbSJames Wright cp *= J_per_kg_K; 3943a8779fbSJames Wright mu *= Pascal * second; 3953a8779fbSJames Wright k *= W_per_m_K; 396493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter; 397493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s; 398b8fb7609SAdeleke O. Bankole reference.pressure *= Pascal; 399b8fb7609SAdeleke O. Bankole for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second; 400b8fb7609SAdeleke O. Bankole reference.temperature *= Kelvin; 4013a8779fbSJames Wright problem->dm_scale = meter; 4023a8779fbSJames Wright 4033a8779fbSJames Wright // -- Solver Settings 4043a8779fbSJames Wright user->phys->implicit = implicit; 4053636f6a4SJames Wright user->phys->state_var = state_var; 4063a8779fbSJames Wright 4073a8779fbSJames Wright // -- QFunction Context 40815a3537eSJed Brown newtonian_ig_ctx->lambda = lambda; 40915a3537eSJed Brown newtonian_ig_ctx->mu = mu; 41015a3537eSJed Brown newtonian_ig_ctx->k = k; 41115a3537eSJed Brown newtonian_ig_ctx->cv = cv; 41215a3537eSJed Brown newtonian_ig_ctx->cp = cp; 41315a3537eSJed Brown newtonian_ig_ctx->c_tau = c_tau; 41415a3537eSJed Brown newtonian_ig_ctx->Ctau_t = Ctau_t; 41515a3537eSJed Brown newtonian_ig_ctx->Ctau_v = Ctau_v; 41615a3537eSJed Brown newtonian_ig_ctx->Ctau_C = Ctau_C; 41715a3537eSJed Brown newtonian_ig_ctx->Ctau_M = Ctau_M; 41815a3537eSJed Brown newtonian_ig_ctx->Ctau_E = Ctau_E; 41915a3537eSJed Brown newtonian_ig_ctx->stabilization = stab; 4208085925cSJames Wright newtonian_ig_ctx->is_implicit = implicit; 4213636f6a4SJames Wright newtonian_ig_ctx->state_var = state_var; 422e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_enable = idl_enable; 423e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * second); 424e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_start = idl_start * meter; 425e7754af5SKenneth E. Jansen newtonian_ig_ctx->idl_length = idl_length * meter; 426a213b8aaSJames Wright newtonian_ig_ctx->idl_pressure = idl_pressure; 4272b916ea7SJeremy L Thompson PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3)); 4283a8779fbSJames Wright 429b8fb7609SAdeleke O. Bankole // -- Setup Context 430b8fb7609SAdeleke O. Bankole setup_context->reference = reference; 431b8fb7609SAdeleke O. Bankole setup_context->gas = *newtonian_ig_ctx; 432b8fb7609SAdeleke O. Bankole setup_context->lx = domain_size[0]; 433b8fb7609SAdeleke O. Bankole setup_context->ly = domain_size[1]; 434b8fb7609SAdeleke O. Bankole setup_context->lz = domain_size[2]; 435b8fb7609SAdeleke O. Bankole setup_context->time = 0; 436b8fb7609SAdeleke O. Bankole 437b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context)); 438b4c37c5cSJames Wright PetscCallCeed(ceed, 439b4c37c5cSJames Wright CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 440b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc)); 441b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, "evaluation time", offsetof(struct SetupContext_, time), 1, 442b4c37c5cSJames Wright "Time of evaluation")); 44315a3537eSJed Brown 444b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context)); 445b4c37c5cSJames Wright PetscCallCeed(ceed, 446b4c37c5cSJames Wright CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx)); 447b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, FreeContextPetsc)); 448b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1, 449b4c37c5cSJames Wright "Size of timestep, delta t")); 450b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift", 451b4c37c5cSJames Wright offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1, 452b4c37c5cSJames Wright "Shift for mass matrix in IJacobian")); 453b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1, 454b4c37c5cSJames Wright "Current solution time")); 455a8bca8acSJames Wright 45615a3537eSJed Brown problem->apply_vol_rhs.qfunction_context = newtonian_ig_context; 457b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ifunction.qfunction_context)); 458b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ijacobian.qfunction_context)); 459b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow.qfunction_context)); 460b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow_jacobian.qfunction_context)); 461f0b65372SJed Brown 4629ed3d70dSJames Wright if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 4639ed3d70dSJames Wright if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference)); 4649ed3d70dSJames Wright if (bc->num_slip > 0) PetscCall(SlipBCSetup(problem, dm, ctx, newtonian_ig_context)); 4659ed3d70dSJames Wright 466f0b65372SJed Brown if (unit_tests) { 467f0b65372SJed Brown PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx)); 468f0b65372SJed Brown } 469d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 4703a8779fbSJames Wright } 4713a8779fbSJames Wright 472991aef52SJames Wright PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx) { 4732d49c0afSJames Wright MPI_Comm comm = user->comm; 474b4c37c5cSJames Wright Ceed ceed = user->ceed; 47515a3537eSJed Brown NewtonianIdealGasContext newtonian_ctx; 47615a3537eSJed Brown 4773a8779fbSJames Wright PetscFunctionBeginUser; 478b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ctx)); 4792b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 48015a3537eSJed Brown " Problem:\n" 48115a3537eSJed Brown " Problem Name : %s\n" 48215a3537eSJed Brown " Stabilization : %s\n", 4832b916ea7SJeremy L Thompson app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization])); 484b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ctx)); 485d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 4863a8779fbSJames Wright } 487