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 "mappedFvPatchBaseBase.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().findIndex(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  thermo_.rho()().name(),
109  thermo_.he().name(),
110  VoF_.U.name()
111  };
112 }
113 
114 
116 {
117  if (curTimeIndex_ == mesh().time().timeIndex())
118  {
119  return;
120  }
121 
122  curTimeIndex_ = mesh().time().timeIndex();
123 
124 
125  const scalar deltaT = mesh().time().deltaTValue();
126 
127  const fvPatch& VoFFilmPatch = mesh().boundary()[filmPatchi_];
128 
129  const scalarField& deltaCoeffs = VoFFilmPatch.deltaCoeffs();
130 
131  const labelList& faceCells = VoFFilmPatch.faceCells();
132 
133 
134  // VoF properties
135 
136  const scalarField& alpha = alpha_.boundaryField()[filmPatchi_];
137 
138 
139  // Film properties
140 
141  const mappedFvPatchBaseBase& VoFFilmPatchMap =
142  refCast<const mappedFvPatchBaseBase>(VoFFilmPatch);
143 
144  const solvers::isothermalFilm& film_
145  (
146  VoFFilmPatchMap.nbrMesh().lookupObject<solvers::isothermalFilm>
147  (
148  solver::typeName
149  )
150  );
151 
152  const label filmVoFPatchi = VoFFilmPatchMap.nbrFvPatch().index();
153 
154  const scalarField delta
155  (
156  VoFFilmPatchMap.fromNeighbour
157  (
158  film_.delta.boundaryField()[filmVoFPatchi]
159  )
160  );
161 
162  transferRate_ = Zero;
163 
164  forAll(faceCells, facei)
165  {
166  const label celli = faceCells[facei];
167 
168  if
169  (
170  alpha[facei] > 0
171  && delta[facei] < 2*deltaFactorToFilm_/deltaCoeffs[facei]
172  && alpha[facei] < alphaToFilm_
173  )
174  {
175  transferRate_[celli] = transferRateCoeff_/deltaT;
176  }
177  }
178 }
179 
180 
181 template<class Type, class TransferRateFunc>
183 inline Foam::fv::VoFFilmTransfer::filmVoFTransferRate
184 (
185  TransferRateFunc transferRateFunc,
186  const dimensionSet& dimProp
187 ) const
188 {
189  const fvPatch& VoFFilmPatch = mesh().boundary()[filmPatchi_];
190 
191  const labelList& faceCells = VoFFilmPatch.faceCells();
192 
193 
194  const mappedFvPatchBaseBase& VoFFilmPatchMap =
195  refCast<const mappedFvPatchBaseBase>(VoFFilmPatch);
196 
197  const Foam::fvModels& fvModels =
198  fvModels::New(VoFFilmPatchMap.nbrMesh());
199 
200  const filmVoFTransfer* filmVoFPtr = nullptr;
201 
202  forAll(fvModels, i)
203  {
204  if (isType<filmVoFTransfer>(fvModels[i]))
205  {
206  filmVoFPtr = &refCast<const filmVoFTransfer>(fvModels[i]);
207  }
208  }
209 
210  if (!filmVoFPtr)
211  {
213  << "Cannot find filmVoFTransfer fvModel for the film region "
214  << VoFFilmPatchMap.nbrMesh().name()
215  << exit(FatalError);
216  }
217 
219  (
221  (
222  "Su",
223  mesh(),
224  dimensioned<Type>(dimProp/dimTime, Zero)
225  )
226  );
227 
228  UIndirectList<Type>(tSu.ref(), faceCells) =
229  VoFFilmPatchMap.fromNeighbour
230  (
231  (filmVoFPtr->*transferRateFunc)()
232  );
233 
234  return tSu/mesh().V();
235 }
236 
237 
239 (
240  const volScalarField& alpha,
241  const volScalarField& rho,
242  fvMatrix<scalar>& eqn
243 ) const
244 {
245  if (debug)
246  {
247  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
248  }
249 
250  if (&rho == &thermo_.rho()())
251  {
252  // Explicit transfer into the VoF
253  eqn +=
254  filmVoFTransferRate<scalar>
255  (
257  dimMass
258  );
259 
260  // Potentially implicit transfer out of the VoF
261  if (&rho == &eqn.psi())
262  {
263  eqn -= fvm::Sp(alpha()*transferRate_, eqn.psi());
264  }
265  else if (&alpha == &eqn.psi())
266  {
267  eqn -= fvm::Sp(rho()*transferRate_, eqn.psi());
268  }
269  else
270  {
271  eqn -= alpha()*rho()*transferRate_;
272  }
273  }
274  else
275  {
277  << "Support for field " << rho.name() << " is not implemented"
278  << exit(FatalError);
279  }
280 }
281 
282 
284 (
285  const volScalarField& alpha,
286  const volScalarField& rho,
287  const volScalarField& he,
288  fvMatrix<scalar>& eqn
289 ) const
290 {
291  if (debug)
292  {
293  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
294  }
295 
296  if (&he == &thermo_.he())
297  {
298  // Explicit transfer into the VoF
299  eqn +=
300  filmVoFTransferRate<scalar>
301  (
303  dimEnergy
304  );
305 
306  // Potentially implicit transfer out of the VoF
307  if (&he == &eqn.psi())
308  {
309  eqn -= fvm::Sp(alpha()*rho()*transferRate_, eqn.psi());
310  }
311  else
312  {
313  eqn -= alpha()*rho()*he*transferRate_;
314  }
315  }
316  else
317  {
319  << "Support for field " << he.name() << " is not implemented"
320  << exit(FatalError);
321  }
322 }
323 
324 
326 (
327  const volScalarField& rho,
328  const volVectorField& U,
329  fvMatrix<vector>& eqn
330 ) const
331 {
332  if (debug)
333  {
334  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
335  }
336 
337  if (&U == &VoF_.U)
338  {
339  // Explicit transfer into the VoF
340  eqn +=
341  filmVoFTransferRate<vector>
342  (
345  );
346 
347  // Potentially implicit transfer out of the VoF
348  if (&U == &eqn.psi())
349  {
350  eqn -= fvm::Sp(alpha_()*rho()*transferRate_, eqn.psi());
351  }
352  else
353  {
354  eqn -= alpha_()*rho()*U*transferRate_;
355  }
356  }
357  else
358  {
360  << "Support for field " << U.name() << " is not implemented"
361  << exit(FatalError);
362  }
363 }
364 
365 
366 template<class Type, class FieldType>
367 inline Foam::tmp<Foam::Field<Type>> Foam::fv::VoFFilmTransfer::TransferRate
368 (
369  const FieldType& f
370 ) const
371 {
372  const fvPatch& VoFFilmPatch = mesh().boundary()[filmPatchi_];
373 
374  const labelList& faceCells = VoFFilmPatch.faceCells();
375 
376 
377  return tmp<Field<Type>>
378  (
379  new Field<Type>
380  (
382  (
383  alpha_()*transferRate_*mesh().V()*f,
384  faceCells
385  )
386  )
387  );
388 }
389 
390 
393 {
394  return TransferRate<scalar>(thermo_.rho()());
395 }
396 
397 
400 {
401  return TransferRate<scalar>(thermo_.rho()()*thermo_.he()());
402 }
403 
404 
407 {
408  return TransferRate<vector>(thermo_.rho()()*VoF_.U());
409 }
410 
411 
413 {
414  transferRate_.setSize(mesh().nCells());
415 }
416 
417 
419 {
420  transferRate_.setSize(mesh().nCells());
421 }
422 
423 
425 {
426  transferRate_.setSize(mesh().nCells());
427 }
428 
429 
431 {
432  return true;
433 }
434 
435 
436 // ************************************************************************* //
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.
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:162
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:125
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:99
const word & name() const
Return reference to name.
Definition: fvMesh.H:432
Finite volume model abstract base class.
Definition: fvModel.H:59
Finite volume models.
Definition: fvModels.H:65
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:156
virtual const scalarField & deltaCoeffs() const
Return the face - cell distance coefficient.
Definition: fvPatch.C:176
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:99
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 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.
virtual void addSup(const volScalarField &alpha, const volScalarField &rho, fvMatrix< scalar > &eqn) const
Add implicit contribution to phase continuity equation.
tmp< scalarField > rhoTransferRate() const
Return the mass transfer rate.
Film<->VoF transfer model.
tmp< scalarField > heTransferRate() const
Return the energy transfer rate.
tmp< vectorField > UTransferRate() const
Return the momentum transfer rate.
tmp< scalarField > rhoTransferRate() const
Return the mass transfer rate.
Base class for fv patches that provide mapping between two fv patches.
const fvPatch & nbrFvPatch() const
Get the patch to map from.
const fvMesh & nbrMesh() const
Get the mesh for the region to map from.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
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.
Abstract base class for run-time selectable region solvers.
Definition: solver.H:56
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Calculate the matrix for implicit and explicit sources.
U
Definition: pEqn.H:72
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:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
messageStream Info
const dimensionSet dimTime
typename VolField< Type >::Internal VolInternalField
Definition: volFieldsFwd.H:61
const dimensionSet dimMass
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
const dimensionSet dimVelocity
error FatalError
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
thermo he()
labelList f(nPoints)
labelList fv(nPoints)
dictionary dict