xref: /honee/problems/newtonian.c (revision 9b103f75867128bb395d4431a2dd4da8eacd1da9)
1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
33a8779fbSJames Wright //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
53a8779fbSJames Wright //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
73a8779fbSJames Wright 
83a8779fbSJames Wright /// @file
93a8779fbSJames Wright /// Utility functions for setting up problems using the Newtonian Qfunction
103a8779fbSJames Wright 
113a8779fbSJames Wright #include "../qfunctions/newtonian.h"
123a8779fbSJames Wright 
13e419654dSJeremy L Thompson #include <ceed.h>
14e419654dSJeremy L Thompson #include <petscdm.h>
15e419654dSJeremy L Thompson 
162b916ea7SJeremy L Thompson #include "../navierstokes.h"
176dd99bedSLeila Ghaffari 
18f7534ca8SJed Brown // For use with PetscOptionsEnum
19c62f7daeSJames Wright static const char *const StateVariables[] = {"conservative", "primitive", "entropy", "StateVariable", "STATEVAR_", NULL};
20f7534ca8SJed Brown 
21c62f7daeSJames Wright static PetscErrorCode CheckQWithTolerance(const CeedScalar Q_s[5], const CeedScalar Q_a[5], const CeedScalar Q_b[5], const char *name,
22c62f7daeSJames Wright                                           PetscReal rtol_0, PetscReal rtol_u, PetscReal rtol_4) {
23c62f7daeSJames Wright   CeedScalar relative_error[5];  // relative error
24c62f7daeSJames Wright   CeedScalar divisor_threshold = 10 * CEED_EPSILON;
2506f41313SJames Wright 
2606f41313SJames Wright   PetscFunctionBeginUser;
27c62f7daeSJames Wright   relative_error[0] = (Q_a[0] - Q_b[0]) / (fabs(Q_s[0]) > divisor_threshold ? Q_s[0] : 1);
28c62f7daeSJames Wright   relative_error[4] = (Q_a[4] - Q_b[4]) / (fabs(Q_s[4]) > divisor_threshold ? Q_s[4] : 1);
29c62f7daeSJames Wright 
30c62f7daeSJames Wright   CeedScalar u_magnitude = sqrt(Square(Q_s[1]) + Square(Q_s[2]) + Square(Q_s[3]));
31c62f7daeSJames Wright   CeedScalar u_divisor   = u_magnitude > divisor_threshold ? u_magnitude : 1;
32c62f7daeSJames Wright   for (int i = 1; i < 4; i++) {
33c62f7daeSJames Wright     relative_error[i] = (Q_a[i] - Q_b[i]) / u_divisor;
342b916ea7SJeremy L Thompson   }
35c62f7daeSJames Wright 
36c62f7daeSJames Wright   if (fabs(relative_error[0]) >= rtol_0) {
37c62f7daeSJames Wright     printf("%s[0] error %g (expected %.10e, got %.10e)\n", name, relative_error[0], Q_s[0], Q_a[0]);
38c62f7daeSJames Wright   }
39c62f7daeSJames Wright   for (int i = 1; i < 4; i++) {
40c62f7daeSJames Wright     if (fabs(relative_error[i]) >= rtol_u) {
41c62f7daeSJames Wright       printf("%s[%d] error %g (expected %.10e, got %.10e)\n", name, i, relative_error[i], Q_s[i], Q_a[i]);
42c62f7daeSJames Wright     }
43c62f7daeSJames Wright   }
44c62f7daeSJames Wright   if (fabs(relative_error[4]) >= rtol_4) {
45c62f7daeSJames Wright     printf("%s[4] error %g (expected %.10e, got %.10e)\n", name, relative_error[4], Q_s[4], Q_a[4]);
46c62f7daeSJames Wright   }
47d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
48f0b65372SJed Brown }
49f0b65372SJed Brown 
50c62f7daeSJames Wright // @brief Verify `StateFromQ` by converting A0 -> B0 -> A0_test, where A0 should equal A0_test
51c62f7daeSJames Wright static PetscErrorCode TestState(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5],
52c62f7daeSJames Wright                                 CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
53c62f7daeSJames Wright   CeedScalar        B0[5], A0_test[5];
54c62f7daeSJames Wright   char              buf[128];
55c62f7daeSJames Wright   const char *const StateVariables_Initial[] = {"U", "Y", "V"};
56c62f7daeSJames Wright 
57c62f7daeSJames Wright   PetscFunctionBeginUser;
58c62f7daeSJames Wright   const char *A_initial = StateVariables_Initial[state_var_A];
59c62f7daeSJames Wright   const char *B_initial = StateVariables_Initial[state_var_B];
60c62f7daeSJames Wright 
61c62f7daeSJames Wright   State state_A0 = StateFromQ(gas, A0, state_var_A);
62c62f7daeSJames Wright   StateToQ(gas, state_A0, B0, state_var_B);
63c62f7daeSJames Wright   State state_B0 = StateFromQ(gas, B0, state_var_B);
64c62f7daeSJames Wright   StateToQ(gas, state_B0, A0_test, state_var_A);
65c62f7daeSJames Wright 
66c62f7daeSJames Wright   snprintf(buf, sizeof buf, "%s->%s->%s: %s", A_initial, B_initial, A_initial, A_initial);
67c62f7daeSJames Wright   PetscCall(CheckQWithTolerance(A0, A0_test, A0, buf, rtol_0, rtol_u, rtol_4));
68c62f7daeSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
69c62f7daeSJames Wright }
70c62f7daeSJames Wright 
71c62f7daeSJames Wright // @brief Verify `StateFromQ_fwd` via a finite difference approximation
72c62f7daeSJames Wright static PetscErrorCode TestState_fwd(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5],
73c62f7daeSJames Wright                                     CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
74c62f7daeSJames Wright   CeedScalar        eps = 4e-7;  // Finite difference step
75c62f7daeSJames Wright   char              buf[128];
76c62f7daeSJames Wright   const char *const StateVariables_Initial[] = {"U", "Y", "V"};
77c62f7daeSJames Wright 
78c62f7daeSJames Wright   PetscFunctionBeginUser;
79c62f7daeSJames Wright   const char *A_initial = StateVariables_Initial[state_var_A];
80c62f7daeSJames Wright   const char *B_initial = StateVariables_Initial[state_var_B];
81c62f7daeSJames Wright   State       state_0   = StateFromQ(gas, A0, state_var_A);
82c62f7daeSJames Wright 
83c62f7daeSJames Wright   for (int i = 0; i < 5; i++) {
84c62f7daeSJames Wright     CeedScalar dB[5] = {0.}, dB_fd[5] = {0.};
85c62f7daeSJames Wright     {  // Calculate dB using State functions
86c62f7daeSJames Wright       CeedScalar dA[5] = {0};
87c62f7daeSJames Wright 
88c62f7daeSJames Wright       dA[i]    = A0[i];
89c62f7daeSJames Wright       State ds = StateFromQ_fwd(gas, state_0, dA, state_var_A);
90c62f7daeSJames Wright       StateToQ(gas, ds, dB, state_var_B);
91c62f7daeSJames Wright     }
92c62f7daeSJames Wright 
93c62f7daeSJames Wright     {  // Calculate dB_fd via finite difference approximation
94c62f7daeSJames Wright       CeedScalar A1[5], B0[5], B1[5];
95c62f7daeSJames Wright 
96c62f7daeSJames Wright       for (int j = 0; j < 5; j++) A1[j] = (1 + eps * (i == j)) * A0[j];
97c62f7daeSJames Wright       State state_1 = StateFromQ(gas, A1, state_var_A);
98c62f7daeSJames Wright       StateToQ(gas, state_0, B0, state_var_B);
99c62f7daeSJames Wright       StateToQ(gas, state_1, B1, state_var_B);
100c62f7daeSJames Wright       for (int j = 0; j < 5; j++) dB_fd[j] = (B1[j] - B0[j]) / eps;
101c62f7daeSJames Wright     }
102c62f7daeSJames Wright 
103c62f7daeSJames 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);
104c62f7daeSJames Wright     PetscCall(CheckQWithTolerance(dB_fd, dB, dB_fd, buf, rtol_0, rtol_u, rtol_4));
105c62f7daeSJames Wright   }
106c62f7daeSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
107c62f7daeSJames Wright }
108c62f7daeSJames Wright 
109c62f7daeSJames Wright // @brief Test the Newtonian State transformation functions, `StateFrom*`
1102b916ea7SJeremy L Thompson static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) {
111f0b65372SJed Brown   Units            units = user->units;
112c62f7daeSJames Wright   const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, K = units->Kelvin;
113c62f7daeSJames Wright 
114f0b65372SJed Brown   PetscFunctionBeginUser;
115c62f7daeSJames Wright   const CeedScalar T          = 200 * K;
116c62f7daeSJames Wright   const CeedScalar rho        = 1.2 * kg / Cube(m);
117c62f7daeSJames Wright   const CeedScalar P          = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T;
118c62f7daeSJames Wright   const CeedScalar u_base     = 40 * m / sec;
119c62f7daeSJames Wright   const CeedScalar u[3]       = {u_base, u_base * 1.1, u_base * 1.2};
120c62f7daeSJames Wright   const CeedScalar e_kinetic  = 0.5 * Dot3(u, u);
121c62f7daeSJames Wright   const CeedScalar e_internal = gas->cv * T;
122c62f7daeSJames Wright   const CeedScalar e_total    = e_kinetic + e_internal;
123c62f7daeSJames Wright   const CeedScalar gamma      = HeatCapacityRatio(gas);
124c62f7daeSJames Wright   const CeedScalar entropy    = log(P) - gamma * log(rho);
125c62f7daeSJames Wright   const CeedScalar rho_div_p  = rho / P;
126c62f7daeSJames Wright   const CeedScalar Y0[5]      = {P, u[0], u[1], u[2], T};
127c62f7daeSJames Wright   const CeedScalar U0[5]      = {rho, rho * u[0], rho * u[1], rho * u[2], rho * e_total};
128c62f7daeSJames 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],
129c62f7daeSJames Wright                                  -rho_div_p};
130c62f7daeSJames Wright 
131c62f7daeSJames Wright   {
132c62f7daeSJames Wright     CeedScalar rtol = 20 * CEED_EPSILON;
133c62f7daeSJames Wright 
134c62f7daeSJames Wright     PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol));
135c62f7daeSJames Wright     PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol));
136c62f7daeSJames Wright     PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol));
137c62f7daeSJames Wright     PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, rtol, rtol, rtol));
138c62f7daeSJames Wright     PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, rtol, rtol, rtol));
139c62f7daeSJames Wright     PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, rtol, rtol, rtol));
140c62f7daeSJames Wright   }
141c62f7daeSJames Wright 
142c62f7daeSJames Wright   {
143c62f7daeSJames Wright     CeedScalar rtol = 5e-6;
144c62f7daeSJames Wright 
145c62f7daeSJames Wright     PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol));
146c62f7daeSJames Wright     PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol));
147c62f7daeSJames Wright     PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol));
148c62f7daeSJames Wright     PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, 10 * rtol, rtol, rtol));
149c62f7daeSJames Wright     PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, 5 * rtol, rtol, rtol));
150c62f7daeSJames Wright     PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, 5 * rtol, 5 * rtol, 5 * rtol));
151f0b65372SJed Brown   }
152d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
153f0b65372SJed Brown }
154f0b65372SJed Brown 
15565dee3d2SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping
15665dee3d2SJames Wright //
15765dee3d2SJames Wright // Only used for SUPG stabilization
15865dee3d2SJames Wright PetscErrorCode CreateKSPMassOperator_NewtonianStabilized(User user, CeedOperator *op_mass) {
15965dee3d2SJames Wright   Ceed                 ceed = user->ceed;
16065dee3d2SJames Wright   CeedInt              num_comp_q, q_data_size;
16165dee3d2SJames Wright   CeedQFunction        qf_mass;
16265dee3d2SJames Wright   CeedElemRestriction  elem_restr_q, elem_restr_qd_i;
16365dee3d2SJames Wright   CeedBasis            basis_q;
16465dee3d2SJames Wright   CeedVector           q_data;
16565dee3d2SJames Wright   CeedQFunctionContext qf_ctx = NULL;
16665dee3d2SJames Wright   PetscInt             dim    = 3;
16765dee3d2SJames Wright 
16865dee3d2SJames Wright   PetscFunctionBeginUser;
16965dee3d2SJames Wright   {  // Get restriction and basis from the RHS function
17065dee3d2SJames Wright     CeedOperator     *sub_ops;
17165dee3d2SJames Wright     CeedOperatorField field;
17265dee3d2SJames Wright     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
17365dee3d2SJames Wright 
17465dee3d2SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops));
17565dee3d2SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
17665dee3d2SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
17765dee3d2SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q));
17865dee3d2SJames Wright 
17965dee3d2SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field));
18065dee3d2SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i));
18165dee3d2SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data));
18265dee3d2SJames Wright 
18365dee3d2SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx));
18465dee3d2SJames Wright   }
18565dee3d2SJames Wright 
18665dee3d2SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
18765dee3d2SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size));
18865dee3d2SJames Wright 
18965dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Newtonian_Conserv, MassFunction_Newtonian_Conserv_loc, &qf_mass));
19065dee3d2SJames Wright 
19165dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx));
19265dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
19365dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
19465dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
19565dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
19665dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
19765dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
19865dee3d2SJames Wright 
19965dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
20065dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
20165dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed));
20265dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data));
20365dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
20465dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
20565dee3d2SJames Wright 
20665dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
20765dee3d2SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
20865dee3d2SJames Wright }
2094c07ec22SJames Wright 
210991aef52SJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
211b7f03f12SJed Brown   SetupContext             setup_context;
2123a8779fbSJames Wright   User                     user   = *(User *)ctx;
2134351d27fSLeila Ghaffari   CeedInt                  degree = user->app_ctx->degree;
2143a8779fbSJames Wright   StabilizationType        stab;
2153636f6a4SJames Wright   StateVariable            state_var;
216b4c37c5cSJames Wright   MPI_Comm                 comm = user->comm;
217b4c37c5cSJames Wright   Ceed                     ceed = user->ceed;
2183a8779fbSJames Wright   PetscBool                implicit;
2190d8cd818SJames Wright   PetscBool                unit_tests;
22015a3537eSJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
22115a3537eSJed Brown   CeedQFunctionContext     newtonian_ig_context;
2223a8779fbSJames Wright 
22315a3537eSJed Brown   PetscFunctionBeginUser;
2242b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &setup_context));
2252b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &newtonian_ig_ctx));
2263a8779fbSJames Wright 
2273a8779fbSJames Wright   // ------------------------------------------------------
2283a8779fbSJames Wright   //           Setup Generic Newtonian IG Problem
2293a8779fbSJames Wright   // ------------------------------------------------------
2303a8779fbSJames Wright   problem->dim               = 3;
231b01ba163SJames Wright   problem->jac_data_size_sur = 11;
2323a8779fbSJames Wright   problem->non_zero_time     = PETSC_FALSE;
233cbe60e31SLeila Ghaffari   problem->print_info        = PRINT_NEWTONIAN;
2347d8a615bSJames Wright   problem->uses_newtonian    = PETSC_TRUE;
2353a8779fbSJames Wright 
2363a8779fbSJames Wright   // ------------------------------------------------------
2373a8779fbSJames Wright   //             Create the libCEED context
2383a8779fbSJames Wright   // ------------------------------------------------------
2393a8779fbSJames Wright   CeedScalar cv         = 717.;          // J/(kg K)
2403a8779fbSJames Wright   CeedScalar cp         = 1004.;         // J/(kg K)
241d9bb1cdbSJames Wright   CeedScalar g[3]       = {0, 0, 0};     // m/s^2
2423a8779fbSJames Wright   CeedScalar lambda     = -2. / 3.;      // -
243bb8a0c61SJames Wright   CeedScalar mu         = 1.8e-5;        // Pa s, dynamic viscosity
2443a8779fbSJames Wright   CeedScalar k          = 0.02638;       // W/(m K)
2454351d27fSLeila Ghaffari   CeedScalar c_tau      = 0.5 / degree;  // -
246bb8a0c61SJames Wright   CeedScalar Ctau_t     = 1.0;           // -
247b5786772SLeila Ghaffari   CeedScalar Cv_func[3] = {36, 60, 128};
248c176988bSLeila Ghaffari   CeedScalar Ctau_v     = Cv_func[(CeedInt)Min(3, degree) - 1];
2492a0eee6eSLeila Ghaffari   CeedScalar Ctau_C     = 0.25 / degree;
2502a0eee6eSLeila Ghaffari   CeedScalar Ctau_M     = 0.25 / degree;
2512a0eee6eSLeila Ghaffari   CeedScalar Ctau_E     = 0.125;
2523a8779fbSJames Wright   PetscReal  domain_min[3], domain_max[3], domain_size[3];
2532b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
254493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
2553a8779fbSJames Wright 
256b8fb7609SAdeleke O. Bankole   StatePrimitive reference      = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15};
257a213b8aaSJames Wright   CeedScalar     idl_decay_time = -1, idl_start = 0, idl_length = 0, idl_pressure = reference.pressure;
258e7754af5SKenneth E. Jansen   PetscBool      idl_enable = PETSC_FALSE;
259b8fb7609SAdeleke O. Bankole 
2603a8779fbSJames Wright   // ------------------------------------------------------
2613a8779fbSJames Wright   //             Create the PETSc context
2623a8779fbSJames Wright   // ------------------------------------------------------
263bb8a0c61SJames Wright   PetscScalar meter    = 1;  // 1 meter in scaled length units
264bb8a0c61SJames Wright   PetscScalar kilogram = 1;  // 1 kilogram in scaled mass units
265bb8a0c61SJames Wright   PetscScalar second   = 1;  // 1 second in scaled time units
2663a8779fbSJames Wright   PetscScalar Kelvin   = 1;  // 1 Kelvin in scaled temperature units
2673a8779fbSJames Wright   PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s;
2683a8779fbSJames Wright 
2693a8779fbSJames Wright   // ------------------------------------------------------
2703a8779fbSJames Wright   //              Command line Options
2713a8779fbSJames Wright   // ------------------------------------------------------
272d9bb1cdbSJames Wright   PetscBool given_option = PETSC_FALSE;
2732b916ea7SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL);
274cbe60e31SLeila Ghaffari   // -- Conservative vs Primitive variables
2752b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE),
2762b916ea7SJeremy L Thompson                              (PetscEnum *)&state_var, NULL));
2773636f6a4SJames Wright 
2783636f6a4SJames Wright   switch (state_var) {
2793636f6a4SJames Wright     case STATEVAR_CONSERVATIVE:
280b8fb7609SAdeleke O. Bankole       problem->ics.qfunction                       = ICsNewtonianIG_Conserv;
281b8fb7609SAdeleke O. Bankole       problem->ics.qfunction_loc                   = ICsNewtonianIG_Conserv_loc;
282cbe60e31SLeila Ghaffari       problem->apply_vol_rhs.qfunction             = RHSFunction_Newtonian;
283cbe60e31SLeila Ghaffari       problem->apply_vol_rhs.qfunction_loc         = RHSFunction_Newtonian_loc;
28476555becSJames Wright       problem->apply_vol_ifunction.qfunction       = IFunction_Newtonian_Conserv;
28576555becSJames Wright       problem->apply_vol_ifunction.qfunction_loc   = IFunction_Newtonian_Conserv_loc;
28676555becSJames Wright       problem->apply_vol_ijacobian.qfunction       = IJacobian_Newtonian_Conserv;
28776555becSJames Wright       problem->apply_vol_ijacobian.qfunction_loc   = IJacobian_Newtonian_Conserv_loc;
288d4559bbeSJames Wright       problem->apply_inflow.qfunction              = BoundaryIntegral_Conserv;
289d4559bbeSJames Wright       problem->apply_inflow.qfunction_loc          = BoundaryIntegral_Conserv_loc;
290d4559bbeSJames Wright       problem->apply_inflow_jacobian.qfunction     = BoundaryIntegral_Jacobian_Conserv;
291d4559bbeSJames Wright       problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc;
2923636f6a4SJames Wright       break;
2933636f6a4SJames Wright 
2943636f6a4SJames Wright     case STATEVAR_PRIMITIVE:
2953636f6a4SJames Wright       problem->ics.qfunction                       = ICsNewtonianIG_Prim;
2963636f6a4SJames Wright       problem->ics.qfunction_loc                   = ICsNewtonianIG_Prim_loc;
2973636f6a4SJames Wright       problem->apply_vol_ifunction.qfunction       = IFunction_Newtonian_Prim;
2983636f6a4SJames Wright       problem->apply_vol_ifunction.qfunction_loc   = IFunction_Newtonian_Prim_loc;
2993636f6a4SJames Wright       problem->apply_vol_ijacobian.qfunction       = IJacobian_Newtonian_Prim;
3003636f6a4SJames Wright       problem->apply_vol_ijacobian.qfunction_loc   = IJacobian_Newtonian_Prim_loc;
3013636f6a4SJames Wright       problem->apply_inflow.qfunction              = BoundaryIntegral_Prim;
3023636f6a4SJames Wright       problem->apply_inflow.qfunction_loc          = BoundaryIntegral_Prim_loc;
3033636f6a4SJames Wright       problem->apply_inflow_jacobian.qfunction     = BoundaryIntegral_Jacobian_Prim;
3043636f6a4SJames Wright       problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc;
3053636f6a4SJames Wright       break;
306c62f7daeSJames Wright     case STATEVAR_ENTROPY:
307*9b103f75SJames Wright       problem->ics.qfunction                       = ICsNewtonianIG_Entropy;
308*9b103f75SJames Wright       problem->ics.qfunction_loc                   = ICsNewtonianIG_Entropy_loc;
309*9b103f75SJames Wright       problem->apply_vol_ifunction.qfunction       = IFunction_Newtonian_Entropy;
310*9b103f75SJames Wright       problem->apply_vol_ifunction.qfunction_loc   = IFunction_Newtonian_Entropy_loc;
311*9b103f75SJames Wright       problem->apply_vol_ijacobian.qfunction       = IJacobian_Newtonian_Entropy;
312*9b103f75SJames Wright       problem->apply_vol_ijacobian.qfunction_loc   = IJacobian_Newtonian_Entropy_loc;
313*9b103f75SJames Wright       problem->apply_inflow.qfunction              = BoundaryIntegral_Entropy;
314*9b103f75SJames Wright       problem->apply_inflow.qfunction_loc          = BoundaryIntegral_Entropy_loc;
315*9b103f75SJames Wright       problem->apply_inflow_jacobian.qfunction     = BoundaryIntegral_Jacobian_Entropy;
316*9b103f75SJames Wright       problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Entropy_loc;
317c62f7daeSJames Wright       break;
318cbe60e31SLeila Ghaffari   }
3191485969bSJeremy L Thompson 
3203a8779fbSJames Wright   // -- Physics
3212b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL));
3222b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL));
3232b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL));
3242b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL));
3252b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL));
3263a8779fbSJames Wright 
327bb8a0c61SJames Wright   PetscInt dim = problem->dim;
328d9bb1cdbSJames Wright   PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL));
329d9bb1cdbSJames Wright   PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option));
330d9bb1cdbSJames Wright   if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim);
331d9bb1cdbSJames Wright 
3322b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
3332b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
3342b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
3352b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL));
3362b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL));
3372b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL));
3382b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL));
3392b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
3402b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL));
3413a8779fbSJames Wright 
342db95602eSJames Wright   dim = 3;
343b8fb7609SAdeleke O. Bankole   PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL));
344b8fb7609SAdeleke O. Bankole   PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL));
345b8fb7609SAdeleke O. Bankole   PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL));
346b8fb7609SAdeleke O. Bankole 
3473a8779fbSJames Wright   // -- Units
3482b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
3493a8779fbSJames Wright   meter = fabs(meter);
3502b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL));
3513a8779fbSJames Wright   kilogram = fabs(kilogram);
3522b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
3533a8779fbSJames Wright   second = fabs(second);
3542b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL));
3553a8779fbSJames Wright   Kelvin = fabs(Kelvin);
3563a8779fbSJames Wright 
3573a8779fbSJames Wright   // -- Warnings
3585d28dccaSJames Wright   PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP,
3595d28dccaSJames Wright              "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n");
360e7754af5SKenneth E. Jansen 
361e7754af5SKenneth E. Jansen   PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point",
362e7754af5SKenneth E. Jansen                                NULL, idl_decay_time, &idl_decay_time, &idl_enable));
3639cbdf780SJames Wright   PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero.");
3649cbdf780SJames Wright   if (idl_decay_time < 0) idl_enable = PETSC_FALSE;
365e7754af5SKenneth E. Jansen   PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL));
366e7754af5SKenneth E. Jansen   PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL));
367a213b8aaSJames Wright   idl_pressure = reference.pressure;
368a213b8aaSJames Wright   PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-reference_pressure`)", NULL, idl_pressure,
369a213b8aaSJames Wright                                &idl_pressure, NULL));
3701485969bSJeremy L Thompson   PetscOptionsEnd();
3713a8779fbSJames Wright 
37265dee3d2SJames Wright   if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized;
37365dee3d2SJames Wright 
3743a8779fbSJames Wright   // ------------------------------------------------------
3753a8779fbSJames Wright   //           Set up the PETSc context
3763a8779fbSJames Wright   // ------------------------------------------------------
3773a8779fbSJames Wright   // -- Define derived units
3783a8779fbSJames Wright   Pascal          = kilogram / (meter * PetscSqr(second));
3793a8779fbSJames Wright   J_per_kg_K      = PetscSqr(meter) / (PetscSqr(second) * Kelvin);
3803a8779fbSJames Wright   m_per_squared_s = meter / PetscSqr(second);
3813a8779fbSJames Wright   W_per_m_K       = kilogram * meter / (pow(second, 3) * Kelvin);
3823a8779fbSJames Wright 
3833a8779fbSJames Wright   user->units->meter           = meter;
3843a8779fbSJames Wright   user->units->kilogram        = kilogram;
3853a8779fbSJames Wright   user->units->second          = second;
3863a8779fbSJames Wright   user->units->Kelvin          = Kelvin;
3873a8779fbSJames Wright   user->units->Pascal          = Pascal;
3883a8779fbSJames Wright   user->units->J_per_kg_K      = J_per_kg_K;
3893a8779fbSJames Wright   user->units->m_per_squared_s = m_per_squared_s;
3903a8779fbSJames Wright   user->units->W_per_m_K       = W_per_m_K;
3913a8779fbSJames Wright 
3923a8779fbSJames Wright   // ------------------------------------------------------
3933a8779fbSJames Wright   //           Set up the libCEED context
3943a8779fbSJames Wright   // ------------------------------------------------------
3953a8779fbSJames Wright   // -- Scale variables to desired units
3963a8779fbSJames Wright   cv *= J_per_kg_K;
3973a8779fbSJames Wright   cp *= J_per_kg_K;
3983a8779fbSJames Wright   mu *= Pascal * second;
3993a8779fbSJames Wright   k *= W_per_m_K;
400493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter;
401493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s;
402b8fb7609SAdeleke O. Bankole   reference.pressure *= Pascal;
403b8fb7609SAdeleke O. Bankole   for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second;
404b8fb7609SAdeleke O. Bankole   reference.temperature *= Kelvin;
4053a8779fbSJames Wright   problem->dm_scale = meter;
4063a8779fbSJames Wright 
4073a8779fbSJames Wright   // -- Solver Settings
4083a8779fbSJames Wright   user->phys->implicit  = implicit;
4093636f6a4SJames Wright   user->phys->state_var = state_var;
4103a8779fbSJames Wright 
4113a8779fbSJames Wright   // -- QFunction Context
41215a3537eSJed Brown   newtonian_ig_ctx->lambda        = lambda;
41315a3537eSJed Brown   newtonian_ig_ctx->mu            = mu;
41415a3537eSJed Brown   newtonian_ig_ctx->k             = k;
41515a3537eSJed Brown   newtonian_ig_ctx->cv            = cv;
41615a3537eSJed Brown   newtonian_ig_ctx->cp            = cp;
41715a3537eSJed Brown   newtonian_ig_ctx->c_tau         = c_tau;
41815a3537eSJed Brown   newtonian_ig_ctx->Ctau_t        = Ctau_t;
41915a3537eSJed Brown   newtonian_ig_ctx->Ctau_v        = Ctau_v;
42015a3537eSJed Brown   newtonian_ig_ctx->Ctau_C        = Ctau_C;
42115a3537eSJed Brown   newtonian_ig_ctx->Ctau_M        = Ctau_M;
42215a3537eSJed Brown   newtonian_ig_ctx->Ctau_E        = Ctau_E;
42315a3537eSJed Brown   newtonian_ig_ctx->stabilization = stab;
4248085925cSJames Wright   newtonian_ig_ctx->is_implicit   = implicit;
4253636f6a4SJames Wright   newtonian_ig_ctx->state_var     = state_var;
426e7754af5SKenneth E. Jansen   newtonian_ig_ctx->idl_enable    = idl_enable;
427e7754af5SKenneth E. Jansen   newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * second);
428e7754af5SKenneth E. Jansen   newtonian_ig_ctx->idl_start     = idl_start * meter;
429e7754af5SKenneth E. Jansen   newtonian_ig_ctx->idl_length    = idl_length * meter;
430a213b8aaSJames Wright   newtonian_ig_ctx->idl_pressure  = idl_pressure;
4312b916ea7SJeremy L Thompson   PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3));
4323a8779fbSJames Wright 
433b8fb7609SAdeleke O. Bankole   // -- Setup Context
434b8fb7609SAdeleke O. Bankole   setup_context->reference = reference;
435b8fb7609SAdeleke O. Bankole   setup_context->gas       = *newtonian_ig_ctx;
436b8fb7609SAdeleke O. Bankole   setup_context->lx        = domain_size[0];
437b8fb7609SAdeleke O. Bankole   setup_context->ly        = domain_size[1];
438b8fb7609SAdeleke O. Bankole   setup_context->lz        = domain_size[2];
439b8fb7609SAdeleke O. Bankole   setup_context->time      = 0;
440b8fb7609SAdeleke O. Bankole 
441b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context));
442b4c37c5cSJames Wright   PetscCallCeed(ceed,
443b4c37c5cSJames Wright                 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
444b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc));
445b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, "evaluation time", offsetof(struct SetupContext_, time), 1,
446b4c37c5cSJames Wright                                                          "Time of evaluation"));
44715a3537eSJed Brown 
448b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context));
449b4c37c5cSJames Wright   PetscCallCeed(ceed,
450b4c37c5cSJames Wright                 CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx));
451b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, FreeContextPetsc));
452b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1,
453b4c37c5cSJames Wright                                                          "Size of timestep, delta t"));
454b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift",
455b4c37c5cSJames Wright                                                          offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1,
456b4c37c5cSJames Wright                                                          "Shift for mass matrix in IJacobian"));
457b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1,
458b4c37c5cSJames Wright                                                          "Current solution time"));
459a8bca8acSJames Wright 
46015a3537eSJed Brown   problem->apply_vol_rhs.qfunction_context = newtonian_ig_context;
461b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ifunction.qfunction_context));
462b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ijacobian.qfunction_context));
463b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow.qfunction_context));
464b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow_jacobian.qfunction_context));
465f0b65372SJed Brown 
4669ed3d70dSJames Wright   if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
4679ed3d70dSJames Wright   if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
4689ed3d70dSJames Wright   if (bc->num_slip > 0) PetscCall(SlipBCSetup(problem, dm, ctx, newtonian_ig_context));
4699ed3d70dSJames Wright 
470f0b65372SJed Brown   if (unit_tests) {
471f0b65372SJed Brown     PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx));
472f0b65372SJed Brown   }
473d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4743a8779fbSJames Wright }
4753a8779fbSJames Wright 
476991aef52SJames Wright PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx) {
4772d49c0afSJames Wright   MPI_Comm                 comm = user->comm;
478b4c37c5cSJames Wright   Ceed                     ceed = user->ceed;
47915a3537eSJed Brown   NewtonianIdealGasContext newtonian_ctx;
48015a3537eSJed Brown 
4813a8779fbSJames Wright   PetscFunctionBeginUser;
482b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ctx));
4832b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(comm,
48415a3537eSJed Brown                         "  Problem:\n"
48515a3537eSJed Brown                         "    Problem Name                       : %s\n"
48615a3537eSJed Brown                         "    Stabilization                      : %s\n",
4872b916ea7SJeremy L Thompson                         app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization]));
488b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ctx));
489d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4903a8779fbSJames Wright }
491