TwoResistanceHeatTransferPhaseSystem.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) 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 \*---------------------------------------------------------------------------*/
25 
28 #include "heatTransferModel.H"
29 #include "fvmSup.H"
30 #include "rhoReactionThermo.H"
31 
32 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
33 
34 template<class BasePhaseSystem>
36 (
37  const phaseSystem::dmdtfTable& dmdtfs,
38  const phaseSystem::dmdtfTable& Tfs,
39  const latentHeatScheme scheme,
40  const latentHeatTransfer transfer,
41  phaseSystem::heatTransferTable& eqns
42 ) const
43 {
44  HeatTransferPhaseSystem<BasePhaseSystem>::addDmdtHefsWithoutL
45  (
46  dmdtfs,
47  Tfs,
48  scheme,
49  eqns
50  );
51 
52  // Loop the pairs
53  forAllConstIter(phaseSystem::dmdtfTable, dmdtfs, dmdtfIter)
54  {
55  const phasePairKey& key = dmdtfIter.key();
56  const phasePair& pair(this->phasePairs_[key]);
57 
58  const volScalarField dmdtf(Pair<word>::compare(pair, key)**dmdtfIter());
59 
60  const volScalarField& Tf = *Tfs[key];
61 
62  const phaseModel& phase1 = pair.phase1();
63  const phaseModel& phase2 = pair.phase2();
64  const rhoThermo& thermo1 = phase1.thermo();
65  const rhoThermo& thermo2 = phase2.thermo();
66 
67  // Transfer coefficients
68  const volScalarField H1(heatTransferModels_[key].first()->K());
69  const volScalarField H2(heatTransferModels_[key].second()->K());
70  const volScalarField H1Fac(H1/(H1 + H2));
71  const volScalarField HEff(H1Fac*H2);
72 
73  // Latent heat contribution
74  switch (transfer)
75  {
76  case latentHeatTransfer::heat:
77  {
78  *eqns[phase1.name()] +=
79  - HEff*(thermo2.T() - thermo1.T()) + H1*(Tf - thermo1.T());
80 
81  *eqns[phase2.name()] +=
82  - HEff*(thermo1.T() - thermo2.T()) + H2*(Tf - thermo2.T());
83 
84  break;
85  }
86  case latentHeatTransfer::mass:
87  {
88  const volScalarField L(this->L(pair, dmdtf, Tf, scheme));
89 
90  *eqns[phase1.name()] += H1Fac*dmdtf*L;
91  *eqns[phase2.name()] += (1 - H1Fac)*dmdtf*L;
92 
93  break;
94  }
95  }
96  }
97 }
98 
99 
100 template<class BasePhaseSystem>
102 (
103  const phaseSystem::dmidtfTable& dmidtfs,
104  const phaseSystem::dmdtfTable& Tfs,
105  const latentHeatScheme scheme,
106  const latentHeatTransfer transfer,
107  phaseSystem::heatTransferTable& eqns
108 ) const
109 {
110  HeatTransferPhaseSystem<BasePhaseSystem>::addDmidtHefsWithoutL
111  (
112  dmidtfs,
113  Tfs,
114  scheme,
115  eqns
116  );
117 
118  // Loop the pairs
119  forAllConstIter(phaseSystem::dmidtfTable, dmidtfs, dmidtfIter)
120  {
121  const phasePairKey& key = dmidtfIter.key();
122  const phasePair& pair(this->phasePairs_[key]);
123 
124  const volScalarField& Tf = *Tfs[key];
125 
126  const phaseModel& phase1 = pair.phase1();
127  const phaseModel& phase2 = pair.phase2();
128  const rhoThermo& thermo1 = phase1.thermo();
129  const rhoThermo& thermo2 = phase2.thermo();
130 
131  // Transfer coefficients
132  const volScalarField H1(heatTransferModels_[key].first()->K());
133  const volScalarField H2(heatTransferModels_[key].second()->K());
134  const volScalarField H1Fac(H1/(H1 + H2));
135  const volScalarField HEff(H1Fac*H2);
136 
137  // Loop the species
138  forAllConstIter(HashPtrTable<volScalarField>, *dmidtfIter(), dmidtfJter)
139  {
140  const word& specie = dmidtfJter.key();
141 
142  // Mass transfer rates
143  const volScalarField dmidtf
144  (
145  Pair<word>::compare(pair, key)**dmidtfJter()
146  );
147 
148  // Latent heat contribution
149  switch (transfer)
150  {
151  case latentHeatTransfer::heat:
152  {
153  // Do nothing. This term is handled outside the specie loop.
154 
155  break;
156  }
157  case latentHeatTransfer::mass:
158  {
159  const volScalarField Li
160  (
161  this->Li(pair, specie, dmidtf, Tf, scheme)
162  );
163 
164  *eqns[phase1.name()] += H1Fac*dmidtf*Li;
165  *eqns[phase2.name()] += (1 - H1Fac)*dmidtf*Li;
166 
167  break;
168  }
169  }
170  }
171 
172  // Latent heat contribution
173  switch (transfer)
174  {
175  case latentHeatTransfer::heat:
176  {
177  *eqns[phase1.name()] +=
178  - HEff*(thermo2.T() - thermo1.T()) + H1*(Tf - thermo1.T());
179 
180  *eqns[phase2.name()] +=
181  - HEff*(thermo1.T() - thermo2.T()) + H2*(Tf - thermo2.T());
182 
183  break;
184  }
185  case latentHeatTransfer::mass:
186  {
187  // Do nothing. This term is handled inside the specie loop.
188 
189  break;
190  }
191  }
192  }
193 }
194 
195 
196 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
197 
198 template<class BasePhaseSystem>
201 (
202  const fvMesh& mesh
203 )
204 :
205  HeatTransferPhaseSystem<BasePhaseSystem>(mesh)
206 {
207  this->generatePairsAndSubModels
208  (
209  "heatTransfer",
210  heatTransferModels_,
211  false
212  );
213 
214  // Check that models have been specified on both sides of the interfaces
216  (
217  heatTransferModelTable,
218  heatTransferModels_,
219  heatTransferModelIter
220  )
221  {
222  const phasePair& pair = this->phasePairs_[heatTransferModelIter.key()];
223 
224  if (!heatTransferModels_[pair].first().valid())
225  {
227  << "A heat transfer model for the " << pair.phase1().name()
228  << " side of the " << pair << " pair is not specified"
229  << exit(FatalError);
230  }
231  if (!heatTransferModels_[pair].second().valid())
232  {
234  << "A heat transfer model for the " << pair.phase2().name()
235  << " side of the " << pair << " pair is not specified"
236  << exit(FatalError);
237  }
238  }
239 }
240 
241 
242 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
243 
244 template<class BasePhaseSystem>
247 {}
248 
249 
250 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
251 
252 template<class BasePhaseSystem>
255 heatTransfer() const
256 {
257  autoPtr<phaseSystem::heatTransferTable> eqnsPtr
258  (
260  );
261 
262  phaseSystem::heatTransferTable& eqns = eqnsPtr();
263 
264  forAll(this->phaseModels_, phasei)
265  {
266  const phaseModel& phase = this->phaseModels_[phasei];
267 
268  eqns.insert
269  (
270  phase.name(),
271  new fvScalarMatrix(phase.thermo().he(), dimEnergy/dimTime)
272  );
273  }
274 
276  (
277  heatTransferModelTable,
278  heatTransferModels_,
279  heatTransferModelIter
280  )
281  {
282  const phasePair& pair = this->phasePairs_[heatTransferModelIter.key()];
283 
284  const phaseModel& phase1 = pair.phase1();
285  const phaseModel& phase2 = pair.phase2();
286  const volScalarField& he1 = phase1.thermo().he();
287  const volScalarField& he2 = phase2.thermo().he();
288  const volScalarField Cpv1(phase1.thermo().Cpv());
289  const volScalarField Cpv2(phase2.thermo().Cpv());
290 
291  const volScalarField H1(heatTransferModelIter().first()->K());
292  const volScalarField H2(heatTransferModelIter().second()->K());
293  const volScalarField HEff(H1*H2/(H1 + H2));
294 
295  *eqns[phase1.name()] +=
296  HEff*(phase2.thermo().T() - phase1.thermo().T())
297  + H1/Cpv1*he1 - fvm::Sp(H1/Cpv1, he1);
298  *eqns[phase2.name()] +=
299  HEff*(phase1.thermo().T() - phase2.thermo().T())
300  + H2/Cpv2*he2 - fvm::Sp(H2/Cpv2, he2);
301  }
302 
303  return eqnsPtr;
304 }
305 
306 
307 template<class BasePhaseSystem>
310 {
311  BasePhaseSystem::correctEnergyTransport();
312 
314 }
315 
316 
317 template<class BasePhaseSystem>
319 {
320  if (BasePhaseSystem::read())
321  {
322  bool readOK = true;
323 
324  // Models ...
325 
326  return readOK;
327  }
328  else
329  {
330  return false;
331  }
332 }
333 
334 
335 // ************************************************************************* //
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
label phasei
Definition: pEqn.H:27
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
static int compare(const Pair< word > &a, const Pair< word > &b)
Compare Pairs.
Definition: Pair.H:153
HashPtrTable< fvScalarMatrix > heatTransferTable
Definition: phaseSystem.H:80
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
CGAL::Exact_predicates_exact_constructions_kernel K
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const dimensionSet dimTime
virtual void correctInterfaceThermo()
Correct the interface thermodynamics.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
const dimensionSet dimEnergy
void addDmdtHefs(const phaseSystem::dmdtfTable &dmdtfs, const phaseSystem::dmdtfTable &Tfs, const latentHeatScheme scheme, const latentHeatTransfer transfer, phaseSystem::heatTransferTable &eqns) const
Add energy transfer terms which result from bulk phase changes.
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
virtual ~TwoResistanceHeatTransferPhaseSystem()
Destructor.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat and Tf.
void addDmidtHefs(const phaseSystem::dmidtfTable &dmidtfs, const phaseSystem::dmdtfTable &Tfs, const latentHeatScheme scheme, const latentHeatTransfer transfer, phaseSystem::heatTransferTable &eqns) const
Add energy transfer terms which result from specie phase changes.
Calculate the matrix for implicit and explicit sources.
TwoResistanceHeatTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual bool read()
Read base phaseProperties dictionary.