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, TopoChangeableMeshObject, 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 a source for an equation
85  template<class Type, class ... AlphaRhoFieldTypes>
86  tmp<fvMatrix<Type>> sourceTerm
87  (
88  const VolField<Type>& eqnField,
89  const dimensionSet& ds,
90  const AlphaRhoFieldTypes& ... alphaRhoFields
91  ) const;
92 
93 
94 
95 protected:
96 
97  friend class DemandDrivenMeshObject
98  <
99  fvMesh,
101  fvModels
102  >;
103 
104  // Protected Constructors
105 
106  //- Construct from components with list of field names
107  explicit fvModels(const fvMesh& mesh);
108 
109 
110 public:
111 
112  //- Runtime type information
113  TypeName("fvModels");
114 
115 
116  // Constructors
117 
118  //- Disallow default bitwise copy construction
119  fvModels(const fvModels&) = delete;
120 
121  //- Inherit the base New method
123  <
124  fvMesh,
126  fvModels
127  >::New;
128 
129 
130  //- Destructor
131  virtual ~fvModels()
132  {}
133 
134 
135  // Member Functions
136 
137  //- Declare fvModels to be a global dictionary
138  virtual bool global() const
139  {
140  return true;
141  }
142 
143 
144  // Checks
145 
146  //- Return true if an fvModel adds a source term to the given
147  // field's transport equation
148  virtual bool addsSupToField(const word& fieldName) const;
149 
150  //- Return the maximum time-step for stable operation
151  virtual scalar maxDeltaT() const;
152 
153 
154  // Sources
155 
156  //- Correct the fvModels
157  // e.g. solve equations, update model, for film, Lagrangian etc.
158  virtual void correct();
159 
160  //- Return source for an equation
161  template<class Type>
163  (
164  const VolField<Type>& eqnField
165  ) const;
166 
167  //- Return source for an equation
168  template<class Type>
170  (
171  const VolField<Type>& field
172  ) const;
173 
174  //- Return source for an equation
175  template<class Type>
177  (
178  const VolField<Type>& field,
179  const VolField<Type>& eqnField
180  ) const;
181 
182  //- Return source for a compressible equation
183  template<class Type>
185  (
186  const volScalarField& rho,
187  const VolField<Type>& field
188  ) const;
189 
190  //- Return source for a compressible equation
191  template<class Type>
193  (
194  const volScalarField& rho,
195  const VolField<Type>& field,
196  const VolField<Type>& eqnField
197  ) const;
198 
199  //- Return source for a phase equation
200  template<class Type>
202  (
203  const volScalarField& alpha,
204  const volScalarField& rho,
205  const VolField<Type>& field
206  ) const;
207 
208  //- Return source for a phase equation
209  template<class Type>
211  (
212  const volScalarField& alpha,
213  const volScalarField& rho,
214  const VolField<Type>& field,
215  const VolField<Type>& eqnField
216  ) const;
217 
218  //- Return source for a phase equation
219  template<class Type>
221  (
222  const volScalarField& alpha,
223  const geometricOneField& rho,
224  const VolField<Type>& field
225  ) const;
226 
227  //- Return source for a phase equation
228  template<class Type>
230  (
231  const geometricOneField& alpha,
232  const volScalarField& rho,
233  const VolField<Type>& field
234  ) const;
235 
236  //- Return source for a phase equation
237  template<class Type>
239  (
240  const geometricOneField& alpha,
241  const geometricOneField& rho,
242  const VolField<Type>& field
243  ) const;
244 
245  //- Return source for an equation with a second time derivative
246  template<class Type>
248  (
249  const VolField<Type>& field
250  ) const;
251 
252 
253  // Mesh changes
254 
255  //- Prepare for mesh update
256  virtual void preUpdateMesh();
257 
258  //- Update for mesh motion
259  virtual bool movePoints();
260 
261  //- Update topology using the given map
262  virtual void topoChange(const polyTopoChangeMap&);
263 
264  //- Update from another mesh using the given map
265  virtual void mapMesh(const polyMeshMap&);
266 
267  //- Redistribute or update using the given distribution map
268  virtual void distribute(const polyDistributionMap&);
269 
270 
271  // IO
272 
273  //- Read the fvModels dictionary if it has changed
274  // and update the models
275  virtual bool read();
276 
277  //- ReadData function which reads the fvModels dictionary
278  // required for regIOobject read operation
279  virtual bool readData(Istream&);
280 
281  //- writeData function required by regIOobject but not used
282  // for this class, writeObject is used instead
283  virtual bool writeData(Ostream& os) const;
284 
285  //- Write the fvModels
286  virtual bool writeObject
287  (
291  const bool write
292  ) const;
293 
294 
295  // Member Operators
296 
297  //- Inherit the PtrListDictionary index operators
299 
300  //- Inherit the PtrListDictionary size function
302 
303  //- Disallow default bitwise assignment
304  void operator=(const fvModels&) = delete;
305 
306 
307  // IOstream Operators
308 
309  //- Ostream operator
310  friend Ostream& operator<<
311  (
312  Ostream& os,
313  const fvModels& models
314  );
315 };
316 
317 
318 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
319 
320 //- Trait for obtaining global status
321 template<>
322 struct typeGlobal<fvModels>
323 {
324  static const bool global = true;
325 };
326 
327 
328 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
329 
330 #ifdef NoRepository
331  #include "fvModelsTemplates.C"
332 #endif
333 
334 
335 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
336 
337 } // End namespace Foam
338 
339 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
340 
341 #endif
342 
343 // ************************************************************************* //
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
Version number type.
Definition: IOstream.H:97
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:87
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:194
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
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
void write(Ostream &, const bool subDict=true) const
Write dictionary, normally with sub-dictionary formatting.
Definition: dictionaryIO.C:211
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:125
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Finite volume models.
Definition: fvModels.H:65
virtual bool movePoints()
Update for mesh motion.
Definition: fvModels.C:271
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool write) const
Write the fvModels.
Definition: fvModels.C:359
tmp< fvMatrix< Type > > d2dt2(const VolField< Type > &field) const
Return source for an equation with a second time derivative.
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:380
virtual bool writeData(Ostream &os) const
writeData function required by regIOobject but not used
Definition: fvModels.C:351
virtual bool global() const
Declare fvModels to be a global dictionary.
Definition: fvModels.H:137
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:130
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
tmp< fvMatrix< Type > > sourceProxy(const VolField< Type > &eqnField) const
Return source for an equation.
virtual bool readData(Istream &)
ReadData function which reads the fvModels dictionary.
Definition: fvModels.C:344
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:319
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 &os, const fvConstraints &constraints)
Trait for obtaining global status.
Definition: IOobject.H:504
static const bool global
Definition: IOobject.H:505