heatTransferModel.H
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) 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 Class
25  Foam::fv::heatTransferModel
26 
27 Description
28  Base class for heat transfer coefficient modelling used in heat transfer
29  fvModels. Area per unit volume [1/m] (AoV) must be provided as a value in
30  the coefficients dictionary or as a field in constant.
31 
32 SourceFiles
33  heatTransferModel.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef heatTransferModel_H
38 #define heatTransferModel_H
39 
40 #include "volFields.H"
41 #include "interRegionModel.H"
42 #include "runTimeSelectionTables.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 namespace fv
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class heatTransferModel Declaration
53 \*---------------------------------------------------------------------------*/
54 
56 {
57  // Private Member Data
58 
59  //- Reference to the mesh
60  const fvMesh& mesh_;
61 
62  //- Dictionary containing source coefficients
63  dictionary coeffs_;
64 
65  //- Area per unit volume [1/m]
66  dimensionedScalar AoV_;
67 
68  //- Area per unit volume [1/m]
70 
71 
72  // Private Member Functions
73 
74  //- Non-virtual read
75  void readCoeffs();
76 
77 
78 public:
79 
80  //- Runtime type information
81  TypeName("heatTransferModel");
82 
83 
84  // Declare runtime construction
85 
87  (
88  autoPtr,
90  mesh,
91  (
92  const dictionary& dict,
93  const fvMesh& mesh
94  ),
95  (dict, mesh)
96  );
97 
99  (
100  autoPtr,
102  model,
103  (
104  const dictionary& dict,
105  const interRegionModel& model
106  ),
107  (dict, model)
108  );
109 
110 
111  // Constructors
112 
113  //- Construct from dictionary and mesh
115  (
116  const word& modelType,
117  const dictionary& dict,
118  const fvMesh& mesh
119  );
120 
121  //- Construct from dictionary and model
123  (
124  const word& modelType,
125  const dictionary& dict,
126  const interRegionModel& model
127  );
128 
129  //- Disallow default bitwise copy construction
130  heatTransferModel(const heatTransferModel&) = delete;
131 
132 
133  // Selectors
134 
135  //- Select from dictionary and mesh
137  (
138  const dictionary& dict,
139  const fvMesh& mesh
140  );
141 
142  //- Select from dictionary and model
144  (
145  const dictionary& dict,
146  const interRegionModel& model
147  );
148 
149 
150  //- Destructor
151  virtual ~heatTransferModel();
152 
153 
154  // Member Functions
155 
156  //- Return reference to the mesh
157  inline const fvMesh& mesh() const
158  {
159  return mesh_;
160  }
161 
162  //- Return coeffs dictionary
163  inline const dictionary& coeffs() const
164  {
165  return coeffs_;
166  }
167 
168  //- Get the area per unit volume
169  tmp<volScalarField> AoV() const;
170 
171  //- Get the heat transfer coefficient
172  virtual tmp<volScalarField> htc() const = 0;
173 
174  //- Correct the heat transfer coefficient
175  virtual void correct() = 0;
176 
177  //- Read dictionary
178  virtual bool read(const dictionary& dict);
179 };
180 
181 
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
183 
184 } // End namespace fv
185 } // End namespace Foam
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 #endif
190 
191 // ************************************************************************* //
tmp< volScalarField > AoV() const
Get the area per unit volume.
dictionary dict
virtual void correct()=0
Correct the heat transfer coefficient.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Base class for heat transfer coefficient modelling used in heat transfer fvModels. Area per unit volume [1/m] (AoV) must be provided as a value in the coefficients dictionary or as a field in constant.
virtual tmp< volScalarField > htc() const =0
Get the heat transfer coefficient.
Base class for inter-region exchange.
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
static autoPtr< heatTransferModel > New(const dictionary &dict, const fvMesh &mesh)
Select from dictionary and mesh.
heatTransferModel(const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from dictionary and mesh.
virtual bool read(const dictionary &dict)
Read dictionary.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
TypeName("heatTransferModel")
Runtime type information.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const fvMesh & mesh() const
Return reference to the mesh.
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: PtrList.H:53
const dictionary & coeffs() const
Return coeffs dictionary.
virtual ~heatTransferModel()
Destructor.
Namespace for OpenFOAM.
declareRunTimeSelectionTable(autoPtr, heatTransferModel, mesh,(const dictionary &dict, const fvMesh &mesh),(dict, mesh))