nutUSpaldingWallFunctionFvPatchScalarFieldDF.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6  This file is modified from OpenFOAM's source code
7  src/TurbulenceModels/turbulenceModels/derivedFvPatchFields/wallFunctions/nutWallFunctions
8  /nutUSpaldingWallFunction/nutUSpaldingWallFunctionFvPatchScalarField.C
9 
10  OpenFOAM: The Open Source CFD Toolbox
11 
12  Copyright (C): 2011-2016 OpenFOAM Foundation
13 
14  OpenFOAM License:
15 
16  OpenFOAM is free software: you can redistribute it and/or modify it
17  under the terms of the GNU General Public License as published by
18  the Free Software Foundation, either version 3 of the License, or
19  (at your option) any later version.
20 
21  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24  for more details.
25 
26  You should have received a copy of the GNU General Public License
27  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28 
29 \*---------------------------------------------------------------------------*/
30 
32 #include "turbulenceModel.H"
33 #include "fvPatchFieldMapper.H"
34 #include "volFields.H"
35 #include "addToRunTimeSelectionTable.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 
42 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
43 
45 {
46  const label patchi = patch().index();
47 
48  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>(
49  IOobject::groupName(
50  turbulenceModel::propertiesName,
51  internalField().group()));
52  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
53  const scalarField magGradU(mag(Uw.snGrad()));
54  const tmp<scalarField> tnuw = turbModel.nu(patchi);
55  const scalarField& nuw = tnuw();
56 
57  // Calculate new viscosity
58  tmp<scalarField> tnutw(
59  max(
60  scalar(0),
61  sqr(calcUTau(magGradU)) / (magGradU + ROOTVSMALL) - nuw));
62 
63  if (tolerance_ != 1.e-14)
64  {
65  // User-specified tolerance. Find out if current nut already satisfies
66  // eqns.
67 
68  // Run ut for one iteration
69  scalarField err;
70  tmp<scalarField> UTau(calcUTau(magGradU, 1, err));
71 
72  // Preserve nutw if the initial error (err) already lower than the
73  // tolerance.
74 
75  scalarField& nutw = tnutw.ref();
76  forAll(err, facei)
77  {
78  if (err[facei] < tolerance_)
79  {
80  nutw[facei] = this->operator[](facei);
81  }
82  }
83  }
84  return tnutw;
85 }
86 
88  const scalarField& magGradU) const
89 {
90  scalarField err;
91  return calcUTau(magGradU, maxIter_, err);
92 }
93 
95  const scalarField& magGradU,
96  const label maxIter,
97  scalarField& err) const
98 {
99  const label patchi = patch().index();
100 
101  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>(
102  IOobject::groupName(
103  turbulenceModel::propertiesName,
104  internalField().group()));
105  const scalarField& y = turbModel.y()[patchi];
106 
107  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
108  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
109 
110  const tmp<scalarField> tnuw = turbModel.nu(patchi);
111  const scalarField& nuw = tnuw();
112 
113  const scalarField& nutw = *this;
114 
115  tmp<scalarField> tuTau(new scalarField(patch().size(), pTraits<scalar>::zero));
116  scalarField& uTau = tuTau.ref();
117 
118  err.setSize(uTau.size());
119  err = 0.0;
120 
121  forAll(uTau, facei)
122  {
123  scalar ut = sqrt((nutw[facei] + nuw[facei]) * magGradU[facei]);
124  // Note: for exact restart seed with laminar viscosity only:
125  //scalar ut = sqrt(nuw[facei]*magGradU[facei]);
126 
127  if (ROOTVSMALL < ut)
128  {
129  int iter = 0;
130 
131  do
132  {
133  scalar kUu = min(kappa_ * magUp[facei] / ut, 50);
134  scalar fkUu = exp(kUu) - 1 - kUu * (1 + 0.5 * kUu);
135 
136  scalar f =
137  -ut * y[facei] / nuw[facei]
138  + magUp[facei] / ut
139  + 1 / E_ * (fkUu - 1.0 / 6.0 * kUu * sqr(kUu));
140 
141  scalar df =
142  y[facei] / nuw[facei]
143  + magUp[facei] / sqr(ut)
144  + 1 / E_ * kUu * fkUu / ut;
145 
146  scalar uTauNew = ut + f / df;
147  err[facei] = mag((ut - uTauNew) / ut);
148  ut = uTauNew;
149 
150  //iterations_++;
151 
152  } while (
153  ut > ROOTVSMALL
154  && err[facei] > tolerance_
155  && ++iter < maxIter);
156 
157  uTau[facei] = max(0.0, ut);
158  }
159  }
160 
161  return tuTau;
162 }
163 
165  Ostream& os) const
166 {
167  nutWallFunctionFvPatchScalarField::writeLocalEntries(os);
168 
169  os.writeEntryIfDifferent<label>("maxIter", 1000, maxIter_);
170  os.writeEntryIfDifferent<scalar>("tolerance", 1.e-14, tolerance_);
171 }
172 
173 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
174 
177  const fvPatch& p,
178  const DimensionedField<scalar, volMesh>& iF)
179  : nutWallFunctionFvPatchScalarField(p, iF),
180  maxIter_(1000),
181  tolerance_(1.e-14)
182 {
183 }
184 
188  const fvPatch& p,
189  const DimensionedField<scalar, volMesh>& iF,
190  const fvPatchFieldMapper& mapper)
191  : nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
192  maxIter_(ptf.maxIter_),
193  tolerance_(ptf.tolerance_)
194 {
195 }
196 
199  const fvPatch& p,
200  const DimensionedField<scalar, volMesh>& iF,
201  const dictionary& dict)
202  : nutWallFunctionFvPatchScalarField(p, iF, dict),
203  maxIter_(dict.lookupOrDefault<label>("maxIter", 1000)),
204  tolerance_(dict.lookupOrDefault<scalar>("tolerance", 1.e-14))
205 {
206 }
207 
211  : nutWallFunctionFvPatchScalarField(wfpsf),
212  maxIter_(wfpsf.maxIter_),
213  tolerance_(wfpsf.tolerance_)
214 {
215 }
216 
220  const DimensionedField<scalar, volMesh>& iF)
221  : nutWallFunctionFvPatchScalarField(wfpsf, iF),
222  maxIter_(wfpsf.maxIter_),
223  tolerance_(wfpsf.tolerance_)
224 {
225 }
226 
227 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
228 
230 {
231  const label patchi = patch().index();
232 
233  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>(
234  IOobject::groupName(
235  turbulenceModel::propertiesName,
236  internalField().group()));
237  const scalarField& y = turbModel.y()[patchi];
238  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
239  const tmp<scalarField> tnuw = turbModel.nu(patchi);
240  const scalarField& nuw = tnuw();
241 
242  return y * calcUTau(mag(Uw.snGrad())) / nuw;
243 }
244 
246 {
247  fvPatchField<scalar>::write(os);
248  writeLocalEntries(os);
249  writeEntry("value", os);
250 }
251 
252 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
253 
255  fvPatchScalarField,
257 
258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259 
260 } // End namespace Foam
261 
262 // ************************************************************************* //
forAll
forAll(pseudoP.boundaryField(), patchI)
Definition: solvePseudoPEqn.H:10
Foam::nutUSpaldingWallFunctionFvPatchScalarFieldDF::calcNut
virtual tmp< scalarField > calcNut() const
Definition: nutUSpaldingWallFunctionFvPatchScalarFieldDF.C:44
Foam::nutUSpaldingWallFunctionFvPatchScalarFieldDF::calcUTau
virtual tmp< scalarField > calcUTau(const scalarField &magGradU) const
Definition: nutUSpaldingWallFunctionFvPatchScalarFieldDF.C:87
Foam::nutUSpaldingWallFunctionFvPatchScalarFieldDF::write
virtual void write(Ostream &os) const
Definition: nutUSpaldingWallFunctionFvPatchScalarFieldDF.C:245
p
volScalarField & p
Definition: createRefsPimple.H:6
Foam::nutUSpaldingWallFunctionFvPatchScalarFieldDF::yPlus
virtual tmp< scalarField > yPlus() const
Definition: nutUSpaldingWallFunctionFvPatchScalarFieldDF.C:229
Foam::nutUSpaldingWallFunctionFvPatchScalarFieldDF::writeLocalEntries
virtual void writeLocalEntries(Ostream &) const
Definition: nutUSpaldingWallFunctionFvPatchScalarFieldDF.C:164
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, multiFreqScalarFvPatchField)
Foam
Definition: multiFreqScalarFvPatchField.C:144
nutUSpaldingWallFunctionFvPatchScalarFieldDF.H
Foam::nutUSpaldingWallFunctionFvPatchScalarFieldDF::maxIter_
const label maxIter_
Definition: nutUSpaldingWallFunctionFvPatchScalarFieldDF.H:126
Foam::nutUSpaldingWallFunctionFvPatchScalarFieldDF::tolerance_
const scalar tolerance_
Definition: nutUSpaldingWallFunctionFvPatchScalarFieldDF.H:129
Foam::nutUSpaldingWallFunctionFvPatchScalarFieldDF::nutUSpaldingWallFunctionFvPatchScalarFieldDF
nutUSpaldingWallFunctionFvPatchScalarFieldDF(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Definition: nutUSpaldingWallFunctionFvPatchScalarFieldDF.C:176
Foam::nutUSpaldingWallFunctionFvPatchScalarFieldDF
Definition: nutUSpaldingWallFunctionFvPatchScalarFieldDF.H:113