oneResistanceHeatTransfer.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) 2024-2025 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 "fvmSup.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
39 }
40 }
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::fv::oneResistanceHeatTransfer::readCoeffs(const dictionary& dict)
45 {
46  // Read and construct the models
47  models_.clear();
48  modelsTable tModels
49  (
50  generateBlendedInterfacialModels<blendedHeatTransferModel>
51  (
52  fluid_,
53  dict,
55  false
56  )
57  );
58  models_.transfer(tModels);
59 
60  // Set the list of affected energy field names
61  wordHashSet heNameSet;
62  forAllConstIter(modelsTable, models_, modelIter)
63  {
64  forAllConstIter(phaseInterface, modelIter()->interface(), phaseIter)
65  {
66  heNameSet.insert(phaseIter().thermo().he().name());
67  }
68  }
69  heNames_ = heNameSet.toc();
70 }
71 
72 
74 Foam::fv::oneResistanceHeatTransfer::Ks() const
75 {
76  if (Ks_.empty())
77  {
78  forAllConstIter(modelsTable, models_, modelIter)
79  {
80  Ks_.insert
81  (
82  modelIter()->interface(),
83  modelIter()->K()
84  );
85  }
86  }
87 
88  return Ks_;
89 }
90 
91 
92 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
93 
95 (
96  const word& name,
97  const word& modelType,
98  const fvMesh& mesh,
99  const dictionary& dict
100 )
101 :
102  fvModel(name, modelType, mesh, dict),
103  fluid_(mesh().lookupObject<phaseSystem>(phaseSystem::propertiesName)),
104  models_(),
105  heNames_(),
106  energyEquationIndex_(-1),
107  Ks_()
108 {
109  readCoeffs(coeffs(dict));
110 }
111 
112 
113 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 
116 (
117  const phaseModel& phase1,
118  const phaseModel& phase2
119 ) const
120 {
121  auto iter = Ks().find(phaseInterface(phase1, phase2));
122 
123  return
124  iter == Ks().end()
126  : tmp<volScalarField>((*iter)());
127 }
128 
129 
131 {
132  return heNames_;
133 }
134 
135 
137 (
138  const volScalarField& alpha,
139  const volScalarField& rho,
140  const volScalarField& he,
141  fvMatrix<scalar>& eqn
142 ) const
143 {
144  // Ensure that the coefficients are up-to date if this is the first call in
145  // the current phase loop
146  if (energyEquationIndex_ % heNames_.size() == 0) Ks_.clear();
147  energyEquationIndex_ ++;
148 
149  // Add sources to the energy equation of this phase
150  const phaseModel& phase = fluid_.phases()[eqn.psi().group()];
151  forAllConstIter(KsTable, Ks(), KIter)
152  {
153  const phaseInterface interface(fluid_, KIter.key());
154 
155  if (!interface.contains(phase)) continue;
156 
157  const phaseModel& otherPhase = interface.otherPhase(phase);
158 
159  const volScalarField& T = phase.thermo().T();
160  const volScalarField& otherT = otherPhase.thermo().T();
161 
162  const volScalarField& K = KIter();
163 
164  const volScalarField stabilisedK
165  (
166  otherPhase/max(otherPhase, otherPhase.residualAlpha())*K
167  );
168 
169  const volScalarField& Cpv = phase.thermo().Cpv();
170 
171  eqn += stabilisedK*(otherT - T + he/Cpv) - fvm::Sp(stabilisedK/Cpv, he);
172  }
173 }
174 
175 
177 {
178  return true;
179 }
180 
181 
183 {}
184 
185 
187 {}
188 
189 
191 {}
192 
193 
195 {
196  energyEquationIndex_ = 0;
197 
198  Ks_.clear();
199 }
200 
201 
203 {
204  if (fvModel::read(dict))
205  {
206  readCoeffs(coeffs(dict));
207  return true;
208  }
209  else
210  {
211  return false;
212  }
213 }
214 
215 
216 // ************************************************************************* //
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
void clear()
Clear all entries from table.
Definition: HashPtrTable.C:105
void transfer(HashTable< T, Key, Hash > &)
Transfer the contents of the argument table into this table.
Definition: HashTable.C:587
virtual const volScalarField & T() const =0
Temperature [K].
virtual const volScalarField & Cpv() const =0
Heat capacity at constant pressure/volume [J/kg/K].
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Finite volume model abstract base class.
Definition: fvModel.H:60
static const dictionary & coeffs(const word &modelType, const dictionary &)
Return the coefficients sub-dictionary for a given model type.
Definition: fvModelI.H:31
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:200
static const wordHashSet keywords
The keywords read by this class.
Definition: fvModel.H:121
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:57
Model for heat transfer between two phases. One heat transfer coefficient is used to calculate the he...
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
void addSup(const volScalarField &alpha, const volScalarField &rho, const volScalarField &he, fvMatrix< scalar > &eqn) const
Add a source term to the phase-compressible energy equation.
tmp< volScalarField > K(const phaseModel &phase1, const phaseModel &phase2) const
Return the heat transfer coefficient between a pair of phases.
virtual void correct()
Correct the fvModel.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read source dictionary.
oneResistanceHeatTransfer(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
const phaseModel & otherPhase(const phaseModel &phase) const
Return the other phase relative to the given phase.
bool contains(const phaseModel &phase) const
Return true if this phaseInterface contains the given phase.
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
Definition: phaseModel.C:169
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
Class to represent a system of phases.
Definition: phaseSystem.H:74
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
Calculate the matrix for implicit and explicit sources.
K
Definition: pEqn.H:75
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
HashSet wordHashSet
A HashSet with word keys.
Definition: HashSet.H:210
thermo he()
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31