1 // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED: http://github.com/ceed
7
8 /// @file
9 /// Utility functions for setting up problems using the Newtonian Qfunction
10
11 #include "../qfunctions/newtonian.h"
12
13 #include <ceed.h>
14 #include <petscdm.h>
15
16 #include "../navierstokes.h"
17
18 // For use with PetscOptionsEnum
19 static const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "ENTROPY", "StateVariable", "STATEVAR_", NULL};
20
CheckQWithTolerance(const CeedScalar Q_s[5],const CeedScalar Q_a[5],const CeedScalar Q_b[5],const char * name,PetscReal rtol_0,PetscReal rtol_u,PetscReal rtol_4)21 static PetscErrorCode CheckQWithTolerance(const CeedScalar Q_s[5], const CeedScalar Q_a[5], const CeedScalar Q_b[5], const char *name,
22 PetscReal rtol_0, PetscReal rtol_u, PetscReal rtol_4) {
23 CeedScalar relative_error[5]; // relative error
24 CeedScalar divisor_threshold = 10 * CEED_EPSILON;
25
26 PetscFunctionBeginUser;
27 relative_error[0] = (Q_a[0] - Q_b[0]) / (fabs(Q_s[0]) > divisor_threshold ? Q_s[0] : 1);
28 relative_error[4] = (Q_a[4] - Q_b[4]) / (fabs(Q_s[4]) > divisor_threshold ? Q_s[4] : 1);
29
30 CeedScalar u_magnitude = sqrt(Square(Q_s[1]) + Square(Q_s[2]) + Square(Q_s[3]));
31 CeedScalar u_divisor = u_magnitude > divisor_threshold ? u_magnitude : 1;
32 for (int i = 1; i < 4; i++) {
33 relative_error[i] = (Q_a[i] - Q_b[i]) / u_divisor;
34 }
35
36 if (fabs(relative_error[0]) >= rtol_0) {
37 printf("%s[0] error %g (expected %.10e, got %.10e)\n", name, relative_error[0], Q_s[0], Q_a[0]);
38 }
39 for (int i = 1; i < 4; i++) {
40 if (fabs(relative_error[i]) >= rtol_u) {
41 printf("%s[%d] error %g (expected %.10e, got %.10e)\n", name, i, relative_error[i], Q_s[i], Q_a[i]);
42 }
43 }
44 if (fabs(relative_error[4]) >= rtol_4) {
45 printf("%s[4] error %g (expected %.10e, got %.10e)\n", name, relative_error[4], Q_s[4], Q_a[4]);
46 }
47 PetscFunctionReturn(PETSC_SUCCESS);
48 }
49
50 // @brief Verify `StateFromQ` by converting A0 -> B0 -> A0_test, where A0 should equal A0_test
TestState(StateVariable state_var_A,StateVariable state_var_B,NewtonianIdealGasContext gas,const CeedScalar A0[5],CeedScalar rtol_0,CeedScalar rtol_u,CeedScalar rtol_4)51 static PetscErrorCode TestState(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5],
52 CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
53 CeedScalar B0[5], A0_test[5];
54 char buf[128];
55 const char *const StateVariables_Initial[] = {"U", "Y", "V"};
56
57 PetscFunctionBeginUser;
58 const char *A_initial = StateVariables_Initial[state_var_A];
59 const char *B_initial = StateVariables_Initial[state_var_B];
60
61 State state_A0 = StateFromQ(gas, A0, state_var_A);
62 StateToQ(gas, state_A0, B0, state_var_B);
63 State state_B0 = StateFromQ(gas, B0, state_var_B);
64 StateToQ(gas, state_B0, A0_test, state_var_A);
65
66 snprintf(buf, sizeof buf, "%s->%s->%s: %s", A_initial, B_initial, A_initial, A_initial);
67 PetscCall(CheckQWithTolerance(A0, A0_test, A0, buf, rtol_0, rtol_u, rtol_4));
68 PetscFunctionReturn(PETSC_SUCCESS);
69 }
70
71 // @brief Verify `StateFromQ_fwd` via a finite difference approximation
TestState_fwd(StateVariable state_var_A,StateVariable state_var_B,NewtonianIdealGasContext gas,const CeedScalar A0[5],CeedScalar rtol_0,CeedScalar rtol_u,CeedScalar rtol_4)72 static PetscErrorCode TestState_fwd(StateVariable state_var_A, StateVariable state_var_B, NewtonianIdealGasContext gas, const CeedScalar A0[5],
73 CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
74 CeedScalar eps = 4e-7; // Finite difference step
75 char buf[128];
76 const char *const StateVariables_Initial[] = {"U", "Y", "V"};
77
78 PetscFunctionBeginUser;
79 const char *A_initial = StateVariables_Initial[state_var_A];
80 const char *B_initial = StateVariables_Initial[state_var_B];
81 State state_0 = StateFromQ(gas, A0, state_var_A);
82
83 for (int i = 0; i < 5; i++) {
84 CeedScalar dB[5] = {0.}, dB_fd[5] = {0.};
85 { // Calculate dB using State functions
86 CeedScalar dA[5] = {0};
87
88 dA[i] = A0[i];
89 State dstate_0 = StateFromQ_fwd(gas, state_0, dA, state_var_A);
90 StateToQ_fwd(gas, state_0, dstate_0, dB, state_var_B);
91 }
92
93 { // Calculate dB_fd via finite difference approximation
94 CeedScalar A1[5], B0[5], B1[5];
95
96 for (int j = 0; j < 5; j++) A1[j] = (1 + eps * (i == j)) * A0[j];
97 State state_1 = StateFromQ(gas, A1, state_var_A);
98 StateToQ(gas, state_0, B0, state_var_B);
99 StateToQ(gas, state_1, B1, state_var_B);
100 for (int j = 0; j < 5; j++) dB_fd[j] = (B1[j] - B0[j]) / eps;
101 }
102
103 snprintf(buf, sizeof buf, "d%s->d%s: StateFrom%s_fwd i=%d: d%s", A_initial, B_initial, A_initial, i, B_initial);
104 PetscCall(CheckQWithTolerance(dB_fd, dB, dB_fd, buf, rtol_0, rtol_u, rtol_4));
105 }
106 PetscFunctionReturn(PETSC_SUCCESS);
107 }
108
109 // @brief Test the Newtonian State transformation functions, `StateFrom*`
UnitTests_Newtonian(User user,NewtonianIdealGasContext gas)110 static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) {
111 Units units = user->units;
112 const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, K = units->Kelvin;
113
114 PetscFunctionBeginUser;
115 const CeedScalar T = 200 * K;
116 const CeedScalar rho = 1.2 * kg / Cube(m);
117 const CeedScalar P = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T;
118 const CeedScalar u_base = 40 * m / sec;
119 const CeedScalar u[3] = {u_base, u_base * 1.1, u_base * 1.2};
120 const CeedScalar e_kinetic = 0.5 * Dot3(u, u);
121 const CeedScalar e_internal = gas->cv * T;
122 const CeedScalar e_total = e_kinetic + e_internal;
123 const CeedScalar gamma = HeatCapacityRatio(gas);
124 const CeedScalar entropy = log(P) - gamma * log(rho);
125 const CeedScalar rho_div_p = rho / P;
126 const CeedScalar Y0[5] = {P, u[0], u[1], u[2], T};
127 const CeedScalar U0[5] = {rho, rho * u[0], rho * u[1], rho * u[2], rho * e_total};
128 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],
129 -rho_div_p};
130
131 {
132 CeedScalar rtol = 40 * CEED_EPSILON;
133
134 PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol));
135 PetscCall(TestState(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol));
136 PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol));
137 PetscCall(TestState(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, rtol, rtol, rtol));
138 PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, rtol, rtol, rtol));
139 PetscCall(TestState(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, rtol, rtol, rtol));
140 }
141
142 {
143 CeedScalar rtol = 5e-6;
144
145 PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_CONSERVATIVE, gas, Y0, rtol, rtol, rtol));
146 PetscCall(TestState_fwd(STATEVAR_PRIMITIVE, STATEVAR_ENTROPY, gas, Y0, rtol, rtol, rtol));
147 PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_PRIMITIVE, gas, U0, rtol, rtol, rtol));
148 PetscCall(TestState_fwd(STATEVAR_CONSERVATIVE, STATEVAR_ENTROPY, gas, U0, 10 * rtol, rtol, rtol));
149 PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_CONSERVATIVE, gas, V0, 5 * rtol, rtol, rtol));
150 PetscCall(TestState_fwd(STATEVAR_ENTROPY, STATEVAR_PRIMITIVE, gas, V0, 5 * rtol, 5 * rtol, 5 * rtol));
151 }
152 PetscFunctionReturn(PETSC_SUCCESS);
153 }
154
155 // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping
156 //
157 // Only used for SUPG stabilization
CreateKSPMassOperator_NewtonianStabilized(User user,CeedOperator * op_mass)158 PetscErrorCode CreateKSPMassOperator_NewtonianStabilized(User user, CeedOperator *op_mass) {
159 Ceed ceed = user->ceed;
160 CeedInt num_comp_q, q_data_size;
161 CeedQFunction qf_mass;
162 CeedElemRestriction elem_restr_q, elem_restr_qd_i;
163 CeedBasis basis_q;
164 CeedVector q_data;
165 CeedQFunctionContext qf_ctx = NULL;
166 PetscInt dim = 3;
167
168 PetscFunctionBeginUser;
169 { // Get restriction and basis from the RHS function
170 CeedOperator *sub_ops;
171 CeedOperatorField field;
172 PetscInt sub_op_index = 0; // will be 0 for the volume op
173
174 PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(user->op_rhs_ctx->op, &sub_ops));
175 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
176 PetscCallCeed(ceed, CeedOperatorFieldGetData(field, NULL, &elem_restr_q, &basis_q, NULL));
177
178 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field));
179 PetscCallCeed(ceed, CeedOperatorFieldGetData(field, NULL, &elem_restr_qd_i, NULL, &q_data));
180
181 PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx));
182 }
183
184 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
185 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size));
186
187 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Newtonian_Conserv, MassFunction_Newtonian_Conserv_loc, &qf_mass));
188
189 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx));
190 PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
191 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
192 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
193 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
194 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
195 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
196
197 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
198 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
199 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed));
200 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data));
201 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
202 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
203
204 PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
205 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
206 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i));
207 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
208 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qf_ctx));
209 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
210 PetscFunctionReturn(PETSC_SUCCESS);
211 }
212
NS_NEWTONIAN_IG(ProblemData problem,DM dm,void * ctx,SimpleBC bc)213 PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
214 SetupContext setup_context;
215 User user = *(User *)ctx;
216 CeedInt degree = user->app_ctx->degree;
217 StabilizationType stab;
218 StateVariable state_var;
219 MPI_Comm comm = user->comm;
220 Ceed ceed = user->ceed;
221 PetscBool implicit;
222 PetscBool unit_tests;
223 NewtonianIdealGasContext newtonian_ig_ctx;
224 CeedQFunctionContext newtonian_ig_context;
225
226 PetscFunctionBeginUser;
227 PetscCall(PetscCalloc1(1, &setup_context));
228 PetscCall(PetscCalloc1(1, &newtonian_ig_ctx));
229
230 // ------------------------------------------------------
231 // Setup Generic Newtonian IG Problem
232 // ------------------------------------------------------
233 problem->dim = 3;
234 problem->jac_data_size_sur = 11;
235 problem->compute_exact_solution_error = PETSC_FALSE;
236 problem->print_info = PRINT_NEWTONIAN;
237 problem->uses_newtonian = PETSC_TRUE;
238
239 // ------------------------------------------------------
240 // Create the libCEED context
241 // ------------------------------------------------------
242 CeedScalar cv = 717.; // J/(kg K)
243 CeedScalar cp = 1004.; // J/(kg K)
244 CeedScalar g[3] = {0, 0, 0}; // m/s^2
245 CeedScalar lambda = -2. / 3.; // -
246 CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity
247 CeedScalar k = 0.02638; // W/(m K)
248 CeedScalar c_tau = 0.5 / degree; // -
249 CeedScalar Ctau_t = 1.0; // -
250 CeedScalar Cv_func[3] = {36, 60, 128};
251 CeedScalar Ctau_v = Cv_func[(CeedInt)Min(3, degree) - 1];
252 CeedScalar Ctau_C = 0.25 / degree;
253 CeedScalar Ctau_M = 0.25 / degree;
254 CeedScalar Ctau_E = 0.125;
255 PetscReal domain_min[3], domain_max[3], domain_size[3];
256 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
257 for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
258
259 StatePrimitive reference = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15};
260 CeedScalar idl_decay_time = -1, idl_start = 0, idl_length = 0, idl_pressure = reference.pressure;
261 PetscBool idl_enable = PETSC_FALSE;
262
263 // ------------------------------------------------------
264 // Create the PETSc context
265 // ------------------------------------------------------
266 PetscScalar meter = 1; // 1 meter in scaled length units
267 PetscScalar kilogram = 1; // 1 kilogram in scaled mass units
268 PetscScalar second = 1; // 1 second in scaled time units
269 PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units
270 PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s;
271
272 // ------------------------------------------------------
273 // Command line Options
274 // ------------------------------------------------------
275 PetscBool given_option = PETSC_FALSE;
276 PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL);
277 // -- Conservative vs Primitive variables
278 PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE),
279 (PetscEnum *)&state_var, NULL));
280
281 switch (state_var) {
282 case STATEVAR_CONSERVATIVE:
283 problem->ics.qfunction = ICsNewtonianIG_Conserv;
284 problem->ics.qfunction_loc = ICsNewtonianIG_Conserv_loc;
285 problem->apply_vol_rhs.qfunction = RHSFunction_Newtonian;
286 problem->apply_vol_rhs.qfunction_loc = RHSFunction_Newtonian_loc;
287 problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Conserv;
288 problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Conserv_loc;
289 problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Conserv;
290 problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Conserv_loc;
291 problem->apply_inflow.qfunction = BoundaryIntegral_Conserv;
292 problem->apply_inflow.qfunction_loc = BoundaryIntegral_Conserv_loc;
293 problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Conserv;
294 problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc;
295 break;
296 case STATEVAR_PRIMITIVE:
297 problem->ics.qfunction = ICsNewtonianIG_Prim;
298 problem->ics.qfunction_loc = ICsNewtonianIG_Prim_loc;
299 problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Prim;
300 problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Prim_loc;
301 problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Prim;
302 problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Prim_loc;
303 problem->apply_inflow.qfunction = BoundaryIntegral_Prim;
304 problem->apply_inflow.qfunction_loc = BoundaryIntegral_Prim_loc;
305 problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Prim;
306 problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc;
307 break;
308 case STATEVAR_ENTROPY:
309 problem->ics.qfunction = ICsNewtonianIG_Entropy;
310 problem->ics.qfunction_loc = ICsNewtonianIG_Entropy_loc;
311 problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Entropy;
312 problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Entropy_loc;
313 problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Entropy;
314 problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Entropy_loc;
315 problem->apply_inflow.qfunction = BoundaryIntegral_Entropy;
316 problem->apply_inflow.qfunction_loc = BoundaryIntegral_Entropy_loc;
317 problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Entropy;
318 problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Entropy_loc;
319 break;
320 }
321
322 // -- Physics
323 PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL));
324 PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL));
325 PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL));
326 PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL));
327 PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL));
328
329 PetscInt dim = problem->dim;
330 PetscCall(PetscOptionsDeprecated("-g", "-gravity", "libCEED 0.11.1", NULL));
331 PetscCall(PetscOptionsRealArray("-gravity", "Gravitational acceleration vector", NULL, g, &dim, &given_option));
332 if (given_option) PetscCheck(dim == 3, comm, PETSC_ERR_ARG_SIZ, "Gravity vector must be size 3, %" PetscInt_FMT " values given", dim);
333
334 PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
335 PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
336 PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
337 PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL));
338 PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL));
339 PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL));
340 PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL));
341 PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
342 PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL));
343
344 dim = 3;
345 PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL));
346 PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL));
347 PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL));
348
349 // -- Units
350 PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
351 meter = fabs(meter);
352 PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL));
353 kilogram = fabs(kilogram);
354 PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
355 second = fabs(second);
356 PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL));
357 Kelvin = fabs(Kelvin);
358
359 // -- Warnings
360 PetscCheck(!(state_var == STATEVAR_PRIMITIVE && !implicit), comm, PETSC_ERR_SUP,
361 "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n");
362
363 PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point",
364 NULL, idl_decay_time, &idl_decay_time, &idl_enable));
365 PetscCheck(!(idl_enable && idl_decay_time == 0), comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero.");
366 if (idl_decay_time < 0) idl_enable = PETSC_FALSE;
367 PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL));
368 PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL));
369 idl_pressure = reference.pressure;
370 PetscCall(PetscOptionsScalar("-idl_pressure", "Pressure IDL uses as reference (default is `-reference_pressure`)", NULL, idl_pressure,
371 &idl_pressure, NULL));
372 PetscOptionsEnd();
373
374 if (stab == STAB_SUPG && !implicit) problem->create_mass_operator = CreateKSPMassOperator_NewtonianStabilized;
375
376 // ------------------------------------------------------
377 // Set up the PETSc context
378 // ------------------------------------------------------
379 // -- Define derived units
380 Pascal = kilogram / (meter * PetscSqr(second));
381 J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin);
382 m_per_squared_s = meter / PetscSqr(second);
383 W_per_m_K = kilogram * meter / (pow(second, 3) * Kelvin);
384
385 user->units->meter = meter;
386 user->units->kilogram = kilogram;
387 user->units->second = second;
388 user->units->Kelvin = Kelvin;
389 user->units->Pascal = Pascal;
390 user->units->J_per_kg_K = J_per_kg_K;
391 user->units->m_per_squared_s = m_per_squared_s;
392 user->units->W_per_m_K = W_per_m_K;
393
394 // ------------------------------------------------------
395 // Set up the libCEED context
396 // ------------------------------------------------------
397 // -- Scale variables to desired units
398 cv *= J_per_kg_K;
399 cp *= J_per_kg_K;
400 mu *= Pascal * second;
401 k *= W_per_m_K;
402 for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter;
403 for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s;
404 reference.pressure *= Pascal;
405 for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second;
406 reference.temperature *= Kelvin;
407 problem->dm_scale = meter;
408
409 // -- Solver Settings
410 user->phys->implicit = implicit;
411 user->phys->state_var = state_var;
412
413 // -- QFunction Context
414 newtonian_ig_ctx->lambda = lambda;
415 newtonian_ig_ctx->mu = mu;
416 newtonian_ig_ctx->k = k;
417 newtonian_ig_ctx->cv = cv;
418 newtonian_ig_ctx->cp = cp;
419 newtonian_ig_ctx->c_tau = c_tau;
420 newtonian_ig_ctx->Ctau_t = Ctau_t;
421 newtonian_ig_ctx->Ctau_v = Ctau_v;
422 newtonian_ig_ctx->Ctau_C = Ctau_C;
423 newtonian_ig_ctx->Ctau_M = Ctau_M;
424 newtonian_ig_ctx->Ctau_E = Ctau_E;
425 newtonian_ig_ctx->stabilization = stab;
426 newtonian_ig_ctx->is_implicit = implicit;
427 newtonian_ig_ctx->state_var = state_var;
428 newtonian_ig_ctx->idl_enable = idl_enable;
429 newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * second);
430 newtonian_ig_ctx->idl_start = idl_start * meter;
431 newtonian_ig_ctx->idl_length = idl_length * meter;
432 newtonian_ig_ctx->idl_pressure = idl_pressure;
433 PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3));
434
435 // -- Setup Context
436 setup_context->reference = reference;
437 setup_context->gas = *newtonian_ig_ctx;
438 setup_context->lx = domain_size[0];
439 setup_context->ly = domain_size[1];
440 setup_context->lz = domain_size[2];
441 setup_context->time = 0;
442
443 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context));
444 PetscCallCeed(ceed,
445 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
446 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc));
447 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, "evaluation time", offsetof(struct SetupContext_, time), 1,
448 "Time of evaluation"));
449
450 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context));
451 PetscCallCeed(ceed,
452 CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx));
453 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, FreeContextPetsc));
454 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1,
455 "Size of timestep, delta t"));
456 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift",
457 offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 1,
458 "Shift for mass matrix in IJacobian"));
459 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(newtonian_ig_context, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1,
460 "Current solution time"));
461
462 problem->apply_vol_rhs.qfunction_context = newtonian_ig_context;
463 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ifunction.qfunction_context));
464 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ijacobian.qfunction_context));
465 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow.qfunction_context));
466 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow_jacobian.qfunction_context));
467
468 if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
469 if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
470 if (bc->num_slip > 0) PetscCall(SlipBCSetup(problem, dm, ctx, newtonian_ig_context));
471
472 if (unit_tests) {
473 PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx));
474 }
475 PetscFunctionReturn(PETSC_SUCCESS);
476 }
477
PRINT_NEWTONIAN(User user,ProblemData problem,AppCtx app_ctx)478 PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx) {
479 MPI_Comm comm = user->comm;
480 Ceed ceed = user->ceed;
481 NewtonianIdealGasContext newtonian_ctx;
482
483 PetscFunctionBeginUser;
484 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ctx));
485 PetscCall(PetscPrintf(comm,
486 " Problem:\n"
487 " Problem Name : %s\n"
488 " Stabilization : %s\n",
489 app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization]));
490 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ctx));
491 PetscFunctionReturn(PETSC_SUCCESS);
492 }
493