All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
HeatTransferPhaseSystem.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) 2015-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 Class
25  Foam::HeatTransferPhaseSystem
26 
27 Description
28  ...
29 
30 SourceFiles
31  HeatTransferPhaseSystem.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef HeatTransferPhaseSystem_H
36 #define HeatTransferPhaseSystem_H
37 
39 #include "phaseSystem.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 /*---------------------------------------------------------------------------*\
47  Class HeatTransferPhaseSystem Declaration
48 \*---------------------------------------------------------------------------*/
49 
50 template<class BasePhaseSystem>
52 :
54  public BasePhaseSystem
55 {
56 protected:
57 
58  // Protected Member Functions
59 
60  //- Add energy transfer terms which result from bulk mass transfers
61  void addDmdtHefs
62  (
63  const phaseSystem::dmdtfTable& dmdtfs,
65  ) const;
66 
67  //- Add energy transfer terms which result from specie mass transfers
68  void addDmidtHefs
69  (
70  const phaseSystem::dmidtfTable& dmidtfs,
72  ) const;
73 
74  //- Add energy transfer terms which result from bulk phase changes,
75  // without the latent heat contribution
77  (
78  const phaseSystem::dmdtfTable& dmdtfs,
79  const phaseSystem::dmdtfTable& Tfs,
82  ) const;
83 
84  //- Add latent heat terms which result from bulk phase changes.
85  // Weight is the proportion of the latent heat contribution applied to
86  // the downwind side of the mass transfer.
87  void addDmdtL
88  (
89  const phaseSystem::dmdtfTable& dmdtfs,
90  const phaseSystem::dmdtfTable& Tfs,
91  const scalar weight,
92  const latentHeatScheme scheme,
94  ) const;
95 
96  //- Add energy transfer terms which result from bulk phase changes.
97  // Weight is the proportion of the latent heat contribution applied to
98  // the downwind side of the mass transfer.
99  void addDmdtHefs
100  (
101  const phaseSystem::dmdtfTable& dmdtfs,
102  const phaseSystem::dmdtfTable& Tfs,
103  const scalar weight,
104  const latentHeatScheme scheme,
106  ) const;
107 
108  //- Add energy transfer terms which result from specie phase changes,
109  // without the latent heat contribution
111  (
112  const phaseSystem::dmidtfTable& dmidtfs,
113  const phaseSystem::dmdtfTable& Tfs,
114  const latentHeatScheme scheme,
116  ) const;
117 
118  //- Add latent heat terms which result from specie phase changes.
119  // Weight is the proportion of the latent heat contribution applied to
120  // the downwind side of the mass transfer.
121  void addDmidtL
122  (
123  const phaseSystem::dmidtfTable& dmidtfs,
124  const phaseSystem::dmdtfTable& Tfs,
125  const scalar weight,
126  const latentHeatScheme scheme,
128  ) const;
129 
130  //- Add energy transfer terms which result from specie phase changes
131  // Weight is the proportion of the latent heat contribution applied to
132  // the downwind side of the mass transfer.
133  void addDmidtHefs
134  (
135  const phaseSystem::dmidtfTable& dmidtfs,
136  const phaseSystem::dmdtfTable& Tfs,
137  const scalar weight,
138  const latentHeatScheme scheme,
140  ) const;
141 
142 
143  // Protected data
144 
145  //- Residual mass fraction used for linearisation of heat transfer
146  // terms that result from mass transfers of individual species
147  scalar residualY_;
148 
149 
150 
151 public:
152 
153  // Constructors
154 
155  //- Construct from fvMesh
157 
158 
159  //- Destructor
160  virtual ~HeatTransferPhaseSystem();
161 
162 
163  // Member Functions
164 
165  //- Return the latent heat for a given pair, mass transfer rate (used
166  // only for it's sign), and interface temperature
167  virtual tmp<volScalarField> L
168  (
169  const phasePair& pair,
170  const volScalarField& dmdtf,
171  const volScalarField& Tf,
172  const latentHeatScheme scheme
173  ) const;
174 
175  //- As above, but for a cell-set
176  virtual tmp<scalarField> L
177  (
178  const phasePair& pair,
179  const scalarField& dmdtf,
180  const scalarField& Tf,
181  const labelUList& cells,
182  const latentHeatScheme scheme
183  ) const;
184 
185  //- Return the latent heat for a given pair, specie, mass transfer rate
186  // (used only for it's sign), and interface temperature
187  virtual tmp<volScalarField> Li
188  (
189  const phasePair& pair,
190  const word& member,
191  const volScalarField& dmidtf,
192  const volScalarField& Tf,
193  const latentHeatScheme scheme
194  ) const;
195 
196  //- As above, but for a cell-set
197  virtual tmp<scalarField> Li
198  (
199  const phasePair& pair,
200  const word& member,
201  const scalarField& dmidtf,
202  const scalarField& Tf,
203  const labelUList& cells,
204  const latentHeatScheme scheme
205  ) const;
206 
207  //- Read base phaseProperties dictionary
208  virtual bool read();
209 };
210 
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 } // End namespace Foam
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 #ifdef NoRepository
219  #include "HeatTransferPhaseSystem.C"
220 #endif
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 #endif
225 
226 // ************************************************************************* //
virtual tmp< volScalarField > L(const phasePair &pair, const volScalarField &dmdtf, const volScalarField &Tf, const latentHeatScheme scheme) const
Return the latent heat for a given pair, mass transfer rate (used.
scalar residualY_
Residual mass fraction used for linearisation of heat transfer.
void addDmdtHefsWithoutL(const phaseSystem::dmdtfTable &dmdtfs, const phaseSystem::dmdtfTable &Tfs, const latentHeatScheme scheme, phaseSystem::heatTransferTable &eqns) const
Add energy transfer terms which result from bulk phase changes,.
static tmp< surfaceInterpolationScheme< Type > > scheme(const surfaceScalarField &faceFlux, Istream &schemeData)
Return weighting factors for scheme given from Istream.
void addDmdtHefs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::heatTransferTable &eqns) const
Add energy transfer terms which result from bulk mass transfers.
void addDmidtHefs(const phaseSystem::dmidtfTable &dmidtfs, phaseSystem::heatTransferTable &eqns) const
Add energy transfer terms which result from specie mass transfers.
const cellShapeList & cells
A class for handling words, derived from string.
Definition: word.H:59
virtual ~HeatTransferPhaseSystem()
Destructor.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
latentHeatScheme
Enumeration for the latent heat scheme.
void addDmidtL(const phaseSystem::dmidtfTable &dmidtfs, const phaseSystem::dmdtfTable &Tfs, const scalar weight, const latentHeatScheme scheme, phaseSystem::heatTransferTable &eqns) const
Add latent heat terms which result from specie phase changes.
void addDmdtL(const phaseSystem::dmdtfTable &dmdtfs, const phaseSystem::dmdtfTable &Tfs, const scalar weight, const latentHeatScheme scheme, phaseSystem::heatTransferTable &eqns) const
Add latent heat terms which result from bulk phase changes.
void addDmidtHefsWithoutL(const phaseSystem::dmidtfTable &dmidtfs, const phaseSystem::dmdtfTable &Tfs, const latentHeatScheme scheme, phaseSystem::heatTransferTable &eqns) const
Add energy transfer terms which result from specie phase changes,.
virtual tmp< volScalarField > Li(const phasePair &pair, const word &member, const volScalarField &dmidtf, const volScalarField &Tf, const latentHeatScheme scheme) const
Return the latent heat for a given pair, specie, mass transfer rate.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A class for managing temporary objects.
Definition: PtrList.H:53
HeatTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual bool read()
Read base phaseProperties dictionary.
Namespace for OpenFOAM.