xref: /libCEED/examples/fluids/problems/newtonian.c (revision dc805cc4a09d29f27b3febd084feb659e74b9d08)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, 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 "../navierstokes.h"
1288b783a1SJames Wright #include "../qfunctions/setupgeo.h"
1388b783a1SJames Wright #include "../qfunctions/newtonian.h"
1488b783a1SJames Wright 
15e334ad8fSJed Brown // Compute relative error |a - b|/|s|
16e334ad8fSJed Brown static PetscErrorCode CheckPrimitiveWithTolerance(StatePrimitive sY,
17e334ad8fSJed Brown     StatePrimitive aY, StatePrimitive bY, const char *name, PetscReal rtol_pressure,
18e334ad8fSJed Brown     PetscReal rtol_velocity, PetscReal rtol_temperature) {
196ef2784eSLeila Ghaffari 
20e334ad8fSJed Brown   PetscFunctionBeginUser;
21e334ad8fSJed Brown   StatePrimitive eY; // relative error
22e334ad8fSJed Brown   eY.pressure = (aY.pressure - bY.pressure) / sY.pressure;
23e334ad8fSJed Brown   PetscScalar u = sqrt(Square(sY.velocity[0]) + Square(sY.velocity[1]) + Square(
24e334ad8fSJed Brown                          sY.velocity[2]));
25e334ad8fSJed Brown   for (int j=0; j<3; j++) eY.velocity[j] = (aY.velocity[j] - bY.velocity[j]) / u;
26e334ad8fSJed Brown   eY.temperature = (aY.temperature - bY.temperature) / sY.temperature;
27a1df05f8SJed Brown   if (fabs(eY.pressure) > rtol_pressure)
28a1df05f8SJed Brown     printf("%s: pressure error %g\n", name, eY.pressure);
29a1df05f8SJed Brown   for (int j=0; j<3; j++)
30a1df05f8SJed Brown     if (fabs(eY.velocity[j]) > rtol_velocity)
31e334ad8fSJed Brown       printf("%s: velocity[%d] error %g\n", name, j, eY.velocity[j]);
32e334ad8fSJed Brown   if (fabs(eY.temperature) > rtol_temperature)
33e334ad8fSJed Brown     printf("%s: temperature error %g\n", name, eY.temperature);
34e334ad8fSJed Brown   PetscFunctionReturn(0);
35e334ad8fSJed Brown }
36e334ad8fSJed Brown 
37e334ad8fSJed Brown static PetscErrorCode UnitTests_Newtonian(User user,
38e334ad8fSJed Brown     NewtonianIdealGasContext gas) {
396ef2784eSLeila Ghaffari 
40e334ad8fSJed Brown   Units units = user->units;
41e334ad8fSJed Brown   const CeedScalar eps    = 1e-6;
42*dc805cc4SLeila Ghaffari   const CeedScalar kg     = units->kilogram,
43*dc805cc4SLeila Ghaffari                    m      = units->meter,
44*dc805cc4SLeila Ghaffari                    sec    = units->second,
45e334ad8fSJed Brown                    Pascal = units->Pascal;
46e334ad8fSJed Brown   PetscFunctionBeginUser;
47*dc805cc4SLeila Ghaffari   const CeedScalar rho = 1.2 * kg / (m*m*m),
48*dc805cc4SLeila Ghaffari                    u   = 40 * m/sec;
49e334ad8fSJed Brown   CeedScalar U[5] = {rho, rho*u, rho *u*1.1, rho *u*1.2, 250e3*Pascal + .5*rho *u*u};
50e334ad8fSJed Brown   const CeedScalar x[3] = {.1, .2, .3};
51e334ad8fSJed Brown   State s = StateFromU(gas, U, x);
52e334ad8fSJed Brown   for (int i=0; i<8; i++) {
53e334ad8fSJed Brown     CeedScalar dU[5] = {0}, dx[3] = {0};
54e334ad8fSJed Brown     if (i < 5) dU[i] = U[i];
55e334ad8fSJed Brown     else dx[i-5] = x[i-5];
56e334ad8fSJed Brown     State ds = StateFromU_fwd(gas, s, dU, x, dx);
57e334ad8fSJed Brown     for (int j=0; j<5; j++) dU[j] = (1 + eps * (i == j)) * U[j];
58e334ad8fSJed Brown     for (int j=0; j<3; j++) dx[j] = (1 + eps * (i == 5 + j)) * x[j];
59e334ad8fSJed Brown     State t = StateFromU(gas, dU, dx);
60e334ad8fSJed Brown     StatePrimitive dY;
61e334ad8fSJed Brown     dY.pressure = (t.Y.pressure - s.Y.pressure) / eps;
62e334ad8fSJed Brown     for (int j=0; j<3; j++)
63e334ad8fSJed Brown       dY.velocity[j] = (t.Y.velocity[j] - s.Y.velocity[j]) / eps;
64e334ad8fSJed Brown     dY.temperature = (t.Y.temperature - s.Y.temperature) / eps;
65e334ad8fSJed Brown     char buf[128];
66e334ad8fSJed Brown     snprintf(buf, sizeof buf, "StateFromU_fwd i=%d", i);
67e334ad8fSJed Brown     PetscCall(CheckPrimitiveWithTolerance(dY, ds.Y, dY, buf, 5e-6, 1e-6, 1e-6));
68e334ad8fSJed Brown   }
69e334ad8fSJed Brown   PetscFunctionReturn(0);
70e334ad8fSJed Brown }
71e334ad8fSJed Brown 
72a0add3c9SJed Brown PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx) {
736ef2784eSLeila Ghaffari 
74a0add3c9SJed Brown   SetupContext      setup_context;
7588b783a1SJames Wright   User              user = *(User *)ctx;
7688b783a1SJames Wright   StabilizationType stab;
7788b783a1SJames Wright   MPI_Comm          comm = PETSC_COMM_WORLD;
7888b783a1SJames Wright   PetscBool         implicit;
79*dc805cc4SLeila Ghaffari   PetscBool         has_curr_time = PETSC_FALSE,
80*dc805cc4SLeila Ghaffari                     prim_var, unit_tests;
8188b783a1SJames Wright   PetscInt          ierr;
82841e4c73SJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
83841e4c73SJed Brown   CeedQFunctionContext newtonian_ig_context;
8488b783a1SJames Wright 
85841e4c73SJed Brown   PetscFunctionBeginUser;
86a0add3c9SJed Brown   ierr = PetscCalloc1(1, &setup_context); CHKERRQ(ierr);
87841e4c73SJed Brown   ierr = PetscCalloc1(1, &newtonian_ig_ctx); CHKERRQ(ierr);
8888b783a1SJames Wright 
8988b783a1SJames Wright   // ------------------------------------------------------
9088b783a1SJames Wright   //           Setup Generic Newtonian IG Problem
9188b783a1SJames Wright   // ------------------------------------------------------
9288b783a1SJames Wright   problem->dim                                  = 3;
9388b783a1SJames Wright   problem->q_data_size_vol                      = 10;
94ba6664aeSJames Wright   problem->q_data_size_sur                      = 10;
950ec2498eSJames Wright   problem->jac_data_size_sur                    = 11;
9691e5af17SJed Brown   problem->setup_vol.qfunction                  = Setup;
9791e5af17SJed Brown   problem->setup_vol.qfunction_loc              = Setup_loc;
9891e5af17SJed Brown   problem->setup_sur.qfunction                  = SetupBoundary;
9991e5af17SJed Brown   problem->setup_sur.qfunction_loc              = SetupBoundary_loc;
100a0add3c9SJed Brown   problem->bc                                   = NULL;
101a0add3c9SJed Brown   problem->bc_ctx                               = setup_context;
10288b783a1SJames Wright   problem->non_zero_time                        = PETSC_FALSE;
103*dc805cc4SLeila Ghaffari   problem->print_info                           = PRINT_NEWTONIAN;
10488b783a1SJames Wright 
10588b783a1SJames Wright   // ------------------------------------------------------
10688b783a1SJames Wright   //             Create the libCEED context
10788b783a1SJames Wright   // ------------------------------------------------------
10888b783a1SJames Wright   CeedScalar cv     = 717.;          // J/(kg K)
10988b783a1SJames Wright   CeedScalar cp     = 1004.;         // J/(kg K)
11088626eedSJames Wright   CeedScalar g[3]   = {0, 0, -9.81}; // m/s^2
11188b783a1SJames Wright   CeedScalar lambda = -2./3.;        // -
11288626eedSJames Wright   CeedScalar mu     = 1.8e-5;        // Pa s, dynamic viscosity
11388b783a1SJames Wright   CeedScalar k      = 0.02638;       // W/(m K)
11488b783a1SJames Wright   CeedScalar c_tau  = 0.5;           // -
11588626eedSJames Wright   CeedScalar Ctau_t = 1.0;           // -
11688626eedSJames Wright   CeedScalar Ctau_v = 36.0;          // TODO make function of degree
11788626eedSJames Wright   CeedScalar Ctau_C = 1.0;           // TODO make function of degree
11888626eedSJames Wright   CeedScalar Ctau_M = 1.0;           // TODO make function of degree
11988626eedSJames Wright   CeedScalar Ctau_E = 1.0;           // TODO make function of degree
12088b783a1SJames Wright   PetscReal domain_min[3], domain_max[3], domain_size[3];
12188b783a1SJames Wright   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
122ba6664aeSJames Wright   for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
12388b783a1SJames Wright 
12488b783a1SJames Wright   // ------------------------------------------------------
12588b783a1SJames Wright   //             Create the PETSc context
12688b783a1SJames Wright   // ------------------------------------------------------
12788626eedSJames Wright   PetscScalar meter    = 1;  // 1 meter in scaled length units
12888626eedSJames Wright   PetscScalar kilogram = 1;  // 1 kilogram in scaled mass units
12988626eedSJames Wright   PetscScalar second   = 1;  // 1 second in scaled time units
13088b783a1SJames Wright   PetscScalar Kelvin   = 1;  // 1 Kelvin in scaled temperature units
13188b783a1SJames Wright   PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s;
13288b783a1SJames Wright 
13388b783a1SJames Wright   // ------------------------------------------------------
13488b783a1SJames Wright   //              Command line Options
13588b783a1SJames Wright   // ------------------------------------------------------
13667490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem",
13767490bc6SJeremy L Thompson                     NULL);
138*dc805cc4SLeila Ghaffari   // -- Conservative vs Primitive variables
139*dc805cc4SLeila Ghaffari   ierr = PetscOptionsBool("-primitive", "Use primitive variables",
140*dc805cc4SLeila Ghaffari                           NULL, prim_var=PETSC_FALSE, &prim_var, NULL); CHKERRQ(ierr);
141*dc805cc4SLeila Ghaffari   if (prim_var) {
142*dc805cc4SLeila Ghaffari     problem->ics.qfunction                        = ICsNewtonianIG_Prim;
143*dc805cc4SLeila Ghaffari     problem->ics.qfunction_loc                    = ICsNewtonianIG_Prim_loc;
144*dc805cc4SLeila Ghaffari     problem->apply_vol_ifunction.qfunction        = IFunction_Newtonian_Prim;
145*dc805cc4SLeila Ghaffari     problem->apply_vol_ifunction.qfunction_loc    = IFunction_Newtonian_Prim_loc;
146*dc805cc4SLeila Ghaffari     problem->apply_vol_ijacobian.qfunction        = IJacobian_Newtonian_Prim;
147*dc805cc4SLeila Ghaffari     problem->apply_vol_ijacobian.qfunction_loc    = IJacobian_Newtonian_Prim_loc;
148*dc805cc4SLeila Ghaffari   } else {
149*dc805cc4SLeila Ghaffari     problem->ics.qfunction                        = ICsNewtonianIG;
150*dc805cc4SLeila Ghaffari     problem->ics.qfunction_loc                    = ICsNewtonianIG_loc;
151*dc805cc4SLeila Ghaffari     problem->apply_vol_rhs.qfunction              = RHSFunction_Newtonian;
152*dc805cc4SLeila Ghaffari     problem->apply_vol_rhs.qfunction_loc          = RHSFunction_Newtonian_loc;
153*dc805cc4SLeila Ghaffari     problem->apply_vol_ifunction.qfunction        = IFunction_Newtonian;
154*dc805cc4SLeila Ghaffari     problem->apply_vol_ifunction.qfunction_loc    = IFunction_Newtonian_loc;
155*dc805cc4SLeila Ghaffari     problem->apply_vol_ijacobian.qfunction        = IJacobian_Newtonian;
156*dc805cc4SLeila Ghaffari     problem->apply_vol_ijacobian.qfunction_loc    = IJacobian_Newtonian_loc;
157*dc805cc4SLeila Ghaffari     problem->apply_inflow.qfunction               = BoundaryIntegral;
158*dc805cc4SLeila Ghaffari     problem->apply_inflow.qfunction_loc           = BoundaryIntegral_loc;
159*dc805cc4SLeila Ghaffari     problem->apply_inflow_jacobian.qfunction      = BoundaryIntegral_Jacobian;
160*dc805cc4SLeila Ghaffari     problem->apply_inflow_jacobian.qfunction_loc  = BoundaryIntegral_Jacobian_loc;
161*dc805cc4SLeila Ghaffari     problem->apply_outflow.qfunction              = PressureOutflow;
162*dc805cc4SLeila Ghaffari     problem->apply_outflow.qfunction_loc          = PressureOutflow_loc;
163*dc805cc4SLeila Ghaffari     problem->apply_outflow_jacobian.qfunction     = PressureOutflow_Jacobian;
164*dc805cc4SLeila Ghaffari     problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_loc;
165*dc805cc4SLeila Ghaffari   }
16667490bc6SJeremy L Thompson 
16788b783a1SJames Wright   // -- Physics
16888b783a1SJames Wright   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
16988b783a1SJames Wright                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
17088b783a1SJames Wright   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
17188b783a1SJames Wright                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
17288b783a1SJames Wright   ierr = PetscOptionsScalar("-lambda",
17388b783a1SJames Wright                             "Stokes hypothesis second viscosity coefficient",
17488b783a1SJames Wright                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
17588b783a1SJames Wright   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
17688b783a1SJames Wright                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
17788b783a1SJames Wright   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
17888b783a1SJames Wright                             NULL, k, &k, NULL); CHKERRQ(ierr);
17988b783a1SJames Wright 
18088626eedSJames Wright   PetscInt dim = problem->dim;
18188626eedSJames Wright   ierr = PetscOptionsRealArray("-g", "Gravitational acceleration",
18288626eedSJames Wright                                NULL, g, &dim, NULL); CHKERRQ(ierr);
18388b783a1SJames Wright   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
18488b783a1SJames Wright                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
18588b783a1SJames Wright                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
18688b783a1SJames Wright   ierr = PetscOptionsScalar("-c_tau", "Stabilization constant",
18788b783a1SJames Wright                             NULL, c_tau, &c_tau, NULL); CHKERRQ(ierr);
18888626eedSJames Wright   ierr = PetscOptionsScalar("-Ctau_t", "Stabilization time constant",
18988626eedSJames Wright                             NULL, Ctau_t, &Ctau_t, NULL); CHKERRQ(ierr);
19088626eedSJames Wright   ierr = PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant",
19188626eedSJames Wright                             NULL, Ctau_v, &Ctau_v, NULL); CHKERRQ(ierr);
19288626eedSJames Wright   ierr = PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant",
19388626eedSJames Wright                             NULL, Ctau_C, &Ctau_C, NULL); CHKERRQ(ierr);
19488626eedSJames Wright   ierr = PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant",
19588626eedSJames Wright                             NULL, Ctau_M, &Ctau_M, NULL); CHKERRQ(ierr);
19688626eedSJames Wright   ierr = PetscOptionsScalar("-Ctau_E", "Stabilization energy constant",
19788626eedSJames Wright                             NULL, Ctau_E, &Ctau_E, NULL); CHKERRQ(ierr);
19888b783a1SJames Wright   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
19988b783a1SJames Wright                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
20088b783a1SJames Wright   CHKERRQ(ierr);
201e334ad8fSJed Brown   ierr = PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests",
202e334ad8fSJed Brown                           NULL, unit_tests=PETSC_FALSE, &unit_tests, NULL);
203e334ad8fSJed Brown   CHKERRQ(ierr);
20488b783a1SJames Wright 
20588b783a1SJames Wright   // -- Units
20688b783a1SJames Wright   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
20788b783a1SJames Wright                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
20888b783a1SJames Wright   meter = fabs(meter);
20988b783a1SJames Wright   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
21088b783a1SJames Wright                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
21188b783a1SJames Wright   kilogram = fabs(kilogram);
21288b783a1SJames Wright   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
21388b783a1SJames Wright                             NULL, second, &second, NULL); CHKERRQ(ierr);
21488b783a1SJames Wright   second = fabs(second);
21588b783a1SJames Wright   ierr = PetscOptionsScalar("-units_Kelvin",
21688b783a1SJames Wright                             "1 Kelvin in scaled temperature units",
21788b783a1SJames Wright                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
21888b783a1SJames Wright   Kelvin = fabs(Kelvin);
21988b783a1SJames Wright 
22088b783a1SJames Wright   // -- Warnings
22188b783a1SJames Wright   if (stab == STAB_SUPG && !implicit) {
22288b783a1SJames Wright     ierr = PetscPrintf(comm,
22388b783a1SJames Wright                        "Warning! Use -stab supg only with -implicit\n");
22488b783a1SJames Wright     CHKERRQ(ierr);
22588b783a1SJames Wright   }
226*dc805cc4SLeila Ghaffari   if (prim_var && !implicit) {
227*dc805cc4SLeila Ghaffari     SETERRQ(comm, PETSC_ERR_ARG_NULL,
228*dc805cc4SLeila Ghaffari             "RHSFunction is not provided for primitive variables (use -primitive only with -implicit)\n");
229*dc805cc4SLeila Ghaffari   }
23067490bc6SJeremy L Thompson   PetscOptionsEnd();
23188b783a1SJames Wright 
23288b783a1SJames Wright   // ------------------------------------------------------
23388b783a1SJames Wright   //           Set up the PETSc context
23488b783a1SJames Wright   // ------------------------------------------------------
23588b783a1SJames Wright   // -- Define derived units
23688b783a1SJames Wright   Pascal          = kilogram / (meter * PetscSqr(second));
23788b783a1SJames Wright   J_per_kg_K      =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
23888b783a1SJames Wright   m_per_squared_s = meter / PetscSqr(second);
23988b783a1SJames Wright   W_per_m_K       = kilogram * meter / (pow(second,3) * Kelvin);
24088b783a1SJames Wright 
24188b783a1SJames Wright   user->units->meter           = meter;
24288b783a1SJames Wright   user->units->kilogram        = kilogram;
24388b783a1SJames Wright   user->units->second          = second;
24488b783a1SJames Wright   user->units->Kelvin          = Kelvin;
24588b783a1SJames Wright   user->units->Pascal          = Pascal;
24688b783a1SJames Wright   user->units->J_per_kg_K      = J_per_kg_K;
24788b783a1SJames Wright   user->units->m_per_squared_s = m_per_squared_s;
24888b783a1SJames Wright   user->units->W_per_m_K       = W_per_m_K;
24988b783a1SJames Wright 
25088b783a1SJames Wright   // ------------------------------------------------------
25188b783a1SJames Wright   //           Set up the libCEED context
25288b783a1SJames Wright   // ------------------------------------------------------
25388b783a1SJames Wright   // -- Scale variables to desired units
25488b783a1SJames Wright   cv     *= J_per_kg_K;
25588b783a1SJames Wright   cp     *= J_per_kg_K;
25688b783a1SJames Wright   mu     *= Pascal * second;
25788b783a1SJames Wright   k      *= W_per_m_K;
258ba6664aeSJames Wright   for (PetscInt i=0; i<3; i++) domain_size[i] *= meter;
259ba6664aeSJames Wright   for (PetscInt i=0; i<3; i++) g[i]           *= m_per_squared_s;
26088b783a1SJames Wright   problem->dm_scale = meter;
26188b783a1SJames Wright 
26288b783a1SJames Wright   // -- Setup Context
26388b783a1SJames Wright   setup_context->cv         = cv;
26488b783a1SJames Wright   setup_context->cp         = cp;
26588b783a1SJames Wright   setup_context->lx         = domain_size[0];
26688b783a1SJames Wright   setup_context->ly         = domain_size[1];
26788b783a1SJames Wright   setup_context->lz         = domain_size[2];
26888b783a1SJames Wright   setup_context->time       = 0;
26988626eedSJames Wright   ierr = PetscArraycpy(setup_context->g, g, 3); CHKERRQ(ierr);
27088b783a1SJames Wright 
27188b783a1SJames Wright   // -- Solver Settings
27288b783a1SJames Wright   user->phys->stab          = stab;
27388b783a1SJames Wright   user->phys->implicit      = implicit;
274*dc805cc4SLeila Ghaffari   user->phys->primitive     = prim_var;
27588b783a1SJames Wright   user->phys->has_curr_time = has_curr_time;
27688b783a1SJames Wright 
27788b783a1SJames Wright   // -- QFunction Context
278841e4c73SJed Brown   newtonian_ig_ctx->lambda        = lambda;
279841e4c73SJed Brown   newtonian_ig_ctx->mu            = mu;
280841e4c73SJed Brown   newtonian_ig_ctx->k             = k;
281841e4c73SJed Brown   newtonian_ig_ctx->cv            = cv;
282841e4c73SJed Brown   newtonian_ig_ctx->cp            = cp;
283841e4c73SJed Brown   newtonian_ig_ctx->c_tau         = c_tau;
284841e4c73SJed Brown   newtonian_ig_ctx->Ctau_t        = Ctau_t;
285841e4c73SJed Brown   newtonian_ig_ctx->Ctau_v        = Ctau_v;
286841e4c73SJed Brown   newtonian_ig_ctx->Ctau_C        = Ctau_C;
287841e4c73SJed Brown   newtonian_ig_ctx->Ctau_M        = Ctau_M;
288841e4c73SJed Brown   newtonian_ig_ctx->Ctau_E        = Ctau_E;
289841e4c73SJed Brown   newtonian_ig_ctx->stabilization = stab;
29065dd5cafSJames Wright   newtonian_ig_ctx->is_implicit   = implicit;
291*dc805cc4SLeila Ghaffari   newtonian_ig_ctx->primitive     = prim_var;
292841e4c73SJed Brown   ierr = PetscArraycpy(newtonian_ig_ctx->g, g, 3); CHKERRQ(ierr);
29388b783a1SJames Wright 
294841e4c73SJed Brown   CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context);
295841e4c73SJed Brown   CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST,
296841e4c73SJed Brown                               CEED_USE_POINTER, sizeof(*setup_context), setup_context);
29717ce10faSJeremy L Thompson   CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context,
29817ce10faSJeremy L Thompson                                      CEED_MEM_HOST,
29917ce10faSJeremy L Thompson                                      FreeContextPetsc);
300841e4c73SJed Brown   CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context,
301841e4c73SJed Brown                                      "evaluation time",
302841e4c73SJed Brown                                      (char *)&setup_context->time - (char *)setup_context, 1, "Time of evaluation");
303841e4c73SJed Brown 
304841e4c73SJed Brown   CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context);
305841e4c73SJed Brown   CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST,
306841e4c73SJed Brown                               CEED_USE_POINTER,
307841e4c73SJed Brown                               sizeof(*newtonian_ig_ctx), newtonian_ig_ctx);
308841e4c73SJed Brown   CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST,
309841e4c73SJed Brown                                      FreeContextPetsc);
310841e4c73SJed Brown   CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size",
311841e4c73SJed Brown                                      offsetof(struct NewtonianIdealGasContext_, dt), 1, "Size of timestep, delta t");
312e334ad8fSJed Brown   CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift",
313e334ad8fSJed Brown                                      offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1,
314e334ad8fSJed Brown                                      "Shift for mass matrix in IJacobian");
315841e4c73SJed Brown   problem->apply_vol_rhs.qfunction_context = newtonian_ig_context;
316841e4c73SJed Brown   CeedQFunctionContextReferenceCopy(newtonian_ig_context,
317841e4c73SJed Brown                                     &problem->apply_vol_ifunction.qfunction_context);
318e334ad8fSJed Brown   CeedQFunctionContextReferenceCopy(newtonian_ig_context,
319e334ad8fSJed Brown                                     &problem->apply_vol_ijacobian.qfunction_context);
32065dd5cafSJames Wright   CeedQFunctionContextReferenceCopy(newtonian_ig_context,
32165dd5cafSJames Wright                                     &problem->apply_inflow.qfunction_context);
32230e9fa81SJames Wright   CeedQFunctionContextReferenceCopy(newtonian_ig_context,
323fc02c281SJames Wright                                     &problem->apply_inflow_jacobian.qfunction_context);
324fc02c281SJames Wright   CeedQFunctionContextReferenceCopy(newtonian_ig_context,
32530e9fa81SJames Wright                                     &problem->apply_outflow.qfunction_context);
32630e9fa81SJames Wright   CeedQFunctionContextReferenceCopy(newtonian_ig_context,
32730e9fa81SJames Wright                                     &problem->apply_outflow_jacobian.qfunction_context);
328e334ad8fSJed Brown 
329e334ad8fSJed Brown   if (unit_tests) {
330e334ad8fSJed Brown     PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx));
331e334ad8fSJed Brown   }
33288b783a1SJames Wright   PetscFunctionReturn(0);
33388b783a1SJames Wright }
33488b783a1SJames Wright 
335*dc805cc4SLeila Ghaffari PetscErrorCode PRINT_NEWTONIAN(ProblemData *problem, AppCtx app_ctx) {
336*dc805cc4SLeila Ghaffari 
337841e4c73SJed Brown   MPI_Comm comm = PETSC_COMM_WORLD;
338841e4c73SJed Brown   PetscErrorCode ierr;
339841e4c73SJed Brown   NewtonianIdealGasContext newtonian_ctx;
340841e4c73SJed Brown 
34188b783a1SJames Wright   PetscFunctionBeginUser;
342841e4c73SJed Brown   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context,
343841e4c73SJed Brown                               CEED_MEM_HOST, &newtonian_ctx);
344841e4c73SJed Brown   ierr = PetscPrintf(comm,
345841e4c73SJed Brown                      "  Problem:\n"
346841e4c73SJed Brown                      "    Problem Name                       : %s\n"
347841e4c73SJed Brown                      "    Stabilization                      : %s\n",
348841e4c73SJed Brown                      app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization]);
349841e4c73SJed Brown   CHKERRQ(ierr);
350841e4c73SJed Brown   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context,
351841e4c73SJed Brown                                   &newtonian_ctx);
35288b783a1SJames Wright   PetscFunctionReturn(0);
35388b783a1SJames Wright }
354