effectivenessHeatExchanger.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) 2013-2025 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::effectivenessHeatExchanger
26 
27 Description
28  Heat exchanger model, based on an effectiveness
29 
30  The total heat exchange source is given by:
31  \f[
32  Q_t = e(\phi, \dot{m}_2) (T_2 - T_1) \phi c_p
33  \f]
34 
35  where:
36  \vartable
37  Q_t | total heat source
38  e(\phi,\dot{m}_2) | effectiveness table
39  \phi | net mass flux entering heat exchanger [kg/s]
40  \dot{m}_2 | secondary mass flow rate [kg/s]
41  T_1 | primary inlet temperature [K]
42  T_2 | secondary inlet temperature [K]
43  c_p | specific heat capacity [J/kg/K]
44  \endvartable
45 
46  The distribution inside the hear exchanger is given by:
47  \f[
48  Q_c = \frac{V_c |U_c| (T_c - T_{ref})}{\sum(V_c |U_c| (T_c - T_{ref}))}
49  \f]
50 
51  where:
52  \vartable
53  Q_c | source for cell
54  V_c | volume of the cell [m^3]
55  U_c | local cell velocity [m/s]
56  T_c | local call temperature [K]
57  T_{ref} | min or max(T) in cell zone depending on the sign of Q_t [K]
58  \endvartable
59 
60 Usage
61  Example usage:
62  \verbatim
63  effectivenessHeatExchanger1
64  {
65  type effectivenessHeatExchanger;
66 
67  cellZone porosity;
68 
69  secondaryMassFlowRate 1.0;
70  secondaryInletT 336;
71  primaryInletT 293;
72 
73  faceZone facesZoneInletOriented;
74 
75  effectiveness <function2>;
76  }
77  \endverbatim
78 
79  Note:
80  - The effectiveness Function2 is described in terms of the primary and
81  secondary mass flow rates and has the same units as the secondary mass
82  flow rate and kg/s for phi.
83  - faceZone is the faces at the inlet of the cellzone, it needs to be
84  created with flip map flags. It is used to integrate the net mass flow
85  rate into the heat exchanger
86 
87 SourceFiles
88  effectivenessHeatExchanger.C
89 
90 \*---------------------------------------------------------------------------*/
91 
92 #ifndef effectivenessHeatExchanger_H
93 #define effectivenessHeatExchanger_H
94 
95 #include "fvModel.H"
96 #include "fvCellZone.H"
97 #include "Function2.H"
98 
99 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
100 
101 namespace Foam
102 {
103 namespace fv
104 {
105 
106 /*---------------------------------------------------------------------------*\
107  Class effectivenessHeatExchanger Declaration
108 \*---------------------------------------------------------------------------*/
109 
110 class effectivenessHeatExchanger
111 :
112  public fvModel
113 {
114  // Private data
115 
116  //- The cellZone the fvConstraint applies to
117  fvCellZone zone_;
118 
119  //- Secondary flow mass rate [kg/s]
120  scalar secondaryMassFlowRate_;
121 
122  //- Inlet secondary temperature [K]
123  scalar secondaryInletT_;
124 
125  //- Primary air temperature at the heat exchanger inlet [K]
126  scalar primaryInletT_;
127 
128  //- 2D function for effectiveness [kg/s]
129  // function of primary and secondary mass flow rates [kg/s]
130  autoPtr<Function2<scalar>> eTable_;
131 
132  //- Name of velocity field; default = U
133  word UName_;
134 
135  //- Name of temperature field; default = T
136  word TName_;
137 
138  //- Name of the flux
139  word phiName_;
140 
141  //- Name of the faceZone at the heat exchange inlet
142  word faceZoneName_;
143 
144  //- Id for the face zone
145  label zoneIndex_;
146 
147  //- Local list of face IDs
148  labelList faceId_;
149 
150  //- Local list of patch ID per face
151  labelList facePatchId_;
152 
153  //- List of +1/-1 representing face flip map (1 use as is, -1 negate)
154  labelList faceSign_;
155 
156  //- Area of the face zone
157  scalar faceZoneArea_;
158 
159 
160 private:
161 
162  // Private Member Functions
163 
164  //- Non-virtual read
165  void readCoeffs(const dictionary& dict);
166 
167  //- Set the zone information
168  void setZone();
169 
170  //- Calculate total area of faceZone across processors
171  void calculateTotalArea(scalar& area) const;
172 
173 
174 public:
175 
176  //- Runtime type information
177  TypeName("effectivenessHeatExchanger");
178 
179 
180  // Constructors
181 
182  //- Construct from components
184  (
185  const word& name,
186  const word& modelType,
187  const fvMesh& mesh,
188  const dictionary& dict
189  );
190 
191  //- Disallow default bitwise copy construction
193  (
195  ) = delete;
196 
197 
198  //- Destructor
200  {}
201 
202 
203  // Member Functions
204 
205  //- Return the list of fields for which the fvModel adds source term
206  // to the transport equation
207  virtual wordList addSupFields() const;
208 
209  //- Explicit and implicit source for compressible equation
210  virtual void addSup
211  (
212  const volScalarField& rho,
213  const volScalarField& he,
214  fvMatrix<scalar>& eqn
215  ) const;
216 
217 
218  // Mesh changes
219 
220  //- Update for mesh motion
221  virtual bool movePoints();
222 
223  //- Update topology using the given map
224  virtual void topoChange(const polyTopoChangeMap&);
225 
226  //- Update from another mesh using the given map
227  virtual void mapMesh(const polyMeshMap&);
228 
229  //- Redistribute or update using the given distribution map
230  virtual void distribute(const polyDistributionMap&);
231 
232 
233  // IO
234 
235  //- Read dictionary
236  virtual bool read(const dictionary& dict);
237 
238 
239  // Member Operators
240 
241  //- Disallow default bitwise assignment
242  void operator=(const effectivenessHeatExchanger&) = delete;
243 };
244 
245 
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 
248 } // End namespace fv
249 } // End namespace Foam
250 
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 
253 #endif
254 
255 // ************************************************************************* //
Generic GeometricField class.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:69
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:57
Heat exchanger model, based on an effectiveness.
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
void operator=(const effectivenessHeatExchanger &)=delete
Disallow default bitwise assignment.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
effectivenessHeatExchanger(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
TypeName("effectivenessHeatExchanger")
Runtime type information.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read dictionary.
virtual void addSup(const volScalarField &rho, const volScalarField &he, fvMatrix< scalar > &eqn) const
Explicit and implicit source for compressible equation.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
thermo he()
labelList fv(nPoints)
dictionary dict