filmVoFTransfer.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::filmVoFTransfer
26 
27 Description
28  Film<->VoF transfer model
29 
30 Usage
31  Example usage:
32  \verbatim
33  filmVoFTransfer
34  {
35  type filmVoFTransfer;
36 
37  deltaFactorToVoF 1.5;
38  alphaToVoF 0.9;
39 
40  transferRateCoeff 0.1;
41  }
42  \endverbatim
43 
44 SourceFiles
45  filmVoFTransfer.C
46 
47 \*---------------------------------------------------------------------------*/
48 
49 #ifndef filmVoFTransfer_H
50 #define filmVoFTransfer_H
51 
52 #include "fvModel.H"
53 #include "isothermalFilm.H"
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 namespace Foam
58 {
59 namespace fv
60 {
61 
62 class VoFFilmTransfer;
63 
64 /*---------------------------------------------------------------------------*\
65  Class filmVoFTransfer Declaration
66 \*---------------------------------------------------------------------------*/
67 
68 class filmVoFTransfer
69 :
70  public fvModel
71 {
72  // Private Data
73 
74  //- The film model
75  const solvers::isothermalFilm& film_;
76 
77  //- Current time index (used for updating)
78  mutable label curTimeIndex_;
79 
80  //- Factor of the cell height above which the film is transferred
81  // to the VoF
82  scalar deltaFactorToVoF_;
83 
84  //- VoF limit above which all of the film is transferred to the VoF
85  scalar alphaToVoF_;
86 
87  //- Transfer rate coefficient
88  scalar transferRateCoeff_;
89 
90  //- Cached transfer rate
91  volScalarField::Internal transferRate_;
92 
93 
94  // Private Member Functions
95 
96  const VoFFilmTransfer& VoFFilm(const fvModels&) const;
97 
98  //- Return the transfer rate from the VoF transferRateFunc
99  template<class Type, class TransferRateFunc>
100  inline tmp<VolInternalField<Type>> VoFToFilmTransferRate
101  (
102  TransferRateFunc transferRateFunc,
103  const dimensionSet& dimProp
104  ) const;
105 
106  //- Return the transfer rate of field f
107  template<class Type, class FieldType>
108  inline tmp<Field<Type>> TransferRate(const FieldType& f) const;
109 
110 
111 public:
112 
113  //- Runtime type information
114  TypeName("filmVoFTransfer");
115 
116 
117  // Constructors
118 
119  //- Construct from explicit source name and mesh
121  (
122  const word& sourceName,
123  const word& modelType,
124  const fvMesh& mesh,
125  const dictionary& dict
126  );
127 
128  //- Disallow default bitwise copy construction
130  (
131  const filmVoFTransfer&
132  ) = delete;
133 
134 
135  // Member Functions
136 
137  // Checks
138 
139  //- Return the list of fields for which the option adds source term
140  // to the transport equation
141  virtual wordList addSupFields() const;
142 
143 
144  // Correct
145 
146  //- Solve the film and update the sources
147  virtual void correct();
148 
149 
150  // Add explicit and implicit contributions to compressible equation
151 
152  //- Add explicit contribution to phase continuity
153  virtual void addSup
154  (
155  const volScalarField& alpha,
156  fvMatrix<scalar>& eqn,
157  const word& fieldName
158  ) const;
159 
160  //- Add explicit contribution to phase energy equation
161  virtual void addSup
162  (
163  const volScalarField& alpha,
164  const volScalarField& rho,
165  fvMatrix<scalar>& eqn,
166  const word& fieldName
167  ) const;
168 
169  //- Add implicit contribution to mixture momentum equation
170  virtual void addSup
171  (
172  const volScalarField& alpha,
173  const volScalarField& rho,
174  fvMatrix<vector>& eqn,
175  const word& fieldName
176  ) const;
177 
178 
179  // Transfer to VoF
180 
181  //- Return the volume transfer rate
183 
184  //- Return the mass transfer rate
186 
187  //- Return the energy transfer rate
189 
190  //- Return the momentum transfer rate
192 
193 
194  // Mesh changes
195 
196  //- Update topology using the given map
197  virtual void topoChange(const polyTopoChangeMap&);
198 
199  //- Update from another mesh using the given map
200  virtual void mapMesh(const polyMeshMap&);
201 
202  //- Redistribute or update using the given distribution map
203  virtual void distribute(const polyDistributionMap&);
204 
205  //- Update for mesh motion
206  virtual bool movePoints();
207 
208 
209  // Member Operators
210 
211  //- Disallow default bitwise assignment
212  void operator=(const filmVoFTransfer&) = delete;
213 };
214 
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 } // End namespace fv
219 } // End namespace Foam
220 
221 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222 
223 #endif
224 
225 // ************************************************************************* //
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
Finite volume models.
Definition: fvModels.H:65
Film<->VoF transfer model.
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.
void operator=(const filmVoFTransfer &)=delete
Disallow default bitwise assignment.
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.
TypeName("filmVoFTransfer")
Runtime type information.
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.
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.
Solver module for flow of compressible isothermal liquid films.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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