VoFFilmTransfer.H
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 Class
25  Foam::fv::VoFFilmTransfer
26 
27 Description
28  Film<->VoF transfer model
29 
30 Usage
31  Example usage:
32  \verbatim
33  VoFFilmTransfer
34  {
35  type VoFFilmTransfer;
36 
37  libs ("libfilmVoFTransfer.so");
38 
39  filmPatch film;
40  phase liquid;
41 
42  deltaFactorToFilm 0.9;
43  alphaToFilm 0.86;
44 
45  transferRateCoeff 0.1;
46  }
47  \endverbatim
48 
49 SourceFiles
50  VoFFilmTransfer.C
51 
52 \*---------------------------------------------------------------------------*/
53 
54 #ifndef VoFFilmTransfer_H
55 #define VoFFilmTransfer_H
56 
57 #include "fvModel.H"
58 #include "compressibleVoF.H"
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 namespace Foam
63 {
64 namespace fv
65 {
66 
67 /*---------------------------------------------------------------------------*\
68  Class VoFFilmTransfer Declaration
69 \*---------------------------------------------------------------------------*/
70 
71 class VoFFilmTransfer
72 :
73  public fvModel
74 {
75  // Private Data
76 
77  //- The VoF model
78  const solvers::compressibleVoF& VoF_;
79 
80  //- Film patch name
81  word filmPatchName_;
82 
83  const label filmPatchi_;
84 
85  //- The name of the transferred phase
86  word phaseName_;
87 
88  //- Reference to the transferred phase thermo
89  const rhoThermo& thermo_;
90 
91  //- Reference to the transferred phase volume fraction
92  const volScalarField& alpha_;
93 
94  //- Current time index (used for updating)
95  mutable label curTimeIndex_;
96 
97  //- Factor of the cell height below which the VoF may be transferred
98  // to the film
99  scalar deltaFactorToFilm_;
100 
101  //- VoF limit below which the VoF may be transferred to the film
102  scalar alphaToFilm_;
103 
104  //- Transfer rate coefficient
105  scalar transferRateCoeff_;
106 
107  volScalarField::Internal transferRate_;
108 
109 
110  // Private Member Functions
111 
112  //- Return the transfer rate from the film transferRateFunc
113  template<class Type, class TransferRateFunc>
114  inline tmp<VolInternalField<Type>> filmVoFTransferRate
115  (
116  TransferRateFunc transferRateFunc,
117  const dimensionSet& dimProp
118  ) const;
119 
120  //- Return the transfer rate of field f
121  template<class Type, class FieldType>
122  inline tmp<Field<Type>> TransferRate(const FieldType& f) const;
123 
124 
125 public:
126 
127  //- Runtime type information
128  TypeName("VoFFilmTransfer");
129 
130 
131  // Constructors
132 
133  //- Construct from explicit source name and mesh
135  (
136  const word& sourceName,
137  const word& modelType,
138  const fvMesh& mesh,
139  const dictionary& dict
140  );
141 
142  //- Disallow default bitwise copy construction
144  (
145  const VoFFilmTransfer&
146  ) = delete;
147 
148 
149  // Member Functions
150 
151  label filmPatchIndex() const
152  {
153  return filmPatchi_;
154  }
155 
156  const volScalarField& alpha() const
157  {
158  return alpha_;
159  }
160 
161 
162  // Checks
163 
164  //- Return the list of fields for which the option adds source term
165  // to the transport equation
166  virtual wordList addSupFields() const;
167 
168 
169  // Correct
170 
171  //- Solve the film and update the sources
172  virtual void correct();
173 
174 
175  // Add explicit and implicit contributions to compressible equation
176 
177  //- Add implicit contribution to phase-fraction equation
178  virtual void addSup
179  (
180  fvMatrix<scalar>& eqn,
181  const word& fieldName
182  ) const;
183 
184  //- Add implicit contribution to phase density equation
185  virtual void addSup
186  (
187  const volScalarField& alpha,
188  fvMatrix<scalar>& eqn,
189  const word& fieldName
190  ) const;
191 
192  //- Add implicit contribution to phase energy equation
193  virtual void addSup
194  (
195  const volScalarField& alpha,
196  const volScalarField& rho,
197  fvMatrix<scalar>& eqn,
198  const word& fieldName
199  ) const;
200 
201  //- Add implicit contribution to mixture momentum equation
202  virtual void addSup
203  (
204  const volScalarField& rho,
205  fvMatrix<vector>& eqn,
206  const word& fieldName
207  ) const;
208 
209 
210  // Transfer to film
211 
212  //- Return the mass transfer rate
214 
215  //- Return the energy transfer rate
217 
218  //- Return the momentum transfer rate
220 
221 
222  // Mesh changes
223 
224  //- Update topology using the given map
225  virtual void topoChange(const polyTopoChangeMap&);
226 
227  //- Update from another mesh using the given map
228  virtual void mapMesh(const polyMeshMap&);
229 
230  //- Redistribute or update using the given distribution map
231  virtual void distribute(const polyDistributionMap&);
232 
233  //- Update for mesh motion
234  virtual bool movePoints();
235 
236 
237  // Member Operators
238 
239  //- Disallow default bitwise assignment
240  void operator=(const VoFFilmTransfer&) = delete;
241 };
242 
243 
244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245 
246 } // End namespace fv
247 } // End namespace Foam
248 
249 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 
251 #endif
252 
253 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Finite volume model abstract base class.
Definition: fvModel.H:59
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
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.
void operator=(const VoFFilmTransfer &)=delete
Disallow default bitwise assignment.
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.
const volScalarField & alpha() const
TypeName("VoFFilmTransfer")
Runtime type information.
tmp< scalarField > rhoTransferRate() const
Return the mass transfer rate.
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.
Base-class for fluid thermodynamic properties based on density.
Definition: rhoThermo.H:55
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
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
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
labelList f(nPoints)
labelList fv(nPoints)
dictionary dict