PhaseTransferPhaseSystem.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) 2018-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 #include "phaseTransferModel.H"
28 #include "fvmSup.H"
29 
30 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
31 
32 template<class BasePhaseSystem>
35 (
36  const phasePairKey& key
37 ) const
38 {
39  if (!rDmdt_.found(key))
40  {
41  return phaseSystem::dmdt(key);
42  }
43 
44  const scalar rDmdtSign(Pair<word>::compare(rDmdt_.find(key).key(), key));
45 
46  return rDmdtSign**rDmdt_[key];
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 template<class BasePhaseSystem>
54 (
55  const fvMesh& mesh
56 )
57 :
58  BasePhaseSystem(mesh)
59 {
60  this->generatePairsAndSubModels
61  (
62  "phaseTransfer",
63  phaseTransferModels_,
64  false
65  );
66 
68  (
69  phaseTransferModelTable,
70  phaseTransferModels_,
71  phaseTransferModelIter
72  )
73  {
74  this->rDmdt_.insert
75  (
76  phaseTransferModelIter.key(),
77  phaseSystem::dmdt(phaseTransferModelIter.key()).ptr()
78  );
79  }
80 }
81 
82 
83 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
84 
85 template<class BasePhaseSystem>
88 {}
89 
90 
91 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
92 
93 template<class BasePhaseSystem>
96 (
97  const phasePairKey& key
98 ) const
99 {
100  return BasePhaseSystem::dmdt(key) + this->rDmdt(key);
101 }
102 
103 
104 template<class BasePhaseSystem>
107 {
108  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
109 
110  forAllConstIter(rDmdtTable, rDmdt_, rDmdtIter)
111  {
112  const phasePair& pair = this->phasePairs_[rDmdtIter.key()];
113  const volScalarField& rDmdt = *rDmdtIter();
114 
115  this->addField(pair.phase1(), "dmdt", rDmdt, dmdts);
116  this->addField(pair.phase2(), "dmdt", - rDmdt, dmdts);
117  }
118 
119  return dmdts;
120 }
121 
122 
123 template<class BasePhaseSystem>
126 {
127  // Create a mass transfer matrix for each species of each phase
128  autoPtr<phaseSystem::massTransferTable> eqnsPtr
129  (
131  );
132 
133  phaseSystem::massTransferTable& eqns = eqnsPtr();
134 
135  forAll(this->phaseModels_, phasei)
136  {
137  const phaseModel& phase = this->phaseModels_[phasei];
138 
139  const PtrList<volScalarField>& Yi = phase.Y();
140 
141  forAll(Yi, i)
142  {
143  eqns.insert
144  (
145  Yi[i].name(),
146  new fvScalarMatrix(Yi[i], dimMass/dimTime)
147  );
148  }
149  }
150 
151  // Mass transfer across the interface
153  (
154  phaseTransferModelTable,
155  phaseTransferModels_,
156  phaseTransferModelIter
157  )
158  {
159  const phasePair& pair(this->phasePairs_[phaseTransferModelIter.key()]);
160 
161  const phaseModel& phase = pair.phase1();
162  const phaseModel& otherPhase = pair.phase2();
163 
164  // Note that the phase YiEqn does not contain a continuity error term,
165  // so these additions represent the entire mass transfer
166 
167  const volScalarField dmdt(this->rDmdt(pair));
168  const volScalarField dmdt12(negPart(dmdt));
169  const volScalarField dmdt21(posPart(dmdt));
170 
171  const PtrList<volScalarField>& Yi = phase.Y();
172 
173  forAll(Yi, i)
174  {
175  const word name
176  (
177  IOobject::groupName(Yi[i].member(), phase.name())
178  );
179 
180  const word otherName
181  (
182  IOobject::groupName(Yi[i].member(), otherPhase.name())
183  );
184 
185  *eqns[name] +=
186  dmdt21*eqns[otherName]->psi()
187  + fvm::Sp(dmdt12, eqns[name]->psi());
188 
189  *eqns[otherName] -=
190  dmdt12*eqns[name]->psi()
191  + fvm::Sp(dmdt21, eqns[otherName]->psi());
192  }
193 
194  }
195 
196  return eqnsPtr;
197 }
198 
199 
200 template<class BasePhaseSystem>
202 {
204 
206  (
207  phaseTransferModelTable,
208  phaseTransferModels_,
209  phaseTransferModelIter
210  )
211  {
212  *rDmdt_[phaseTransferModelIter.key()] =
214  }
215 
217  (
218  phaseTransferModelTable,
219  phaseTransferModels_,
220  phaseTransferModelIter
221  )
222  {
223  *rDmdt_[phaseTransferModelIter.key()] +=
224  phaseTransferModelIter()->dmdt();
225  }
226 }
227 
228 
229 template<class BasePhaseSystem>
231 {
232  if (BasePhaseSystem::read())
233  {
234  bool readOK = true;
235 
236  // Models ...
237 
238  return readOK;
239  }
240  else
241  {
242  return false;
243  }
244 }
245 
246 
247 // ************************************************************************* //
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
PhaseTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
label phasei
Definition: pEqn.H:27
HashPtrTable< fvScalarMatrix > massTransferTable
Definition: phaseSystem.H:79
static int compare(const Pair< word > &a, const Pair< word > &b)
Compare Pairs.
Definition: Pair.H:142
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the mass transfer rate for a pair.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
dimensionedScalar posPart(const dimensionedScalar &ds)
virtual tmp< volScalarField > rDmdt(const phasePairKey &key) const
Return the representation mass transfer rate.
virtual void correct()
Correct the mass transfer rates.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual bool read()
Read base phaseProperties dictionary.
static word groupName(Name name, const word &group)
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-7/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimDensity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const volScalarField & psi
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: PtrList.H:53
virtual autoPtr< phaseSystem::massTransferTable > massTransfer() const
Return the mass transfer matrices.
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the mass transfer rate for a pair.
Calculate the matrix for implicit and explicit sources.
dimensionedScalar negPart(const dimensionedScalar &ds)
virtual ~PhaseTransferPhaseSystem()
Destructor.