PhaseTransferPhaseSystem.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) 2018 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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#include "phaseTransferModel.H"
30#include "fvmSup.H"
31
32// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
33
34template<class BasePhaseSystem>
37(
38 const phasePairKey& key
39) const
40{
41 if (!rDmdt_.found(key))
42 {
43 return phaseSystem::dmdt(key);
44 }
45
46 const scalar rDmdtSign(Pair<word>::compare(rDmdt_.find(key).key(), key));
47
48 return rDmdtSign**rDmdt_[key];
49}
50
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
54template<class BasePhaseSystem>
56(
57 const fvMesh& mesh
58)
59:
60 BasePhaseSystem(mesh)
61{
62 this->generatePairsAndSubModels
63 (
64 "phaseTransfer",
66 false
67 );
68
70 (
73 phaseTransferModelIter
74 )
75 {
76 this->rDmdt_.set
77 (
78 phaseTransferModelIter.key(),
79 phaseSystem::dmdt(phaseTransferModelIter.key()).ptr()
80 );
81 }
82}
83
84
85// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
86
87template<class BasePhaseSystem>
90{}
91
92
93// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
94
95template<class BasePhaseSystem>
98(
99 const phasePairKey& key
100) const
101{
102 return BasePhaseSystem::dmdt(key) + this->rDmdt(key);
103}
104
105
106template<class BasePhaseSystem>
109{
110 PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
111
112 forAllConstIter(rDmdtTable, rDmdt_, rDmdtIter)
113 {
114 const phasePair& pair = this->phasePairs_[rDmdtIter.key()];
115 const volScalarField& rDmdt = *rDmdtIter();
116
117 this->addField(pair.phase1(), "dmdt", rDmdt, dmdts);
118 this->addField(pair.phase2(), "dmdt", - rDmdt, dmdts);
119 }
120
121 return dmdts;
122}
123
124
125template<class BasePhaseSystem>
128{
129 // Create a mass transfer matrix for each species of each phase
131 (
133 );
134
135 phaseSystem::massTransferTable& eqns = eqnsPtr();
136
137 forAll(this->phaseModels_, phasei)
138 {
139 const phaseModel& phase = this->phaseModels_[phasei];
140
141 const PtrList<volScalarField>& Yi = phase.Y();
142
143 forAll(Yi, i)
144 {
145 eqns.set
146 (
147 Yi[i].name(),
148 new fvScalarMatrix(Yi[i], dimMass/dimTime)
149 );
150 }
151 }
152
153 // Mass transfer across the interface
155 (
157 phaseTransferModels_,
158 phaseTransferModelIter
159 )
160 {
161 const phasePair& pair(this->phasePairs_[phaseTransferModelIter.key()]);
162
163 const phaseModel& phase = pair.phase1();
164 const phaseModel& otherPhase = pair.phase2();
165
166 // Note that the phase YiEqn does not contain a continuity error term,
167 // so these additions represent the entire mass transfer
168
169 const volScalarField dmdt(this->rDmdt(pair));
170 const volScalarField dmdt12(negPart(dmdt));
171 const volScalarField dmdt21(posPart(dmdt));
172
173 const PtrList<volScalarField>& Yi = phase.Y();
174
175 forAll(Yi, i)
176 {
177 const word name
178 (
179 IOobject::groupName(Yi[i].member(), phase.name())
180 );
181
182 const word otherName
183 (
184 IOobject::groupName(Yi[i].member(), otherPhase.name())
185 );
186
187 *eqns[name] +=
188 dmdt21*eqns[otherName]->psi()
189 + fvm::Sp(dmdt12, eqns[name]->psi());
190
191 *eqns[otherName] -=
192 dmdt12*eqns[name]->psi()
193 + fvm::Sp(dmdt21, eqns[otherName]->psi());
194 }
195
196 }
197
198 return eqnsPtr;
199}
200
201
202template<class BasePhaseSystem>
204{
205 BasePhaseSystem::correct();
206
208 (
210 phaseTransferModels_,
211 phaseTransferModelIter
212 )
213 {
214 *rDmdt_[phaseTransferModelIter.key()] =
216 }
217
219 (
221 phaseTransferModels_,
222 phaseTransferModelIter
223 )
224 {
225 *rDmdt_[phaseTransferModelIter.key()] +=
226 phaseTransferModelIter()->dmdt();
227 }
228}
229
230
231template<class BasePhaseSystem>
233{
234 if (BasePhaseSystem::read())
235 {
236 bool readOK = true;
237
238 // Models ...
239
240 return readOK;
241 }
242 else
243 {
244 return false;
245 }
246}
247
248
249// ************************************************************************* //
tmp< volScalarField > dmdt() const
Return the blended mass transfer rate.
bool set(const Key &key, T *ptr)
Assign a new entry, overwriting existing entries.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:69
Class which models non-thermally-coupled mass transfers; i.e., representation changes,...
virtual void correct()
Correct the mass transfer rates.
rDmdtTable rDmdt_
Mass transfer rates.
virtual autoPtr< phaseSystem::massTransferTable > massTransfer() const
Return the mass transfer matrices.
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
virtual ~PhaseTransferPhaseSystem()
Destructor.
phaseTransferModelTable phaseTransferModels_
Mass transfer models.
virtual tmp< volScalarField > rDmdt(const phasePairKey &key) const
Return the representation mass transfer rate.
virtual bool read()
Read base phaseProperties dictionary.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
const word & name() const
Definition: phaseModel.H:144
An ordered or unorder pair of phase names. Typically specified as follows.
Definition: phasePairKey.H:68
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:56
const multiphaseInter::phaseModel & phase1() const
Definition: phasePairI.H:30
const multiphaseInter::phaseModel & phase2() const
Definition: phasePairI.H:36
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:57
const word & name() const
Definition: phase.H:111
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
const volScalarField & psi
dynamicFvMesh & mesh
Calculate the finiteVolume matrix for implicit and explicit sources.
label phasei
Definition: pEqn.H:27
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:45
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
dimensionedScalar negPart(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
const dimensionSet dimDensity
dimensionedScalar posPart(const dimensionedScalar &ds)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:381