xref: /libCEED/examples/fluids/problems/newtonian.c (revision 49aac155e7a09736f56fb3abac0f57dab29f7cbf)
1 // Copyright (c) 2017-2022, 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 #include "../qfunctions/setupgeo.h"
18 
19 // For use with PetscOptionsEnum
20 static const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "StateVariable", "STATEVAR_", NULL};
21 
22 // Compute relative error |a - b|/|s|
23 static PetscErrorCode CheckPrimitiveWithTolerance(StatePrimitive sY, StatePrimitive aY, StatePrimitive bY, const char *name, PetscReal rtol_pressure,
24                                                   PetscReal rtol_velocity, PetscReal rtol_temperature) {
25   PetscFunctionBeginUser;
26   StatePrimitive eY;  // relative error
27   eY.pressure   = (aY.pressure - bY.pressure) / sY.pressure;
28   PetscScalar u = sqrt(Square(sY.velocity[0]) + Square(sY.velocity[1]) + Square(sY.velocity[2]));
29   for (int j = 0; j < 3; j++) eY.velocity[j] = (aY.velocity[j] - bY.velocity[j]) / u;
30   eY.temperature = (aY.temperature - bY.temperature) / sY.temperature;
31   if (fabs(eY.pressure) > rtol_pressure) printf("%s: pressure error %g\n", name, eY.pressure);
32   for (int j = 0; j < 3; j++) {
33     if (fabs(eY.velocity[j]) > rtol_velocity) printf("%s: velocity[%d] error %g\n", name, j, eY.velocity[j]);
34   }
35   if (fabs(eY.temperature) > rtol_temperature) printf("%s: temperature error %g\n", name, eY.temperature);
36   PetscFunctionReturn(0);
37 }
38 
39 static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) {
40   Units            units = user->units;
41   const CeedScalar eps   = 1e-6;
42   const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, Pascal = units->Pascal;
43   PetscFunctionBeginUser;
44   const CeedScalar rho = 1.2 * kg / (m * m * m), u = 40 * m / sec;
45   CeedScalar       U[5] = {rho, rho * u, rho * u * 1.1, rho * u * 1.2, 250e3 * Pascal + .5 * rho * u * u};
46   const CeedScalar x[3] = {.1, .2, .3};
47   State            s    = StateFromU(gas, U, x);
48   for (int i = 0; i < 8; i++) {
49     CeedScalar dU[5] = {0}, dx[3] = {0};
50     if (i < 5) dU[i] = U[i];
51     else dx[i - 5] = x[i - 5];
52     State ds = StateFromU_fwd(gas, s, dU, x, dx);
53     for (int j = 0; j < 5; j++) dU[j] = (1 + eps * (i == j)) * U[j];
54     for (int j = 0; j < 3; j++) dx[j] = (1 + eps * (i == 5 + j)) * x[j];
55     State          t = StateFromU(gas, dU, dx);
56     StatePrimitive dY;
57     dY.pressure = (t.Y.pressure - s.Y.pressure) / eps;
58     for (int j = 0; j < 3; j++) dY.velocity[j] = (t.Y.velocity[j] - s.Y.velocity[j]) / eps;
59     dY.temperature = (t.Y.temperature - s.Y.temperature) / eps;
60     char buf[128];
61     snprintf(buf, sizeof buf, "StateFromU_fwd i=%d", i);
62     PetscCall(CheckPrimitiveWithTolerance(dY, ds.Y, dY, buf, 5e-6, 1e-6, 1e-6));
63   }
64   PetscFunctionReturn(0);
65 }
66 
67 PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
68   SetupContext             setup_context;
69   User                     user   = *(User *)ctx;
70   CeedInt                  degree = user->app_ctx->degree;
71   StabilizationType        stab;
72   StateVariable            state_var;
73   MPI_Comm                 comm = PETSC_COMM_WORLD;
74   PetscBool                implicit;
75   PetscBool                has_curr_time = PETSC_FALSE, unit_tests;
76   NewtonianIdealGasContext newtonian_ig_ctx;
77   CeedQFunctionContext     newtonian_ig_context;
78 
79   PetscFunctionBeginUser;
80   PetscCall(PetscCalloc1(1, &setup_context));
81   PetscCall(PetscCalloc1(1, &newtonian_ig_ctx));
82 
83   // ------------------------------------------------------
84   //           Setup Generic Newtonian IG Problem
85   // ------------------------------------------------------
86   problem->dim                     = 3;
87   problem->q_data_size_vol         = 10;
88   problem->q_data_size_sur         = 10;
89   problem->jac_data_size_sur       = 11;
90   problem->setup_vol.qfunction     = Setup;
91   problem->setup_vol.qfunction_loc = Setup_loc;
92   problem->setup_sur.qfunction     = SetupBoundary;
93   problem->setup_sur.qfunction_loc = SetupBoundary_loc;
94   problem->bc                      = NULL;
95   problem->bc_ctx                  = NULL;
96   problem->non_zero_time           = PETSC_FALSE;
97   problem->print_info              = PRINT_NEWTONIAN;
98 
99   // ------------------------------------------------------
100   //             Create the libCEED context
101   // ------------------------------------------------------
102   CeedScalar cv         = 717.;           // J/(kg K)
103   CeedScalar cp         = 1004.;          // J/(kg K)
104   CeedScalar g[3]       = {0, 0, -9.81};  // m/s^2
105   CeedScalar lambda     = -2. / 3.;       // -
106   CeedScalar mu         = 1.8e-5;         // Pa s, dynamic viscosity
107   CeedScalar k          = 0.02638;        // W/(m K)
108   CeedScalar c_tau      = 0.5 / degree;   // -
109   CeedScalar Ctau_t     = 1.0;            // -
110   CeedScalar Cv_func[3] = {36, 60, 128};
111   CeedScalar Ctau_v     = Cv_func[(CeedInt)Min(3, degree) - 1];
112   CeedScalar Ctau_C     = 0.25 / degree;
113   CeedScalar Ctau_M     = 0.25 / degree;
114   CeedScalar Ctau_E     = 0.125;
115   PetscReal  domain_min[3], domain_max[3], domain_size[3];
116   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
117   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
118 
119   StatePrimitive reference      = {.pressure = 1.01e5, .velocity = {0}, .temperature = 288.15};
120   CeedScalar     idl_decay_time = -1, idl_start = 0, idl_length = 0;
121   PetscBool      idl_enable = PETSC_FALSE;
122 
123   // ------------------------------------------------------
124   //             Create the PETSc context
125   // ------------------------------------------------------
126   PetscScalar meter    = 1;  // 1 meter in scaled length units
127   PetscScalar kilogram = 1;  // 1 kilogram in scaled mass units
128   PetscScalar second   = 1;  // 1 second in scaled time units
129   PetscScalar Kelvin   = 1;  // 1 Kelvin in scaled temperature units
130   PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s;
131 
132   // ------------------------------------------------------
133   //              Command line Options
134   // ------------------------------------------------------
135   PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL);
136   // -- Conservative vs Primitive variables
137   PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE),
138                              (PetscEnum *)&state_var, NULL));
139 
140   switch (state_var) {
141     case STATEVAR_CONSERVATIVE:
142       problem->ics.qfunction                       = ICsNewtonianIG_Conserv;
143       problem->ics.qfunction_loc                   = ICsNewtonianIG_Conserv_loc;
144       problem->apply_vol_rhs.qfunction             = RHSFunction_Newtonian;
145       problem->apply_vol_rhs.qfunction_loc         = RHSFunction_Newtonian_loc;
146       problem->apply_vol_ifunction.qfunction       = IFunction_Newtonian_Conserv;
147       problem->apply_vol_ifunction.qfunction_loc   = IFunction_Newtonian_Conserv_loc;
148       problem->apply_vol_ijacobian.qfunction       = IJacobian_Newtonian_Conserv;
149       problem->apply_vol_ijacobian.qfunction_loc   = IJacobian_Newtonian_Conserv_loc;
150       problem->apply_inflow.qfunction              = BoundaryIntegral_Conserv;
151       problem->apply_inflow.qfunction_loc          = BoundaryIntegral_Conserv_loc;
152       problem->apply_inflow_jacobian.qfunction     = BoundaryIntegral_Jacobian_Conserv;
153       problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc;
154       break;
155 
156     case STATEVAR_PRIMITIVE:
157       problem->ics.qfunction                       = ICsNewtonianIG_Prim;
158       problem->ics.qfunction_loc                   = ICsNewtonianIG_Prim_loc;
159       problem->apply_vol_ifunction.qfunction       = IFunction_Newtonian_Prim;
160       problem->apply_vol_ifunction.qfunction_loc   = IFunction_Newtonian_Prim_loc;
161       problem->apply_vol_ijacobian.qfunction       = IJacobian_Newtonian_Prim;
162       problem->apply_vol_ijacobian.qfunction_loc   = IJacobian_Newtonian_Prim_loc;
163       problem->apply_inflow.qfunction              = BoundaryIntegral_Prim;
164       problem->apply_inflow.qfunction_loc          = BoundaryIntegral_Prim_loc;
165       problem->apply_inflow_jacobian.qfunction     = BoundaryIntegral_Jacobian_Prim;
166       problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc;
167       break;
168   }
169 
170   // -- Physics
171   PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL));
172   PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL));
173   PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL));
174   PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL));
175   PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL));
176 
177   PetscInt dim = problem->dim;
178   PetscCall(PetscOptionsRealArray("-g", "Gravitational acceleration", NULL, g, &dim, NULL));
179   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
180   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
181   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
182   PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL));
183   PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL));
184   PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL));
185   PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL));
186   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
187   PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL));
188 
189   PetscCall(PetscOptionsScalar("-reference_pressure", "Reference/initial pressure", NULL, reference.pressure, &reference.pressure, NULL));
190   PetscCall(PetscOptionsScalarArray("-reference_velocity", "Reference/initial velocity", NULL, reference.velocity, &dim, NULL));
191   PetscCall(PetscOptionsScalar("-reference_temperature", "Reference/initial temperature", NULL, reference.temperature, &reference.temperature, NULL));
192 
193   // -- Units
194   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
195   meter = fabs(meter);
196   PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL));
197   kilogram = fabs(kilogram);
198   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
199   second = fabs(second);
200   PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL));
201   Kelvin = fabs(Kelvin);
202 
203   // -- Warnings
204   if (stab == STAB_SUPG && !implicit) {
205     PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n"));
206   }
207   if (state_var == STATEVAR_PRIMITIVE && !implicit) {
208     SETERRQ(comm, PETSC_ERR_SUP, "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n");
209   }
210 
211   PetscCall(PetscOptionsScalar("-idl_decay_time", "Characteristic timescale of the pressure deviance decay. The timestep is good starting point",
212                                NULL, idl_decay_time, &idl_decay_time, &idl_enable));
213   if (idl_enable && idl_decay_time == 0) SETERRQ(comm, PETSC_ERR_SUP, "idl_decay_time may not be equal to zero.");
214   else if (idl_decay_time < 0) idl_enable = PETSC_FALSE;
215   PetscCall(PetscOptionsScalar("-idl_start", "Start of IDL in the x direction", NULL, idl_start, &idl_start, NULL));
216   PetscCall(PetscOptionsScalar("-idl_length", "Length of IDL in the positive x direction", NULL, idl_length, &idl_length, NULL));
217   PetscOptionsEnd();
218 
219   // ------------------------------------------------------
220   //           Set up the PETSc context
221   // ------------------------------------------------------
222   // -- Define derived units
223   Pascal          = kilogram / (meter * PetscSqr(second));
224   J_per_kg_K      = PetscSqr(meter) / (PetscSqr(second) * Kelvin);
225   m_per_squared_s = meter / PetscSqr(second);
226   W_per_m_K       = kilogram * meter / (pow(second, 3) * Kelvin);
227 
228   user->units->meter           = meter;
229   user->units->kilogram        = kilogram;
230   user->units->second          = second;
231   user->units->Kelvin          = Kelvin;
232   user->units->Pascal          = Pascal;
233   user->units->J_per_kg_K      = J_per_kg_K;
234   user->units->m_per_squared_s = m_per_squared_s;
235   user->units->W_per_m_K       = W_per_m_K;
236 
237   // ------------------------------------------------------
238   //           Set up the libCEED context
239   // ------------------------------------------------------
240   // -- Scale variables to desired units
241   cv *= J_per_kg_K;
242   cp *= J_per_kg_K;
243   mu *= Pascal * second;
244   k *= W_per_m_K;
245   for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter;
246   for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s;
247   reference.pressure *= Pascal;
248   for (PetscInt i = 0; i < 3; i++) reference.velocity[i] *= meter / second;
249   reference.temperature *= Kelvin;
250   problem->dm_scale = meter;
251 
252   // -- Solver Settings
253   user->phys->stab          = stab;
254   user->phys->implicit      = implicit;
255   user->phys->state_var     = state_var;
256   user->phys->has_curr_time = has_curr_time;
257 
258   // -- QFunction Context
259   newtonian_ig_ctx->lambda        = lambda;
260   newtonian_ig_ctx->mu            = mu;
261   newtonian_ig_ctx->k             = k;
262   newtonian_ig_ctx->cv            = cv;
263   newtonian_ig_ctx->cp            = cp;
264   newtonian_ig_ctx->c_tau         = c_tau;
265   newtonian_ig_ctx->Ctau_t        = Ctau_t;
266   newtonian_ig_ctx->Ctau_v        = Ctau_v;
267   newtonian_ig_ctx->Ctau_C        = Ctau_C;
268   newtonian_ig_ctx->Ctau_M        = Ctau_M;
269   newtonian_ig_ctx->Ctau_E        = Ctau_E;
270   newtonian_ig_ctx->P0            = reference.pressure;
271   newtonian_ig_ctx->stabilization = stab;
272   newtonian_ig_ctx->P0            = reference.pressure;
273   newtonian_ig_ctx->is_implicit   = implicit;
274   newtonian_ig_ctx->state_var     = state_var;
275   newtonian_ig_ctx->idl_enable    = idl_enable;
276   newtonian_ig_ctx->idl_amplitude = 1 / (idl_decay_time * second);
277   newtonian_ig_ctx->idl_start     = idl_start * meter;
278   newtonian_ig_ctx->idl_length    = idl_length * meter;
279   PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3));
280 
281   // -- Setup Context
282   setup_context->reference = reference;
283   setup_context->gas       = *newtonian_ig_ctx;
284   setup_context->lx        = domain_size[0];
285   setup_context->ly        = domain_size[1];
286   setup_context->lz        = domain_size[2];
287   setup_context->time      = 0;
288 
289   if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
290   if (bc->num_outflow > 0) PetscCall(OutflowBCSetup(problem, dm, ctx, newtonian_ig_ctx, &reference));
291 
292   CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context);
293   CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context);
294   CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc);
295   CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, "evaluation time", offsetof(struct SetupContext_, time), 1,
296                                      "Time of evaluation");
297 
298   CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context);
299   CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx);
300   CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, FreeContextPetsc);
301   CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1,
302                                      "Size of timestep, delta t");
303   CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift", offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift),
304                                      1, "Shift for mass matrix in IJacobian");
305   CeedQFunctionContextRegisterDouble(newtonian_ig_context, "solution time", offsetof(struct NewtonianIdealGasContext_, time), 1,
306                                      "Current solution time");
307 
308   problem->apply_vol_rhs.qfunction_context = newtonian_ig_context;
309   CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ifunction.qfunction_context);
310   CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ijacobian.qfunction_context);
311   CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow.qfunction_context);
312   CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow_jacobian.qfunction_context);
313 
314   if (unit_tests) {
315     PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx));
316   }
317   PetscFunctionReturn(0);
318 }
319 
320 PetscErrorCode PRINT_NEWTONIAN(ProblemData *problem, AppCtx app_ctx) {
321   MPI_Comm                 comm = PETSC_COMM_WORLD;
322   NewtonianIdealGasContext newtonian_ctx;
323 
324   PetscFunctionBeginUser;
325   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ctx);
326   PetscCall(PetscPrintf(comm,
327                         "  Problem:\n"
328                         "    Problem Name                       : %s\n"
329                         "    Stabilization                      : %s\n",
330                         app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization]));
331   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ctx);
332   PetscFunctionReturn(0);
333 }
334