3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
freeflow/rans/twoeq/sst/volumevariables.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_SST_VOLUME_VARIABLES_HH
26#define DUMUX_SST_VOLUME_VARIABLES_HH
27
32
33namespace Dumux {
34
39template <class Traits, class NSVolumeVariables>
41: public RANSVolumeVariables<Traits, NSVolumeVariables>
42{
44
45 using Scalar = typename Traits::PrimaryVariables::value_type;
46 using DimVector = Dune::FieldVector<Scalar, Traits::ModelTraits::dim()>;
47
48public:
50 using Indices = typename Traits::ModelTraits::Indices;
51
60 template<class ElementSolution, class Problem, class Element, class SubControlVolume>
61 void update(const ElementSolution &elemSol,
62 const Problem &problem,
63 const Element &element,
64 const SubControlVolume& scv)
65 {
66 RANSParentType::updateNavierStokesVolVars(elemSol, problem, element, scv);
67 updateRANSProperties(elemSol, problem, element, scv);
68 }
69
80 template<class ElementSolution, class Problem, class Element, class SubControlVolume>
81 void updateRANSProperties(const ElementSolution &elemSol,
82 const Problem &problem,
83 const Element &element,
84 const SubControlVolume& scv)
85 {
86 RANSParentType::updateRANSProperties(elemSol, problem, element, scv);
87
88 turbulentKineticEnergy_ = elemSol[0][Indices::turbulentKineticEnergyIdx];
89 dissipation_ = elemSol[0][Indices::dissipationIdx];
90
91 storedDissipation_ = problem.storedDissipation(RANSParentType::elementIdx());
92 storedTurbulentKineticEnergy_ = problem.storedTurbulentKineticEnergy(RANSParentType::elementIdx());
93
94 storedDissipationGradient_ = problem.storedDissipationGradient(RANSParentType::elementIdx());
95 storedTurbulentKineticEnergyGradient_ = problem.storedTurbulentKineticEnergyGradient(RANSParentType::elementIdx());
96
97 stressTensorScalarProduct_ = problem.stressTensorScalarProduct(RANSParentType::elementIdx());
98 vorticityTensorScalarProduct_ = problem.vorticityTensorScalarProduct(RANSParentType::elementIdx());
99 wallDistance_ = problem.wallDistance(RANSParentType::elementIdx());
100 kinematicViscosity_ = problem.kinematicViscosity(RANSParentType::elementIdx());
101
102 if (problem.useStoredEddyViscosity())
103 RANSParentType::setDynamicEddyViscosity_(problem.storedDynamicEddyViscosity(RANSParentType::elementIdx()));
104 else
106
109 }
110
114 template<class Problem>
115 Scalar calculateEddyViscosity(const Problem& problem)
116 {
117 using std::sqrt;
118 using std::max;
119
120 if(problem.sstModelVersion() == SSTModel::BSL)
122 else if(problem.sstModelVersion() == SSTModel::SST)
123 {
124 const Scalar dividend = a1SST() * turbulentKineticEnergy();
125 const Scalar possibleMax1 = a1SST()*dissipation();
126 const Scalar possibleMax2 = absoluteValueVorticity() * F2();
127 const Scalar divisor = std::max(possibleMax1,possibleMax2);
128 return dividend / divisor * RANSParentType::density();
129 }
130 else
131 DUNE_THROW(Dune::NotImplemented, "\nThis SST Model is not implemented.\n");
132 }
133
134
137 { return turbulentKineticEnergy_; }
138
140 Scalar dissipation() const
141 { return dissipation_; }
142
146
148 Scalar storedDissipation() const
149 { return storedDissipation_; }
150
154
158
162
166
168 Scalar kinematicViscosity() const
169 { return kinematicViscosity_; }
170
171 Scalar wallDistance() const
172 { return wallDistance_; }
173
175 const Scalar betaOmega() const
176 { return 0.0708; }
177
179 const Scalar alpha() const
180 { return 0.52; }
181
183 const Scalar sigmaK() const
184 { return 0.6; }
185
187 const Scalar sigmaOmega() const
188 { return 0.5; }
189
191 const Scalar betaK() const
192 { return 0.09; }
193
195 const Scalar sigmaK1BSL() const
196 { return 0.5; }
197
199 const Scalar sigmaK1SST() const
200 { return 0.85; }
201
203 const Scalar sigmaK2() const
204 { return 1.0; }
205
207 const Scalar sigmaOmega1BSL() const
208 { return 0.5; }
209
211 const Scalar sigmaOmega1SST() const
212 { return 0.5; }
213
215 const Scalar sigmaOmega2() const
216 { return 0.856; }
217
219 const Scalar beta1BSL() const
220 { return 0.0750; }
221
223 const Scalar beta1SST() const
224 { return 0.0750; }
225
227 const Scalar beta2() const
228 { return 0.0820; }
229
231 const Scalar betaStar1BSL() const
232 { return 0.09; }
233
235 const Scalar betaStar1SST() const
236 { return 0.09; }
237
239 const Scalar betaStar2() const
240 { return 0.09; }
241
243 const Scalar kappa1BSL() const
244 { return 0.41; }
245
247 const Scalar kappa1SST() const
248 { return 0.41; }
249
251 const Scalar kappa2() const
252 { return 0.41; }
253
255 const Scalar gamma1BSL() const
256 { return (beta1BSL()/betaStar1BSL()) - ((sigmaOmega1BSL()*kappa1BSL()*kappa1BSL())/(std::sqrt(betaStar1BSL()))); }
257
259 const Scalar gamma1SST() const
260 { return (beta1SST()/betaStar1SST()) - ((sigmaOmega1SST()*kappa1SST()*kappa1SST())/(std::sqrt(betaStar1SST()))); }
261
263 const Scalar gamma2() const
264 { return (beta2()/betaStar2()) - ((sigmaOmega2()*kappa2()*kappa2())/(std::sqrt(betaStar2()))); }
265
267 const Scalar a1SST() const
268 { return 0.31; }
269
272 { return std::sqrt(2.0 * vorticityTensorScalarProduct()); }
273
275 Scalar F1() const
276 {
277 Scalar gradientProduct = 0.0;
278 for (unsigned int i = 0; i < Traits::ModelTraits::dim(); ++i){
280
281 Scalar positiveCrossDiffusion = 2.0 * RANSParentType::density() * sigmaOmega2() / storedDissipation() * gradientProduct;
282 positiveCrossDiffusion = std::max(positiveCrossDiffusion, 1e-20);
283
284 const Scalar possibleMin2 = (4.0 * RANSParentType::density() * sigmaOmega2() * storedTurbulentKineticEnergy())
285 / (positiveCrossDiffusion * wallDistance() * wallDistance());
286
287 const Scalar possibleMax1 = (std::sqrt(storedTurbulentKineticEnergy())) / (0.09 * storedDissipation() * wallDistance());
288 const Scalar possibleMax2 = (500.00 * kinematicViscosity()) / (wallDistance() * wallDistance() * storedDissipation());
289 const Scalar possibleMin1 = std::max(possibleMax1, possibleMax2);
290
291 const Scalar argument = std::min(possibleMin1,possibleMin2);
292
293 return tanh(argument * argument * argument * argument);
294 }
295
296
298 Scalar F2() const
299 {
300 const Scalar possibleMax1 = 2.0 * (std::sqrt(storedTurbulentKineticEnergy())) / (0.09 * storedDissipation() * wallDistance());
301 const Scalar possibleMax2 = (500.0 * kinematicViscosity()) / (wallDistance() * wallDistance() * storedDissipation());
302 const Scalar argument = std::max(possibleMax1, possibleMax2);
303
304 return tanh(argument * argument);
305 }
306
307 //------------------------Constants for SST--------------------------------
309 const Scalar sigmaKSST() const
310 { return F1()*sigmaK1SST() + (1-F1())*sigmaK2(); }
311
313 const Scalar sigmaOmegaSST() const
314 { return F1()*sigmaOmega1SST() + (1-F1())*sigmaOmega2(); }
315
317 const Scalar betaSST() const
318 { return F1()*beta1SST() + (1-F1())*beta2(); }
319
321 const Scalar betaStarSST() const
322 { return F1()*betaStar1SST() + (1-F1())*betaStar2(); }
323
325 const Scalar kappaSST() const
326 { return F1()*kappa1SST() + (1-F1())*kappa2(); }
327
329 const Scalar gammaSST() const
330 { return F1()*gamma1SST() + (1-F1())*gamma2(); }
331
332
333
334 //------------------------Constants for BSL--------------------------------
336 const Scalar sigmaKBSL() const
337 { return F1()*sigmaK1BSL() + (1-F1())*sigmaK2(); }
338
340 const Scalar sigmaOmegaBSL() const
341 { return F1()*sigmaOmega1BSL() + (1-F1())*sigmaOmega2(); }
342
344 const Scalar betaBSL() const
345 { return F1()*beta1BSL() + (1-F1())*beta2(); }
346
348 const Scalar betaStarBSL() const
349 { return F1()*betaStar1BSL() + (1-F1())*betaStar2(); }
350
352 const Scalar kappaBSL() const
353 { return F1()*kappa1BSL() + (1-F1())*kappa2(); }
354
356 const Scalar gammaBSL() const
357 { return F1()*gamma1BSL() + (1-F1())*gamma2(); }
358
359
360protected:
361 Scalar betaOmega_ = 0.0;
362
363 Scalar dissipation_ = 0.0;
364 Scalar storedDissipation_ = 0.0;
365 DimVector storedDissipationGradient_ = DimVector(0.0);
368 DimVector storedTurbulentKineticEnergyGradient_ = DimVector(0.0);
369
372 Scalar wallDistance_ = 0.0;
374};
375}
376
377#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available free flow turbulence models in Dumux.
Definition: adapt.hh:29
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Volume variables for the isothermal single-phase SST 2-Eq model.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:42
const Scalar kappa1SST() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:247
Scalar storedDissipation_
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:364
const Scalar sigmaK1SST() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:199
Scalar storedTurbulentKineticEnergy() const
Returns the turbulent kinetic energy .
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:144
const Scalar gamma1SST() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:259
const Scalar sigmaOmega1SST() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:211
const Scalar gamma1BSL() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:255
const Scalar gammaBSL() const
Returns the constant for the BSL-model.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:356
const Scalar sigmaOmega() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:187
Scalar wallDistance() const
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:171
void updateRANSProperties(const ElementSolution &elemSol, const Problem &problem, const Element &element, const SubControlVolume &scv)
Update all turbulent quantities for a given control volume.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:81
const Scalar beta2() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:227
Scalar absoluteValueVorticity() const
Returns the absolute value of the vorticity .
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:271
typename Traits::ModelTraits::Indices Indices
export the indices type
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:50
Scalar vorticityTensorScalarProduct_
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:371
const Scalar sigmaOmegaSST() const
Returns the constant for the SST-model.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:313
const Scalar sigmaK2() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:203
const Scalar betaStar1BSL() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:231
const Scalar betaStarSST() const
Returns the constant for the SST-model. is the same for all models.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:321
Scalar F1() const
Returns the transformation function for the constants of the BSL- and SST-model.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:275
const Scalar kappa1BSL() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:243
Scalar dissipation_
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:363
const Scalar betaK() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:191
DimVector storedTurbulentKineticEnergyGradient() const
Returns the gradient of the turbulent kinetic energy .
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:152
const Scalar betaBSL() const
Returns the constant for the BSL-model.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:344
const Scalar kappa2() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:251
const Scalar betaSST() const
Returns the constant for the SST-model.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:317
const Scalar betaStar1SST() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:235
Scalar F2() const
Returns the transformation function for the eddy viscosity of the SST-model.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:298
const Scalar betaStar2() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:239
Scalar storedTurbulentKineticEnergy_
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:367
DimVector storedDissipationGradient() const
Returns the gradient of the effective dissipation .
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:156
Scalar calculateEddyViscosity(const Problem &problem)
Returns the dynamic eddy viscosity for the SST-model.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:115
const Scalar kappaSST() const
Returns the constant for the SST-model. is the same for all models.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:325
Scalar turbulentKineticEnergy() const
Returns the turbulent kinetic energy .
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:136
const Scalar betaOmega() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:175
Scalar vorticityTensorScalarProduct() const
Returns the scalar product of the vorticity tensor.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:164
const Scalar sigmaK() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:183
const Scalar beta1BSL() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:219
Scalar betaOmega_
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:361
const Scalar sigmaKBSL() const
Returns the constant for the BSL-model.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:336
Scalar kinematicViscosity() const
Returns the kinematic viscosity.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:168
Scalar dissipation() const
Returns an effective dissipation .
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:140
Scalar kinematicViscosity_
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:373
Scalar stressTensorScalarProduct() const
Returns the scalar product of the stress tensor.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:160
Scalar turbulentKineticEnergy_
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:366
const Scalar a1SST() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:267
const Scalar sigmaK1BSL() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:195
Scalar storedDissipation() const
Returns an effective dissipation .
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:148
const Scalar betaStarBSL() const
Returns the constant for the BSL-model. is the same for all models.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:348
Scalar stressTensorScalarProduct_
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:370
const Scalar sigmaOmegaBSL() const
Returns the constant for the BSL-model.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:340
void update(const ElementSolution &elemSol, const Problem &problem, const Element &element, const SubControlVolume &scv)
Update all quantities for a given control volume.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:61
const Scalar sigmaOmega1BSL() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:207
const Scalar gammaSST() const
Returns the constant for the SST-model.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:329
DimVector storedTurbulentKineticEnergyGradient_
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:368
const Scalar gamma2() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:263
const Scalar sigmaOmega2() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:215
const Scalar sigmaKSST() const
Returns the constant for the SST-model.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:309
const Scalar alpha() const
Returns the value.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:179
const Scalar kappaBSL() const
Returns the constant for the BSL-model. is the same for all models.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:352
const Scalar beta1SST() const
Returns the constant.
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:223
DimVector storedDissipationGradient_
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:365
Scalar wallDistance_
Definition: freeflow/rans/twoeq/sst/volumevariables.hh:372
Volume variables for the isothermal single-phase Reynolds-Averaged Navier-Stokes models.
Definition: freeflow/rans/volumevariables.hh:42
void calculateEddyThermalConductivity(const Problem &problem)
Calculates the eddy thermal conductivity based on the kinematic eddy viscosity and the turbulent Pra...
Definition: freeflow/rans/volumevariables.hh:221
void updateNavierStokesVolVars(const ElementSolution &elemSol, const Problem &problem, const Element &element, const SubControlVolume &scv)
Update all quantities for a given control volume.
Definition: freeflow/rans/volumevariables.hh:65
Scalar setDynamicEddyViscosity_(Scalar value)
Sets the dynamic eddy viscosity .
Definition: freeflow/rans/volumevariables.hh:266
void updateRANSProperties(const ElementSolution &elemSol, const Problem &problem, const Element &element, const SubControlVolume &scv)
Update all turbulent quantities for a given control volume.
Definition: freeflow/rans/volumevariables.hh:85
void calculateEddyDiffusivity(const Problem &problem)
Calculates the eddy diffusivity based on the kinematic eddy viscosity and the turbulent Schmidt numb...
Definition: freeflow/rans/volumevariables.hh:210
unsigned int elementIdx() const
Return the element Idx of the control volume.
Definition: freeflow/rans/volumevariables.hh:120
Declares all properties used in Dumux.