fvModels.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) 2021-2022 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::fvModels
26 
27 Description
28  Finite volume models
29 
30 SourceFiles
31  fvModels.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef fvModels_H
36 #define fvModels_H
37 
38 #include "fvModel.H"
39 #include "PtrListDictionary.H"
40 #include "MeshObject.H"
41 #include "HashSet.H"
42 #include "volFields.H"
43 #include "geometricOneField.H"
44 #include "fvMesh.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 // Forward declaration of friend functions and operators
52 class fvModels;
53 
54 Ostream& operator<<(Ostream& os, const fvModels& models);
55 
56 
57 /*---------------------------------------------------------------------------*\
58  Class fvModels Declaration
59 \*---------------------------------------------------------------------------*/
60 
61 class fvModels
62 :
63  public MeshObject<fvMesh, UpdateableMeshObject, fvModels>,
64  public dictionary,
65  public PtrListDictionary<fvModel>
66 {
67  // Private Member Data
68 
69  //- Time index to check that all defined sources have been applied
70  mutable label checkTimeIndex_;
71 
72  //- Sets of the fields that have had sources added by the fvModels
73  mutable PtrList<wordHashSet> addSupFields_;
74 
75 
76  // Private Member Functions
77 
78  //- Create IO object if dictionary is present
79  IOobject createIOobject(const fvMesh& mesh) const;
80 
81  //- Check that all fvModels have been applied
82  void checkApplied() const;
83 
84  //- Return source for equation with specified name and dimensions
85  template<class Type, class ... AlphaRhoFieldTypes>
86  tmp<fvMatrix<Type>> source
87  (
89  const word& fieldName,
90  const dimensionSet& ds,
91  const AlphaRhoFieldTypes& ... alphaRhos
92  ) const;
93 
94 
95 public:
96 
97  //- Runtime type information
98  TypeName("fvModels");
99 
100 
101  // Constructors
102 
103  //- Construct from components with list of field names
104  fvModels(const fvMesh& mesh);
105 
106  //- Disallow default bitwise copy construction
107  fvModels(const fvModels&) = delete;
108 
109  //- Inherit the base New method
111 
112 
113  //- Destructor
114  virtual ~fvModels()
115  {}
116 
117 
118  // Member Functions
119 
120  //- Declare fvModels to be a global dictionary
121  virtual bool global() const
122  {
123  return true;
124  }
125 
126 
127  // Checks
128 
129  //- Return true if an fvModel adds a source term to the given
130  // field's transport equation
131  virtual bool addsSupToField(const word& fieldName) const;
132 
133  //- Return the maximum time-step for stable operation
134  virtual scalar maxDeltaT() const;
135 
136 
137  // Sources
138 
139  //- Correct the fvModels
140  // e.g. solve equations, update model, for film, Lagrangian etc.
141  virtual void correct();
142 
143  //- Return source for an equation
144  template<class Type>
145  tmp<fvMatrix<Type>> source
146  (
148  ) const;
149 
150  //- Return source for an equation with a specified name
151  template<class Type>
152  tmp<fvMatrix<Type>> source
153  (
155  const word& fieldName
156  ) const;
157 
158  //- Return source for a compressible equation
159  template<class Type>
160  tmp<fvMatrix<Type>> source
161  (
162  const volScalarField& rho,
164  ) const;
165 
166  //- Return source for a compressible equation with a specified name
167  template<class Type>
168  tmp<fvMatrix<Type>> source
169  (
170  const volScalarField& rho,
172  const word& fieldName
173  ) const;
174 
175  //- Return source for a phase equation
176  template<class Type>
177  tmp<fvMatrix<Type>> source
178  (
179  const volScalarField& alpha,
180  const volScalarField& rho,
182  ) const;
183 
184  //- Return source for a phase equation with a specified name
185  template<class Type>
186  tmp<fvMatrix<Type>> source
187  (
188  const volScalarField& alpha,
189  const volScalarField& rho,
191  const word& fieldName
192  ) const;
193 
194  //- Return source for a phase equation
195  template<class Type>
196  tmp<fvMatrix<Type>> source
197  (
198  const volScalarField& alpha,
199  const geometricOneField& rho,
201  ) const;
202 
203  //- Return source for a phase equation
204  template<class Type>
205  tmp<fvMatrix<Type>> source
206  (
207  const geometricOneField& alpha,
208  const volScalarField& rho,
210  ) const;
211 
212  //- Return source for a phase equation
213  template<class Type>
214  tmp<fvMatrix<Type>> source
215  (
216  const geometricOneField& alpha,
217  const geometricOneField& rho,
219  ) const;
220 
221  //- Return source for an equation with a second time derivative
222  template<class Type>
224  (
226  ) const;
227 
228  //- Return source for an equation with a second time derivative
229  template<class Type>
231  (
233  const word& fieldName
234  ) const;
235 
236 
237  // Mesh changes
238 
239  //- Prepare for mesh update
240  virtual void preUpdateMesh();
241 
242  //- Update for mesh motion
243  virtual bool movePoints();
244 
245  //- Update topology using the given map
246  virtual void topoChange(const polyTopoChangeMap&);
247 
248  //- Update from another mesh using the given map
249  virtual void mapMesh(const polyMeshMap&);
250 
251  //- Redistribute or update using the given distribution map
252  virtual void distribute(const polyDistributionMap&);
253 
254 
255  // IO
256 
257  //- ReadData function which reads the fvModels dictionary
258  // required for regIOobject read operation
259  virtual bool readData(Istream&);
260 
261  //- Write data to Ostream
262  virtual bool writeData(Ostream& os) const;
263 
264  //- Read the fvModels dictionary if it has changed
265  // and update the models
266  virtual bool read();
267 
268 
269  // Member Operators
270 
271  //- Inherit the PtrListDictionary index operators
273 
274  //- Disallow default bitwise assignment
275  void operator=(const fvModels&) = delete;
276 
277 
278  // IOstream Operators
279 
280  //- Ostream operator
281  friend Ostream& operator<<
282  (
283  Ostream& os,
284  const fvModels& models
285  );
286 };
287 
288 
289 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
290 
291 //- Template function for obtaining global status
292 template<>
293 inline bool typeGlobal<fvModels>()
294 {
295  return true;
296 }
297 
298 
299 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
300 
301 #ifdef NoRepository
302  #include "fvModelsTemplates.C"
303 #endif
304 
305 
306 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
307 
308 } // End namespace Foam
309 
310 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
311 
312 #endif
313 
314 // ************************************************************************* //
virtual bool global() const
Declare fvModels to be a global dictionary.
Definition: fvModels.H:120
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
void operator=(const fvModels &)=delete
Disallow default bitwise assignment.
virtual scalar maxDeltaT() const
Return the maximum time-step for stable operation.
Definition: fvModels.C:240
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:353
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModels.C:255
Generic GeometricField class.
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:87
Dimension set for the base types.
Definition: dimensionSet.H:121
Template dictionary class which manages the storage associated with it.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
A class for handling words, derived from string.
Definition: word.H:59
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: fvModels.C:303
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: fvModels.C:281
Finite volume models.
Definition: fvModels.H:60
TypeName("fvModels")
Runtime type information.
tmp< fvMatrix< Type > > d2dt2(const GeometricField< Type, fvPatchField, volMesh > &field) const
Return source for an equation with a second time derivative.
bool typeGlobal< fvModels >()
Template function for obtaining global status.
Definition: fvModels.H:292
virtual bool readData(Istream &)
ReadData function which reads the fvModels dictionary.
Definition: fvModels.C:314
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
virtual ~fvModels()
Destructor.
Definition: fvModels.H:113
virtual bool writeData(Ostream &os) const
Write data to Ostream.
Definition: fvModels.C:321
Foam::fvModels & fvModels
virtual bool movePoints()
Update for mesh motion.
Definition: fvModels.C:266
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
Ostream & operator<<(Ostream &, const ensightPart &)
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: fvModels.C:292
rDeltaTY field()
A class for managing temporary objects.
Definition: PtrList.H:53
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
virtual bool addsSupToField(const word &fieldName) const
Return true if an fvModel adds a source term to the given.
Definition: fvModels.C:224
fvModels(const fvMesh &mesh)
Construct from components with list of field names.
Definition: fvModels.C:139
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
virtual bool read()
Read the fvModels dictionary if it has changed.
Definition: fvModels.C:328
Namespace for OpenFOAM.