All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
constantHeatTransfer.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-2018 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 "constantHeatTransfer.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace fv
34 {
35  defineTypeNameAndDebug(constantHeatTransfer, 0);
37  (
38  option,
39  constantHeatTransfer,
40  dictionary
41  );
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const word& name,
51  const word& modelType,
52  const dictionary& dict,
53  const fvMesh& mesh
54 )
55 :
56  interRegionHeatTransferModel(name, modelType, dict, mesh),
57  htcConst_(),
58  AoV_()
59 {
60  if (active() && master_)
61  {
62  htcConst_.reset
63  (
64  new volScalarField
65  (
66  IOobject
67  (
68  "htcConst",
69  mesh_.time().timeName(),
70  mesh_,
73  ),
74  mesh_
75  )
76  );
77 
78  AoV_.reset
79  (
80  new volScalarField
81  (
82  IOobject
83  (
84  "AoV",
85  mesh_.time().timeName(),
86  mesh_,
89  ),
90  mesh_
91  )
92  );
93 
94  htc_ = htcConst_()*AoV_();
95  }
96 }
97 
98 
99 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
100 
102 {}
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
108 {}
109 
110 
112 {
114  {
115  return true;
116  }
117  else
118  {
119  return false;
120  }
121 }
122 
123 
124 // ************************************************************************* //
defineTypeNameAndDebug(option, 0)
Base class for inter region heat exchange. The derived classes must provide the heat transfer coeffic...
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
Macros for easy insertion into run-time selection tables.
virtual ~constantHeatTransfer()
Destructor.
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
virtual bool read(const dictionary &dict)
Read dictionary.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual void calculateHtc()
Calculate the heat transfer coefficient.
virtual bool read(const dictionary &dict)
Read dictionary.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.
constantHeatTransfer(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from dictionary.