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-2018 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 "fvmSup.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class BasePhaseSystem>
34 (
35  const fvMesh& mesh
36 )
37 :
38  BasePhaseSystem(mesh)
39 {
40  this->generatePairsAndSubModels
41  (
42  "heatTransfer",
43  heatTransferModels_,
44  false
45  );
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
50 
51 template<class BasePhaseSystem>
54 {}
55 
56 
57 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
58 
59 template<class BasePhaseSystem>
62 heatTransfer() const
63 {
64  autoPtr<phaseSystem::heatTransferTable> eqnsPtr
65  (
67  );
68 
69  phaseSystem::heatTransferTable& eqns = eqnsPtr();
70 
71  forAll(this->phaseModels_, phasei)
72  {
73  const phaseModel& phase = this->phaseModels_[phasei];
74 
75  eqns.insert
76  (
77  phase.name(),
78  new fvScalarMatrix(phase.thermo().he(), dimEnergy/dimTime)
79  );
80  }
81 
82  // Heat transfer across the interface
84  (
85  heatTransferModelTable,
86  heatTransferModels_,
87  heatTransferModelIter
88  )
89  {
90  const volScalarField K(heatTransferModelIter()->K());
91 
92  const phasePair& pair(this->phasePairs_[heatTransferModelIter.key()]);
93 
94  forAllConstIter(phasePair, pair, iter)
95  {
96  const phaseModel& phase = iter();
97  const phaseModel& otherPhase = iter.otherPhase();
98 
99  const volScalarField& he(phase.thermo().he());
100  volScalarField Cpv(phase.thermo().Cpv());
101 
102  *eqns[phase.name()] +=
103  K*(otherPhase.thermo().T() - phase.thermo().T() + he/Cpv)
104  - fvm::Sp(K/Cpv, he);
105  }
106  }
107 
108  // Source term due to mass transfer
110  (
112  this->phasePairs_,
113  phasePairIter
114  )
115  {
116  const phasePair& pair(phasePairIter());
117 
118  if (pair.ordered())
119  {
120  continue;
121  }
122 
123  const phaseModel& phase1 = pair.phase1();
124  const phaseModel& phase2 = pair.phase2();
125 
126  const volScalarField& he1(phase1.thermo().he());
127  const volScalarField& he2(phase2.thermo().he());
128 
129  const volScalarField K1(phase1.K());
130  const volScalarField K2(phase2.K());
131 
132  // Note that the phase heEqn contains a continuity error term, which
133  // implicitly adds a mass transfer term of fvm::Sp(dmdt, he). These
134  // additions do not include this term.
135 
136  const volScalarField dmdt(this->dmdt(pair));
137  const volScalarField dmdt21(posPart(dmdt));
138  const volScalarField dmdt12(negPart(dmdt));
139 
140  *eqns[phase1.name()] +=
141  dmdt21*he2 - fvm::Sp(dmdt21, he1) + dmdt21*(K2 - K1);
142 
143  *eqns[phase2.name()] -=
144  dmdt12*he1 - fvm::Sp(dmdt12, he2) + dmdt12*(K1 - K2);
145  }
146 
147  return eqnsPtr;
148 }
149 
150 
151 template<class BasePhaseSystem>
153 {
154  if (BasePhaseSystem::read())
155  {
156  bool readOK = true;
157 
158  // Models ...
159 
160  return readOK;
161  }
162  else
163  {
164  return false;
165  }
166 }
167 
168 
169 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:50
OneResistanceHeatTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
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
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
Definition: phaseSystem.H:87
#define K2
Definition: SHA1.C:168
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)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
volScalarField & he2
Definition: EEqns.H:3
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.