volumeFractionSource.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) 2019-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::volumeFractionSource
26 
27 Description
28  This fvModel adds transport terms into the equations to account for the
29  presence of a constant volume fraction. The volume fraction is read from
30  constant/alpha.<volumePhase>, where <volumePhase> is given as a parameter
31  to the fvModel. Both advective and diffusive terms are added, and the
32  resulting solution is time-accurate. The flux and velocity are treated as
33  superficial.
34 
35  This can be used to represent the effect of porous media that are caused
36  purely by the reduction in volume of the fluid phase; i.e., additional
37  blockage, and changes to transport and diffusion rates. It does not
38  represent losses or transfers with the porous media. That requires separate
39  sub-modelling.
40 
41 Usage
42  \table
43  Property | Description | Req'd? | Default
44  phi | Name of the flux field | no | phi
45  rho | Name of the density field | no | rho
46  U | Name of the velocity field | no | U
47  volumePhase | Name of the phase associated with the volume fraction \\
48  | yes |
49  \endtable
50 
51  Example specification:
52  \verbatim
53  volumeFractionSource1
54  {
55  type volumeFractionSource;
56 
57  phi phi;
58  rho rho;
59  U U;
60 
61  volumePhase solid;
62  }
63  \endverbatim
64 
65 SourceFiles
66  volumeFractionSource.C
67 
68 \*---------------------------------------------------------------------------*/
69 
70 #ifndef volumeFractionSource_H
71 #define volumeFractionSource_H
72 
73 #include "fvModel.H"
74 #include "surfaceFields.H"
75 #include "volFields.H"
76 
77 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78 
79 namespace Foam
80 {
81 namespace fv
82 {
83 
84 /*---------------------------------------------------------------------------*\
85  Class volumeFractionSource Declaration
86 \*---------------------------------------------------------------------------*/
87 
88 class volumeFractionSource
89 :
90  public fvModel
91 {
92  // Private Member Data
93 
94  //- The name of the flux field
95  word phiName_;
96 
97  //- The name of the density field
98  word rhoName_;
99 
100  //- The name of the velocity field
101  word UName_;
102 
103  //- The name of the volume fraction phase
104  word volumePhaseName_;
105 
106 
107  // Private Member Functions
108 
109  //- Non-virtual read
110  void readCoeffs();
111 
112  //- Get the volume fraction field
113  const volScalarField& volumeAlpha() const;
114 
115  //- Get the diffusivity for a given field
116  tmp<volScalarField> D(const word& fieldName) const;
117 
118 
119  // Sources
120 
121  //- Add source terms to an equation
122  template<class Type, class AlphaFieldType>
123  void addGeneralSup
124  (
125  const AlphaFieldType& alpha,
127  const word& fieldName
128  ) const;
129 
130  //- Add a source term to an equation
131  template<class Type, class AlphaFieldType>
132  void addAlphaSupType
133  (
134  const AlphaFieldType& alpha,
135  fvMatrix<Type>& eqn,
136  const word& fieldName
137  ) const;
138 
139  //- Add a source term to a scalar equation
140  template<class AlphaFieldType>
141  void addAlphaSupType
142  (
143  const AlphaFieldType& alpha,
144  fvMatrix<scalar>& eqn,
145  const word& fieldName
146  ) const;
147 
148  //- Add a source term to a vector equation
149  template<class AlphaFieldType>
150  void addAlphaSupType
151  (
152  const AlphaFieldType& alpha,
153  fvMatrix<vector>& eqn,
154  const word& fieldName
155  ) const;
156 
157  //- Add a source term to an equation
158  template<class Type>
159  void addSupType(fvMatrix<Type>& eqn, const word& fieldName) const;
160 
161  //- Add a source term to a compressible equation
162  template<class Type>
163  void addSupType
164  (
165  const volScalarField& rho,
166  fvMatrix<Type>& eqn,
167  const word& fieldName
168  ) const;
169 
170  //- Add a source term to a phase equation
171  template<class Type>
172  void addSupType
173  (
174  const volScalarField& alpha,
175  const volScalarField& rho,
176  fvMatrix<Type>& eqn,
177  const word& fieldName
178  ) const;
179 
180 
181 
182 
183 public:
184 
185  //- Runtime type information
186  TypeName("volumeFractionSource");
187 
188 
189  // Constructors
190 
191  //- Construct from components
193  (
194  const word& name,
195  const word& modelType,
196  const fvMesh& mesh,
197  const dictionary& dict
198  );
199 
200  //- Disallow default bitwise copy construction
202 
203 
204  //- Destructor
205  virtual ~volumeFractionSource();
206 
207 
208  // Member Functions
209 
210  // Checks
211 
212  //- Return true if the fvModel adds a source term to the given
213  // field's transport equation
214  virtual bool addsSupToField(const word& fieldName) const;
215 
216  //- Return the list of fields for which the fvModel adds source term
217  // to the transport equation
218  virtual wordList addSupFields() const;
219 
220 
221  // Sources
222 
223  //- Add a source term to an equation
225 
226  //- Add a source term to a compressible equation
228 
229  //- Add a source term to a phase equation
231 
232 
233  // Mesh changes
234 
235  //- Update for mesh motion
236  virtual bool movePoints();
237 
238  //- Update topology using the given map
239  virtual void topoChange(const polyTopoChangeMap&);
240 
241  //- Update from another mesh using the given map
242  virtual void mapMesh(const polyMeshMap&);
243 
244  //- Redistribute or update using the given distribution map
245  virtual void distribute(const polyDistributionMap&);
246 
247 
248  // IO
249 
250  //- Read dictionary
251  virtual bool read(const dictionary& dict);
252 
253 
254  // Member Operators
255 
256  //- Disallow default bitwise assignment
257  void operator=(const volumeFractionSource&) = delete;
258 };
259 
260 
261 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
262 
263 } // End namespace fv
264 } // End namespace Foam
265 
266 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
267 
268 #endif
269 
270 // ************************************************************************* //
Generic GeometricField class.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
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
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:28
This fvModel adds transport terms into the equations to account for the presence of a constant volume...
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
virtual ~volumeFractionSource()
Destructor.
TypeName("volumeFractionSource")
Runtime type information.
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 bool read(const dictionary &dict)
Read dictionary.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
volumeFractionSource(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
FOR_ALL_FIELD_TYPES(DEFINE_FV_MODEL_ADD_SUP)
Add a source term to an equation.
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
void operator=(const volumeFractionSource &)=delete
Disallow default bitwise assignment.
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.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define DEFINE_FV_MODEL_ADD_SUP(Type, nullArg)
Definition: fvModelM.H:26
#define DEFINE_FV_MODEL_ADD_RHO_SUP(Type, nullArg)
Definition: fvModelM.H:43
#define DEFINE_FV_MODEL_ADD_ALPHA_RHO_SUP(Type, nullArg)
Definition: fvModelM.H:62
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
labelList fv(nPoints)
dictionary dict
Foam::surfaceFields.