InterfaceCompositionPhaseChangePhaseSystem.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 | Copyright (C) 2015-2016 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 
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class BasePhaseSystem>
35 (
36  const fvMesh& mesh
37 )
38 :
39  HeatAndMassTransferPhaseSystem<BasePhaseSystem>(mesh)
40 {
41  this->generatePairsAndSubModels
42  (
43  "interfaceComposition",
44  interfaceCompositionModels_
45  );
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
50 
51 template<class BasePhaseSystem>
54 {}
55 
56 
57 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
58 
59 template<class BasePhaseSystem>
62 massTransfer() const
63 {
64  // Create a mass transfer matrix for each species of each phase
65  autoPtr<phaseSystem::massTransferTable> eqnsPtr
66  (
68  );
69 
70  phaseSystem::massTransferTable& eqns = eqnsPtr();
71 
72  forAll(this->phaseModels_, phasei)
73  {
74  const phaseModel& phase = this->phaseModels_[phasei];
75 
76  const PtrList<volScalarField>& Yi = phase.Y();
77 
78  forAll(Yi, i)
79  {
80  eqns.insert
81  (
82  Yi[i].name(),
83  new fvScalarMatrix(Yi[i], dimMass/dimTime)
84  );
85  }
86  }
87 
88  // Reset the interfacial mass flow rates
90  (
92  this->phasePairs_,
93  phasePairIter
94  )
95  {
96  const phasePair& pair(phasePairIter());
97 
98  if (pair.ordered())
99  {
100  continue;
101  }
102 
103  *this->dmdt_[pair] =
104  *this->dmdtExplicit_[pair];
105 
106  *this->dmdtExplicit_[pair] =
108  }
109 
110  // Sum up the contribution from each interface composition model
112  (
113  interfaceCompositionModelTable,
114  interfaceCompositionModels_,
115  interfaceCompositionModelIter
116  )
117  {
118  const interfaceCompositionModel& compositionModel
119  (
120  interfaceCompositionModelIter()
121  );
122 
123  const phasePair& pair
124  (
125  this->phasePairs_[interfaceCompositionModelIter.key()]
126  );
127  const phaseModel& phase = pair.phase1();
128  const phaseModel& otherPhase = pair.phase2();
129  const phasePairKey key(phase.name(), otherPhase.name());
130 
131  const volScalarField& Tf(*this->Tf_[key]);
132 
133  volScalarField& dmdtExplicit(*this->dmdtExplicit_[key]);
134  volScalarField& dmdt(*this->dmdt_[key]);
135 
136  scalar dmdtSign(Pair<word>::compare(this->dmdt_.find(key).key(), key));
137 
138  const volScalarField K
139  (
140  this->massTransferModels_[key][phase.name()]->K()
141  );
142 
144  (
145  hashedWordList,
146  compositionModel.species(),
147  memberIter
148  )
149  {
150  const word& member = *memberIter;
151 
152  const word name
153  (
154  IOobject::groupName(member, phase.name())
155  );
156 
157  const word otherName
158  (
159  IOobject::groupName(member, otherPhase.name())
160  );
161 
162  const volScalarField KD
163  (
164  K*compositionModel.D(member)
165  );
166 
167  const volScalarField Yf
168  (
169  compositionModel.Yf(member, Tf)
170  );
171 
172  // Implicit transport through the phase
173  *eqns[name] +=
174  phase.rho()*KD*Yf
175  - fvm::Sp(phase.rho()*KD, eqns[name]->psi());
176 
177  // Sum the mass transfer rate
178  dmdtExplicit += dmdtSign*phase.rho()*KD*Yf;
179  dmdt -= dmdtSign*phase.rho()*KD*eqns[name]->psi();
180 
181  // Explicit transport out of the other phase
182  if (eqns.found(otherName))
183  {
184  *eqns[otherName] -=
185  otherPhase.rho()*KD*compositionModel.dY(member, Tf);
186  }
187  }
188  }
189 
190  return eqnsPtr;
191 }
192 
193 
194 template<class BasePhaseSystem>
197 {
199 
200  // This loop solves for the interface temperatures, Tf, and updates the
201  // interface composition models.
202  //
203  // The rate of heat transfer to the interface must equal the latent heat
204  // consumed at the interface, i.e.:
205  //
206  // H1*(T1 - Tf) + H2*(T2 - Tf) == mDotL
207  // == K*rho*(Yfi - Yi)*Li
208  //
209  // Yfi is likely to be a strong non-linear (typically exponential) function
210  // of Tf, so the solution for the temperature is newton-accelerated
211 
213  (
215  this->phasePairs_,
216  phasePairIter
217  )
218  {
219  const phasePair& pair(phasePairIter());
220 
221  if (pair.ordered())
222  {
223  continue;
224  }
225 
226  const phasePairKey key12(pair.first(), pair.second(), true);
227  const phasePairKey key21(pair.second(), pair.first(), true);
228 
229  volScalarField H1(this->heatTransferModels_[pair][pair.first()]->K());
230  volScalarField H2(this->heatTransferModels_[pair][pair.second()]->K());
231  dimensionedScalar HSmall("small", heatTransferModel::dimK, SMALL);
232 
233  volScalarField mDotL
234  (
235  IOobject
236  (
237  "mDotL",
238  this->mesh().time().timeName(),
239  this->mesh()
240  ),
241  this->mesh(),
243  );
244  volScalarField mDotLPrime
245  (
246  IOobject
247  (
248  "mDotLPrime",
249  this->mesh().time().timeName(),
250  this->mesh()
251  ),
252  this->mesh(),
253  dimensionedScalar("zero", mDotL.dimensions()/dimTemperature, 0)
254  );
255 
256  volScalarField& Tf = *this->Tf_[pair];
257 
258  // Add latent heats from forward and backward models
259  if (this->interfaceCompositionModels_.found(key12))
260  {
261  this->interfaceCompositionModels_[key12]->addMDotL
262  (
263  this->massTransferModels_[pair][pair.first()]->K(),
264  Tf,
265  mDotL,
266  mDotLPrime
267  );
268  }
269  if (this->interfaceCompositionModels_.found(key21))
270  {
271  this->interfaceCompositionModels_[key21]->addMDotL
272  (
273  this->massTransferModels_[pair][pair.second()]->K(),
274  Tf,
275  mDotL,
276  mDotLPrime
277  );
278  }
279 
280  // Update the interface temperature by applying one step of newton's
281  // method to the interface relation
282  Tf -=
283  (
284  H1*(Tf - pair.phase1().thermo().T())
285  + H2*(Tf - pair.phase2().thermo().T())
286  + mDotL
287  )
288  /(
289  max(H1 + H2 + mDotLPrime, HSmall)
290  );
291 
293 
294  Info<< "Tf." << pair.name()
295  << ": min = " << min(Tf.primitiveField())
296  << ", mean = " << average(Tf.primitiveField())
297  << ", max = " << max(Tf.primitiveField())
298  << endl;
299 
300  // Update the interface compositions
301  if (this->interfaceCompositionModels_.found(key12))
302  {
303  this->interfaceCompositionModels_[key12]->update(Tf);
304  }
305  if (this->interfaceCompositionModels_.found(key21))
306  {
307  this->interfaceCompositionModels_[key21]->update(Tf);
308  }
309  }
310 }
311 
312 
313 template<class BasePhaseSystem>
315 {
316  if (BasePhaseSystem::read())
317  {
318  bool readOK = true;
319 
320  // Models ...
321 
322  return readOK;
323  }
324  else
325  {
326  return false;
327  }
328 }
329 
330 
331 // ************************************************************************* //
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual ~InterfaceCompositionPhaseChangePhaseSystem()
Destructor.
label phasei
Definition: pEqn.H:27
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
Definition: phaseSystem.H:133
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
InterfaceCompositionPhaseChangePhaseSystem(const fvMesh &)
Construct from fvMesh.
static int compare(const Pair< word > &a, const Pair< word > &b)
Compare Pairs.
Definition: Pair.H:141
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the interfacial mass flow rate.
CGAL::Exact_predicates_exact_constructions_kernel K
static const dimensionSet dimK
Coefficient dimensions.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
dynamicFvMesh & mesh
heatTransferModelTable heatTransferModels_
Heat transfer models.
static word groupName(Name name, const word &group)
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
word timeName
Definition: getTimeIndex.H:3
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimEnergy
const dimensionSet dimDensity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
void correctBoundaryConditions()
Correct boundary field.
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
HashPtrTable< fvScalarMatrix, word, string::hash > massTransferTable
Definition: phaseSystem.H:118
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
virtual autoPtr< phaseSystem::massTransferTable > massTransfer() const
Return the mass transfer matrices.
virtual bool read()
Read base phaseProperties dictionary.
mixture correctThermo()
virtual void correctThermo()
Correct the thermodynamics.