xref: /honee/problems/newtonian.c (revision 991aef52b7f76ea745826c507d0f3e4bba35bb34)
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"
172b916ea7SJeremy L Thompson #include "../qfunctions/setupgeo.h"
186dd99bedSLeila Ghaffari 
19f7534ca8SJed Brown // For use with PetscOptionsEnum
20f7534ca8SJed Brown static const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "StateVariable", "STATEVAR_", NULL};
21f7534ca8SJed Brown 
222b916ea7SJeremy L Thompson // Compute relative error |a - b|/|s|
232b916ea7SJeremy L Thompson static PetscErrorCode CheckPrimitiveWithTolerance(StatePrimitive sY, StatePrimitive aY, StatePrimitive bY, const char *name, PetscReal rtol_pressure,
242b916ea7SJeremy L Thompson                                                   PetscReal rtol_velocity, PetscReal rtol_temperature) {
25f0b65372SJed Brown   StatePrimitive eY;  // relative error
2606f41313SJames Wright 
2706f41313SJames Wright   PetscFunctionBeginUser;
28f0b65372SJed Brown   eY.pressure   = (aY.pressure - bY.pressure) / sY.pressure;
292b916ea7SJeremy L Thompson   PetscScalar u = sqrt(Square(sY.velocity[0]) + Square(sY.velocity[1]) + Square(sY.velocity[2]));
30f0b65372SJed Brown   for (int j = 0; j < 3; j++) eY.velocity[j] = (aY.velocity[j] - bY.velocity[j]) / u;
31f0b65372SJed Brown   eY.temperature = (aY.temperature - bY.temperature) / sY.temperature;
322b916ea7SJeremy L Thompson   if (fabs(eY.pressure) > rtol_pressure) printf("%s: pressure error %g\n", name, eY.pressure);
332b916ea7SJeremy L Thompson   for (int j = 0; j < 3; j++) {
342b916ea7SJeremy L Thompson     if (fabs(eY.velocity[j]) > rtol_velocity) printf("%s: velocity[%d] error %g\n", name, j, eY.velocity[j]);
352b916ea7SJeremy L Thompson   }
362b916ea7SJeremy L Thompson   if (fabs(eY.temperature) > rtol_temperature) printf("%s: temperature error %g\n", name, eY.temperature);
37d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
38f0b65372SJed Brown }
39f0b65372SJed Brown 
402b916ea7SJeremy L Thompson static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) {
41f0b65372SJed Brown   Units            units = user->units;
42f0b65372SJed Brown   const CeedScalar eps   = 1e-6;
432b916ea7SJeremy L Thompson   const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, Pascal = units->Pascal;
44f0b65372SJed Brown   PetscFunctionBeginUser;
452b916ea7SJeremy L Thompson   const CeedScalar rho = 1.2 * kg / (m * m * m), u = 40 * m / sec;
46f0b65372SJed Brown   CeedScalar       U[5] = {rho, rho * u, rho * u * 1.1, rho * u * 1.2, 250e3 * Pascal + .5 * rho * u * u};
47edcfef1bSKenneth E. Jansen   State            s    = StateFromU(gas, U);
48f0b65372SJed Brown   for (int i = 0; i < 8; i++) {
49c4f7305dSKenneth E. Jansen     CeedScalar dU[5] = {0};
50f0b65372SJed Brown     if (i < 5) dU[i] = U[i];
51edcfef1bSKenneth E. Jansen     State ds = StateFromU_fwd(gas, s, dU);
52f0b65372SJed Brown     for (int j = 0; j < 5; j++) dU[j] = (1 + eps * (i == j)) * U[j];
53edcfef1bSKenneth E. Jansen     State          t = StateFromU(gas, dU);
54f0b65372SJed Brown     StatePrimitive dY;
55f0b65372SJed Brown     dY.pressure = (t.Y.pressure - s.Y.pressure) / eps;
562b916ea7SJeremy L Thompson     for (int j = 0; j < 3; j++) dY.velocity[j] = (t.Y.velocity[j] - s.Y.velocity[j]) / eps;
57f0b65372SJed Brown     dY.temperature = (t.Y.temperature - s.Y.temperature) / eps;
58f0b65372SJed Brown     char buf[128];
59f0b65372SJed Brown     snprintf(buf, sizeof buf, "StateFromU_fwd i=%d", i);
60f0b65372SJed Brown     PetscCall(CheckPrimitiveWithTolerance(dY, ds.Y, dY, buf, 5e-6, 1e-6, 1e-6));
61f0b65372SJed Brown   }
62d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
63f0b65372SJed Brown }
64f0b65372SJed Brown 
6565dee3d2SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping
6665dee3d2SJames Wright //
6765dee3d2SJames Wright // Only used for SUPG stabilization
6865dee3d2SJames Wright PetscErrorCode CreateKSPMassOperator_NewtonianStabilized(User user, CeedOperator *op_mass) {
6965dee3d2SJames Wright   Ceed                 ceed = user->ceed;
7065dee3d2SJames Wright   CeedInt              num_comp_q, q_data_size;
7165dee3d2SJames Wright   CeedQFunction        qf_mass;
7265dee3d2SJames Wright   CeedElemRestriction  elem_restr_q, elem_restr_qd_i;
7365dee3d2SJames Wright   CeedBasis            basis_q;
7465dee3d2SJames Wright   CeedVector           q_data;
7565dee3d2SJames Wright   CeedQFunctionContext qf_ctx = NULL;
7665dee3d2SJames Wright   PetscInt             dim    = 3;
7765dee3d2SJames Wright 
7865dee3d2SJames Wright   PetscFunctionBeginUser;
7965dee3d2SJames Wright   {  // Get restriction and basis from the RHS function
8065dee3d2SJames Wright     CeedOperator     *sub_ops;
8165dee3d2SJames Wright     CeedOperatorField field;
8265dee3d2SJames Wright     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
8365dee3d2SJames Wright 
8465dee3d2SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops));
8565dee3d2SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
8665dee3d2SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
8765dee3d2SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q));
8865dee3d2SJames Wright 
8965dee3d2SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field));
9065dee3d2SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i));
9165dee3d2SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data));
9265dee3d2SJames Wright 
9365dee3d2SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx));
9465dee3d2SJames Wright   }
9565dee3d2SJames Wright 
9665dee3d2SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
9765dee3d2SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size));
9865dee3d2SJames Wright 
9965dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Newtonian_Conserv, MassFunction_Newtonian_Conserv_loc, &qf_mass));
10065dee3d2SJames Wright 
10165dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx));
10265dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
10365dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
10465dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
10565dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
10665dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
10765dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
10865dee3d2SJames Wright 
10965dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
11065dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
11165dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed));
11265dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data));
11365dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
11465dee3d2SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
11565dee3d2SJames Wright 
11665dee3d2SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
11765dee3d2SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
11865dee3d2SJames Wright }
119*991aef52SJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
120b7f03f12SJed Brown   SetupContext             setup_context;
1213a8779fbSJames Wright   User                     user   = *(User *)ctx;
1224351d27fSLeila Ghaffari   CeedInt                  degree = user->app_ctx->degree;
1233a8779fbSJames Wright   StabilizationType        stab;
1243636f6a4SJames Wright   StateVariable            state_var;
125b4c37c5cSJames Wright   MPI_Comm                 comm = user->comm;
126b4c37c5cSJames Wright   Ceed                     ceed = user->ceed;
1273a8779fbSJames Wright   PetscBool                implicit;
1280d8cd818SJames Wright   PetscBool                unit_tests;
12915a3537eSJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
13015a3537eSJed Brown   CeedQFunctionContext     newtonian_ig_context;
1313a8779fbSJames Wright 
13215a3537eSJed Brown   PetscFunctionBeginUser;
1332b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &setup_context));
1342b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &newtonian_ig_ctx));
1353a8779fbSJames Wright 
1363a8779fbSJames Wright   // ------------------------------------------------------
1373a8779fbSJames Wright   //           Setup Generic Newtonian IG Problem
1383a8779fbSJames Wright   // ------------------------------------------------------
1393a8779fbSJames Wright   problem->dim                     = 3;
1403a8779fbSJames Wright   problem->q_data_size_vol         = 10;
141493642f1SJames Wright   problem->q_data_size_sur         = 10;
142b01ba163SJames Wright   problem->jac_data_size_sur       = 11;
1439785fe93SJed Brown   problem->setup_vol.qfunction     = Setup;
1449785fe93SJed Brown   problem->setup_vol.qfunction_loc = Setup_loc;
1459785fe93SJed Brown   problem->setup_sur.qfunction     = SetupBoundary;
1469785fe93SJed Brown   problem->setup_sur.qfunction_loc = SetupBoundary_loc;
1473a8779fbSJames Wright   problem->non_zero_time           = PETSC_FALSE;
148cbe60e31SLeila Ghaffari   problem->print_info              = PRINT_NEWTONIAN;
1497d8a615bSJames Wright   problem->uses_newtonian          = PETSC_TRUE;
1503a8779fbSJames Wright 
1513a8779fbSJames Wright   // ------------------------------------------------------
1523a8779fbSJames Wright   //             Create the libCEED context
1533a8779fbSJames Wright   // ------------------------------------------------------
1543a8779fbSJames Wright   CeedScalar cv         = 717.;          // J/(kg K)
1553a8779fbSJames Wright   CeedScalar cp         = 1004.;         // J/(kg K)
156d9bb1cdbSJames Wright   CeedScalar g[3]       = {0, 0, 0};     // m/s^2
1573a8779fbSJames Wright   CeedScalar lambda     = -2. / 3.;      // -
158bb8a0c61SJames Wright   CeedScalar mu         = 1.8e-5;        // Pa s, dynamic viscosity
1593a8779fbSJames Wright   CeedScalar k          = 0.02638;       // W/(m K)
1604351d27fSLeila Ghaffari   CeedScalar c_tau      = 0.5 / degree;  // -
161bb8a0c61SJames Wright   CeedScalar Ctau_t     = 1.0;           // -
162b5786772SLeila Ghaffari   CeedScalar Cv_func[3] = {36, 60, 128};
163c176988bSLeila Ghaffari   CeedScalar Ctau_v     = Cv_func[(CeedInt)Min(3, degree) - 1];
1642a0eee6eSLeila Ghaffari   CeedScalar Ctau_C     = 0.25 / degree;
1652a0eee6eSLeila Ghaffari   CeedScalar Ctau_M     = 0.25 / degree;
1662a0eee6eSLeila Ghaffari   CeedScalar Ctau_E     = 0.125;
1673a8779fbSJames Wright   PetscReal  domain_min[3], domain_max[3], domain_size[3];
1682b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
169493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
1703a8779fbSJames Wright 
171b8fb7609SAdeleke O. Bankole   StatePrimitive reference      = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15};
172e7754af5SKenneth E. Jansen   CeedScalar     idl_decay_time = -1, idl_start = 0, idl_length = 0;
173e7754af5SKenneth E. Jansen   PetscBool      idl_enable = PETSC_FALSE;
174b8fb7609SAdeleke O. Bankole 
1753a8779fbSJames Wright   // ------------------------------------------------------
1763a8779fbSJames Wright   //             Create the PETSc context
1773a8779fbSJames Wright   // ------------------------------------------------------
178bb8a0c61SJames Wright   PetscScalar meter    = 1;  // 1 meter in scaled length units
179bb8a0c61SJames Wright   PetscScalar kilogram = 1;  // 1 kilogram in scaled mass units
180bb8a0c61SJames Wright   PetscScalar second   = 1;  // 1 second in scaled time units
1813a8779fbSJames Wright   PetscScalar Kelvin   = 1;  // 1 Kelvin in scaled temperature units
1823a8779fbSJames Wright   PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s;
1833a8779fbSJames Wright 
1843a8779fbSJames Wright   // ------------------------------------------------------
1853a8779fbSJames Wright   //              Command line Options
1863a8779fbSJames Wright   // ------------------------------------------------------
187d9bb1cdbSJames Wright   PetscBool given_option = PETSC_FALSE;
1882b916ea7SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL);
189cbe60e31SLeila Ghaffari   // -- Conservative vs Primitive variables
1902b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE),
1912b916ea7SJeremy L Thompson                              (PetscEnum *)&state_var, NULL));
1923636f6a4SJames Wright 
1933636f6a4SJames Wright   switch (state_var) {
1943636f6a4SJames Wright     case STATEVAR_CONSERVATIVE:
195b8fb7609SAdeleke O. Bankole       problem->ics.qfunction                       = ICsNewtonianIG_Conserv;
196b8fb7609SAdeleke O. Bankole       problem->ics.qfunction_loc                   = ICsNewtonianIG_Conserv_loc;
197cbe60e31SLeila Ghaffari       problem->apply_vol_rhs.qfunction             = RHSFunction_Newtonian;
198cbe60e31SLeila Ghaffari       problem->apply_vol_rhs.qfunction_loc         = RHSFunction_Newtonian_loc;
19976555becSJames Wright       problem->apply_vol_ifunction.qfunction       = IFunction_Newtonian_Conserv;
20076555becSJames Wright       problem->apply_vol_ifunction.qfunction_loc   = IFunction_Newtonian_Conserv_loc;
20176555becSJames Wright       problem->apply_vol_ijacobian.qfunction       = IJacobian_Newtonian_Conserv;
20276555becSJames Wright       problem->apply_vol_ijacobian.qfunction_loc   = IJacobian_Newtonian_Conserv_loc;
203d4559bbeSJames Wright       problem->apply_inflow.qfunction              = BoundaryIntegral_Conserv;
204d4559bbeSJames Wright       problem->apply_inflow.qfunction_loc          = BoundaryIntegral_Conserv_loc;
205d4559bbeSJames Wright       problem->apply_inflow_jacobian.qfunction     = BoundaryIntegral_Jacobian_Conserv;
206d4559bbeSJames Wright       problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc;
2073636f6a4SJames Wright       break;
2083636f6a4SJames Wright 
2093636f6a4SJames Wright     case STATEVAR_PRIMITIVE:
2103636f6a4SJames Wright       problem->ics.qfunction                       = ICsNewtonianIG_Prim;
2113636f6a4SJames Wright       problem->ics.qfunction_loc                   = ICsNewtonianIG_Prim_loc;
2123636f6a4SJames Wright       problem->apply_vol_ifunction.qfunction       = IFunction_Newtonian_Prim;
2133636f6a4SJames Wright       problem->apply_vol_ifunction.qfunction_loc   = IFunction_Newtonian_Prim_loc;
2143636f6a4SJames Wright       problem->apply_vol_ijacobian.qfunction       = IJacobian_Newtonian_Prim;
2153636f6a4SJames Wright       problem->apply_vol_ijacobian.qfunction_loc   = IJacobian_Newtonian_Prim_loc;
2163636f6a4SJames Wright       problem->apply_inflow.qfunction              = BoundaryIntegral_Prim;
2173636f6a4SJames Wright       problem->apply_inflow.qfunction_loc          = BoundaryIntegral_Prim_loc;
2183636f6a4SJames Wright       problem->apply_inflow_jacobian.qfunction     = BoundaryIntegral_Jacobian_Prim;
2193636f6a4SJames Wright       problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc;
2203636f6a4SJames Wright       break;
221cbe60e31SLeila Ghaffari   }
2221485969bSJeremy L Thompson 
2233a8779fbSJames Wright   // -- Physics
2242b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL));
2252b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL));
2262b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL));
2272b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL));
2282b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL));
2293a8779fbSJames Wright 
230bb8a0c61SJames Wright   PetscInt dim = problem->dim;
231d9bb1cdbSJames Wright   PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL));
232d9bb1cdbSJames Wright   PetscCall(PetscOptionsRealArray("-g", "Gravitational acceleration vector", NULL, g, &dim, &given_option));
233d9bb1cdbSJames Wright   dim = problem->dim;
234d9bb1cdbSJames Wright   PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option));
235d9bb1cdbSJames Wright   if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim);
236d9bb1cdbSJames Wright 
2372b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
2382b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
2392b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
2402b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL));
2412b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL));
2422b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL));
2432b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL));
2442b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
2452b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL));
2463a8779fbSJames Wright 
247db95602eSJames Wright   dim = 3;
248b8fb7609SAdeleke O. Bankole   PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL));
249b8fb7609SAdeleke O. Bankole   PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL));
250b8fb7609SAdeleke O. Bankole   PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL));
251b8fb7609SAdeleke O. Bankole 
2523a8779fbSJames Wright   // -- Units
2532b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
2543a8779fbSJames Wright   meter = fabs(meter);
2552b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL));
2563a8779fbSJames Wright   kilogram = fabs(kilogram);
2572b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
2583a8779fbSJames Wright   second = fabs(second);
2592b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL));
2603a8779fbSJames Wright   Kelvin = fabs(Kelvin);
2613a8779fbSJames Wright 
2623a8779fbSJames Wright   // -- Warnings
2635d28dccaSJames Wright   PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP,
2645d28dccaSJames Wright              "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n");
265e7754af5SKenneth E. Jansen 
266e7754af5SKenneth E. Jansen   PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point",
267e7754af5SKenneth E. Jansen                                NULL, idl_decay_time, &idl_decay_time, &idl_enable));
2689cbdf780SJames Wright   PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero.");
2699cbdf780SJames Wright   if (idl_decay_time < 0) idl_enable = PETSC_FALSE;
270e7754af5SKenneth E. Jansen   PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL));
271e7754af5SKenneth E. Jansen   PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL));
2721485969bSJeremy L Thompson   PetscOptionsEnd();
2733a8779fbSJames Wright 
27465dee3d2SJames Wright   if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized;
27565dee3d2SJames Wright 
2763a8779fbSJames Wright   // ------------------------------------------------------
2773a8779fbSJames Wright   //           Set up the PETSc context
2783a8779fbSJames Wright   // ------------------------------------------------------
2793a8779fbSJames Wright   // -- Define derived units
2803a8779fbSJames Wright   Pascal          = kilogram / (meter * PetscSqr(second));
2813a8779fbSJames Wright   J_per_kg_K      = PetscSqr(meter) / (PetscSqr(second) * Kelvin);
2823a8779fbSJames Wright   m_per_squared_s = meter / PetscSqr(second);
2833a8779fbSJames Wright   W_per_m_K       = kilogram * meter / (pow(second, 3) * Kelvin);
2843a8779fbSJames Wright 
2853a8779fbSJames Wright   user->units->meter           = meter;
2863a8779fbSJames Wright   user->units->kilogram        = kilogram;
2873a8779fbSJames Wright   user->units->second          = second;
2883a8779fbSJames Wright   user->units->Kelvin          = Kelvin;
2893a8779fbSJames Wright   user->units->Pascal          = Pascal;
2903a8779fbSJames Wright   user->units->J_per_kg_K      = J_per_kg_K;
2913a8779fbSJames Wright   user->units->m_per_squared_s = m_per_squared_s;
2923a8779fbSJames Wright   user->units->W_per_m_K       = W_per_m_K;
2933a8779fbSJames Wright 
2943a8779fbSJames Wright   // ------------------------------------------------------
2953a8779fbSJames Wright   //           Set up the libCEED context
2963a8779fbSJames Wright   // ------------------------------------------------------
2973a8779fbSJames Wright   // -- Scale variables to desired units
2983a8779fbSJames Wright   cv *= J_per_kg_K;
2993a8779fbSJames Wright   cp *= J_per_kg_K;
3003a8779fbSJames Wright   mu *= Pascal * second;
3013a8779fbSJames Wright   k *= W_per_m_K;
302493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter;
303493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s;
304b8fb7609SAdeleke O. Bankole   reference.pressure *= Pascal;
305b8fb7609SAdeleke O. Bankole   for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second;
306b8fb7609SAdeleke O. Bankole   reference.temperature *= Kelvin;
3073a8779fbSJames Wright   problem->dm_scale = meter;
3083a8779fbSJames Wright 
3093a8779fbSJames Wright   // -- Solver Settings
3103a8779fbSJames Wright   user->phys->implicit  = implicit;
3113636f6a4SJames Wright   user->phys->state_var = state_var;
3123a8779fbSJames Wright 
3133a8779fbSJames Wright   // -- QFunction Context
31415a3537eSJed Brown   newtonian_ig_ctx->lambda        = lambda;
31515a3537eSJed Brown   newtonian_ig_ctx->mu            = mu;
31615a3537eSJed Brown   newtonian_ig_ctx->k             = k;
31715a3537eSJed Brown   newtonian_ig_ctx->cv            = cv;
31815a3537eSJed Brown   newtonian_ig_ctx->cp            = cp;
31915a3537eSJed Brown   newtonian_ig_ctx->c_tau         = c_tau;
32015a3537eSJed Brown   newtonian_ig_ctx->Ctau_t        = Ctau_t;
32115a3537eSJed Brown   newtonian_ig_ctx->Ctau_v        = Ctau_v;
32215a3537eSJed Brown   newtonian_ig_ctx->Ctau_C        = Ctau_C;
32315a3537eSJed Brown   newtonian_ig_ctx->Ctau_M        = Ctau_M;
32415a3537eSJed Brown   newtonian_ig_ctx->Ctau_E        = Ctau_E;
325e7754af5SKenneth E. Jansen   newtonian_ig_ctx->P0            = reference.pressure;
32615a3537eSJed Brown   newtonian_ig_ctx->stabilization = stab;
327c8c30d87SJed Brown   newtonian_ig_ctx->P0            = reference.pressure;
3288085925cSJames Wright   newtonian_ig_ctx->is_implicit   = implicit;
3293636f6a4SJames Wright   newtonian_ig_ctx->state_var     = state_var;
330e7754af5SKenneth E. Jansen   newtonian_ig_ctx->idl_enable    = idl_enable;
331e7754af5SKenneth E. Jansen   newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * second);
332e7754af5SKenneth E. Jansen   newtonian_ig_ctx->idl_start     = idl_start * meter;
333e7754af5SKenneth E. Jansen   newtonian_ig_ctx->idl_length    = idl_length * meter;
3342b916ea7SJeremy L Thompson   PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3));
3353a8779fbSJames Wright 
336b8fb7609SAdeleke O. Bankole   // -- Setup Context
337b8fb7609SAdeleke O. Bankole   setup_context->reference = reference;
338b8fb7609SAdeleke O. Bankole   setup_context->gas       = *newtonian_ig_ctx;
339b8fb7609SAdeleke O. Bankole   setup_context->lx        = domain_size[0];
340b8fb7609SAdeleke O. Bankole   setup_context->ly        = domain_size[1];
341b8fb7609SAdeleke O. Bankole   setup_context->lz        = domain_size[2];
342b8fb7609SAdeleke O. Bankole   setup_context->time      = 0;
343b8fb7609SAdeleke O. Bankole 
344b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context));
345b4c37c5cSJames Wright   PetscCallCeed(ceed,
346b4c37c5cSJames Wright                 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
347b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc));
348b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, "evaluation time", offsetof(struct SetupContext_, time), 1,
349b4c37c5cSJames Wright                                                          "Time of evaluation"));
35015a3537eSJed Brown 
351b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context));
352b4c37c5cSJames Wright   PetscCallCeed(ceed,
353b4c37c5cSJames Wright                 CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx));
354b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, FreeContextPetsc));
355b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1,
356b4c37c5cSJames Wright                                                          "Size of timestep, delta t"));
357b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift",
358b4c37c5cSJames Wright                                                          offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1,
359b4c37c5cSJames Wright                                                          "Shift for mass matrix in IJacobian"));
360b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1,
361b4c37c5cSJames Wright                                                          "Current solution time"));
362a8bca8acSJames Wright 
36315a3537eSJed Brown   problem->apply_vol_rhs.qfunction_context = newtonian_ig_context;
364b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ifunction.qfunction_context));
365b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ijacobian.qfunction_context));
366b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow.qfunction_context));
367b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow_jacobian.qfunction_context));
368f0b65372SJed Brown 
3699ed3d70dSJames Wright   if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
3709ed3d70dSJames Wright   if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
3719ed3d70dSJames Wright   if (bc->num_slip > 0) PetscCall(SlipBCSetup(problem, dm, ctx, newtonian_ig_context));
3729ed3d70dSJames Wright 
373f0b65372SJed Brown   if (unit_tests) {
374f0b65372SJed Brown     PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx));
375f0b65372SJed Brown   }
376d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3773a8779fbSJames Wright }
3783a8779fbSJames Wright 
379*991aef52SJames Wright PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx) {
3802d49c0afSJames Wright   MPI_Comm                 comm = user->comm;
381b4c37c5cSJames Wright   Ceed                     ceed = user->ceed;
38215a3537eSJed Brown   NewtonianIdealGasContext newtonian_ctx;
38315a3537eSJed Brown 
3843a8779fbSJames Wright   PetscFunctionBeginUser;
385b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ctx));
3862b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(comm,
38715a3537eSJed Brown                         "  Problem:\n"
38815a3537eSJed Brown                         "    Problem Name                       : %s\n"
38915a3537eSJed Brown                         "    Stabilization                      : %s\n",
3902b916ea7SJeremy L Thompson                         app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization]));
391b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ctx));
392d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3933a8779fbSJames Wright }
394