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