Loading...
Searching...
No Matches
kOmegaSSTLM.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2016-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "kOmegaSSTLM.H"
30#include "fvOptions.H"
31#include "bound.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
37namespace RASModels
38{
39
40// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41
42template<class BasicTurbulenceModel>
44(
45 const volScalarField& CDkOmega
46) const
47{
48 const volScalarField Ry(this->y_*sqrt(this->k_)/this->nu());
49 const volScalarField F3(exp(-pow(Ry/120.0, 8)));
50
52}
53
54
55template<class BasicTurbulenceModel>
57(
59) const
60{
62}
63
64
65template<class BasicTurbulenceModel>
67(
68 const volScalarField& F1,
69 const volTensorField& gradU
70) const
71{
72 return
75}
76
77
78template<class BasicTurbulenceModel>
80(
82 const volScalarField::Internal& Omega,
84) const
85{
86 const volScalarField::Internal& omega = this->omega_();
87 const volScalarField::Internal& y = this->y_();
88
89 dimensionedScalar deltaMin("deltaMin", dimLength, SMALL);
91 (
92 max(375*Omega*nu*ReThetat_()*y/sqr(Us), deltaMin)
93 );
94
95 const volScalarField::Internal ReOmega(sqr(y)*omega/nu);
96 const volScalarField::Internal Fwake(exp(-sqr(ReOmega/1e5)));
97
99 (
100 IOobject::groupName("Fthetat", this->alphaRhoPhi_.group()),
102 (
103 min
104 (
105 max
106 (
107 Fwake*exp(-pow4((y/delta))),
108 (1 - sqr((gammaInt_() - 1.0/ce2_)/(1 - 1.0/ce2_)))
109 ),
110 scalar(1)
111 )
113 );
114}
115
116
117template<class BasicTurbulenceModel>
120{
121 auto tReThetac = volScalarField::Internal::New
122 (
123 IOobject::groupName("ReThetac", this->alphaRhoPhi_.group()),
125 this->mesh_,
126 dimless
127 );
128 auto& ReThetac = tReThetac.ref();
129
130 forAll(ReThetac, celli)
131 {
132 const scalar ReThetat = ReThetat_[celli];
133
134 ReThetac[celli] =
135 ReThetat <= 1870
136 ?
137 ReThetat
138 - 396.035e-2
139 + 120.656e-4*ReThetat
140 - 868.230e-6*sqr(ReThetat)
141 + 696.506e-9*pow3(ReThetat)
142 - 174.105e-12*pow4(ReThetat)
143 :
144 ReThetat - 593.11 - 0.482*(ReThetat - 1870);
146
147 return tReThetac;
148}
149
150
151template<class BasicTurbulenceModel>
153(
155) const
156{
157 auto tFlength = volScalarField::Internal::New
158 (
159 IOobject::groupName("Flength", this->alphaRhoPhi_.group()),
161 this->mesh_,
162 dimless
163 );
164 auto& Flength = tFlength.ref();
165
166 const volScalarField::Internal& omega = this->omega_();
167 const volScalarField::Internal& y = this->y_();
168
169 forAll(ReThetat_, celli)
170 {
171 const scalar ReThetat = ReThetat_[celli];
172
173 if (ReThetat < 400)
174 {
175 Flength[celli] =
176 398.189e-1
177 - 119.270e-4*ReThetat
178 - 132.567e-6*sqr(ReThetat);
179 }
180 else if (ReThetat < 596)
181 {
182 Flength[celli] =
183 263.404
184 - 123.939e-2*ReThetat
185 + 194.548e-5*sqr(ReThetat)
186 - 101.695e-8*pow3(ReThetat);
187 }
188 else if (ReThetat < 1200)
189 {
190 Flength[celli] = 0.5 - 3e-4*(ReThetat - 596);
191 }
192 else
193 {
194 Flength[celli] = 0.3188;
195 }
196
197 const scalar Fsublayer =
198 exp(-sqr(sqr(y[celli])*omega[celli]/(200*nu[celli])));
199
200 Flength[celli] = Flength[celli]*(1 - Fsublayer) + 40*Fsublayer;
202
203 return tFlength;
204}
205
206
207template<class BasicTurbulenceModel>
209(
211 const volScalarField::Internal& dUsds,
213) const
214{
215 auto tReThetat0 = volScalarField::Internal::New
216 (
217 IOobject::groupName("ReThetat0", this->alphaRhoPhi_.group()),
219 this->mesh_,
220 dimless
221 );
222 auto& ReThetat0 = tReThetat0.ref();
223
224 const volScalarField& k = this->k_;
225
226 label maxIter = 0;
227
228 forAll(ReThetat0, celli)
229 {
230 const scalar Tu
231 (
232 max(100*sqrt((2.0/3.0)*k[celli])/Us[celli], scalar(0.027))
233 );
234
235 // Initialize lambda to zero.
236 // If lambda were cached between time-steps convergence would be faster
237 // starting from the previous time-step value.
238 scalar lambda = 0;
239
240 scalar lambdaErr;
241 scalar thetat;
242 label iter = 0;
243
244 do
245 {
246 // Previous iteration lambda for convergence test
247 const scalar lambda0 = lambda;
248
249 if (Tu <= 1.3)
250 {
251 const scalar Flambda =
252 dUsds[celli] <= 0
253 ?
254 1
255 - (
256 - 12.986*lambda
257 - 123.66*sqr(lambda)
258 - 405.689*pow3(lambda)
259 )*exp(-pow(Tu/1.5, 1.5))
260 :
261 1
262 + 0.275*(1 - exp(-35*lambda))
263 *exp(-Tu/0.5);
264
265 thetat =
266 (1173.51 - 589.428*Tu + 0.2196/sqr(Tu))
267 *Flambda*nu[celli]
268 /Us[celli];
269 }
270 else
271 {
272 const scalar Flambda =
273 dUsds[celli] <= 0
274 ?
275 1
276 - (
277 -12.986*lambda
278 -123.66*sqr(lambda)
279 -405.689*pow3(lambda)
280 )*exp(-pow(Tu/1.5, 1.5))
281 :
282 1
283 + 0.275*(1 - exp(-35*lambda))
284 *exp(-2*Tu);
285
286 thetat =
287 331.50*pow((Tu - 0.5658), -0.671)
288 *Flambda*nu[celli]/Us[celli];
289 }
290
291 lambda = sqr(thetat)/nu[celli]*dUsds[celli];
292 lambda = clamp(lambda, scalar(-0.1), scalar(0.1));
293
294 lambdaErr = mag(lambda - lambda0);
295
296 maxIter = max(maxIter, ++iter);
297
298 } while (lambdaErr > lambdaErr_);
299
300 ReThetat0[celli] = max(thetat*Us[celli]/nu[celli], scalar(20));
301 }
302
303 if (maxIter > maxLambdaIter_)
304 {
306 << "Number of lambda iterations exceeds maxLambdaIter("
307 << maxLambdaIter_ << ')'<< endl;
309
310 return tReThetat0;
311}
312
313
314template<class BasicTurbulenceModel>
316(
317 const volScalarField::Internal& Rev,
318 const volScalarField::Internal& ReThetac,
320) const
321{
322 const volScalarField::Internal Fonset1(Rev/(2.193*ReThetac));
323
324 const volScalarField::Internal Fonset2
325 (
326 min(max(Fonset1, pow4(Fonset1)), scalar(2))
327 );
328
329 const volScalarField::Internal Fonset3(max(1 - pow3(RT/2.5), scalar(0)));
330
332 (
333 IOobject::groupName("Fonset", this->alphaRhoPhi_.group()),
335 max(Fonset2 - Fonset3, scalar(0))
336 );
337}
338
339
340// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
341
342template<class BasicTurbulenceModel>
343kOmegaSSTLM<BasicTurbulenceModel>::kOmegaSSTLM
344(
345 const alphaField& alpha,
346 const rhoField& rho,
347 const volVectorField& U,
348 const surfaceScalarField& alphaRhoPhi,
349 const surfaceScalarField& phi,
350 const transportModel& transport,
351 const word& propertiesName,
352 const word& type
353)
354:
355 kOmegaSST<BasicTurbulenceModel>
356 (
357 alpha,
358 rho,
359 U,
360 alphaRhoPhi,
361 phi,
362 transport,
363 propertiesName,
365 ),
366
367 ca1_
368 (
369 dimensionedScalar::getOrAddToDict
370 (
371 "ca1",
372 this->coeffDict_,
373 2
374 )
375 ),
376 ca2_
377 (
378 dimensionedScalar::getOrAddToDict
379 (
380 "ca2",
381 this->coeffDict_,
382 0.06
383 )
384 ),
385 ce1_
386 (
387 dimensionedScalar::getOrAddToDict
388 (
389 "ce1",
390 this->coeffDict_,
391 1
392 )
393 ),
394 ce2_
395 (
396 dimensionedScalar::getOrAddToDict
397 (
398 "ce2",
399 this->coeffDict_,
400 50
401 )
402 ),
403 cThetat_
404 (
405 dimensionedScalar::getOrAddToDict
406 (
407 "cThetat",
408 this->coeffDict_,
409 0.03
410 )
411 ),
412 sigmaThetat_
413 (
414 dimensionedScalar::getOrAddToDict
415 (
416 "sigmaThetat",
417 this->coeffDict_,
418 2
419 )
420 ),
421 lambdaErr_
422 (
423 this->coeffDict_.template getOrDefault<scalar>("lambdaErr", 1e-6)
424 ),
425 maxLambdaIter_
426 (
427 this->coeffDict_.template getOrDefault<label>("maxLambdaIter", 10)
428 ),
429 deltaU_("deltaU", dimVelocity, SMALL),
430
431 ReThetat_
432 (
434 (
435 IOobject::groupName("ReThetat", alphaRhoPhi.group()),
436 this->runTime_.timeName(),
437 this->mesh_,
438 IOobject::MUST_READ,
439 IOobject::AUTO_WRITE,
440 IOobject::REGISTER
441 ),
442 this->mesh_
443 ),
444
445 gammaInt_
446 (
448 (
449 IOobject::groupName("gammaInt", alphaRhoPhi.group()),
450 this->runTime_.timeName(),
451 this->mesh_,
452 IOobject::MUST_READ,
453 IOobject::AUTO_WRITE,
454 IOobject::REGISTER
455 ),
456 this->mesh_
457 ),
458
459 gammaIntEff_
460 (
462 (
463 IOobject::groupName("gammaIntEff", alphaRhoPhi.group()),
464 this->runTime_.timeName(),
465 this->mesh_
466 ),
467 this->mesh_,
469 )
470{
471 if (type == typeName)
472 {
473 this->printCoeffs(type);
475}
476
477
478// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
479
480template<class BasicTurbulenceModel>
482{
484 {
485 ca1_.readIfPresent(this->coeffDict());
486 ca2_.readIfPresent(this->coeffDict());
487 ce1_.readIfPresent(this->coeffDict());
488 ce2_.readIfPresent(this->coeffDict());
489 sigmaThetat_.readIfPresent(this->coeffDict());
490 cThetat_.readIfPresent(this->coeffDict());
491 this->coeffDict().readIfPresent("lambdaErr", lambdaErr_);
492 this->coeffDict().readIfPresent("maxLambdaIter", maxLambdaIter_);
493
494 return true;
496
497 return false;
498}
499
500
501template<class BasicTurbulenceModel>
503{
504 // Local references
505 const alphaField& alpha = this->alpha_;
506 const rhoField& rho = this->rho_;
507 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
508 const volVectorField& U = this->U_;
509 const volScalarField& k = this->k_;
510 const volScalarField& omega = this->omega_;
511 const tmp<volScalarField> tnu = this->nu();
512 const volScalarField::Internal& nu = tnu()();
513 const volScalarField::Internal& y = this->y_();
515
516 // Fields derived from the velocity gradient
518 const volScalarField::Internal Omega(sqrt(2*magSqr(skew(tgradU()()))));
519 const volScalarField::Internal S(sqrt(2*magSqr(symm(tgradU()()))));
520 const volScalarField::Internal Us(max(mag(U()), deltaU_));
521 const volScalarField::Internal dUsds((U() & (U() & tgradU()()))/sqr(Us));
522 tgradU.clear();
523
524 const volScalarField::Internal Fthetat(this->Fthetat(Us, Omega, nu));
525
526 {
527 const volScalarField::Internal t(500*nu/sqr(Us));
528 const volScalarField::Internal Pthetat
529 (
530 alpha()*rho()*(cThetat_/t)*(1 - Fthetat)
531 );
532
533 // Transition onset momentum-thickness Reynolds number equation
534 tmp<fvScalarMatrix> ReThetatEqn
535 (
536 fvm::ddt(alpha, rho, ReThetat_)
537 + fvm::div(alphaRhoPhi, ReThetat_)
538 - fvm::laplacian(alpha*rho*DReThetatEff(), ReThetat_)
539 ==
540 Pthetat*ReThetat0(Us, dUsds, nu) - fvm::Sp(Pthetat, ReThetat_)
541 + fvOptions(alpha, rho, ReThetat_)
542 );
543
544 ReThetatEqn.ref().relax();
545 fvOptions.constrain(ReThetatEqn.ref());
546 solve(ReThetatEqn);
547 fvOptions.correct(ReThetat_);
548 bound(ReThetat_, 0);
549 }
550
551 const volScalarField::Internal ReThetac(this->ReThetac());
552 const volScalarField::Internal Rev(sqr(y)*S/nu);
553 const volScalarField::Internal RT(k()/(nu*omega()));
554
555 {
556 const volScalarField::Internal Pgamma
557 (
558 alpha()*rho()
559 *ca1_*Flength(nu)*S*sqrt(gammaInt_()*Fonset(Rev, ReThetac, RT))
560 );
561
562 const volScalarField::Internal Fturb(exp(-pow4(0.25*RT)));
563
564 const volScalarField::Internal Egamma
565 (
566 alpha()*rho()*ca2_*Omega*Fturb*gammaInt_()
567 );
568
569 // Intermittency equation
570 tmp<fvScalarMatrix> gammaIntEqn
571 (
572 fvm::ddt(alpha, rho, gammaInt_)
573 + fvm::div(alphaRhoPhi, gammaInt_)
574 - fvm::laplacian(alpha*rho*DgammaIntEff(), gammaInt_)
575 ==
576 Pgamma - fvm::Sp(ce1_*Pgamma, gammaInt_)
577 + Egamma - fvm::Sp(ce2_*Egamma, gammaInt_)
578 + fvOptions(alpha, rho, gammaInt_)
579 );
580
581 gammaIntEqn.ref().relax();
582 fvOptions.constrain(gammaIntEqn.ref());
583 solve(gammaIntEqn);
584 fvOptions.correct(gammaInt_);
585 bound(gammaInt_, 0);
586 }
587
588 const volScalarField::Internal Freattach(exp(-pow4(RT/20.0)));
589 const volScalarField::Internal gammaSep
590 (
591 min(2*max(Rev/(3.235*ReThetac) - 1, scalar(0))*Freattach, scalar(2))
592 *Fthetat
593 );
594
595 gammaIntEff_ = max(gammaInt_(), gammaSep);
596}
597
598
599template<class BasicTurbulenceModel>
601{
602 if (!this->turbulence_)
603 {
604 return;
605 }
606
607 // Correct k and omega
609
610 // Correct ReThetat and gammaInt
611 correctReThetatGammaInt();
612}
613
614
615// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
616
617} // End namespace RASModels
618} // End namespace Foam
619
620// ************************************************************************* //
scalar y
label k
scalar delta
#define F1(B, C, D)
Definition SHA1.C:149
#define F3(B, C, D)
Definition SHA1.C:151
Bound the given scalar field if it has gone unbounded.
fv::options & fvOptions
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
DimensionedField< scalar, volMesh > Internal
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
label maxLambdaIter_
Maximum number of iterations to converge the lambda/thetat loop.
BasicTurbulenceModel::alphaField alphaField
tmp< volScalarField > DReThetatEff() const
Return the effective diffusivity for transition onset.
volScalarField gammaInt_
Intermittency.
BasicTurbulenceModel::rhoField rhoField
tmp< volScalarField::Internal > Flength(const volScalarField::Internal &nu) const
Empirical correlation that controls the length of the.
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
Modified form of the k-omega SST F1 function.
Definition kOmegaSSTLM.C:37
volScalarField::Internal gammaIntEff_
Effective intermittency.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField::Internal > Fthetat(const volScalarField::Internal &Us, const volScalarField::Internal &Omega, const volScalarField::Internal &nu) const
Freestream blending-function.
Definition kOmegaSSTLM.C:73
void correctReThetatGammaInt()
Solve the turbulence equations and correct the turbulence viscosity.
const dimensionedScalar deltaU_
Stabilization for division by the magnitude of the velocity.
tmp< volScalarField > DgammaIntEff() const
Return the effective diffusivity for intermittency.
tmp< volScalarField::Internal > ReThetac() const
Empirical correlation for critical Reynolds number where the.
BasicTurbulenceModel::transportModel transportModel
volScalarField ReThetat_
Transition onset momentum-thickness Reynolds number.
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField &F1, const volTensorField &gradU) const
Modified form of the k-omega SST epsilon/k.
Definition kOmegaSSTLM.C:60
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Modified form of the k-omega SST k production rate.
Definition kOmegaSSTLM.C:50
scalar lambdaErr_
Convergence criterion for the lambda/thetat loop.
tmp< volScalarField::Internal > Fonset(const volScalarField::Internal &Rev, const volScalarField::Internal &ReThetac, const volScalarField::Internal &RT) const
Transition onset location control function.
dimensionedScalar cThetat_
dimensionedScalar sigmaThetat_
tmp< volScalarField::Internal > ReThetat0(const volScalarField::Internal &Us, const volScalarField::Internal &dUsds, const volScalarField::Internal &nu) const
Return the transition onset momentum-thickness Reynolds number.
virtual bool read()
Re-read model coefficients if they have changed.
const volScalarField & ReThetat() const
Access function transition onset momentum-thickness Reynolds number.
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows.
Definition kOmegaSST.H:128
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
Definition fvOptions.H:69
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Definition fvOptions.C:116
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField &F1, const volTensorField &gradU) const
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
virtual bool read()=0
Re-read model coefficients if they have changed.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
word timeName
Definition getTimeIndex.H:3
#define WarningInFunction
Report a warning using Foam::Warning.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition bound.C:29
GeometricField< scalar, fvPatchField, volMesh > volScalarField
tensor Ry(const scalar omega)
Rotational transformation tensor about the y-axis by omega radians.
Definition transform.H:99
const dimensionSet dimVelocity
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
dimensionedScalar pow4(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedTensor skew(const dimensionedTensor &dt)
volScalarField & nu
volScalarField & alpha
CEqn solve()
volScalarField & e
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299