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-2022 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 "heatTransferModel.H"
28 #include "fvmSup.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class BasePhaseSystem>
35 (
36  const fvMesh& mesh
37 )
38 :
39  HeatTransferPhaseSystem<BasePhaseSystem>(mesh)
40 {
41  this->generateInterfacialModels(heatTransferModels_);
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
46 
47 template<class BasePhaseSystem>
50 {}
51 
52 
53 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
54 
55 template<class BasePhaseSystem>
58 heatTransfer() const
59 {
60  autoPtr<phaseSystem::heatTransferTable> eqnsPtr
61  (
63  );
64 
65  phaseSystem::heatTransferTable& eqns = eqnsPtr();
66 
67  forAll(this->phaseModels_, phasei)
68  {
69  const phaseModel& phase = this->phaseModels_[phasei];
70 
71  eqns.insert
72  (
73  phase.name(),
74  new fvScalarMatrix(phase.thermo().he(), dimEnergy/dimTime)
75  );
76  }
77 
78  // Heat transfer across the interface
80  (
81  heatTransferModelTable,
82  heatTransferModels_,
83  heatTransferModelIter
84  )
85  {
86  const phaseInterface interface(*this, heatTransferModelIter.key());
87 
88  const volScalarField K(heatTransferModelIter()->K());
89 
90  forAllConstIter(phaseInterface, interface, iter)
91  {
92  const phaseModel& phase = iter();
93  const phaseModel& otherPhase = iter.otherPhase();
94 
95  const volScalarField& he = phase.thermo().he();
96  const volScalarField Cpv(phase.thermo().Cpv());
97 
98  *eqns[phase.name()] +=
99  K*(otherPhase.thermo().T() - phase.thermo().T() + he/Cpv)
100  - fvm::Sp(K/Cpv, he);
101  }
102  }
103 
104  return eqnsPtr;
105 }
106 
107 
108 template<class BasePhaseSystem>
110 {
112  {
113  bool readOK = true;
114 
115  // Models ...
116 
117  return readOK;
118  }
119  else
120  {
121  return false;
122  }
123 }
124 
125 
126 // ************************************************************************* //
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 forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
HashPtrTable< fvScalarMatrix > heatTransferTable
Definition: phaseSystem.H:78
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvMesh & mesh
CGAL::Exact_predicates_exact_constructions_kernel K
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const dimensionSet dimTime
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.