xref: /libCEED/examples/fluids/problems/newtonian.c (revision 013a555157cefb861605c4e4ea6ce4e03bb7c3d7)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
388b783a1SJames Wright //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
588b783a1SJames Wright //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
788b783a1SJames Wright 
888b783a1SJames Wright /// @file
988b783a1SJames Wright /// Utility functions for setting up problems using the Newtonian Qfunction
1088b783a1SJames Wright 
1188b783a1SJames Wright #include "../qfunctions/newtonian.h"
1288b783a1SJames Wright 
1349aac155SJeremy L Thompson #include <ceed.h>
1449aac155SJeremy L Thompson #include <petscdm.h>
1549aac155SJeremy L Thompson 
162b730f8bSJeremy L Thompson #include "../navierstokes.h"
176ef2784eSLeila Ghaffari 
18d64246eaSJed Brown // For use with PetscOptionsEnum
19d64246eaSJed Brown static const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "StateVariable", "STATEVAR_", NULL};
20d64246eaSJed Brown 
212b730f8bSJeremy L Thompson // Compute relative error |a - b|/|s|
222b730f8bSJeremy L Thompson static PetscErrorCode CheckPrimitiveWithTolerance(StatePrimitive sY, StatePrimitive aY, StatePrimitive bY, const char *name, PetscReal rtol_pressure,
232b730f8bSJeremy L Thompson                                                   PetscReal rtol_velocity, PetscReal rtol_temperature) {
24e334ad8fSJed Brown   StatePrimitive eY;  // relative error
25f17d818dSJames Wright 
26f17d818dSJames Wright   PetscFunctionBeginUser;
27e334ad8fSJed Brown   eY.pressure   = (aY.pressure - bY.pressure) / sY.pressure;
282b730f8bSJeremy L Thompson   PetscScalar u = sqrt(Square(sY.velocity[0]) + Square(sY.velocity[1]) + Square(sY.velocity[2]));
29e334ad8fSJed Brown   for (int j = 0; j < 3; j++) eY.velocity[j] = (aY.velocity[j] - bY.velocity[j]) / u;
30e334ad8fSJed Brown   eY.temperature = (aY.temperature - bY.temperature) / sY.temperature;
312b730f8bSJeremy L Thompson   if (fabs(eY.pressure) > rtol_pressure) printf("%s: pressure error %g\n", name, eY.pressure);
322b730f8bSJeremy L Thompson   for (int j = 0; j < 3; j++) {
332b730f8bSJeremy L Thompson     if (fabs(eY.velocity[j]) > rtol_velocity) printf("%s: velocity[%d] error %g\n", name, j, eY.velocity[j]);
342b730f8bSJeremy L Thompson   }
352b730f8bSJeremy L Thompson   if (fabs(eY.temperature) > rtol_temperature) printf("%s: temperature error %g\n", name, eY.temperature);
36ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
37e334ad8fSJed Brown }
38e334ad8fSJed Brown 
392b730f8bSJeremy L Thompson static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) {
40e334ad8fSJed Brown   Units            units = user->units;
41e334ad8fSJed Brown   const CeedScalar eps   = 1e-6;
422b730f8bSJeremy L Thompson   const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, Pascal = units->Pascal;
43e334ad8fSJed Brown   PetscFunctionBeginUser;
442b730f8bSJeremy L Thompson   const CeedScalar rho = 1.2 * kg / (m * m * m), u = 40 * m / sec;
45e334ad8fSJed Brown   CeedScalar       U[5] = {rho, rho * u, rho * u * 1.1, rho * u * 1.2, 250e3 * Pascal + .5 * rho * u * u};
463bd61617SKenneth E. Jansen   State            s    = StateFromU(gas, U);
47e334ad8fSJed Brown   for (int i = 0; i < 8; i++) {
48242a7fdaSKenneth E. Jansen     CeedScalar dU[5] = {0};
49e334ad8fSJed Brown     if (i < 5) dU[i] = U[i];
503bd61617SKenneth E. Jansen     State ds = StateFromU_fwd(gas, s, dU);
51e334ad8fSJed Brown     for (int j = 0; j < 5; j++) dU[j] = (1 + eps * (i == j)) * U[j];
523bd61617SKenneth E. Jansen     State          t = StateFromU(gas, dU);
53e334ad8fSJed Brown     StatePrimitive dY;
54e334ad8fSJed Brown     dY.pressure = (t.Y.pressure - s.Y.pressure) / eps;
552b730f8bSJeremy L Thompson     for (int j = 0; j < 3; j++) dY.velocity[j] = (t.Y.velocity[j] - s.Y.velocity[j]) / eps;
56e334ad8fSJed Brown     dY.temperature = (t.Y.temperature - s.Y.temperature) / eps;
57e334ad8fSJed Brown     char buf[128];
58e334ad8fSJed Brown     snprintf(buf, sizeof buf, "StateFromU_fwd i=%d", i);
59e334ad8fSJed Brown     PetscCall(CheckPrimitiveWithTolerance(dY, ds.Y, dY, buf, 5e-6, 1e-6, 1e-6));
60e334ad8fSJed Brown   }
61ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
62e334ad8fSJed Brown }
63e334ad8fSJed Brown 
640fcbc436SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping
650fcbc436SJames Wright //
660fcbc436SJames Wright // Only used for SUPG stabilization
670fcbc436SJames Wright PetscErrorCode CreateKSPMassOperator_NewtonianStabilized(User user, CeedOperator *op_mass) {
680fcbc436SJames Wright   Ceed                 ceed = user->ceed;
690fcbc436SJames Wright   CeedInt              num_comp_q, q_data_size;
700fcbc436SJames Wright   CeedQFunction        qf_mass;
710fcbc436SJames Wright   CeedElemRestriction  elem_restr_q, elem_restr_qd_i;
720fcbc436SJames Wright   CeedBasis            basis_q;
730fcbc436SJames Wright   CeedVector           q_data;
740fcbc436SJames Wright   CeedQFunctionContext qf_ctx = NULL;
750fcbc436SJames Wright   PetscInt             dim    = 3;
760fcbc436SJames Wright 
770fcbc436SJames Wright   PetscFunctionBeginUser;
780fcbc436SJames Wright   {  // Get restriction and basis from the RHS function
790fcbc436SJames Wright     CeedOperator     *sub_ops;
800fcbc436SJames Wright     CeedOperatorField field;
810fcbc436SJames Wright     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
820fcbc436SJames Wright 
830fcbc436SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops));
840fcbc436SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
850fcbc436SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
860fcbc436SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q));
870fcbc436SJames Wright 
880fcbc436SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field));
890fcbc436SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i));
900fcbc436SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data));
910fcbc436SJames Wright 
920fcbc436SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx));
930fcbc436SJames Wright   }
940fcbc436SJames Wright 
950fcbc436SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
960fcbc436SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size));
970fcbc436SJames Wright 
980fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Newtonian_Conserv, MassFunction_Newtonian_Conserv_loc, &qf_mass));
990fcbc436SJames Wright 
1000fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx));
1010fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
1020fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
1030fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
1040fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
1050fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
1060fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
1070fcbc436SJames Wright 
1080fcbc436SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
1090fcbc436SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
1100fcbc436SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed));
1110fcbc436SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data));
1120fcbc436SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
1130fcbc436SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
1140fcbc436SJames Wright 
1150fcbc436SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
1160fcbc436SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1170fcbc436SJames Wright }
118*013a5551SJames Wright 
119731c13d7SJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
120a0add3c9SJed Brown   SetupContext             setup_context;
12188b783a1SJames Wright   User                     user   = *(User *)ctx;
122ebd2ea64SLeila Ghaffari   CeedInt                  degree = user->app_ctx->degree;
12388b783a1SJames Wright   StabilizationType        stab;
12497baf651SJames Wright   StateVariable            state_var;
125a424bcd0SJames Wright   MPI_Comm                 comm = user->comm;
126a424bcd0SJames Wright   Ceed                     ceed = user->ceed;
12788b783a1SJames Wright   PetscBool                implicit;
128a2f6637eSJames Wright   PetscBool                unit_tests;
129841e4c73SJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
130841e4c73SJed Brown   CeedQFunctionContext     newtonian_ig_context;
13188b783a1SJames Wright 
132841e4c73SJed Brown   PetscFunctionBeginUser;
1332b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &setup_context));
1342b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &newtonian_ig_ctx));
13588b783a1SJames Wright 
13688b783a1SJames Wright   // ------------------------------------------------------
13788b783a1SJames Wright   //           Setup Generic Newtonian IG Problem
13888b783a1SJames Wright   // ------------------------------------------------------
13988b783a1SJames Wright   problem->dim               = 3;
1400ec2498eSJames Wright   problem->jac_data_size_sur = 11;
14188b783a1SJames Wright   problem->non_zero_time     = PETSC_FALSE;
142dc805cc4SLeila Ghaffari   problem->print_info        = PRINT_NEWTONIAN;
143752a08a7SJames Wright   problem->uses_newtonian    = PETSC_TRUE;
14488b783a1SJames Wright 
14588b783a1SJames Wright   // ------------------------------------------------------
14688b783a1SJames Wright   //             Create the libCEED context
14788b783a1SJames Wright   // ------------------------------------------------------
14888b783a1SJames Wright   CeedScalar cv         = 717.;          // J/(kg K)
14988b783a1SJames Wright   CeedScalar cp         = 1004.;         // J/(kg K)
150a2726bdbSJames Wright   CeedScalar g[3]       = {0, 0, 0};     // m/s^2
15188b783a1SJames Wright   CeedScalar lambda     = -2. / 3.;      // -
15288626eedSJames Wright   CeedScalar mu         = 1.8e-5;        // Pa s, dynamic viscosity
15388b783a1SJames Wright   CeedScalar k          = 0.02638;       // W/(m K)
154ebd2ea64SLeila Ghaffari   CeedScalar c_tau      = 0.5 / degree;  // -
15588626eedSJames Wright   CeedScalar Ctau_t     = 1.0;           // -
15694c01735SLeila Ghaffari   CeedScalar Cv_func[3] = {36, 60, 128};
157c5dde687SLeila Ghaffari   CeedScalar Ctau_v     = Cv_func[(CeedInt)Min(3, degree) - 1];
15875538696SLeila Ghaffari   CeedScalar Ctau_C     = 0.25 / degree;
15975538696SLeila Ghaffari   CeedScalar Ctau_M     = 0.25 / degree;
16075538696SLeila Ghaffari   CeedScalar Ctau_E     = 0.125;
16188b783a1SJames Wright   PetscReal  domain_min[3], domain_max[3], domain_size[3];
1622b730f8bSJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
163ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
16488b783a1SJames Wright 
165d310b3d3SAdeleke O. Bankole   StatePrimitive reference      = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15};
1662249ac91SJames Wright   CeedScalar     idl_decay_time = -1, idl_start = 0, idl_length = 0, idl_pressure = reference.pressure;
167530ad8c4SKenneth E. Jansen   PetscBool      idl_enable = PETSC_FALSE;
168d310b3d3SAdeleke O. Bankole 
16988b783a1SJames Wright   // ------------------------------------------------------
17088b783a1SJames Wright   //             Create the PETSc context
17188b783a1SJames Wright   // ------------------------------------------------------
17288626eedSJames Wright   PetscScalar meter    = 1;  // 1 meter in scaled length units
17388626eedSJames Wright   PetscScalar kilogram = 1;  // 1 kilogram in scaled mass units
17488626eedSJames Wright   PetscScalar second   = 1;  // 1 second in scaled time units
17588b783a1SJames Wright   PetscScalar Kelvin   = 1;  // 1 Kelvin in scaled temperature units
17688b783a1SJames Wright   PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s;
17788b783a1SJames Wright 
17888b783a1SJames Wright   // ------------------------------------------------------
17988b783a1SJames Wright   //              Command line Options
18088b783a1SJames Wright   // ------------------------------------------------------
181a2726bdbSJames Wright   PetscBool given_option = PETSC_FALSE;
1822b730f8bSJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL);
183dc805cc4SLeila Ghaffari   // -- Conservative vs Primitive variables
1842b730f8bSJeremy L Thompson   PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE),
1852b730f8bSJeremy L Thompson                              (PetscEnum *)&state_var, NULL));
18697baf651SJames Wright 
18797baf651SJames Wright   switch (state_var) {
18897baf651SJames Wright     case STATEVAR_CONSERVATIVE:
189d310b3d3SAdeleke O. Bankole       problem->ics.qfunction                       = ICsNewtonianIG_Conserv;
190d310b3d3SAdeleke O. Bankole       problem->ics.qfunction_loc                   = ICsNewtonianIG_Conserv_loc;
191dc805cc4SLeila Ghaffari       problem->apply_vol_rhs.qfunction             = RHSFunction_Newtonian;
192dc805cc4SLeila Ghaffari       problem->apply_vol_rhs.qfunction_loc         = RHSFunction_Newtonian_loc;
1933d02368aSJames Wright       problem->apply_vol_ifunction.qfunction       = IFunction_Newtonian_Conserv;
1943d02368aSJames Wright       problem->apply_vol_ifunction.qfunction_loc   = IFunction_Newtonian_Conserv_loc;
1953d02368aSJames Wright       problem->apply_vol_ijacobian.qfunction       = IJacobian_Newtonian_Conserv;
1963d02368aSJames Wright       problem->apply_vol_ijacobian.qfunction_loc   = IJacobian_Newtonian_Conserv_loc;
19720840d50SJames Wright       problem->apply_inflow.qfunction              = BoundaryIntegral_Conserv;
19820840d50SJames Wright       problem->apply_inflow.qfunction_loc          = BoundaryIntegral_Conserv_loc;
19920840d50SJames Wright       problem->apply_inflow_jacobian.qfunction     = BoundaryIntegral_Jacobian_Conserv;
20020840d50SJames Wright       problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc;
20197baf651SJames Wright       break;
20297baf651SJames Wright 
20397baf651SJames Wright     case STATEVAR_PRIMITIVE:
20497baf651SJames Wright       problem->ics.qfunction                       = ICsNewtonianIG_Prim;
20597baf651SJames Wright       problem->ics.qfunction_loc                   = ICsNewtonianIG_Prim_loc;
20697baf651SJames Wright       problem->apply_vol_ifunction.qfunction       = IFunction_Newtonian_Prim;
20797baf651SJames Wright       problem->apply_vol_ifunction.qfunction_loc   = IFunction_Newtonian_Prim_loc;
20897baf651SJames Wright       problem->apply_vol_ijacobian.qfunction       = IJacobian_Newtonian_Prim;
20997baf651SJames Wright       problem->apply_vol_ijacobian.qfunction_loc   = IJacobian_Newtonian_Prim_loc;
21097baf651SJames Wright       problem->apply_inflow.qfunction              = BoundaryIntegral_Prim;
21197baf651SJames Wright       problem->apply_inflow.qfunction_loc          = BoundaryIntegral_Prim_loc;
21297baf651SJames Wright       problem->apply_inflow_jacobian.qfunction     = BoundaryIntegral_Jacobian_Prim;
21397baf651SJames Wright       problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc;
21497baf651SJames Wright       break;
215dc805cc4SLeila Ghaffari   }
21667490bc6SJeremy L Thompson 
21788b783a1SJames Wright   // -- Physics
2182b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL));
2192b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL));
2202b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL));
2212b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL));
2222b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL));
22388b783a1SJames Wright 
22488626eedSJames Wright   PetscInt dim = problem->dim;
225a2726bdbSJames Wright   PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL));
226a2726bdbSJames Wright   PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option));
227a2726bdbSJames Wright   if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim);
228a2726bdbSJames Wright 
2292b730f8bSJeremy L Thompson   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
2302b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
2312b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
2322b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL));
2332b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL));
2342b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL));
2352b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL));
2362b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
2372b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL));
23888b783a1SJames Wright 
239e30587cbSJames Wright   dim = 3;
240d310b3d3SAdeleke O. Bankole   PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL));
241d310b3d3SAdeleke O. Bankole   PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL));
242d310b3d3SAdeleke O. Bankole   PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL));
243d310b3d3SAdeleke O. Bankole 
24488b783a1SJames Wright   // -- Units
2452b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
24688b783a1SJames Wright   meter = fabs(meter);
2472b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL));
24888b783a1SJames Wright   kilogram = fabs(kilogram);
2492b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
25088b783a1SJames Wright   second = fabs(second);
2512b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL));
25288b783a1SJames Wright   Kelvin = fabs(Kelvin);
25388b783a1SJames Wright 
25488b783a1SJames Wright   // -- Warnings
2550e654f56SJames Wright   PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP,
2560e654f56SJames Wright              "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n");
257530ad8c4SKenneth E. Jansen 
258530ad8c4SKenneth E. Jansen   PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point",
259530ad8c4SKenneth E. Jansen                                NULL, idl_decay_time, &idl_decay_time, &idl_enable));
260fe1e732eSJames Wright   PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero.");
261fe1e732eSJames Wright   if (idl_decay_time < 0) idl_enable = PETSC_FALSE;
262530ad8c4SKenneth E. Jansen   PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL));
263530ad8c4SKenneth E. Jansen   PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL));
2642249ac91SJames Wright   idl_pressure = reference.pressure;
2652249ac91SJames Wright   PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-reference_pressure`)", NULL, idl_pressure,
2662249ac91SJames Wright                                &idl_pressure, NULL));
26767490bc6SJeremy L Thompson   PetscOptionsEnd();
26888b783a1SJames Wright 
2690fcbc436SJames Wright   if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized;
2700fcbc436SJames Wright 
27188b783a1SJames Wright   // ------------------------------------------------------
27288b783a1SJames Wright   //           Set up the PETSc context
27388b783a1SJames Wright   // ------------------------------------------------------
27488b783a1SJames Wright   // -- Define derived units
27588b783a1SJames Wright   Pascal          = kilogram / (meter * PetscSqr(second));
27688b783a1SJames Wright   J_per_kg_K      = PetscSqr(meter) / (PetscSqr(second) * Kelvin);
27788b783a1SJames Wright   m_per_squared_s = meter / PetscSqr(second);
27888b783a1SJames Wright   W_per_m_K       = kilogram * meter / (pow(second, 3) * Kelvin);
27988b783a1SJames Wright 
28088b783a1SJames Wright   user->units->meter           = meter;
28188b783a1SJames Wright   user->units->kilogram        = kilogram;
28288b783a1SJames Wright   user->units->second          = second;
28388b783a1SJames Wright   user->units->Kelvin          = Kelvin;
28488b783a1SJames Wright   user->units->Pascal          = Pascal;
28588b783a1SJames Wright   user->units->J_per_kg_K      = J_per_kg_K;
28688b783a1SJames Wright   user->units->m_per_squared_s = m_per_squared_s;
28788b783a1SJames Wright   user->units->W_per_m_K       = W_per_m_K;
28888b783a1SJames Wright 
28988b783a1SJames Wright   // ------------------------------------------------------
29088b783a1SJames Wright   //           Set up the libCEED context
29188b783a1SJames Wright   // ------------------------------------------------------
29288b783a1SJames Wright   // -- Scale variables to desired units
29388b783a1SJames Wright   cv *= J_per_kg_K;
29488b783a1SJames Wright   cp *= J_per_kg_K;
29588b783a1SJames Wright   mu *= Pascal * second;
29688b783a1SJames Wright   k *= W_per_m_K;
297ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter;
298ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s;
299d310b3d3SAdeleke O. Bankole   reference.pressure *= Pascal;
300d310b3d3SAdeleke O. Bankole   for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second;
301d310b3d3SAdeleke O. Bankole   reference.temperature *= Kelvin;
30288b783a1SJames Wright   problem->dm_scale = meter;
30388b783a1SJames Wright 
30488b783a1SJames Wright   // -- Solver Settings
30588b783a1SJames Wright   user->phys->implicit  = implicit;
30697baf651SJames Wright   user->phys->state_var = state_var;
30788b783a1SJames Wright 
30888b783a1SJames Wright   // -- QFunction Context
309841e4c73SJed Brown   newtonian_ig_ctx->lambda        = lambda;
310841e4c73SJed Brown   newtonian_ig_ctx->mu            = mu;
311841e4c73SJed Brown   newtonian_ig_ctx->k             = k;
312841e4c73SJed Brown   newtonian_ig_ctx->cv            = cv;
313841e4c73SJed Brown   newtonian_ig_ctx->cp            = cp;
314841e4c73SJed Brown   newtonian_ig_ctx->c_tau         = c_tau;
315841e4c73SJed Brown   newtonian_ig_ctx->Ctau_t        = Ctau_t;
316841e4c73SJed Brown   newtonian_ig_ctx->Ctau_v        = Ctau_v;
317841e4c73SJed Brown   newtonian_ig_ctx->Ctau_C        = Ctau_C;
318841e4c73SJed Brown   newtonian_ig_ctx->Ctau_M        = Ctau_M;
319841e4c73SJed Brown   newtonian_ig_ctx->Ctau_E        = Ctau_E;
320841e4c73SJed Brown   newtonian_ig_ctx->stabilization = stab;
32165dd5cafSJames Wright   newtonian_ig_ctx->is_implicit   = implicit;
32297baf651SJames Wright   newtonian_ig_ctx->state_var     = state_var;
323530ad8c4SKenneth E. Jansen   newtonian_ig_ctx->idl_enable    = idl_enable;
324530ad8c4SKenneth E. Jansen   newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * second);
325530ad8c4SKenneth E. Jansen   newtonian_ig_ctx->idl_start     = idl_start * meter;
326530ad8c4SKenneth E. Jansen   newtonian_ig_ctx->idl_length    = idl_length * meter;
3272249ac91SJames Wright   newtonian_ig_ctx->idl_pressure  = idl_pressure;
3282b730f8bSJeremy L Thompson   PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3));
32988b783a1SJames Wright 
330d310b3d3SAdeleke O. Bankole   // -- Setup Context
331d310b3d3SAdeleke O. Bankole   setup_context->reference = reference;
332d310b3d3SAdeleke O. Bankole   setup_context->gas       = *newtonian_ig_ctx;
333d310b3d3SAdeleke O. Bankole   setup_context->lx        = domain_size[0];
334d310b3d3SAdeleke O. Bankole   setup_context->ly        = domain_size[1];
335d310b3d3SAdeleke O. Bankole   setup_context->lz        = domain_size[2];
336d310b3d3SAdeleke O. Bankole   setup_context->time      = 0;
337d310b3d3SAdeleke O. Bankole 
338a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context));
339a424bcd0SJames Wright   PetscCallCeed(ceed,
340a424bcd0SJames Wright                 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
341a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc));
342a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, "evaluation time", offsetof(struct SetupContext_, time), 1,
343a424bcd0SJames Wright                                                          "Time of evaluation"));
344841e4c73SJed Brown 
345a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context));
346a424bcd0SJames Wright   PetscCallCeed(ceed,
347a424bcd0SJames Wright                 CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx));
348a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, FreeContextPetsc));
349a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1,
350a424bcd0SJames Wright                                                          "Size of timestep, delta t"));
351a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift",
352a424bcd0SJames Wright                                                          offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1,
353a424bcd0SJames Wright                                                          "Shift for mass matrix in IJacobian"));
354a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1,
355a424bcd0SJames Wright                                                          "Current solution time"));
356f179a332SJames Wright 
357841e4c73SJed Brown   problem->apply_vol_rhs.qfunction_context = newtonian_ig_context;
358a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ifunction.qfunction_context));
359a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ijacobian.qfunction_context));
360a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow.qfunction_context));
361a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow_jacobian.qfunction_context));
362e334ad8fSJed Brown 
3639f844368SJames Wright   if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
3649f844368SJames Wright   if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
3659f844368SJames Wright   if (bc->num_slip > 0) PetscCall(SlipBCSetup(problem, dm, ctx, newtonian_ig_context));
3669f844368SJames Wright 
367e334ad8fSJed Brown   if (unit_tests) {
368e334ad8fSJed Brown     PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx));
369e334ad8fSJed Brown   }
370ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
37188b783a1SJames Wright }
37288b783a1SJames Wright 
373731c13d7SJames Wright PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx) {
374367c849eSJames Wright   MPI_Comm                 comm = user->comm;
375a424bcd0SJames Wright   Ceed                     ceed = user->ceed;
376841e4c73SJed Brown   NewtonianIdealGasContext newtonian_ctx;
377841e4c73SJed Brown 
37888b783a1SJames Wright   PetscFunctionBeginUser;
379a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ctx));
3802b730f8bSJeremy L Thompson   PetscCall(PetscPrintf(comm,
381841e4c73SJed Brown                         "  Problem:\n"
382841e4c73SJed Brown                         "    Problem Name                       : %s\n"
383841e4c73SJed Brown                         "    Stabilization                      : %s\n",
3842b730f8bSJeremy L Thompson                         app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization]));
385a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ctx));
386ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
38788b783a1SJames Wright }
388