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-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::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 "DemandDrivenMeshObject.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 DemandDrivenMeshObject<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  (
88  const VolField<Type>& field,
89  const word& fieldName,
90  const dimensionSet& ds,
91  const AlphaRhoFieldTypes& ... alphaRhos
92  ) const;
93 
94 
95 
96 protected:
97 
99 
100  // Protected Constructors
101 
102  //- Construct from components with list of field names
103  explicit fvModels(const fvMesh& mesh);
104 
105 
106 public:
107 
108  //- Runtime type information
109  TypeName("fvModels");
110 
111 
112  // Constructors
113 
114  //- Disallow default bitwise copy construction
115  fvModels(const fvModels&) = delete;
116 
117  //- Inherit the base New method
119  <
120  fvMesh,
122  fvModels
123  >::New;
124 
125 
126  //- Destructor
127  virtual ~fvModels()
128  {}
129 
130 
131  // Member Functions
132 
133  //- Declare fvModels to be a global dictionary
134  virtual bool global() const
135  {
136  return true;
137  }
138 
139 
140  // Checks
141 
142  //- Return true if an fvModel adds a source term to the given
143  // field's transport equation
144  virtual bool addsSupToField(const word& fieldName) const;
145 
146  //- Return the maximum time-step for stable operation
147  virtual scalar maxDeltaT() const;
148 
149 
150  // Sources
151 
152  //- Correct the fvModels
153  // e.g. solve equations, update model, for film, Lagrangian etc.
154  virtual void correct();
155 
156  //- Return source for an equation
157  template<class Type>
158  tmp<fvMatrix<Type>> source
159  (
160  const VolField<Type>& field
161  ) const;
162 
163  //- Return source for an equation with a specified name
164  template<class Type>
165  tmp<fvMatrix<Type>> source
166  (
167  const VolField<Type>& field,
168  const word& fieldName
169  ) const;
170 
171  //- Return source for a compressible equation
172  template<class Type>
173  tmp<fvMatrix<Type>> source
174  (
175  const volScalarField& rho,
176  const VolField<Type>& field
177  ) const;
178 
179  //- Return source for a compressible equation with a specified name
180  template<class Type>
181  tmp<fvMatrix<Type>> source
182  (
183  const volScalarField& rho,
184  const VolField<Type>& field,
185  const word& fieldName
186  ) const;
187 
188  //- Return source for a phase equation
189  template<class Type>
190  tmp<fvMatrix<Type>> source
191  (
192  const volScalarField& alpha,
193  const volScalarField& rho,
194  const VolField<Type>& field
195  ) const;
196 
197  //- Return source for a phase equation with a specified name
198  template<class Type>
199  tmp<fvMatrix<Type>> source
200  (
201  const volScalarField& alpha,
202  const volScalarField& rho,
203  const VolField<Type>& field,
204  const word& fieldName
205  ) const;
206 
207  //- Return source for a phase equation
208  template<class Type>
209  tmp<fvMatrix<Type>> source
210  (
211  const volScalarField& alpha,
212  const geometricOneField& rho,
213  const VolField<Type>& field
214  ) const;
215 
216  //- Return source for a phase equation
217  template<class Type>
218  tmp<fvMatrix<Type>> source
219  (
220  const geometricOneField& alpha,
221  const volScalarField& rho,
222  const VolField<Type>& field
223  ) const;
224 
225  //- Return source for a phase equation
226  template<class Type>
227  tmp<fvMatrix<Type>> source
228  (
229  const geometricOneField& alpha,
230  const geometricOneField& rho,
231  const VolField<Type>& field
232  ) const;
233 
234  //- Return source for an equation with a second time derivative
235  template<class Type>
237  (
238  const VolField<Type>& field
239  ) const;
240 
241  //- Return source for an equation with a second time derivative
242  template<class Type>
244  (
245  const VolField<Type>& field,
246  const word& fieldName
247  ) const;
248 
249 
250  // Mesh changes
251 
252  //- Prepare for mesh update
253  virtual void preUpdateMesh();
254 
255  //- Update for mesh motion
256  virtual bool movePoints();
257 
258  //- Update topology using the given map
259  virtual void topoChange(const polyTopoChangeMap&);
260 
261  //- Update from another mesh using the given map
262  virtual void mapMesh(const polyMeshMap&);
263 
264  //- Redistribute or update using the given distribution map
265  virtual void distribute(const polyDistributionMap&);
266 
267 
268  // IO
269 
270  //- ReadData function which reads the fvModels dictionary
271  // required for regIOobject read operation
272  virtual bool readData(Istream&);
273 
274  //- Write data to Ostream
275  virtual bool writeData(Ostream& os) const;
276 
277  //- Read the fvModels dictionary if it has changed
278  // and update the models
279  virtual bool read();
280 
281 
282  // Member Operators
283 
284  //- Inherit the PtrListDictionary index operators
286 
287  //- Inherit the PtrListDictionary size function
289 
290  //- Disallow default bitwise assignment
291  void operator=(const fvModels&) = delete;
292 
293 
294  // IOstream Operators
295 
296  //- Ostream operator
297  friend Ostream& operator<<
298  (
299  Ostream& os,
300  const fvModels& models
301  );
302 };
303 
304 
305 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
306 
307 //- Trait for obtaining global status
308 template<>
309 struct typeGlobal<fvModels>
310 {
311  static const bool global = true;
312 };
313 
314 
315 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
316 
317 #ifdef NoRepository
318  #include "fvModelsTemplates.C"
319 #endif
320 
321 
322 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
323 
324 } // End namespace Foam
325 
326 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
327 
328 #endif
329 
330 // ************************************************************************* //
Templated abstract base-class for demand-driven mesh objects used to automate their allocation to the...
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Template dictionary class which manages the storage associated with it.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
UpdateableMeshObject(regIOobject &io, const fvMesh &mesh)
Definition: MeshObjects.H:150
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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Finite volume models.
Definition: fvModels.H:65
virtual bool movePoints()
Update for mesh motion.
Definition: fvModels.C:271
tmp< fvMatrix< Type > > d2dt2(const VolField< Type > &field) const
Return source for an equation with a second time derivative.
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:358
virtual bool writeData(Ostream &os) const
Write data to Ostream.
Definition: fvModels.C:326
virtual bool global() const
Declare fvModels to be a global dictionary.
Definition: fvModels.H:133
TypeName("fvModels")
Runtime type information.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: fvModels.C:286
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: fvModels.C:308
virtual scalar maxDeltaT() const
Return the maximum time-step for stable operation.
Definition: fvModels.C:245
virtual ~fvModels()
Destructor.
Definition: fvModels.H:126
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModels.C:260
void operator=(const fvModels &)=delete
Disallow default bitwise assignment.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: fvModels.C:297
virtual bool readData(Istream &)
ReadData function which reads the fvModels dictionary.
Definition: fvModels.C:319
virtual bool addsSupToField(const word &fieldName) const
Return true if an fvModel adds a source term to the given.
Definition: fvModels.C:229
fvModels(const fvMesh &mesh)
Construct from components with list of field names.
Definition: fvModels.C:144
virtual bool read()
Read the fvModels dictionary if it has changed.
Definition: fvModels.C:333
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
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
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
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
Ostream & operator<<(Ostream &, const ensightPart &)
Trait for obtaining global status.
Definition: IOobject.H:504
static const bool global
Definition: IOobject.H:505