kappaEff.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) 2020-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 
26 #include "kappaEff.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace wallHeatTransferCoeffModels
35 {
38  (
40  kappaEff,
41  word
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const word& name,
52  const fvMesh& mesh,
53  const dictionary& dict
54 )
55 :
56  wallHeatTransferCoeffModel(name, mesh, dict),
57  mesh_(mesh),
58  Pr_
59  (
60  "Pr",
61  dimless,
62  dict.lookupBackwardsCompatible({"Pr", "Prl"})
63  ),
64  Prt_("Prt", dimless, dict),
65  Lchar_("Lchar", dimLength, Zero),
66  isCharLength_(false)
67 {
68  read(dict);
69 }
70 
71 
72 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
73 
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
81 {
82  Pr_.read(dict);
83  Prt_.read(dict);
84  isCharLength_ = dict.found("Lchar");
85 
86  if (isCharLength_)
87  {
88  Lchar_.read(dict);
89  }
90 
91  return true;
92 }
93 
94 
97 (
98  const momentumTransportModel& mmtm,
99  const labelHashSet& patches
100 ) const
101 {
102  tmp<volScalarField> thtcByRhoCp
103  (
105  (
106  type(),
107  mesh_,
108  isCharLength_
111  )
112  );
113 
114  volScalarField::Boundary& thtcByRhoCpBf =
115  thtcByRhoCp.ref().boundaryFieldRef();
116 
117  forAllConstIter(labelHashSet, patches, iter)
118  {
119  label patchi = iter.key();
120 
121  if (!thtcByRhoCpBf[patchi].coupled())
122  {
123  thtcByRhoCpBf[patchi] =
124  mmtm.nu(patchi)/Pr_.value() + mmtm.nut(patchi)/Prt_.value();
125 
126  if (isCharLength_)
127  {
128  thtcByRhoCpBf[patchi] /= Lchar_.value();
129  }
130  }
131  }
132  return thtcByRhoCp;
133 
134 }
135 
136 // ************************************************************************* //
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
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
const dimensionSet dimArea
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Abstract base class for run time selection of heat transfer coefficient models.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
const dimensionSet dimless
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
const dimensionSet dimTime
bool read(const char *, int32_t &)
Definition: int32IO.C:85
A class for handling words, derived from string.
Definition: word.H:59
virtual tmp< volScalarField > htcByRhoCp(const momentumTransportModel &mmtm, const labelHashSet &patches) const
Calculate the wall heat transfer coefficient.
Definition: kappaEff.C:97
static const zero Zero
Definition: zero.H:97
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
addToRunTimeSelectionTable(wallHeatTransferCoeffModel, kappaEff, word)
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
Abstract base class for turbulence models (RAS, LES and laminar).
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual bool read(const dictionary &)
Read the kappaEff data.
Definition: kappaEff.C:80
A class for managing temporary objects.
Definition: PtrList.H:53
Calculates the estimated flow heat transfer coefficient at wall patches as the volScalarField field &#39;...
Definition: kappaEff.H:194
ITstream & lookupBackwardsCompatible(const wordList &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream, trying a list of keywords.
Definition: dictionary.C:855
Namespace for OpenFOAM.
kappaEff(const word &name, const fvMesh &mesh, const dictionary &)
Construct from name, mesh and dict.
Definition: kappaEff.C:50