alphatWallFunctionFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 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 
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace compressible
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
40 (
41  const fvPatch& p,
43 )
44 :
45  fixedValueFvPatchScalarField(p, iF),
46  Prt_(0.85)
47 {}
48 
49 
51 (
52  const fvPatch& p,
54  const dictionary& dict
55 )
56 :
57  fixedValueFvPatchScalarField(p, iF, dict),
58  Prt_(dict.lookupOrDefault<scalar>("Prt", 0.85))
59 {}
60 
61 
63 (
65  const fvPatch& p,
67  const fvPatchFieldMapper& mapper
68 )
69 :
70  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
71  Prt_(ptf.Prt_)
72 {}
73 
74 
76 (
79 )
80 :
81  fixedValueFvPatchScalarField(awfpsf, iF),
82  Prt_(awfpsf.Prt_)
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
89 {
90  if (updated())
91  {
92  return;
93  }
94 
95  const label patchi = patch().index();
96 
97  const thermophysicalTransportModel& ttm =
98  db().lookupObject<thermophysicalTransportModel>
99  (
101  (
102  thermophysicalTransportModel::typeName,
103  internalField().group()
104  )
105  );
106 
107  const compressibleMomentumTransportModel& turbModel =
108  ttm.momentumTransport();
109 
110  const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
111 
112  const tmp<scalarField> tnutw = turbModel.nut(patchi);
113 
114  operator==(rhow*tnutw/Prt_);
115 
116  fixedValueFvPatchScalarField::updateCoeffs();
117 }
118 
119 
121 {
123  writeEntry(os, "Prt", Prt_);
124  writeEntry(os, "value", *this);
125 }
126 
127 
128 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
129 
131 (
134 );
135 
136 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137 
138 } // End namespace compressible
139 } // End namespace Foam
140 
141 // ************************************************************************* //
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:156
const compressibleMomentumTransportModel & momentumTransport() const
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
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:260
Macros for easy insertion into run-time selection tables.
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:54
alphatWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
This boundary condition provides a turbulent thermal diffusivity condition when using wall functions...
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Abstract base class for thermophysical transport models (RAS, LES and laminar).
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.
A class for managing temporary objects.
Definition: PtrList.H:53
Abstract base class for turbulence models (RAS, LES and laminar).
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
Namespace for OpenFOAM.