alphatWallFunctionFvPatchScalarField.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 | Copyright (C) 2011-2016 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 "fvPatchFieldMapper.H"
29 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace compressible
37 {
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42 (
43  const fvPatch& p,
45 )
46 :
47  fixedValueFvPatchScalarField(p, iF),
48  Prt_(0.85)
49 {}
50 
51 
53 (
55  const fvPatch& p,
57  const fvPatchFieldMapper& mapper
58 )
59 :
60  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
61  Prt_(ptf.Prt_)
62 {}
63 
64 
66 (
67  const fvPatch& p,
69  const dictionary& dict
70 )
71 :
72  fixedValueFvPatchScalarField(p, iF, dict),
73  Prt_(dict.lookupOrDefault<scalar>("Prt", 0.85))
74 {}
75 
76 
78 (
80 )
81 :
82  fixedValueFvPatchScalarField(awfpsf),
83  Prt_(awfpsf.Prt_)
84 {}
85 
86 
88 (
91 )
92 :
93  fixedValueFvPatchScalarField(awfpsf, iF),
94  Prt_(awfpsf.Prt_)
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
101 {
102  if (updated())
103  {
104  return;
105  }
106 
107  const label patchi = patch().index();
108 
109  // Retrieve turbulence properties from model
110  const compressibleTurbulenceModel& turbModel =
111  db().lookupObject<compressibleTurbulenceModel>
112  (
114  (
116  internalField().group()
117  )
118  );
119 
120  const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
121  const tmp<scalarField> tnutw = turbModel.nut(patchi);
122 
123  operator==(rhow*tnutw/Prt_);
124 
125  fixedValueFvPatchScalarField::updateCoeffs();
126 }
127 
128 
130 {
132  os.writeKeyword("Prt") << Prt_ << token::END_STATEMENT << nl;
133  writeEntry("value", os);
134 }
135 
136 
137 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
138 
140 (
143 );
144 
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146 
147 } // End namespace compressible
148 } // End namespace Foam
149 
150 // ************************************************************************* //
const char *const group
Group name for atomic constants.
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:362
Macros for easy insertion into run-time selection tables.
static const word propertiesName
Default name of the turbulence properties dictionary.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
const volScalarField & rho() const
Return the density field.
static const char nl
Definition: Ostream.H:262
Abstract base class for turbulence models (RAS, LES and laminar).
alphatWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
This boundary condition provides a turbulent thermal diffusivity conditon when using wall functions...
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
A class for managing temporary objects.
Definition: PtrList.H:53
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
Namespace for OpenFOAM.
bool compressible
Definition: pEqn.H:30