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