epsilonLowReWallFunctionFvPatchScalarField.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 | Copyright (C) 2012-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "turbulenceModel.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
31 
32 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
33 
35 (
36  const scalar kappa,
37  const scalar E
38 )
39 {
40  scalar ypl = 11.0;
41 
42  for (int i=0; i<10; i++)
43  {
44  ypl = log(max(E*ypl, 1))/kappa;
45  }
46 
47  return ypl;
48 }
49 
50 
52 (
53  const turbulenceModel& turbModel,
54  const List<scalar>& cornerWeights,
55  const fvPatch& patch,
56  scalarField& G0,
58 )
59 {
60  const label patchi = patch.index();
61 
62  const scalarField& y = turbModel.y()[patchi];
63 
64  const scalar Cmu25 = pow025(Cmu_);
65  const scalar Cmu75 = pow(Cmu_, 0.75);
66 
67  const tmp<volScalarField> tk = turbModel.k();
68  const volScalarField& k = tk();
69 
70  const tmp<scalarField> tnuw = turbModel.nu(patchi);
71  const scalarField& nuw = tnuw();
72 
73  const tmp<scalarField> tnutw = turbModel.nut(patchi);
74  const scalarField& nutw = tnutw();
75 
76  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
77 
78  const scalarField magGradUw(mag(Uw.snGrad()));
79 
81  db().lookupObject<DimensionedField<scalar, volMesh> >
82  (
83  turbModel.GName()
84  );
85 
86  // Set epsilon and G
87  forAll(nutw, facei)
88  {
89  label celli = patch.faceCells()[facei];
90 
91  scalar yPlus = Cmu25*sqrt(k[celli])*y[facei]/nuw[facei];
92 
93  scalar w = cornerWeights[facei];
94 
95  if (yPlus > yPlusLam_)
96  {
97  epsilon0[celli] += w*Cmu75*pow(k[celli], 1.5)/(kappa_*y[facei]);
98 
99  G0[celli] +=
100  w
101  *(nutw[facei] + nuw[facei])
102  *magGradUw[facei]
103  *Cmu25*sqrt(k[celli])
104  /(kappa_*y[facei]);
105  }
106  else
107  {
108  epsilon0[celli] += w*2.0*k[celli]*nuw[facei]/sqr(y[facei]);
109  G0[celli] += G[celli];
110  }
111  }
112 }
113 
114 
115 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
116 
119 (
120  const fvPatch& p,
122 )
123 :
125  yPlusLam_(yPlusLam(kappa_, E_))
126 {}
127 
128 
131 (
133  const fvPatch& p,
135  const fvPatchFieldMapper& mapper
136 )
137 :
138  epsilonWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
139  yPlusLam_(ptf.yPlusLam_)
140 {}
141 
142 
145 (
146  const fvPatch& p,
148  const dictionary& dict
149 )
150 :
152  yPlusLam_(yPlusLam(kappa_, E_))
153 {}
154 
155 
158 (
160 )
161 :
163  yPlusLam_(ewfpsf.yPlusLam_)
164 {}
165 
166 
169 (
172 )
173 :
175  yPlusLam_(ewfpsf.yPlusLam_)
176 {}
177 
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 namespace Foam
182 {
184  (
187  );
188 }
189 
190 
191 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:200
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
w
volScalarField w(IOobject("w", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("w", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0))
p
p
Definition: pEqn.H:62
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::epsilonLowReWallFunctionFvPatchScalarField::epsilonLowReWallFunctionFvPatchScalarField
epsilonLowReWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: epsilonLowReWallFunctionFvPatchScalarField.C:119
Foam::epsilonLowReWallFunctionFvPatchScalarField
This boundary condition provides a turbulence dissipation wall function condition for low- and high-R...
Definition: epsilonLowReWallFunctionFvPatchScalarField.H:99
Foam::turbulenceModel::nut
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::epsilonLowReWallFunctionFvPatchScalarField::yPlusLam_
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
Definition: epsilonLowReWallFunctionFvPatchScalarField.H:109
fvPatchFieldMapper.H
Foam::turbulenceModel::y
const nearWallDist & y() const
Return the near wall distances.
Definition: turbulenceModel.H:157
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::constant::electromagnetic::kappa
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Foam::epsilonLowReWallFunctionFvPatchScalarField::yPlusLam
scalar yPlusLam(const scalar kappa, const scalar E)
Calculate the Y+ at the edge of the laminar sublayer.
Definition: epsilonLowReWallFunctionFvPatchScalarField.C:35
Foam::constant::electromagnetic::epsilon0
const dimensionedScalar epsilon0
Electric constant: default SI units: [F/m].
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:131
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::constant::electromagnetic::G0
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::epsilonWallFunctionFvPatchScalarField
This boundary condition provides a turbulence dissipation wall function condition for high Reynolds n...
Definition: epsilonWallFunctionFvPatchScalarField.H:114
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
Foam::turbulenceModel::GName
word GName() const
Helper function to return the name of the turbulence G field.
Definition: turbulenceModel.H:136
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:60
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::fvPatch::index
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:179
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:253
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::turbulenceModel::k
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::makePatchTypeField
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::List< scalar >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::yPlus
This function object evaluates and outputs turbulence y+ for turbulence models. The field is stored o...
Definition: yPlus.H:110
patchi
label patchi
Definition: getPatchFieldScalar.H:1
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::fvPatch::faceCells
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::epsilonLowReWallFunctionFvPatchScalarField::calculate
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &epsilon)
Calculate the epsilon and G.
Definition: epsilonLowReWallFunctionFvPatchScalarField.C:52
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
epsilonLowReWallFunctionFvPatchScalarField.H
turbulenceModel.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
Foam::turbulenceModel::U
const volVectorField & U() const
Access function to velocity field.
Definition: turbulenceModel.H:142
y
scalar y
Definition: LISASMDCalcMethod1.H:14