heatTransferLimitedPhaseChange.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 
27 #include "fvmSup.H"
28 #include "multiphaseEuler.H"
31 
32 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace fv
37 {
40  (
41  fvModel,
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void Foam::fv::heatTransferLimitedPhaseChange::readCoeffs
52 (
53  const dictionary& dict
54 )
55 {
56  saturationModelPtr_.reset
57  (
59  (
60  "saturationTemperature",
61  dict
62  ).ptr()
63  );
64 
65  pressureImplicit_ = dict.lookup<bool>("pressureImplicit");
66 
67  if (pressureImplicit_ && !dmDotdpPtr_.valid())
68  {
69  dmDotdpPtr_.set
70  (
72  (
73  IOobject
74  (
75  name() + ":dmDotdp",
76  mesh().time().name(),
77  mesh(),
80  ),
81  mesh(),
83  )
84  );
85  }
86 
87  if (!pressureImplicit_ && dmDotdpPtr_.valid())
88  {
89  dmDotdpPtr_.clear();
90  }
91 }
92 
93 
94 void Foam::fv::heatTransferLimitedPhaseChange::correctMDot() const
95 {
96  Info<< type() << ": " << name() << endl << incrIndent;
97 
98  const volScalarField::Internal& p = this->p();
99  const volScalarField::Internal& T1 = phase1_.thermo().T();
100  const volScalarField::Internal& T2 = phase2_.thermo().T();
101 
102  // Saturation temperature
103  const volScalarField::Internal Tsat(saturationModelPtr_->Tsat(p));
104 
105  Info<< indent << "min/mean/max Tsat"
106  << " = " << gMin(Tsat.primitiveField())
107  << '/' << gAverage(Tsat.primitiveField())
108  << '/' << gMax(Tsat.primitiveField())
109  << endl;
110 
111  // Latent heat
112  const volScalarField::Internal L(this->L(Tsat));
113 
114  // Heat transfer coefficients
115  const Pair<tmp<volScalarField>> Hs =
116  solver_.heatTransfer.Hs(phase1_, phase2_, scalar(0));
117  const volScalarField::Internal& H1 = Hs.first();
118  const volScalarField::Internal& H2 = Hs.second();
119 
120  // Relaxation factor
121  const scalar f = mesh().solution().fieldRelaxationFactor(mDot_.member());
122 
123  // Recalculate the phase change rate
124  mDot_ = (1 - f)*mDot_ + f*(H1*(T1 - Tsat) + H2*(T2 - Tsat))/L;
125 
126  Info<< indent << "min/mean/max mDot"
127  << " = " << gMin(mDot_.primitiveField())
128  << '/' << gAverage(mDot_.primitiveField())
129  << '/' << gMax(mDot_.primitiveField())
130  << endl;
131 
132  // Recalculate the phase change rate derivative w.r.t. pressure
133  if (pressureImplicit_)
134  {
135  // Saturation temperature derivative w.r.t. pressure
136  const volScalarField::Internal TsatPrime
137  (
138  saturationModelPtr_->TsatPrime(p)
139  );
140 
141  dmDotdpPtr_() = (1 - f)*dmDotdpPtr_() - f*(H1 + H2)*TsatPrime/L;
142  }
143 
144  Info<< decrIndent;
145 }
146 
147 
148 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
149 
151 (
152  const word& name,
153  const word& modelType,
154  const fvMesh& mesh,
155  const dictionary& dict
156 )
157 :
158  phaseChange(name, modelType, mesh, dict, wordList()),
159  solver_(mesh().lookupObject<solvers::multiphaseEuler>(solver::typeName)),
160  fluid_(solver_.fluid),
161  phase1_(fluid_.phases()[phaseNames().first()]),
162  phase2_(fluid_.phases()[phaseNames().second()]),
163  saturationModelPtr_(nullptr),
164  pressureImplicit_(false),
165  pressureEquationIndex_(-1),
166  mDot_
167  (
168  IOobject
169  (
170  name + ":mDot",
171  mesh.time().name(),
172  mesh,
173  IOobject::READ_IF_PRESENT,
174  IOobject::AUTO_WRITE
175  ),
176  mesh,
178  ),
179  dmDotdpPtr_(nullptr)
180 {
181  readCoeffs(coeffs(dict));
182 }
183 
184 
185 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
186 
189 {
190  // Heat transfer coefficients
191  const Pair<tmp<volScalarField>> Hs =
192  solver_.heatTransfer.Hs(phase1_, phase2_);
193  const volScalarField::Internal& H1 = Hs.first();
194  const volScalarField::Internal& H2 = Hs.second();
195 
196  return H2/(H1 + H2);
197 }
198 
199 
202 {
203  return mDot_;
204 }
205 
206 
208 (
209  const volScalarField& alpha,
210  const volScalarField& rho,
211  fvMatrix<scalar>& eqn
212 ) const
213 {
214  // Pressure equation (i.e., continuity, linearised in the pressure)
215  if
216  (
217  (&alpha == &phase1_ || &alpha == &phase2_)
218  && (&rho == &phase1_.rho() || &rho == &phase2_.rho())
219  && &eqn.psi() == &solver_.p_rgh
220  )
221  {
222  // Ensure that the source is up-to date if this is the first call in
223  // the current phase loop
224  if (pressureEquationIndex_ % 2 == 0) correctMDot();
225  pressureEquationIndex_ ++;
226 
227  // Add linearisation if necessary
228  if (pressureImplicit_)
229  {
230  const label s = &alpha == &phase1_ ? -1 : +1;
231 
232  eqn += correction(fvm::Sp(s*dmDotdpPtr_(), eqn.psi()));
233  }
234  }
235 
236  // Let the base class add the actual source
238 }
239 
240 
242 {
243  // Reset the p_rgh equation solution counter
244  pressureEquationIndex_ = 0;
245 
246  // Correct the total phase change rate
247  correctMDot();
248 }
249 
250 
252 {
253  if (phaseChange::read(dict))
254  {
255  readCoeffs(coeffs(dict));
256  return true;
257  }
258  else
259  {
260  return false;
261  }
262 }
263 
264 
265 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
DimensionedField< Type, GeoMesh, PrimitiveField > Internal
Type of the internal field from which this GeometricField is derived.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
An ordered pair of two objects of type <Type> with first() and second() elements.
Definition: Pair.H:66
const Type & second() const
Return second.
Definition: PairI.H:121
const Type & first() const
Return first.
Definition: PairI.H:107
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
const fvSolution & solution() const
Return the fvSolution.
Definition: fvMesh.C:1810
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
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:69
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:57
Model for heat transfer rate limited phase change between two phases.
virtual tmp< DimensionedField< scalar, volMesh > > Lfraction() const
Return the fraction of the latent heat that is transferred into.
virtual tmp< DimensionedField< scalar, volMesh > > mDot() const
Return the mass transfer rate.
virtual bool read(const dictionary &dict)
Read source dictionary.
void addSup(const volScalarField &alpha, const volScalarField &rho, const volScalarField &heOrYi, fvMatrix< scalar > &eqn) const
Use phaseChange's source functions.
Definition: phaseChange.C:551
heatTransferLimitedPhaseChange(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
virtual void addSup(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
Definition: massTransfer.C:225
Base class for phase change models.
Definition: phaseChange.H:61
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: phaseChange.C:624
static autoPtr< saturationTemperatureModel > New(const dictionary &dict)
Select with dictionary.
scalar fieldRelaxationFactor(const word &name) const
Return the relaxation factor for the given field.
Definition: solution.C:271
Abstract base class for run-time selectable region solvers.
Definition: solver.H:56
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.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
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.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:242
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const dimensionSet dimPressure
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:235
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
const dimensionSet dimTime
const dimensionSet dimDensity
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
Type gAverage(const FieldField< Field, Type > &f)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
Type gMin(const FieldField< Field, Type > &f)
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:228
Type gMax(const FieldField< Field, Type > &f)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
labelList f(nPoints)
labelList fv(nPoints)
dictionary dict
volScalarField & p