ThermalPhaseChangePhaseSystem.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 #include "fvcVolumeIntegrate.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class BasePhaseSystem>
35 (
36  const fvMesh& mesh
37 )
38 :
39  HeatAndMassTransferPhaseSystem<BasePhaseSystem>(mesh),
40  volatile_(this->lookup("volatile")),
41  saturationModel_(saturationModel::New(this->subDict("saturationModel"))),
42  massTransfer_(this->lookup("massTransfer"))
43 {
44 
46  (
48  this->phasePairs_,
49  phasePairIter
50  )
51  {
52  const phasePair& pair(phasePairIter());
53 
54  if (pair.ordered())
55  {
56  continue;
57  }
58 
59  // Initially assume no mass transfer
60  iDmdt_.insert
61  (
62  pair,
63  new volScalarField
64  (
65  IOobject
66  (
67  IOobject::groupName("iDmdt", pair.name()),
68  this->mesh().time().timeName(),
69  this->mesh(),
72  ),
73  this->mesh(),
75  )
76  );
77  }
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
82 
83 template<class BasePhaseSystem>
86 {}
87 
88 
89 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
90 
91 template<class BasePhaseSystem>
94 {
95  return saturationModel_();
96 }
97 
98 
99 template<class BasePhaseSystem>
102 {
103  typedef compressible::alphatPhaseChangeWallFunctionFvPatchScalarField
104  alphatPhaseChangeWallFunction;
105 
106  autoPtr<phaseSystem::heatTransferTable> eqnsPtr =
108 
109  phaseSystem::heatTransferTable& eqns = eqnsPtr();
110 
111  // Accumulate mDotL contributions from boundaries
113  (
115  this->phasePairs_,
116  phasePairIter
117  )
118  {
119  const phasePair& pair(phasePairIter());
120 
121  if (pair.ordered())
122  {
123  continue;
124  }
125 
126  const phaseModel& phase = pair.phase1();
127  const phaseModel& otherPhase = pair.phase2();
128 
129  volScalarField mDotL
130  (
131  IOobject
132  (
133  "mDotL",
134  phase.mesh().time().timeName(),
135  phase.mesh(),
138  false
139  ),
140  phase.mesh(),
141  dimensionedScalar("",dimensionSet(1,-1,-3,0,0),0.0)
142  );
143 
144  if
145  (
146  otherPhase.mesh().foundObject<volScalarField>
147  (
148  "alphat." + otherPhase.name()
149  )
150  )
151  {
152  const volScalarField& alphat =
153  otherPhase.mesh().lookupObject<volScalarField>
154  (
155  "alphat." + otherPhase.name()
156  );
157 
158  const fvPatchList& patches = this->mesh().boundary();
159  forAll(patches, patchi)
160  {
161  const fvPatch& currPatch = patches[patchi];
162 
163  if
164  (
165  isA<alphatPhaseChangeWallFunction>
166  (
167  alphat.boundaryField()[patchi]
168  )
169  )
170  {
171  const scalarField& patchMDotL =
172  refCast<const alphatPhaseChangeWallFunction>
173  (
174  alphat.boundaryField()[patchi]
175  ).mDotL();
176 
177  forAll(patchMDotL,facei)
178  {
179  label faceCelli = currPatch.faceCells()[facei];
180  mDotL[faceCelli] = patchMDotL[facei];
181  }
182  }
183  }
184  }
185 
186  *eqns[otherPhase.name()] -= mDotL;
187 
188  }
189 
190  return eqnsPtr;
191 }
192 
193 
194 template<class BasePhaseSystem>
197 {
198  // Create a mass transfer matrix for each species of each phase
199  autoPtr<phaseSystem::massTransferTable> eqnsPtr
200  (
202  );
203 
204  phaseSystem::massTransferTable& eqns = eqnsPtr();
205 
206  forAll(this->phaseModels_, phasei)
207  {
208  const phaseModel& phase = this->phaseModels_[phasei];
209 
210  const PtrList<volScalarField>& Yi = phase.Y();
211 
212  forAll(Yi, i)
213  {
214  eqns.insert
215  (
216  Yi[i].name(),
217  new fvScalarMatrix(Yi[i], dimMass/dimTime)
218  );
219  }
220  }
221 
223  (
225  this->phasePairs_,
226  phasePairIter
227  )
228  {
229  const phasePair& pair(phasePairIter());
230 
231  if (pair.ordered())
232  {
233  continue;
234  }
235  const phaseModel& phase = pair.phase1();
236  const phaseModel& otherPhase = pair.phase2();
237 
238  const word name
239  (
240  IOobject::groupName(volatile_, phase.name())
241  );
242 
243  const word otherName
244  (
245  IOobject::groupName(volatile_, otherPhase.name())
246  );
247 
248  const volScalarField dmdt(this->dmdt(pair));
249  const volScalarField dmdt12(posPart(dmdt));
250  const volScalarField dmdt21(negPart(dmdt));
251 
252  *eqns[name] += fvm::Sp(dmdt21, eqns[name]->psi()) - dmdt21;
253  *eqns[otherName] += dmdt12 - fvm::Sp(dmdt12, eqns[otherName]->psi());
254  }
255 
256  return eqnsPtr;
257 }
258 
259 
260 template<class BasePhaseSystem>
262 {
263  typedef compressible::alphatPhaseChangeWallFunctionFvPatchScalarField
264  alphatPhaseChangeWallFunction;
265 
266  BasePhaseSystem::correctThermo();
267 
268 
269 
271  (
273  this->phasePairs_,
274  phasePairIter
275  )
276  {
277  const phasePair& pair(phasePairIter());
278 
279  if (pair.ordered())
280  {
281  continue;
282  }
283 
284  const phaseModel& phase1 = pair.phase1();
285  const phaseModel& phase2 = pair.phase2();
286 
287  Info<< phase1.name() << " min/max T "
288  << min(phase1.thermo().T()).value()
289  << " - "
290  << max(phase1.thermo().T()).value()
291  << endl;
292 
293  Info<< phase2.name() << " min/max T "
294  << min(phase2.thermo().T()).value()
295  << " - "
296  << max(phase2.thermo().T()).value()
297  << endl;
298 
299  const volScalarField& T1(phase1.thermo().T());
300  const volScalarField& T2(phase2.thermo().T());
301 
302  const volScalarField& he1(phase1.thermo().he());
303  const volScalarField& he2(phase2.thermo().he());
304 
305  volScalarField& dmdt(*this->dmdt_[pair]);
306  volScalarField& iDmdt(*this->iDmdt_[pair]);
307 
308  volScalarField& Tf = *this->Tf_[pair];
309 
310  volScalarField hef1(phase1.thermo().he(phase1.thermo().p(), Tf));
311  volScalarField hef2(phase2.thermo().he(phase2.thermo().p(), Tf));
312 
313  if (massTransfer_ )
314  {
315  volScalarField H1
316  (
317  this->heatTransferModels_[pair][pair.first()]->K(0)
318  );
319 
320  volScalarField H2
321  (
322  this->heatTransferModels_[pair][pair.second()]->K(0)
323  );
324 
325  Tf = saturationModel_->Tsat(phase1.thermo().p());
326 
327  scalar iDmdtRelax(this->mesh().fieldRelaxationFactor("iDmdt"));
328 
329  iDmdt =
330  (1 - iDmdtRelax)*iDmdt
331  + iDmdtRelax*(H1*(Tf - T1) + H2*(Tf - T2))
332  /min
333  (
334  (pos(iDmdt)*he2 + neg(iDmdt)*hef2)
335  - (neg(iDmdt)*he1 + pos(iDmdt)*hef1),
336  0.3*mag(hef2 - hef1)
337  );
338 
339  Info<< "iDmdt." << pair.name()
340  << ": min = " << min(iDmdt.primitiveField())
341  << ", mean = " << average(iDmdt.primitiveField())
342  << ", max = " << max(iDmdt.primitiveField())
343  << ", integral = " << fvc::domainIntegrate(iDmdt).value()
344  << endl;
345  }
346  else
347  {
348  iDmdt == dimensionedScalar("0", dmdt.dimensions(), 0);
349  }
350 
351  volScalarField H1(this->heatTransferModels_[pair][pair.first()]->K());
352  volScalarField H2(this->heatTransferModels_[pair][pair.second()]->K());
353 
354  // Limit the H[12] boundary field to avoid /0
355  const scalar HLimit = 1e-4;
356  H1.boundaryFieldRef() =
357  max(H1.boundaryField(), phase1.boundaryField()*HLimit);
358  H2.boundaryFieldRef() =
359  max(H2.boundaryField(), phase2.boundaryField()*HLimit);
360 
361  volScalarField mDotL
362  (
363  iDmdt*
364  (
365  (pos(iDmdt)*he2 + neg(iDmdt)*hef2)
366  - (neg(iDmdt)*he1 + pos(iDmdt)*hef1)
367  )
368  );
369 
370  Tf = (H1*T1 + H2*T2 + mDotL)/(H1 + H2);
371 
372  Info<< "Tf." << pair.name()
373  << ": min = " << min(Tf.primitiveField())
374  << ", mean = " << average(Tf.primitiveField())
375  << ", max = " << max(Tf.primitiveField())
376  << endl;
377 
378  // Accumulate dmdt contributions from boundaries
379  volScalarField wDmdt
380  (
381  IOobject
382  (
383  IOobject::groupName("wDmdt", pair.name()),
384  this->mesh().time().timeName(),
385  this->mesh(),
388  false
389  ),
390  this->mesh(),
392  );
393 
394  if
395  (
396  phase2.mesh().foundObject<volScalarField>
397  (
398  "alphat." + phase2.name()
399  )
400  )
401  {
402  const volScalarField& alphat =
403  phase2.mesh().lookupObject<volScalarField>
404  (
405  "alphat." + phase2.name()
406  );
407 
408  const fvPatchList& patches = this->mesh().boundary();
409  forAll(patches, patchi)
410  {
411  const fvPatch& currPatch = patches[patchi];
412 
413  if
414  (
415  isA<alphatPhaseChangeWallFunction>
416  (
417  alphat.boundaryField()[patchi]
418  )
419  )
420  {
421  const scalarField& patchDmdt =
422  refCast<const alphatPhaseChangeWallFunction>
423  (
424  alphat.boundaryField()[patchi]
425  ).dmdt();
426 
427  forAll(patchDmdt,facei)
428  {
429  label faceCelli = currPatch.faceCells()[facei];
430  wDmdt[faceCelli] += patchDmdt[facei];
431  }
432  }
433  }
434 
435  Info<< "wDmdt." << pair.name()
436  << ": min = " << min(wDmdt.primitiveField())
437  << ", mean = " << average(wDmdt.primitiveField())
438  << ", max = " << max(wDmdt.primitiveField())
439  << ", integral = " << fvc::domainIntegrate(wDmdt).value()
440  << endl;
441  }
442 
443  dmdt = wDmdt + iDmdt;
444 
445  Info<< "dmdt." << pair.name()
446  << ": min = " << min(dmdt.primitiveField())
447  << ", mean = " << average(dmdt.primitiveField())
448  << ", max = " << max(dmdt.primitiveField())
449  << ", integral = " << fvc::domainIntegrate(dmdt).value()
450  << endl;
451  }
452 }
453 
454 
455 template<class BasePhaseSystem>
457 {
458  if (BasePhaseSystem::read())
459  {
460  bool readOK = true;
461 
462  // Models ...
463 
464  return readOK;
465  }
466  else
467  {
468  return false;
469  }
470 }
471 
472 
473 // ************************************************************************* //
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
label phasei
Definition: pEqn.H:27
const double e
Elementary charge.
Definition: doubleFloat.H:78
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
Definition: phaseSystem.H:133
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
ThermalPhaseChangePhaseSystem(const fvMesh &)
Construct from fvMesh.
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
virtual void correctThermo()
Correct the thermodynamics.
const saturationModel & saturation() const
Return the saturationModel.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
virtual autoPtr< phaseSystem::massTransferTable > massTransfer() const
Return the mass transfer matrices.
dimensionedScalar posPart(const dimensionedScalar &ds)
dimensionedScalar neg(const dimensionedScalar &ds)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
HashPtrTable< fvScalarMatrix, word, string::hash > heatTransferTable
Definition: phaseSystem.H:109
dimensionedScalar pos(const dimensionedScalar &ds)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
Volume integrate volField creating a volField.
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
volScalarField & he2
Definition: EEqns.H:3
virtual bool read()
Read base phaseProperties dictionary.
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 dimDensity
PtrList< fvPatch > fvPatchList
container classes for fvPatch
Definition: fvPatchList.H:45
label patchi
virtual ~ThermalPhaseChangePhaseSystem()
Destructor.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Mesh & mesh() const
Return mesh.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
const volScalarField & psi
HashPtrTable< fvScalarMatrix, word, string::hash > massTransferTable
Definition: phaseSystem.H:118
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
dimensionedScalar negPart(const dimensionedScalar &ds)