Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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"
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38namespace 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),
106:
107 fixedValueFvPatchScalarField(ptf, p, iF, mapper),
121:
122 fixedValueFvPatchScalarField(p, iF, dict),
134:
135 fixedValueFvPatchScalarField(wfpsf),
148:
149 fixedValueFvPatchScalarField(wfpsf, iF),
150 UName_(wfpsf.UName_),
151 wallCoeffs_(wfpsf.wallCoeffs_)
152{
154}
155
156
157// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
158
161(
162 const turbulenceModel& turbModel,
163 const label patchi
164)
165{
166 return
169 turbModel.nut()().boundaryField()[patchi],
170 patchi
171 );
172}
173
174
176{
177 if (updated())
178 {
179 return;
180 }
183
184 fixedValueFvPatchScalarField::updateCoeffs();
185}
186
187
189(
190 Ostream& os
191) const
192{
194 writeLocalEntries(os);
195}
196
197
198// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
compressible::turbulenceModel & turb
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
The class nutWallFunction is an abstract base class that hosts calculation methods and common functi...
virtual const volVectorField & U(const turbulenceModel &turb) const
Helper to return the velocity field either from the turbulence model (default) or the mesh database.
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
void writeLocalEntries(Ostream &) const
Write local wall function variables.
virtual tmp< scalarField > calcNut() const =0
Calculate the turbulent viscosity.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
virtual void checkType()
Check the type of the patch.
nutWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
GeometricField< vector, fvPatchField, volMesh > volVectorField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
errorManip< error > abort(error &err)
Definition errorManip.H:139
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict