OneResistanceHeatTransferPhaseSystem.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) 2015-2020 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 "heatTransferModel.H"
29 #include "fvmSup.H"
30 #include "rhoReactionThermo.H"
31 
32 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
33 
34 template<class BasePhaseSystem>
36 (
37  const phaseSystem::dmdtfTable& dmdtfs,
38  phaseSystem::heatTransferTable& eqns
39 ) const
40 {
41  forAllConstIter(phaseSystem::dmdtfTable, dmdtfs, dmdtfIter)
42  {
43  const phasePairKey& key = dmdtfIter.key();
44  const phasePair& pair(this->phasePairs_[key]);
45  const volScalarField dmdtf(Pair<word>::compare(pair, key)**dmdtfIter());
46  const volScalarField dmdtf21(posPart(dmdtf));
47  const volScalarField dmdtf12(negPart(dmdtf));
48 
49  const phaseModel& phase1 = pair.phase1();
50  const phaseModel& phase2 = pair.phase2();
51  const rhoThermo& thermo1 = phase1.thermo();
52  const rhoThermo& thermo2 = phase2.thermo();
53  const volScalarField& he1(thermo1.he());
54  const volScalarField& he2(thermo2.he());
55  const volScalarField K1(phase1.K());
56  const volScalarField K2(phase2.K());
57 
58  // Transfer of energy from bulk to bulk
59  *eqns[phase1.name()] += dmdtf21*he2 + fvm::Sp(dmdtf12, he1);
60  *eqns[phase2.name()] -= dmdtf12*he1 + fvm::Sp(dmdtf21, he2);
61 
62  // Transfer of kinetic energy
63  *eqns[phase1.name()] += dmdtf21*K2 + dmdtf12*K1;
64  *eqns[phase2.name()] -= dmdtf12*K1 + dmdtf21*K2;
65  }
66 }
67 
68 
69 template<class BasePhaseSystem>
71 (
72  const phaseSystem::dmidtfTable& dmidtfs,
73  phaseSystem::heatTransferTable& eqns
74 ) const
75 {
76  forAllConstIter(phaseSystem::dmidtfTable, dmidtfs, dmidtfIter)
77  {
78  const phasePairKey& key = dmidtfIter.key();
79  const phasePair& pair(this->phasePairs_[key]);
80 
81  const phaseModel& phase1 = pair.phase1();
82  const phaseModel& phase2 = pair.phase2();
83  const rhoThermo& thermo1 = phase1.thermo();
84  const rhoThermo& thermo2 = phase2.thermo();
85  const volScalarField& he1(thermo1.he());
86  const volScalarField& he2(thermo2.he());
87  const volScalarField K1(phase1.K());
88  const volScalarField K2(phase2.K());
89 
90  forAllConstIter(HashPtrTable<volScalarField>, *dmidtfIter(), dmidtfJter)
91  {
92  const word& member = dmidtfJter.key();
93 
94  const volScalarField dmidtf
95  (
96  Pair<word>::compare(pair, key)**dmidtfJter()
97  );
98  const volScalarField dmidtf21(posPart(dmidtf));
99  const volScalarField dmidtf12(negPart(dmidtf));
100 
101  // Create the energies for the transferring specie
102  volScalarField hei1(he1);
103  if (isA<rhoReactionThermo>(thermo1))
104  {
105  const basicSpecieMixture& composition1 =
106  refCast<const rhoReactionThermo>(thermo1).composition();
107  hei1 =
108  composition1.HE
109  (
110  composition1.species()[member],
111  thermo1.p(),
112  thermo1.T()
113  );
114  }
115  volScalarField hei2(he2);
116  if (isA<rhoReactionThermo>(thermo2))
117  {
118  const basicSpecieMixture& composition2 =
119  refCast<const rhoReactionThermo>(thermo2).composition();
120  hei2 =
121  composition2.HE
122  (
123  composition2.species()[member],
124  thermo2.p(),
125  thermo2.T()
126  );
127  }
128 
129  // Transfer of energy from bulk to bulk
130  *eqns[phase1.name()] += dmidtf21*hei2 + fvm::Sp(dmidtf12, he1);
131  *eqns[phase2.name()] -= dmidtf12*hei1 + fvm::Sp(dmidtf21, he2);
132 
133  // Transfer of kinetic energy
134  *eqns[phase1.name()] += dmidtf21*K2 + dmidtf12*K1;
135  *eqns[phase2.name()] -= dmidtf12*K1 + dmidtf21*K2;
136  }
137  }
138 }
139 
140 
141 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
142 
143 template<class BasePhaseSystem>
146 (
147  const fvMesh& mesh
148 )
149 :
150  BasePhaseSystem(mesh)
151 {
152  this->generatePairsAndSubModels
153  (
154  "heatTransfer",
155  heatTransferModels_,
156  false
157  );
158 }
159 
160 
161 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
162 
163 template<class BasePhaseSystem>
166 {}
167 
168 
169 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
170 
171 template<class BasePhaseSystem>
174 heatTransfer() const
175 {
176  autoPtr<phaseSystem::heatTransferTable> eqnsPtr
177  (
179  );
180 
181  phaseSystem::heatTransferTable& eqns = eqnsPtr();
182 
183  forAll(this->phaseModels_, phasei)
184  {
185  const phaseModel& phase = this->phaseModels_[phasei];
186 
187  eqns.insert
188  (
189  phase.name(),
190  new fvScalarMatrix(phase.thermo().he(), dimEnergy/dimTime)
191  );
192  }
193 
194  // Heat transfer across the interface
196  (
197  heatTransferModelTable,
198  heatTransferModels_,
199  heatTransferModelIter
200  )
201  {
202  const volScalarField K(heatTransferModelIter()->K());
203 
204  const phasePair& pair(this->phasePairs_[heatTransferModelIter.key()]);
205 
206  forAllConstIter(phasePair, pair, iter)
207  {
208  const phaseModel& phase = iter();
209  const phaseModel& otherPhase = iter.otherPhase();
210 
211  const volScalarField& he(phase.thermo().he());
212  volScalarField Cpv(phase.thermo().Cpv());
213 
214  *eqns[phase.name()] +=
215  K*(otherPhase.thermo().T() - phase.thermo().T() + he/Cpv)
216  - fvm::Sp(K/Cpv, he);
217  }
218  }
219 
220  return eqnsPtr;
221 }
222 
223 
224 template<class BasePhaseSystem>
226 {
227  if (BasePhaseSystem::read())
228  {
229  bool readOK = true;
230 
231  // Models ...
232 
233  return readOK;
234  }
235  else
236  {
237  return false;
238  }
239 }
240 
241 
242 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:50
OneResistanceHeatTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
void addDmidtHef(const phaseSystem::dmidtfTable &dmidtfs, phaseSystem::heatTransferTable &eqns) const
Add energy transfer terms which result from specie mass transfers.
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
label phasei
Definition: pEqn.H:27
#define K1
Definition: SHA1.C:167
basicSpecieMixture & composition
#define K2
Definition: SHA1.C:168
static int compare(const Pair< word > &a, const Pair< word > &b)
Compare Pairs.
Definition: Pair.H:142
HashPtrTable< fvScalarMatrix > heatTransferTable
Definition: phaseSystem.H:77
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
CGAL::Exact_predicates_exact_constructions_kernel K
dimensionedScalar posPart(const dimensionedScalar &ds)
void addDmdtHefs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::heatTransferTable &eqns) const
Add energy transfer terms which result from bulk mass transfers.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dynamicFvMesh & mesh
mixture thermo2().correctRho(psi2 *(p_rgh - p_rgh_0))
mixture thermo1().correctRho(psi1 *(p_rgh - p_rgh_0))
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
const dimensionSet dimEnergy
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
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
Calculate the matrix for implicit and explicit sources.
virtual bool read()
Read base phaseProperties dictionary.
dimensionedScalar negPart(const dimensionedScalar &ds)
virtual ~OneResistanceHeatTransferPhaseSystem()
Destructor.
zeroField Sp
Definition: alphaSuSp.H:2