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-2019 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 
27 
28 #include "BlendedInterfacialModel.H"
29 #include "heatTransferModel.H"
30 
31 #include "HashPtrTable.H"
32 
33 #include "fvcDiv.H"
34 #include "fvmSup.H"
35 #include "fvMatrix.H"
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
40 template<class BasePhaseSystem>
43 (
44  const fvMesh& mesh
45 )
46 :
47  BasePhaseSystem(mesh)
48 {
49  this->generatePairsAndSubModels
50  (
51  "heatTransfer",
52  heatTransferModels_,
53  false
54  );
55 
56  // Check that models have been specified on both sides of the interfaces
58  (
59  heatTransferModelTable,
60  heatTransferModels_,
61  heatTransferModelIter
62  )
63  {
64  const phasePair& pair = this->phasePairs_[heatTransferModelIter.key()];
65 
66  if (!heatTransferModels_[pair].first().valid())
67  {
69  << "A heat transfer model for the " << pair.phase1().name()
70  << " side of the " << pair << " pair is not specified"
71  << exit(FatalError);
72  }
73  if (!heatTransferModels_[pair].second().valid())
74  {
76  << "A heat transfer model for the " << pair.phase2().name()
77  << " side of the " << pair << " pair is not specified"
78  << exit(FatalError);
79  }
80  }
81 
82  // Calculate initial Tf-s as if there is no mass transfer
84  (
85  heatTransferModelTable,
86  heatTransferModels_,
87  heatTransferModelIter
88  )
89  {
90  const phasePair& pair = this->phasePairs_[heatTransferModelIter.key()];
91 
92  const phaseModel& phase1 = pair.phase1();
93  const phaseModel& phase2 = pair.phase2();
94 
95  const volScalarField& T1(phase1.thermo().T());
96  const volScalarField& T2(phase2.thermo().T());
97 
98  volScalarField H1(heatTransferModels_[pair].first()->K());
99  volScalarField H2(heatTransferModels_[pair].second()->K());
100  dimensionedScalar HSmall("small", heatTransferModel::dimK, small);
101 
102  Tf_.insert
103  (
104  pair,
105  new volScalarField
106  (
107  IOobject
108  (
109  IOobject::groupName("Tf", pair.name()),
110  this->mesh().time().timeName(),
111  this->mesh(),
114  ),
115  (H1*T1 + H2*T2)/max(H1 + H2, HSmall)
116  )
117  );
118  }
119 }
120 
121 
122 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
123 
124 template<class BasePhaseSystem>
127 {}
128 
129 
130 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
131 
132 template<class BasePhaseSystem>
135 heatTransfer() const
136 {
137  autoPtr<phaseSystem::heatTransferTable> eqnsPtr
138  (
140  );
141 
142  phaseSystem::heatTransferTable& eqns = eqnsPtr();
143 
144  forAll(this->phaseModels_, phasei)
145  {
146  const phaseModel& phase = this->phaseModels_[phasei];
147 
148  eqns.insert
149  (
150  phase.name(),
151  new fvScalarMatrix(phase.thermo().he(), dimEnergy/dimTime)
152  );
153  }
154 
155  // Heat transfer with the interface
157  (
158  heatTransferModelTable,
159  heatTransferModels_,
160  heatTransferModelIter
161  )
162  {
163  const phasePair& pair
164  (
165  this->phasePairs_[heatTransferModelIter.key()]
166  );
167 
168  const volScalarField& Tf(*Tf_[pair]);
169 
170  const Pair<tmp<volScalarField>> Ks
171  (
172  heatTransferModelIter().first()->K(),
173  heatTransferModelIter().second()->K()
174  );
175 
176  forAllConstIter(phasePair, pair, iter)
177  {
178  const phaseModel& phase = iter();
179 
180  const volScalarField& he(phase.thermo().he());
181  const volScalarField Cpv(phase.thermo().Cpv());
182  const volScalarField& K(Ks[iter.index()]);
183 
184  *eqns[phase.name()] +=
185  K*(Tf - phase.thermo().T())
186  + K/Cpv*he - fvm::Sp(K/Cpv, he);
187  }
188  }
189 
190  // Source term due to mass transfer
192  (
194  this->phasePairs_,
195  phasePairIter
196  )
197  {
198  const phasePair& pair(phasePairIter());
199 
200  if (pair.ordered())
201  {
202  continue;
203  }
204 
205  const phaseModel& phase1 = pair.phase1();
206  const phaseModel& phase2 = pair.phase2();
207 
208  const volScalarField& he1(phase1.thermo().he());
209  const volScalarField& he2(phase2.thermo().he());
210 
211  const volScalarField K1(phase1.K());
212  const volScalarField K2(phase2.K());
213 
214  const volScalarField dmdt(this->dmdt(pair));
215  const volScalarField dmdt21(posPart(dmdt));
216  const volScalarField dmdt12(negPart(dmdt));
217 
218  *eqns[phase1.name()] += - fvm::Sp(dmdt21, he1) + dmdt21*(K2 - K1);
219 
220  *eqns[phase2.name()] -= - fvm::Sp(dmdt12, he2) + dmdt12*(K1 - K2);
221 
222  if (this->heatTransferModels_.found(phasePairIter.key()))
223  {
224  const volScalarField& Tf(*Tf_[pair]);
225 
226  *eqns[phase1.name()] +=
227  dmdt21*phase1.thermo().he(phase1.thermo().p(), Tf);
228 
229  *eqns[phase2.name()] -=
230  dmdt12*phase2.thermo().he(phase2.thermo().p(), Tf);
231  }
232  else
233  {
234  *eqns[phase1.name()] += dmdt21*he2;
235 
236  *eqns[phase2.name()] -= dmdt12*he1;
237  }
238  }
239 
240  return eqnsPtr;
241 }
242 
243 
244 template<class BasePhaseSystem>
247 {
248  BasePhaseSystem::correctEnergyTransport();
249 
251 }
252 
253 
254 template<class BasePhaseSystem>
257 {
259  (
260  heatTransferModelTable,
261  heatTransferModels_,
262  heatTransferModelIter
263  )
264  {
265  const phasePair& pair
266  (
267  this->phasePairs_[heatTransferModelIter.key()]
268  );
269 
270  const phaseModel& phase1 = pair.phase1();
271  const phaseModel& phase2 = pair.phase2();
272 
273  const volScalarField& p(phase1.thermo().p());
274 
275  const volScalarField& T1(phase1.thermo().T());
276  const volScalarField& T2(phase2.thermo().T());
277 
278  volScalarField& Tf(*this->Tf_[pair]);
279 
280  const volScalarField L
281  (
282  phase1.thermo().he(p, Tf) - phase2.thermo().he(p, Tf)
283  );
284 
285  const volScalarField dmdt(this->dmdt(pair));
286 
287  volScalarField H1
288  (
289  this->heatTransferModels_[pair].first()->K()
290  );
291 
292  volScalarField H2
293  (
294  this->heatTransferModels_[pair].second()->K()
295  );
296 
297  // Limit the H[12] to avoid /0
298  H1.max(small);
299  H2.max(small);
300 
301  Tf = (H1*T1 + H2*T2 + dmdt*L)/(H1 + H2);
302 
303  Info<< "Tf." << pair.name()
304  << ": min = " << min(Tf.primitiveField())
305  << ", mean = " << average(Tf.primitiveField())
306  << ", max = " << max(Tf.primitiveField())
307  << endl;
308  }
309 }
310 
311 
312 template<class BasePhaseSystem>
314 {
315  if (BasePhaseSystem::read())
316  {
317  bool readOK = true;
318 
319  // Models ...
320 
321  return readOK;
322  }
323  else
324  {
325  return false;
326  }
327 }
328 
329 
330 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:50
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
#define K1
Definition: SHA1.C:167
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
Definition: phaseSystem.H:87
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
#define K2
Definition: SHA1.C:168
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the mass transfer rate for a pair.
HashPtrTable< fvScalarMatrix > heatTransferTable
Definition: phaseSystem.H:77
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
dimensionedScalar posPart(const dimensionedScalar &ds)
static const dimensionSet dimK
Coefficient dimensions.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual void correctInterfaceThermo()
Correct the interface thermodynamics.
dynamicFvMesh & mesh
static word groupName(Name name, const word &group)
virtual void correctInterfaceThermo()
Correct the interface thermodynamics.
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Calculate the divergence of the given field.
volScalarField & he2
Definition: EEqns.H:3
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const dimensionSet dimEnergy
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual ~TwoResistanceHeatTransferPhaseSystem()
Destructor.
messageStream Info
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
volScalarField & p
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat and Tf.
Calculate the matrix for implicit and explicit sources.
dimensionedScalar negPart(const dimensionedScalar &ds)
TwoResistanceHeatTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual bool read()
Read base phaseProperties dictionary.