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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2020-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 
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace compressible
40 {
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const fvPatch& p,
48 )
49 :
50  fixedValueFvPatchScalarField(p, iF),
51  Prt_(0.85)
52 {}
53 
54 
56 (
58  const fvPatch& p,
60  const fvPatchFieldMapper& mapper
61 )
62 :
63  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
64  Prt_(ptf.Prt_)
65 {}
66 
67 
69 (
70  const fvPatch& p,
72  const dictionary& dict
73 )
74 :
75  fixedValueFvPatchScalarField(p, iF, dict),
76  Prt_(dict.getOrDefault<scalar>("Prt", 0.85))
77 {}
78 
79 
81 (
83 )
84 :
85  fixedValueFvPatchScalarField(awfpsf),
86  Prt_(awfpsf.Prt_)
87 {}
88 
89 
91 (
94 )
95 :
96  fixedValueFvPatchScalarField(awfpsf, iF),
97  Prt_(awfpsf.Prt_)
98 {}
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
104 {
105  if (updated())
106  {
107  return;
108  }
109 
110  const label patchi = patch().index();
111 
112  // Retrieve turbulence properties from model
113  const auto& turbModel = db().lookupObject<compressibleTurbulenceModel>
114  (
116  (
118  internalField().group()
119  )
120  );
121 
122  const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
123  const tmp<scalarField> tnutw = turbModel.nut(patchi);
125  operator==(rhow*tnutw/Prt_);
126 
127  fixedValueFvPatchScalarField::updateCoeffs();
128 }
129 
130 
132 {
134  os.writeEntryIfDifferent<scalar>("Prt", 0.85, Prt_);
136 }
137 
138 
139 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
140 
142 (
144  alphatWallFunctionFvPatchScalarField
145 );
146 
147 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
148 
149 } // End namespace compressible
150 } // End namespace Foam
151 
152 // ************************************************************************* //
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
Macros for easy insertion into run-time selection tables.
constexpr const char *const group
Group name for atomic constants.
makePatchTypeField(fvPatchScalarField, alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField)
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
static const word propertiesName
Default name of the turbulence properties dictionary.
A FieldMapper for finite-volume patch fields.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
bool compressible
Definition: pEqn.H:2
Abstract base class for turbulence models (RAS, LES and laminar).
alphatWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
OBJstream os(runTime.globalPath()/outputName)
This boundary condition provides a turbulent thermal diffusivity condition when using wall functions...
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.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Namespace for OpenFOAM.