epsilonWallFunctionFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 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 
28 #include "momentumTransportModel.H"
29 #include "fvMatrix.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
35 Foam::epsilonWallFunctionFvPatchScalarField::calculate
36 (
37  const momentumTransportModel& mtm
38 ) const
39 {
40  const label patchi = patch().index();
41 
42  const nutWallFunctionFvPatchScalarField& nutw =
44 
45  const scalarField& y = mtm.yb()[patchi];
46 
47  const tmp<scalarField> tnuw = mtm.nu(patchi);
48  const scalarField& nuw = tnuw();
49 
50  const tmp<volScalarField> tk = mtm.k();
51  const volScalarField& k = tk();
52 
53  const fvPatchVectorField& Uw = mtm.U().boundaryField()[patchi];
54 
55  const scalarField magGradUw(mag(Uw.snGrad()));
56 
57  const scalar Cmu25 = pow025(nutw.Cmu());
58  const scalar Cmu75 = pow(nutw.Cmu(), 0.75);
59 
60  // Initialise the returned field pair
61  Pair<tmp<scalarField>> GandEpsilon
62  (
63  tmp<scalarField>(new scalarField(patch().size())),
64  tmp<scalarField>(new scalarField(patch().size()))
65  );
66  scalarField& G = GandEpsilon.first().ref();
67  scalarField& epsilon = GandEpsilon.second().ref();
68 
69  // Calculate G and epsilon
70  forAll(nutw, facei)
71  {
72  const label celli = patch().faceCells()[facei];
73 
74  const scalar yPlus = Cmu25*y[facei]*sqrt(k[celli])/nuw[facei];
75 
76  if (yPlus > nutw.yPlusLam())
77  {
78  G[facei] =
79  (nutw[facei] + nuw[facei])
80  *magGradUw[facei]
81  *Cmu25*sqrt(k[celli])
82  /(nutw.kappa()*y[facei]);
83 
84  epsilon[facei] =
85  Cmu75*k[celli]*sqrt(k[celli])/(nutw.kappa()*y[facei]);
86  }
87  else
88  {
89  G[facei] = 0;
90 
91  epsilon[facei] = 2*k[celli]*nuw[facei]/sqr(y[facei]);
92  }
93  }
94 
95  return GandEpsilon;
96 }
97 
98 
99 void Foam::epsilonWallFunctionFvPatchScalarField::updateCoeffsMaster()
100 {
101  if (patch().index() != masterPatchIndex())
102  {
103  return;
104  }
105 
106  // Lookup the momentum transport model
107  const momentumTransportModel& mtm =
108  db().lookupType<momentumTransportModel>(internalField().group());
109 
110  // Get mutable references to the turbulence fields
111  VolInternalField<scalar>& G =
112  db().lookupObjectRef<VolInternalField<scalar>>(mtm.GName());
114  const_cast<volScalarField&>
115  (
116  static_cast<const volScalarField&>
117  (
118  internalField()
119  )
120  );
121 
122  const volScalarField::Boundary& epsilonBf = epsilon.boundaryField();
123 
124  // Make all processors build the near wall distances
125  mtm.yb();
126 
127  // Evaluate all the wall functions
128  PtrList<scalarField> Gpfs(epsilonBf.size());
129  PtrList<scalarField> epsilonPfs(epsilonBf.size());
130  forAll(epsilonBf, patchi)
131  {
132  if (isA<epsilonWallFunctionFvPatchScalarField>(epsilonBf[patchi]))
133  {
134  const epsilonWallFunctionFvPatchScalarField& epsilonPf =
135  refCast<const epsilonWallFunctionFvPatchScalarField>
136  (epsilonBf[patchi]);
137 
138  Pair<tmp<scalarField>> GandEpsilon = epsilonPf.calculate(mtm);
139 
140  Gpfs.set(patchi, GandEpsilon.first().ptr());
141  epsilonPfs.set(patchi, GandEpsilon.second().ptr());
142  }
143  }
144 
145  // Average the values into the wall-adjacent cells and store
146  wallCellGPtr_.reset(patchFieldsToWallCellField(Gpfs).ptr());
147  wallCellEpsilonPtr_.reset(patchFieldsToWallCellField(epsilonPfs).ptr());
148 
149  // Set the fractional proportion of the value in the wall cells
150  UIndirectList<scalar>(G, wallCells()) =
151  (1 - wallCellFraction())*scalarField(G, wallCells())
152  + wallCellFraction()*wallCellGPtr_();
153  UIndirectList<scalar>(epsilon, wallCells()) =
154  (1 - wallCellFraction())*scalarField(epsilon, wallCells())
155  + wallCellFraction()*wallCellEpsilonPtr_();
156 }
157 
158 
159 void Foam::epsilonWallFunctionFvPatchScalarField::manipulateMatrixMaster
160 (
161  fvMatrix<scalar>& matrix
162 )
163 {
164  if (patch().index() != masterPatchIndex())
165  {
166  return;
167  }
168 
169  matrix.setValues(wallCells(), wallCellEpsilonPtr_(), wallCellFraction());
170 }
171 
172 
173 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
174 
177 (
178  const fvPatch& p,
180  const dictionary& dict
181 )
182 :
184  wallCellGPtr_(nullptr),
185  wallCellEpsilonPtr_(nullptr)
186 {}
187 
188 
191 (
193  const fvPatch& p,
195  const fieldMapper& mapper
196 )
197 :
198  wallCellWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
199  wallCellGPtr_(nullptr),
200  wallCellEpsilonPtr_(nullptr)
201 {}
202 
203 
206 (
209 )
210 :
212  wallCellGPtr_(nullptr),
213  wallCellEpsilonPtr_(nullptr)
214 {}
215 
216 
217 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
218 
220 (
221  const fvPatchScalarField& ptf,
222  const fieldMapper& mapper
223 )
224 {
226  wallCellGPtr_.clear();
227  wallCellEpsilonPtr_.clear();
228 }
229 
230 
232 (
233  const fvPatchScalarField& ptf
234 )
235 {
237  wallCellGPtr_.clear();
238  wallCellEpsilonPtr_.clear();
239 }
240 
241 
243 {
244  if (updated())
245  {
246  return;
247  }
248 
249  initMaster();
250 
251  updateCoeffsMaster();
252 
253  operator==(patchInternalField());
254 
256 }
257 
258 
260 (
261  fvMatrix<scalar>& matrix
262 )
263 {
264  if (manipulatedMatrix())
265  {
266  return;
267  }
268 
269  if (masterPatchIndex() == -1)
270  {
272  << "updateCoeffs must be called before manipulateMatrix"
273  << exit(FatalError);
274  }
275 
276  manipulateMatrixMaster(matrix);
277 
279 }
280 
281 
282 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
283 
284 namespace Foam
285 {
287  (
290  );
291 }
292 
293 
294 // ************************************************************************* //
scalar y
label k
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
label size() const
Return the number of elements in the UList.
Definition: ListI.H:171
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:65
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
This boundary condition provides a turbulence dissipation wall constraint for low- and high-Reynolds ...
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
epsilonWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
Abstract base class for field mapping.
Definition: fieldMapper.H:48
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:202
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:222
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:356
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:156
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:99
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
Base class for wall functions that modify cell values.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const scalar yPlus
const scalar epsilon
label patchi
compressibleMomentumTransportModel momentumTransportModel
const dimensionedScalar G
Newtonian constant of gravitation.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
dimensioned< scalar > mag(const dimensioned< Type > &)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
fvPatchField< vector > fvPatchVectorField
dimensionedScalar pow025(const dimensionedScalar &ds)
dictionary dict
volScalarField & p