filmVoFTransfer.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 "filmVoFTransfer.H"
27 #include "VoFFilmTransfer.H"
28 #include "mappedPatchBase.H"
29 #include "compressibleVoF.H"
30 #include "fvmSup.H"
32 
33 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  namespace fv
38  {
40 
42  (
43  fvModel,
46  );
47  }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const word& sourceName,
56  const word& modelType,
57  const fvMesh& mesh,
58  const dictionary& dict
59 )
60 :
61  fvModel(sourceName, modelType, mesh, dict),
62  film_(mesh.lookupObject<solvers::isothermalFilm>(solver::typeName)),
63  curTimeIndex_(-1),
64  deltaFactorToVoF_
65  (
66  dict.lookupOrDefault<scalar>("deltaFactorToVoF", 1.0)
67  ),
68  alphaToVoF_
69  (
70  dict.lookupOrDefault<scalar>("alphaToVoF", 0.5)
71  ),
72  transferRateCoeff_
73  (
74  dict.lookupOrDefault<scalar>("transferRateCoeff", 0.1)
75  ),
76  transferRate_
77  (
78  volScalarField::Internal::New
79  (
80  "transferRate",
81  mesh,
83  )
84  )
85 {}
86 
87 
88 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89 
91 {
92  return wordList
93  {
94  film_.alpha.name(),
95  film_.thermo.he().name(),
96  film_.U.name()
97  };
98 }
99 
100 
102 {
103  if (curTimeIndex_ == mesh().time().timeIndex())
104  {
105  return;
106  }
107 
108  curTimeIndex_ = mesh().time().timeIndex();
109 
110  const scalar deltaT = mesh().time().deltaTValue();
111 
112 
113  // Film properties
114 
115  const labelList& faceCells = film_.surfacePatch().faceCells();
116  const scalarField& delta = film_.delta();
117 
118 
119  // VoF properties
120 
121  const solvers::compressibleVoF& VoF_
122  (
123  film_.surfacePatchMap().nbrMesh().lookupObject<solvers::compressibleVoF>
124  (
125  solver::typeName
126  )
127  );
128 
129  const label patchiVoF = film_.surfacePatchMap().nbrPolyPatch().index();
130 
131  const VoFFilmTransfer& VoFFilm(this->VoFFilm(VoF_.fvModels()));
132 
133  const scalarField alphaVoF
134  (
135  film_.surfacePatchMap().fromNeighbour
136  (
137  VoFFilm.alpha().boundaryField()[patchiVoF]
138  )
139  );
140 
141  const scalarField deltaCoeffsVoF
142  (
143  film_.surfacePatchMap().fromNeighbour
144  (
145  VoF_.mesh.boundary()[patchiVoF].deltaCoeffs()
146  )
147  );
148 
149  // Reset the transfer rate
150  transferRate_ = Zero;
151 
152  forAll(faceCells, facei)
153  {
154  const label celli = faceCells[facei];
155 
156  if
157  (
158  delta[celli] > 2*deltaFactorToVoF_/deltaCoeffsVoF[facei]
159  || alphaVoF[facei] > alphaToVoF_
160  )
161  {
162  transferRate_[celli] = transferRateCoeff_/deltaT;
163  }
164  }
165 }
166 
167 
168 const Foam::fv::VoFFilmTransfer& Foam::fv::filmVoFTransfer::VoFFilm
169 (
170  const Foam::fvModels& fvModels
171 ) const
172 {
173  const VoFFilmTransfer* VoFFilmPtr = nullptr;
174 
175  forAll(fvModels, i)
176  {
177  if (isType<VoFFilmTransfer>(fvModels[i]))
178  {
179  const VoFFilmTransfer& VoFFilm
180  (
181  refCast<const VoFFilmTransfer>(fvModels[i])
182  );
183 
184  if
185  (
186  VoFFilm.filmPatchIndex()
187  == film_.surfacePatchMap().nbrPolyPatch().index()
188  )
189  {
190  VoFFilmPtr = &VoFFilm;
191  }
192  }
193  }
194 
195  if (!VoFFilmPtr)
196  {
198  << "Cannot find VoFFilmTransfer fvModel for this film "
199  "in VoF region " << film_.surfacePatchMap().nbrMesh().name()
200  << exit(FatalError);
201  }
202 
203  return *VoFFilmPtr;
204 }
205 
206 
207 template<class Type, class TransferRateFunc>
209 inline Foam::fv::filmVoFTransfer::VoFToFilmTransferRate
210 (
211  TransferRateFunc transferRateFunc,
212  const dimensionSet& dimProp
213 ) const
214 {
215  const Foam::fvModels& fvModels
216  (
218  (
219  refCast<const fvMesh>(film_.surfacePatchMap().nbrMesh())
220  )
221  );
222 
224  (
226  (
227  "Su",
228  mesh(),
229  dimensioned<Type>(dimProp/dimTime, Zero)
230  )
231  );
232 
233  UIndirectList<Type>(tSu.ref(), film_.surfacePatch().faceCells()) =
234  film_.surfacePatchMap().fromNeighbour
235  (
236  (VoFFilm(fvModels).*transferRateFunc)()
237  );
238 
239  return tSu/mesh().V();
240 }
241 
242 
244 (
245  const volScalarField& rho,
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 == film_.alpha.name())
256  {
257  eqn +=
258  VoFToFilmTransferRate<scalar>
259  (
261  dimMass
262  )
263  - fvm::Sp(transferRate_*rho(), 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  const volScalarField& rho,
278  fvMatrix<scalar>& eqn,
279  const word& fieldName
280 ) const
281 {
282  if (debug)
283  {
284  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
285  }
286 
287  if (fieldName == film_.thermo.he().name())
288  {
289  eqn +=
290  VoFToFilmTransferRate<scalar>
291  (
293  dimEnergy
294  )
295  - fvm::Sp(alpha()*rho()*transferRate_, eqn.psi());
296  }
297  else
298  {
300  << "Support for field " << fieldName << " is not implemented"
301  << exit(FatalError);
302  }
303 }
304 
305 
307 (
308  const volScalarField& alpha,
309  const volScalarField& rho,
310  fvMatrix<vector>& eqn,
311  const word& fieldName
312 ) const
313 {
314  if (debug)
315  {
316  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
317  }
318 
319  eqn +=
320  VoFToFilmTransferRate<vector>
321  (
324  )
325  - fvm::Sp(alpha()*rho()*transferRate_, eqn.psi());
326 }
327 
328 
329 template<class Type, class FieldType>
330 inline Foam::tmp<Foam::Field<Type>> Foam::fv::filmVoFTransfer::TransferRate
331 (
332  const FieldType& f
333 ) const
334 {
335  const labelList& faceCells = film_.surfacePatch().faceCells();
336 
337  return tmp<Field<Type>>
338  (
339  new Field<Type>
340  (
342  (
343  film_.alpha()*transferRate_*mesh().V()*f,
344  faceCells
345  )
346  )
347  );
348 }
349 
350 
353 {
354  return TransferRate<scalar>(oneField());
355 }
356 
357 
360 {
361  return TransferRate<scalar>(film_.thermo.rho()());
362 }
363 
364 
367 {
368  return TransferRate<scalar>(film_.thermo.rho()()*film_.thermo.he()());
369 }
370 
371 
374 {
375  return TransferRate<vector>(film_.thermo.rho()()*film_.U());
376 }
377 
378 
380 {
381  transferRate_.setSize(mesh().nCells());
382 }
383 
384 
386 {
387  transferRate_.setSize(mesh().nCells());
388 }
389 
390 
392 {
393  transferRate_.setSize(mesh().nCells());
394 }
395 
396 
398 {
399  return true;
400 }
401 
402 
403 // ************************************************************************* //
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: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
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:893
Finite volume model abstract base class.
Definition: fvModel.H:59
Finite volume models.
Definition: fvModels.H:65
Film<->VoF transfer model.
tmp< scalarField > heTransferRate() const
Return the energy transfer rate.
tmp< vectorField > UTransferRate() const
Return the momentum transfer rate.
const volScalarField & alpha() const
tmp< scalarField > rhoTransferRate() const
Return the mass transfer rate.
Film<->VoF transfer model.
filmVoFTransfer(const word &sourceName, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
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.
virtual void correct()
Solve the film and update the sources.
tmp< scalarField > transferRate() const
Return the volume transfer rate.
virtual void addSup(const volScalarField &alpha, fvMatrix< scalar > &eqn, const word &fieldName) const
Add explicit contribution to phase continuity.
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.
A class representing the concept of a field of 1 used to avoid unnecessary manipulations for objects ...
Definition: oneField.H:53
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:55
Foam::fvModels & fvModels() const
Return the fvModels that are created on demand.
Definition: solver.C:85
const fvMesh & mesh
Region mesh.
Definition: solver.H:94
Solver module for for 2 compressible, non-isothermal immiscible fluids using a VOF (volume of fluid) ...
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 dimMomentum
const dimensionSet dimTime
typename VolField< Type >::Internal VolInternalField
Definition: volFieldsFwd.H:58
const dimensionSet dimMass
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
labelList f(nPoints)
labelList fv(nPoints)
dictionary dict