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) 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 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class BasePhaseSystem>
36 (
37  const fvMesh& mesh
38 )
39 :
40  HeatTransferPhaseSystem<BasePhaseSystem>(mesh)
41 {
42  this->generatePairsAndSubModels
43  (
44  "heatTransfer",
45  heatTransferModels_,
46  false
47  );
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
52 
53 template<class BasePhaseSystem>
56 {}
57 
58 
59 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
60 
61 template<class BasePhaseSystem>
64 heatTransfer() const
65 {
66  autoPtr<phaseSystem::heatTransferTable> eqnsPtr
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.insert
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  const 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  return eqnsPtr;
111 }
112 
113 
114 template<class BasePhaseSystem>
116 {
118  {
119  bool readOK = true;
120 
121  // Models ...
122 
123  return readOK;
124  }
125  else
126  {
127  return false;
128  }
129 }
130 
131 
132 // ************************************************************************* //
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
HashPtrTable< fvScalarMatrix > heatTransferTable
Definition: phaseSystem.H:80
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
CGAL::Exact_predicates_exact_constructions_kernel K
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const dimensionSet dimTime
dynamicFvMesh & mesh
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
const dimensionSet dimEnergy
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
virtual bool read()
Read base phaseProperties dictionary.
Calculate the matrix for implicit and explicit sources.
virtual bool read()
Read base phaseProperties dictionary.
virtual ~OneResistanceHeatTransferPhaseSystem()
Destructor.