nutWallFunctionFvPatchScalarField.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) 2011-2016, 2019 OpenFOAM Foundation
9  Copyright (C) 2019-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
32 #include "wallFvPatch.H"
33 #include "turbulenceModel.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
41 }
42 
43 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
44 
46 {
47  if (!isA<wallFvPatch>(patch()))
48  {
50  << "Invalid wall function specification" << nl
51  << " Patch type for patch " << patch().name()
52  << " must be wall" << nl
53  << " Current patch type is " << patch().type() << nl << endl
54  << abort(FatalError);
55  }
56 }
57 
58 
60 (
61  const turbulenceModel& turb
62 ) const
63 {
64  if (UName_.empty())
65  {
66  return turb.U();
67  }
68 
69  return db().lookupObject<volVectorField>(UName_);
70 }
71 
72 
74 (
75  Ostream& os
76 ) const
77 {
78  os.writeEntryIfDifferent<word>("U", word::null, UName_);
79  wallCoeffs_.writeEntries(os);
80 }
81 
82 
83 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84 
86 (
87  const fvPatch& p,
89 )
90 :
91  fixedValueFvPatchScalarField(p, iF),
92  UName_(),
93  wallCoeffs_()
94 {
95  checkType();
96 }
97 
98 
100 (
102  const fvPatch& p,
104  const fvPatchFieldMapper& mapper
105 )
106 :
107  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
108  UName_(ptf.UName_),
109  wallCoeffs_(ptf.wallCoeffs_)
110 {
111  checkType();
112 }
113 
114 
116 (
117  const fvPatch& p,
119  const dictionary& dict
120 )
121 :
122  fixedValueFvPatchScalarField(p, iF, dict),
123  UName_(dict.getOrDefault<word>("U", word::null)),
124  wallCoeffs_(dict)
125 {
126  checkType();
127 }
128 
129 
131 (
133 )
134 :
135  fixedValueFvPatchScalarField(wfpsf),
136  UName_(wfpsf.UName_),
137  wallCoeffs_(wfpsf.wallCoeffs_)
138 {
139  checkType();
140 }
141 
142 
144 (
147 )
148 :
149  fixedValueFvPatchScalarField(wfpsf, iF),
150  UName_(wfpsf.UName_),
151  wallCoeffs_(wfpsf.wallCoeffs_)
152 {
153  checkType();
154 }
155 
156 
157 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
158 
161 (
162  const turbulenceModel& turbModel,
163  const label patchi
164 )
165 {
166  return
167  refCast<const nutWallFunctionFvPatchScalarField>
168  (
169  turbModel.nut()().boundaryField()[patchi],
170  patchi
171  );
172 }
173 
174 
176 {
177  if (updated())
178  {
179  return;
180  }
181 
182  operator==(calcNut());
183 
184  fixedValueFvPatchScalarField::updateCoeffs();
185 }
186 
187 
189 (
190  Ostream& os
191 ) const
192 {
194  writeLocalEntries(os);
195 }
196 
197 
198 // ************************************************************************* //
dictionary dict
virtual const volVectorField & U(const turbulenceModel &turb) const
Helper to return the velocity field either from the turbulence model (default) or the mesh database...
void writeLocalEntries(Ostream &) const
Write local wall function variables.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
compressible::turbulenceModel & turb
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
The class nutWallFunction is an abstract base class that hosts calculation methods and common functi...
A class for handling words, derived from Foam::string.
Definition: word.H:63
A FieldMapper for finite-volume patch fields.
static const word null
An empty word.
Definition: word.H:84
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:336
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
virtual void checkType()
Check the type of the patch.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
nutWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
volScalarField & p
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
Namespace for OpenFOAM.