nutUSpaldingWallFunctionFvPatchScalarFieldDF.H
Go to the documentation of this file.
1 
2 /*---------------------------------------------------------------------------*\
3 
4  DAFoam : Discrete Adjoint with OpenFOAM
5  Version : v3
6 
7  This file is modified from OpenFOAM's source code
8  src/TurbulenceModels/turbulenceModels/derivedFvPatchFields/wallFunctions/nutWallFunctions
9  /nutUSpaldingWallFunction/nutUSpaldingWallFunctionFvPatchScalarField.H
10 
11  OpenFOAM: The Open Source CFD Toolbox
12 
13  Copyright (C): 2011-2016 OpenFOAM Foundation
14 
15  OpenFOAM License:
16 
17  OpenFOAM is free software: you can redistribute it and/or modify it
18  under the terms of the GNU General Public License as published by
19  the Free Software Foundation, either version 3 of the License, or
20  (at your option) any later version.
21 
22  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
23  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
24  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
25  for more details.
26 
27  You should have received a copy of the GNU General Public License
28  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
29 
30  Description:
31 
32  NOTE: It's the nutUSpaldingWallFunction implementation in OpenFOAM-v2006,
33  which supports two additional inputs: maxIter and tolerance. In DAFoam,
34  we need to make sure this boundary condition provides as accurate as
35  possible boundary values, so the default values are maxIter = 1000
36  and tolerance = 1e-14. If this causes very slow simulations,
37  decrease maxIter or increate tolerance.
38 
39  This boundary condition provides a wall constraint on the turbulent
40  viscosity, i.e. \c nut, based on velocity, i.e. \c U. Using Spalding's
41  law gives a continuous \c nut profile to the wall.
42 
43  \f[
44  y^+ = u^+ + \frac{1}{E} \left[exp(\kappa u^+) - 1 - \kappa u^+\,
45  - 0.5 (\kappa u^+)^2 - \frac{1}{6} (\kappa u^+)^3\right]
46  \f]
47 
48  where
49  \vartable
50  y^+ | Non-dimensional position
51  u^+ | Non-dimensional velocity
52  \kappa | von Kármán constant
53  \endvartable
54 
55 
56  Usage
57  Example of the boundary condition specification:
58  \verbatim
59  <patchName>
60  {
61  // Mandatory entries (unmodifiable)
62  type nutUSpaldingWallFunction;
63 
64  // Optional entries (unmodifiable)
65  maxIter 100;
66  tolerance 1e-8;
67 
68  // Optional (inherited) entries
69  ...
70  }
71  \endverbatim
72 
73  where the entries mean:
74  \table
75  Property | Description | Type | Req'd | Dflt
76  type | Type name: nutUBlendedWallFunction | word | yes | -
77  maxIter | Number of Newton-Raphson iterations | label | no | 1000
78  tolerance | Convergence tolerance | scalar | no | 1e-14
79  \endtable
80 
81  The inherited entries are elaborated in:
82  - \link nutWallFunctionFvPatchScalarField.H \endlink
83 
84  Note
85  - Suffers from non-exact restart since \c correctNut() (called through
86  \c turbulence->validate) returns a slightly different value every time
87  it is called. This is since the seed for the Newton-Raphson iteration
88  uses the current value of \c *this (\c =nut ).
89  - This can be avoided by overriding the tolerance. This also switches on
90  a pre-detection whether the current nut already satisfies the turbulence
91  conditions and if so does not change it at all. This means that the nut
92  only changes if it 'has really changed'. This probably should be used with
93  a tight tolerance, to make sure to kick every iteration, e.g.
94  maxIter 100;
95  tolerance 1e-7;
96 
97 \*---------------------------------------------------------------------------*/
98 
99 #ifndef nutUSpaldingWallFunctionFvPatchScalarFieldDF_H
100 #define nutUSpaldingWallFunctionFvPatchScalarFieldDF_H
101 
102 #include "nutWallFunctionFvPatchScalarField.H"
103 
104 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
105 
106 namespace Foam
107 {
108 
109 /*---------------------------------------------------------------------------*\
110  Class nutUSpaldingWallFunctionFvPatchScalarFieldDF Declaration
111 \*---------------------------------------------------------------------------*/
112 
114  : public nutWallFunctionFvPatchScalarField
115 {
116 protected:
117  // Protected Member Functions
118 
119  //- Calculate the turbulence viscosity
120  virtual tmp<scalarField> calcNut() const;
121 
122  //- Calculate the friction velocity
123  virtual tmp<scalarField> calcUTau(const scalarField& magGradU) const;
124 
125  //- Max iterations in calcNut
126  const label maxIter_;
127 
128  //- Convergence tolerance
129  const scalar tolerance_;
130 
131  //- Calculate the friction velocity and number of iterations for
132  //- convergence
133  virtual tmp<scalarField> calcUTau(
134  const scalarField& magGradU,
135  const label maxIter,
136  scalarField& err) const;
137 
138  //- Write local wall function variables
139  virtual void writeLocalEntries(Ostream&) const;
140 
141 public:
142  //- Runtime type information
143  TypeName("nutUSpaldingWallFunctionDF");
144 
145  // Constructors
146 
147  //- Construct from patch and internal field
149  const fvPatch&,
150  const DimensionedField<scalar, volMesh>&);
151 
152  //- Construct from patch, internal field and dictionary
154  const fvPatch&,
155  const DimensionedField<scalar, volMesh>&,
156  const dictionary&);
157 
158  //- Construct by mapping given
159  // nutUSpaldingWallFunctionFvPatchScalarFieldDF
160  // onto a new patch
163  const fvPatch&,
164  const DimensionedField<scalar, volMesh>&,
165  const fvPatchFieldMapper&);
166 
167  //- Construct as copy
170 
171  //- Construct and return a clone
172  virtual tmp<fvPatchScalarField> clone() const
173  {
174  return tmp<fvPatchScalarField>(
176  }
177 
178  //- Construct as copy setting internal field reference
181  const DimensionedField<scalar, volMesh>&);
182 
183  //- Construct and return a clone setting internal field reference
184  virtual tmp<fvPatchScalarField> clone(
185  const DimensionedField<scalar, volMesh>& iF) const
186  {
187  return tmp<fvPatchScalarField>(
189  }
190 
191  //- Destructor
193  {}
194 
195  // Member functions
196 
197  // Evaluation functions
198 
199  //- Calculate and return the yPlus at the boundary
200  virtual tmp<scalarField> yPlus() const;
201 
202  // I-O
203 
204  //- Write
205  virtual void write(Ostream& os) const;
206 };
207 
208 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209 
210 } // End namespace Foam
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 #endif
215 
216 // ************************************************************************* //
Foam::nutUSpaldingWallFunctionFvPatchScalarFieldDF::calcNut
virtual tmp< scalarField > calcNut() const
Definition: nutUSpaldingWallFunctionFvPatchScalarFieldDF.C:44
Foam::nutUSpaldingWallFunctionFvPatchScalarFieldDF::TypeName
TypeName("nutUSpaldingWallFunctionDF")
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
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
Definition: multiFreqScalarFvPatchField.C:144
Foam::nutUSpaldingWallFunctionFvPatchScalarFieldDF::clone
virtual tmp< fvPatchScalarField > clone(const DimensionedField< scalar, volMesh > &iF) const
Definition: nutUSpaldingWallFunctionFvPatchScalarFieldDF.H:184
Foam::nutUSpaldingWallFunctionFvPatchScalarFieldDF::~nutUSpaldingWallFunctionFvPatchScalarFieldDF
virtual ~nutUSpaldingWallFunctionFvPatchScalarFieldDF()
Definition: nutUSpaldingWallFunctionFvPatchScalarFieldDF.H:192
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::clone
virtual tmp< fvPatchScalarField > clone() const
Definition: nutUSpaldingWallFunctionFvPatchScalarFieldDF.H:172
Foam::nutUSpaldingWallFunctionFvPatchScalarFieldDF
Definition: nutUSpaldingWallFunctionFvPatchScalarFieldDF.H:113