VoFFilmTransfer.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) 2023 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 
26 #include "VoFFilmTransfer.H"
27 #include "filmVoFTransfer.H"
28 #include "mappedPatchBase.H"
29 #include "fvmSup.H"
31 
32 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace fv
37  {
39 
41  (
42  fvModel,
45  );
46  }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const word& sourceName,
55  const word& modelType,
56  const fvMesh& mesh,
57  const dictionary& dict
58 )
59 :
60  fvModel(sourceName, modelType, mesh, dict),
61  VoF_(mesh.lookupObject<solvers::compressibleVoF>(solver::typeName)),
62  filmPatchName_(dict.lookup("filmPatch")),
63  filmPatchi_(mesh.boundaryMesh().findPatchID(filmPatchName_)),
64  phaseName_(dict.lookup("phase")),
65  thermo_
66  (
67  phaseName_ == VoF_.mixture.phase1Name()
68  ? VoF_.mixture.thermo1()
69  : VoF_.mixture.thermo2()
70  ),
71  alpha_
72  (
73  phaseName_ == VoF_.mixture.phase1Name()
74  ? VoF_.mixture.alpha1()
75  : VoF_.mixture.alpha2()
76  ),
77  curTimeIndex_(-1),
78  deltaFactorToFilm_
79  (
80  dict.lookupOrDefault<scalar>("deltaFactorToFilm", 0.5)
81  ),
82  alphaToFilm_
83  (
84  dict.lookupOrDefault<scalar>("alphaToFilm", 0.1)
85  ),
86  transferRateCoeff_
87  (
88  dict.lookupOrDefault<scalar>("transferRateCoeff", 0.1)
89  ),
90  transferRate_
91  (
92  volScalarField::Internal::New
93  (
94  "transferRate",
95  mesh,
97  )
98  )
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
105 {
106  return wordList
107  {
108  alpha_.name(),
109  thermo_.rho()().name(),
110  thermo_.he().name(),
111  VoF_.U.name()
112  };
113 }
114 
115 
117 {
118  if (curTimeIndex_ == mesh().time().timeIndex())
119  {
120  return;
121  }
122 
123  curTimeIndex_ = mesh().time().timeIndex();
124 
125 
126  const scalar deltaT = mesh().time().deltaTValue();
127 
128  const polyPatch& VoFFilmPatch = mesh().boundaryMesh()[filmPatchi_];
129 
130 
131  // VoF properties
132 
133  const scalarField& alpha = alpha_.boundaryField()[filmPatchi_];
134 
135  const scalarField& deltaCoeffs =
136  mesh().boundary()[filmPatchi_].deltaCoeffs();
137 
138  const labelList& faceCells = mesh().boundary()[filmPatchi_].faceCells();
139 
140 
141  // Film properties
142 
143  const mappedPatchBase& VoFFilmPatchMap = refCast<const mappedPatchBase>
144  (
145  VoFFilmPatch
146  );
147 
148  const solvers::isothermalFilm& film_
149  (
150  VoFFilmPatchMap.nbrMesh().lookupObject<solvers::isothermalFilm>
151  (
152  solver::typeName
153  )
154  );
155 
156  const label filmVoFPatchi = VoFFilmPatchMap.nbrPolyPatch().index();
157 
158  const scalarField delta
159  (
160  VoFFilmPatchMap.fromNeighbour
161  (
162  film_.delta.boundaryField()[filmVoFPatchi]
163  )
164  );
165 
166  transferRate_ = Zero;
167 
168  forAll(faceCells, facei)
169  {
170  const label celli = faceCells[facei];
171 
172  if
173  (
174  alpha[facei] > 0
175  && delta[facei] < 2*deltaFactorToFilm_/deltaCoeffs[facei]
176  && alpha[facei] < alphaToFilm_
177  )
178  {
179  transferRate_[celli] = transferRateCoeff_/deltaT;
180  }
181  }
182 }
183 
184 
185 template<class Type, class TransferRateFunc>
187 inline Foam::fv::VoFFilmTransfer::filmVoFTransferRate
188 (
189  TransferRateFunc transferRateFunc,
190  const dimensionSet& dimProp
191 ) const
192 {
193  const mappedPatchBase& VoFFilmPatchMap = refCast<const mappedPatchBase>
194  (
195  mesh().boundaryMesh()[filmPatchi_]
196  );
197 
198  const Foam::fvModels& fvModels
199  (
201  (
202  refCast<const fvMesh>(VoFFilmPatchMap.nbrMesh())
203  )
204  );
205 
206  const filmVoFTransfer* filmVoFPtr = nullptr;
207 
208  forAll(fvModels, i)
209  {
210  if (isType<filmVoFTransfer>(fvModels[i]))
211  {
212  filmVoFPtr = &refCast<const filmVoFTransfer>(fvModels[i]);
213  }
214  }
215 
216  if (!filmVoFPtr)
217  {
219  << "Cannot find filmVoFTransfer fvModel for the film region "
220  << VoFFilmPatchMap.nbrMesh().name()
221  << exit(FatalError);
222  }
223 
225  (
227  (
228  "Su",
229  mesh(),
230  dimensioned<Type>(dimProp/dimTime, Zero)
231  )
232  );
233 
234  UIndirectList<Type>(tSu.ref(), mesh().boundary()[filmPatchi_].faceCells()) =
235  VoFFilmPatchMap.fromNeighbour
236  (
237  (filmVoFPtr->*transferRateFunc)()
238  );
239 
240  return tSu/mesh().V();
241 }
242 
243 
245 (
246  fvMatrix<scalar>& eqn,
247  const word& fieldName
248 ) const
249 {
250  if (debug)
251  {
252  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
253  }
254 
255  if (fieldName == alpha_.name())
256  {
257  eqn +=
258  filmVoFTransferRate<scalar>
259  (
261  dimVolume
262  )
263  - fvm::Sp(transferRate_, eqn.psi());
264  }
265  else
266  {
268  << "Support for field " << fieldName << " is not implemented"
269  << exit(FatalError);
270  }
271 }
272 
273 
275 (
276  const volScalarField& alpha,
277  fvMatrix<scalar>& eqn,
278  const word& fieldName
279 ) const
280 {
281  if (debug)
282  {
283  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
284  }
285 
286  if (fieldName == thermo_.rho()().name())
287  {
288  eqn +=
289  filmVoFTransferRate<scalar>
290  (
292  dimMass
293  )
294  - fvm::Sp(alpha()*transferRate_, eqn.psi());
295  }
296  else
297  {
299  << "Support for field " << fieldName << " is not implemented"
300  << exit(FatalError);
301  }
302 }
303 
304 
306 (
307  const volScalarField& alpha,
308  const volScalarField& rho,
309  fvMatrix<scalar>& eqn,
310  const word& fieldName
311 ) const
312 {
313  if (debug)
314  {
315  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
316  }
317 
318  if (fieldName == thermo_.he().name())
319  {
320  eqn +=
321  filmVoFTransferRate<scalar>
322  (
324  dimEnergy
325  )
326  - fvm::Sp(alpha()*rho()*transferRate_, eqn.psi());
327  }
328  else
329  {
331  << "Support for field " << fieldName << " is not implemented"
332  << exit(FatalError);
333  }
334 }
335 
336 
338 (
339  const volScalarField& rho,
340  fvMatrix<vector>& eqn,
341  const word& fieldName
342 ) const
343 {
344  if (debug)
345  {
346  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
347  }
348 
349  eqn +=
350  filmVoFTransferRate<vector>
351  (
354  )
355  - fvm::Sp(alpha_()*thermo_.rho()()*transferRate_, eqn.psi());
356 }
357 
358 
359 template<class Type, class FieldType>
360 inline Foam::tmp<Foam::Field<Type>> Foam::fv::VoFFilmTransfer::TransferRate
361 (
362  const FieldType& f
363 ) const
364 {
365  const labelList& faceCells = mesh().boundary()[filmPatchi_].faceCells();
366 
367  return tmp<Field<Type>>
368  (
369  new Field<Type>
370  (
372  (
373  alpha_()*transferRate_*mesh().V()*f,
374  faceCells
375  )
376  )
377  );
378 }
379 
380 
383 {
384  return TransferRate<scalar>(thermo_.rho()());
385 }
386 
387 
390 {
391  return TransferRate<scalar>(thermo_.rho()()*thermo_.he()());
392 }
393 
394 
397 {
398  return TransferRate<vector>(thermo_.rho()()*VoF_.U());
399 }
400 
401 
403 {
404  transferRate_.setSize(mesh().nCells());
405 }
406 
407 
409 {
410  transferRate_.setSize(mesh().nCells());
411 }
412 
413 
415 {
416  transferRate_.setSize(mesh().nCells());
417 }
418 
419 
421 {
422  return true;
423 }
424 
425 
426 // ************************************************************************* //
scalar delta
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const word & name() const
Return name.
Definition: IOobject.H:310
A List with indirect addressing.
Definition: UIndirectList.H:60
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
Dimension set for the base types.
Definition: dimensionSet.H:122
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:101
Finite volume model abstract base class.
Definition: fvModel.H:59
Finite volume models.
Definition: fvModels.H:65
Film<->VoF transfer model.
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
tmp< scalarField > heTransferRate() const
Return the energy transfer rate.
VoFFilmTransfer(const word &sourceName, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
virtual void correct()
Solve the film and update the sources.
virtual void addSup(fvMatrix< scalar > &eqn, const word &fieldName) const
Add implicit contribution to phase-fraction equation.
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 void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
tmp< vectorField > UTransferRate() const
Return the momentum transfer rate.
tmp< scalarField > rhoTransferRate() const
Return the mass transfer rate.
Film<->VoF transfer model.
tmp< scalarField > heTransferRate() const
Return the energy transfer rate.
tmp< scalarField > transferRate() const
Return the volume transfer rate.
tmp< vectorField > UTransferRate() const
Return the momentum transfer rate.
tmp< scalarField > rhoTransferRate() const
Return the mass transfer rate.
Engine which provides mapping between two patches.
const polyPatch & nbrPolyPatch() const
Get the patch to map from.
const polyMesh & nbrMesh() const
Get the mesh for the region to map from.
tmp< Field< Type > > fromNeighbour(const Field< Type > &nbrFld) const
Map/interpolate the neighbour patch field to this patch.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
label index() const
Return the index of this patch in the boundaryMesh.
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
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:270
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Abstract base class for run-time selectable region solvers.
Definition: solver.H:55
Solver module for flow of compressible isothermal liquid films.
const volScalarField & delta
Film thickness.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Calculate the matrix for implicit and explicit sources.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
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 dimEnergy
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const dimensionSet dimless
messageStream Info
const dimensionSet dimTime
const dimensionSet dimVolume
typename VolField< Type >::Internal VolInternalField
Definition: volFieldsFwd.H:58
const dimensionSet dimMass
const dimensionSet dimVelocity
error FatalError
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
label timeIndex
Definition: getTimeIndex.H:4
faceListList boundary(nPatches)
labelList f(nPoints)
labelList fv(nPoints)
dictionary dict