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