twoResistanceHeatTransfer.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::twoResistanceHeatTransfer::readCoeffs(const dictionary& dict)
45 {
46  // Read and construct the models
47  models_.clear();
48  modelsTable tModels
49  (
50  generateBlendedInterfacialModels<blendedSidedHeatTransferModel>
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::twoResistanceHeatTransfer::Ks() const
75 {
76  if (Ks_.empty())
77  {
78  forAllConstIter(modelsTable, models_, modelIter)
79  {
80  const phaseInterface& interface = modelIter()->interface();
81 
82  Pair<tmp<volScalarField>> interfaceKs
83  (
84  modelIter()->KinThe(interface.phase1()),
85  modelIter()->KinThe(interface.phase2())
86  );
87 
88  Ks_.insert(interface, interfaceKs);
89  }
90  }
91 
92  return Ks_;
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
97 
99 (
100  const word& name,
101  const word& modelType,
102  const fvMesh& mesh,
103  const dictionary& dict
104 )
105 :
106  fvModel(name, modelType, mesh, dict),
107  fluid_(mesh().lookupObject<phaseSystem>(phaseSystem::propertiesName)),
108  models_(),
109  heNames_(),
110  energyEquationIndex_(-1),
111  Ks_()
112 {
113  readCoeffs(coeffs(dict));
114 }
115 
116 
117 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 
120 Foam::fv::twoResistanceHeatTransfer::Ks
121 (
122  const phaseModel& phase1,
123  const phaseModel& phase2
124 ) const
125 {
126  const phaseInterface interface(phase1, phase2);
127 
128  auto iter = Ks().find(interface);
129 
130  if (iter == Ks().end()) return Pair<tmp<volScalarField>>();
131 
132  Pair<tmp<volScalarField>> Ks(iter().first()(), iter().second()());
133 
134  return interface.index(phase2) ? Ks : reverse(Ks);
135 }
136 
137 
139 Foam::fv::twoResistanceHeatTransfer::Ks
140 (
141  const phaseModel& phase1,
142  const phaseModel& phase2,
143  const scalar residualAlpha
144 ) const
145 {
146  const phaseInterface interface(phase1, phase2);
147 
148  auto iter = models_.find(interface);
149 
150  if (iter == models_.end()) return Pair<tmp<volScalarField>>();
151 
152  return
154  (
155  iter()->KinThe(phase1, residualAlpha),
156  iter()->KinThe(phase2, residualAlpha)
157  );
158 }
159 
160 
162 {
163  return heNames_;
164 }
165 
166 
168 (
169  const volScalarField& alpha,
170  const volScalarField& rho,
171  const volScalarField& he,
172  fvMatrix<scalar>& eqn
173 ) const
174 {
175  // Ensure that the coefficients are up-to date if this is the first call in
176  // the current phase loop
177  if (energyEquationIndex_ % heNames_.size() == 0) Ks_.clear();
178  energyEquationIndex_ ++;
179 
180  // Add sources to the energy equation of this phase
181  const phaseModel& phase = fluid_.phases()[eqn.psi().group()];
182  forAllConstIter(KsTable, Ks(), KIter)
183  {
184  const phaseInterface interface(fluid_, KIter.key());
185 
186  if (!interface.contains(phase)) continue;
187 
188  const phaseModel& otherPhase = interface.otherPhase(phase);
189 
190  const volScalarField& T = phase.thermo().T();
191  const volScalarField& otherT = otherPhase.thermo().T();
192 
193  const volScalarField& K = KIter()[interface.index(phase)];
194  const volScalarField& otherK = KIter()[interface.index(otherPhase)];
195 
196  const volScalarField KEff(K*otherK/(K + otherK));
197 
198  const volScalarField& Cpv = phase.thermo().Cpv();
199 
200  eqn += KEff*(otherT - T) + K/Cpv*he - fvm::Sp(K/Cpv, he);
201  }
202 }
203 
204 
206 {
207  return true;
208 }
209 
210 
212 {}
213 
214 
216 {}
217 
218 
220 {}
221 
222 
224 {
225  energyEquationIndex_ = 0;
226 
227  Ks_.clear();
228 }
229 
230 
232 {
233  if (fvModel::read(dict))
234  {
235  readCoeffs(coeffs(dict));
236  return true;
237  }
238  else
239  {
240  return false;
241  }
242 }
243 
244 
245 // ************************************************************************* //
#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
An ordered pair of two objects of type <Type> with first() and second() elements.
Definition: Pair.H:66
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. Two heat transfer coefficients are used to calculate the ...
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.
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.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
twoResistanceHeatTransfer(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
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.
label index(const phaseModel &phase) const
Return the index of the given phase. Generates a FatalError if.
bool contains(const phaseModel &phase) const
Return true if this phaseInterface contains the given phase.
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 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)
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
void reverse(UList< T > &, const label n)
Definition: UListI.H:334
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
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