Loading...
Searching...
No Matches
OneResistanceHeatTransferPhaseSystem.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) 2015-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\*---------------------------------------------------------------------------*/
29#include "fvmSup.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class BasePhaseSystem>
36(
37 const fvMesh& mesh
38)
39:
40 BasePhaseSystem(mesh)
41{
42 this->generatePairsAndSubModels
43 (
44 "heatTransfer",
46 false
47 );
48}
49
50
51// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
52
53template<class BasePhaseSystem>
56{}
57
58
59// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
60
61template<class BasePhaseSystem>
64heatTransfer() const
65{
67 (
69 );
70
71 phaseSystem::heatTransferTable& eqns = eqnsPtr();
72
73 forAll(this->phaseModels_, phasei)
74 {
75 const phaseModel& phase = this->phaseModels_[phasei];
76
77 eqns.set
78 (
79 phase.name(),
80 new fvScalarMatrix(phase.thermo().he(), dimEnergy/dimTime)
81 );
82 }
83
84 // Heat transfer across the interface
86 (
87 heatTransferModelTable,
88 heatTransferModels_,
89 heatTransferModelIter
90 )
91 {
92 const volScalarField K(heatTransferModelIter()->K());
93
94 const phasePair& pair(this->phasePairs_[heatTransferModelIter.key()]);
95
96 forAllConstIter(phasePair, pair, iter)
97 {
98 const phaseModel& phase = iter();
99 const phaseModel& otherPhase = iter.otherPhase();
100
101 const volScalarField& he(phase.thermo().he());
102 volScalarField Cpv(phase.thermo().Cpv());
103
104 *eqns[phase.name()] +=
105 K*(otherPhase.thermo().T() - phase.thermo().T() + he/Cpv)
106 - fvm::Sp(K/Cpv, he);
107 }
108 }
109
110 // Source term due to mass transfer
112 (
114 this->phasePairs_,
115 phasePairIter
116 )
117 {
118 const phasePair& pair(phasePairIter());
119
120 if (pair.ordered())
121 {
122 continue;
123 }
124
125 const phaseModel& phase1 = pair.phase1();
126 const phaseModel& phase2 = pair.phase2();
127
128 const volScalarField& he1(phase1.thermo().he());
129 const volScalarField& he2(phase2.thermo().he());
130
131 const volScalarField K1(phase1.K());
132 const volScalarField K2(phase2.K());
133
134 // Note that the phase heEqn contains a continuity error term, which
135 // implicitly adds a mass transfer term of fvm::Sp(dmdt, he). These
136 // additions do not include this term.
137
138 const volScalarField dmdt(this->dmdt(pair));
139 const volScalarField dmdt21(posPart(dmdt));
140 const volScalarField dmdt12(negPart(dmdt));
141
142 *eqns[phase1.name()] +=
143 dmdt21*he2 - fvm::Sp(dmdt21, he1) + dmdt21*(K2 - K1);
144
145 *eqns[phase2.name()] -=
146 dmdt12*he1 - fvm::Sp(dmdt12, he2) + dmdt12*(K1 - K2);
148
149 return eqnsPtr;
150}
151
152
153template<class BasePhaseSystem>
155{
156 if (BasePhaseSystem::read())
157 {
158 bool readOK = true;
159
160 // Models ...
161
162 return readOK;
163 }
164 else
165 {
166 return false;
167 }
168}
169
170
171// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
#define K1
Definition SHA1.C:143
#define K2
Definition SHA1.C:144
volScalarField & he
Definition YEEqn.H:52
phaseModel & phase1
phaseModel & phase2
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field).
bool set(const Key &key, T *ptr)
Assign a new entry, overwrites existing.
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
HashTable< autoPtr< BlendedInterfacialModel< heatTransferModel > >, phasePairKey, phasePairKey::hash > heatTransferModelTable
heatTransferModelTable heatTransferModels_
Heat transfer models.
OneResistanceHeatTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual bool read()
Read base phaseProperties dictionary.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
virtual tmp< volScalarField > K() const =0
Return the phase kinetic energy.
const phaseModel & otherPhase() const
Return the other phase in this two-phase system.
Definition phaseModel.C:211
const word & name() const
Definition phaseModel.H:166
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
bool ordered() const noexcept
Return the ordered flag.
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition phasePair.H:52
const multiphaseInter::phaseModel & phase1() const
Definition phasePairI.H:23
const multiphaseInter::phaseModel & phase2() const
Definition phasePairI.H:29
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
Definition phaseSystem.H:89
HashPtrTable< fvScalarMatrix > heatTransferTable
Definition phaseSystem.H:79
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phase.H:53
const word & name() const
Definition phase.H:114
const volScalarField & T
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.
const dimensionSet dimEnergy
fvMatrix< scalar > fvScalarMatrix
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar negPart(const dimensionedScalar &ds)
dimensionedScalar posPart(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition stdFoam.H:356
volScalarField & he2
Definition EEqns.H:3