HeatTransferModel.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-2020 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 "HeatTransferModel.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class CloudType>
32 :
34  BirdCorrection_(false)
35 {}
36 
37 
38 template<class CloudType>
40 (
41  const dictionary& dict,
42  CloudType& owner,
43  const word& type
44 )
45 :
46  CloudSubModelBase<CloudType>(owner, dict, typeName, type),
47  BirdCorrection_(this->coeffDict().lookup("BirdCorrection"))
48 {}
49 
50 
51 template<class CloudType>
53 (
55 )
56 :
58  BirdCorrection_(htm.BirdCorrection_)
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 
64 template<class CloudType>
66 {}
67 
68 
69 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70 
71 template<class CloudType>
73 {
74  return BirdCorrection_;
75 }
76 
77 
78 template<class CloudType>
80 (
81  const scalar dp,
82  const scalar Re,
83  const scalar Pr,
84  const scalar kappa,
85  const scalar NCpW
86 ) const
87 {
88  const scalar Nu = this->Nu(Re, Pr);
89 
90  scalar htc = Nu*kappa/dp;
91 
92  if (BirdCorrection_ && (mag(htc) > rootVSmall) && (mag(NCpW) > rootVSmall))
93  {
94  const scalar phit = min(NCpW/htc, 50);
95  if (phit > 0.001)
96  {
97  htc *= phit/(exp(phit) - 1.0);
98  }
99  }
100 
101  return htc;
102 }
103 
104 
105 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
106 
107 #include "HeatTransferModelNew.C"
108 
109 // ************************************************************************* //
Base class for cloud sub-models.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:79
Templated heat transfer model class.
virtual ~HeatTransferModel()
Destructor.
virtual scalar htc(const scalar dp, const scalar Re, const scalar Pr, const scalar kappa, const scalar NCpW) const
Return heat transfer coefficient.
HeatTransferModel(CloudType &owner)
Construct null from owner.
const Switch & BirdCorrection() const
Return the Bird htc correction flag.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
A class for handling words, derived from string.
Definition: word.H:62
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
dimensionedScalar exp(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< scalar > mag(const dimensioned< Type > &)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
dictionary dict