TwoResistanceHeatTransferPhaseSystem.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2019 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 
30 #include "BlendedInterfacialModel.H"
31 #include "heatTransferModel.H"
32 
33 #include "HashPtrTable.H"
34 
35 #include "fvcDiv.H"
36 #include "fvmSup.H"
37 #include "fvMatrix.H"
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class BasePhaseSystem>
45 (
46  const fvMesh& mesh
47 )
48 :
49  BasePhaseSystem(mesh)
50 {
51  this->generatePairsAndSubModels
52  (
53  "heatTransfer",
54  heatTransferModels_,
55  false
56  );
57 
58  // Check that models have been specified on both sides of the interfaces
60  (
62  heatTransferModels_,
63  heatTransferModelIter
64  )
65  {
66  const phasePair& pair = this->phasePairs_[heatTransferModelIter.key()];
67 
68  if (!heatTransferModels_[pair].first())
69  {
71  << "A heat transfer model for the " << pair.phase1().name()
72  << " side of the " << pair << " pair is not specified"
73  << exit(FatalError);
74  }
75  if (!heatTransferModels_[pair].second())
76  {
78  << "A heat transfer model for the " << pair.phase2().name()
79  << " side of the " << pair << " pair is not specified"
80  << exit(FatalError);
81  }
82  }
83 
84  // Calculate initial Tf-s as if there is no mass transfer
86  (
88  heatTransferModels_,
89  heatTransferModelIter
90  )
91  {
92  const phasePair& pair = this->phasePairs_[heatTransferModelIter.key()];
93 
94  const phaseModel& phase1 = pair.phase1();
95  const phaseModel& phase2 = pair.phase2();
96 
97  const volScalarField& T1(phase1.thermo().T());
98  const volScalarField& T2(phase2.thermo().T());
99 
100  volScalarField H1(heatTransferModels_[pair].first()->K());
101  volScalarField H2(heatTransferModels_[pair].second()->K());
102  dimensionedScalar HSmall("small", heatTransferModel::dimK, SMALL);
103 
104  Tf_.set
105  (
106  pair,
107  new volScalarField
108  (
109  IOobject
110  (
111  IOobject::groupName("Tf", pair.name()),
112  this->mesh().time().timeName(),
113  this->mesh(),
114  IOobject::NO_READ,
115  IOobject::AUTO_WRITE
116  ),
117  (H1*T1 + H2*T2)/max(H1 + H2, HSmall)
118  )
119  );
120  }
121 }
122 
123 
124 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
125 
126 template<class BasePhaseSystem>
129 {}
130 
131 
132 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
133 
134 template<class BasePhaseSystem>
138 {
140  (
142  );
143 
144  phaseSystem::heatTransferTable& eqns = eqnsPtr();
145 
146  forAll(this->phaseModels_, phasei)
147  {
148  const phaseModel& phase = this->phaseModels_[phasei];
149 
150  eqns.set
151  (
152  phase.name(),
153  new fvScalarMatrix(phase.thermo().he(), dimEnergy/dimTime)
154  );
155  }
156 
157  // Heat transfer with the interface
159  (
161  heatTransferModels_,
162  heatTransferModelIter
163  )
164  {
165  const phasePair& pair
166  (
167  this->phasePairs_[heatTransferModelIter.key()]
168  );
169 
170  const volScalarField& Tf(*Tf_[pair]);
171 
172  const Pair<tmp<volScalarField>> Ks
173  (
174  heatTransferModelIter().first()->K(),
175  heatTransferModelIter().second()->K()
176  );
177 
178  forAllConstIter(phasePair, pair, iter)
179  {
180  const phaseModel& phase = iter();
181 
182  const volScalarField& he(phase.thermo().he());
183  const volScalarField Cpv(phase.thermo().Cpv());
184  const volScalarField& K(Ks[iter.index()]);
185 
186  *eqns[phase.name()] +=
187  K*(Tf - phase.thermo().T())
188  + K/Cpv*he - fvm::Sp(K/Cpv, he);
189  }
190  }
191 
192  // Source term due to mass transfer
194  (
196  this->phasePairs_,
197  phasePairIter
198  )
199  {
200  const phasePair& pair(phasePairIter());
201 
202  if (pair.ordered())
203  {
204  continue;
205  }
206 
207  const phaseModel& phase1 = pair.phase1();
208  const phaseModel& phase2 = pair.phase2();
209 
210  const volScalarField& he1(phase1.thermo().he());
211  const volScalarField& he2(phase2.thermo().he());
212 
213  const volScalarField K1(phase1.K());
214  const volScalarField K2(phase2.K());
215 
216  const volScalarField dmdt(this->dmdt(pair));
217  const volScalarField dmdt21(posPart(dmdt));
218  const volScalarField dmdt12(negPart(dmdt));
219 
220  *eqns[phase1.name()] += - fvm::Sp(dmdt21, he1) + dmdt21*(K2 - K1);
221 
222  *eqns[phase2.name()] -= - fvm::Sp(dmdt12, he2) + dmdt12*(K1 - K2);
223 
224  if (this->heatTransferModels_.found(phasePairIter.key()))
225  {
226  const volScalarField& Tf(*Tf_[pair]);
227 
228  *eqns[phase1.name()] +=
229  dmdt21*phase1.thermo().he(phase1.thermo().p(), Tf);
230 
231  *eqns[phase2.name()] -=
232  dmdt12*phase2.thermo().he(phase2.thermo().p(), Tf);
233  }
234  else
235  {
236  *eqns[phase1.name()] += dmdt21*he2;
237 
238  *eqns[phase2.name()] -= dmdt12*he1;
239  }
240  }
241 
242  return eqnsPtr;
243 }
244 
245 
246 template<class BasePhaseSystem>
249 {
250  BasePhaseSystem::correctEnergyTransport();
251 
252  correctInterfaceThermo();
253 }
254 
255 
256 template<class BasePhaseSystem>
259 {
261  (
263  heatTransferModels_,
264  heatTransferModelIter
265  )
266  {
267  const phasePair& pair
268  (
269  this->phasePairs_[heatTransferModelIter.key()]
270  );
271 
272  const phaseModel& phase1 = pair.phase1();
273  const phaseModel& phase2 = pair.phase2();
274 
275  const volScalarField& p(phase1.thermo().p());
276 
277  const volScalarField& T1(phase1.thermo().T());
278  const volScalarField& T2(phase2.thermo().T());
279 
280  volScalarField& Tf(*this->Tf_[pair]);
281 
282  const volScalarField L
283  (
284  phase1.thermo().he(p, Tf) - phase2.thermo().he(p, Tf)
285  );
286 
287  const volScalarField dmdt(this->dmdt(pair));
288 
289  volScalarField H1
290  (
291  this->heatTransferModels_[pair].first()->K()
292  );
293 
294  volScalarField H2
295  (
296  this->heatTransferModels_[pair].second()->K()
297  );
298 
299  // Limit the H[12] to avoid /0
300  H1.max(SMALL);
301  H2.max(SMALL);
302 
303  Tf = (H1*T1 + H2*T2 + dmdt*L)/(H1 + H2);
304 
305  Info<< "Tf." << pair.name()
306  << ": min = " << min(Tf.primitiveField())
307  << ", mean = " << average(Tf.primitiveField())
308  << ", max = " << max(Tf.primitiveField())
309  << endl;
310  }
311 }
312 
313 
314 template<class BasePhaseSystem>
316 {
317  if (BasePhaseSystem::read())
318  {
319  bool readOK = true;
320 
321  // Models ...
322 
323  return readOK;
324  }
325  else
326  {
327  return false;
328  }
329 }
330 
331 
332 // ************************************************************************* //
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
L
const vector L(dict.get< vector >("L"))
Foam::phasePair
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:53
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::TwoResistanceHeatTransferPhaseSystem::TwoResistanceHeatTransferPhaseSystem
TwoResistanceHeatTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
Definition: TwoResistanceHeatTransferPhaseSystem.C:45
K1
#define K1
Definition: SHA1.C:144
Foam::phase
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:54
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::posPart
dimensionedScalar posPart(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:221
Sp
zeroField Sp
Definition: alphaSuSp.H:2
fvcDiv.H
Calculate the divergence of the given field.
phasei
label phasei
Definition: pEqn.H:27
fvMatrix.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
he2
volScalarField & he2
Definition: EEqns.H:3
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAllConstIter
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:344
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:58
TwoResistanceHeatTransferPhaseSystem.H
Foam::blockMeshTools::read
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:57
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::phasePairKey::ordered
bool ordered() const noexcept
Return the ordered flag.
Definition: phasePairKey.H:98
Foam::HashPtrTable::set
bool set(const Key &key, T *ptr)
Assign a new entry, overwriting existing entries.
Definition: HashPtrTableI.H:135
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:44
Foam::fac::average
tmp< GeometricField< Type, faPatchField, areaMesh > > average(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Area-weighted average a edgeField creating a areaField.
Definition: facAverage.C:47
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::TwoResistanceHeatTransferPhaseSystem::correctEnergyTransport
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat and Tf.
Definition: TwoResistanceHeatTransferPhaseSystem.C:248
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
K2
#define K2
Definition: SHA1.C:145
phase2
phaseModel & phase2
Definition: setRegionFluidFields.H:6
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned< scalar >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::TwoResistanceHeatTransferPhaseSystem::heatTransfer
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
Definition: TwoResistanceHeatTransferPhaseSystem.C:137
Foam::TwoResistanceHeatTransferPhaseSystem::~TwoResistanceHeatTransferPhaseSystem
virtual ~TwoResistanceHeatTransferPhaseSystem()
Destructor.
Definition: TwoResistanceHeatTransferPhaseSystem.C:128
fvmSup.H
Calculate the matrix for implicit and explicit sources.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::phase::name
const word & name() const
Definition: phase.H:111
Foam::GeometricField::T
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
Definition: GeometricField.C:1046
Foam::HashTable< Pair< autoPtr< BlendedInterfacialModel< heatTransferModel > > >, phasePairKey, phasePairKey::hash >
Foam::phaseModel::name
const word & name() const
Definition: phaseModel.H:140
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::negPart
dimensionedScalar negPart(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:232
he
volScalarField & he
Definition: YEEqn.H:52
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:54
Foam::HashPtrTable
A HashTable of pointers to objects of type <T>, with deallocation management of the pointers.
Definition: HashPtrTable.H:54
Foam::TwoResistanceHeatTransferPhaseSystem::read
virtual bool read()
Read base phaseProperties dictionary.
Definition: TwoResistanceHeatTransferPhaseSystem.C:315
Foam::TwoResistanceHeatTransferPhaseSystem::correctInterfaceThermo
virtual void correctInterfaceThermo()
Correct the interface thermodynamics.
Definition: TwoResistanceHeatTransferPhaseSystem.C:258
Foam::phasePair::name
virtual word name() const
Pair name.
Definition: phasePair.C:69
HashPtrTable.H
Foam::phasePair::phase2
const phaseModel & phase2() const
Definition: phasePairI.H:36
Foam::GeometricField::max
void max(const dimensioned< Type > &dt)
Use the maximum of the field and specified value.
Definition: GeometricField.C:1143
phase1
phaseModel & phase1
Definition: setRegionFluidFields.H:5
Foam::phasePair::phase1
const phaseModel & phase1() const
Definition: phasePairI.H:30
Foam::GeometricField< scalar, fvPatchField, volMesh >
zeroGradientFvPatchFields.H