fvModel.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-2025 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::fvModel
26 
27 Description
28  Finite volume model abstract base class.
29 
30 SourceFiles
31  fvModel.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef fvModel_H
36 #define fvModel_H
37 
38 #include "fvMatricesFwd.H"
39 #include "volFieldsFwd.H"
40 #include "dictionary.H"
41 #include "dimensionSet.H"
42 #include "HashSet.H"
43 #include "fvModelM.H"
44 #include "geometricOneField.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 class fvMesh;
52 class polyTopoChangeMap;
53 class polyMeshMap;
54 class polyDistributionMap;
55 
56 /*---------------------------------------------------------------------------*\
57  Class fvModel Declaration
58 \*---------------------------------------------------------------------------*/
59 
60 class fvModel
61 {
62  // Private Member Data
63 
64  //- Model name
65  const word name_;
66 
67  //- Model type
68  const word modelType_;
69 
70  //- Reference to the mesh database
71  const fvMesh& mesh_;
72 
73 
74 protected:
75 
76  // Protected Member Functions
77 
78  //- Add a source term to an equation
79  template<class Type>
80  void addSupType
81  (
82  const VolField<Type>& field,
83  fvMatrix<Type>& eqn
84  ) const;
85 
86  //- Add a source term to a compressible equation
87  template<class Type>
88  void addSupType
89  (
90  const volScalarField& rho,
91  const VolField<Type>& field,
92  fvMatrix<Type>& eqn
93  ) const;
94 
95  //- Add a source term to a phase equation
96  template<class Type>
97  void addSupType
98  (
99  const volScalarField& alpha,
100  const volScalarField& rho,
101  const VolField<Type>& field,
102  fvMatrix<Type>& eqn
103  ) const;
104 
105  //- Return a source for an equation
106  template<class Type, class ... AlphaRhoFieldTypes>
108  (
109  const VolField<Type>& eqnField,
110  const dimensionSet& ds,
111  const AlphaRhoFieldTypes& ... alphaRhoFields
112  ) const;
113 
114 
115 public:
116 
117  //- Runtime type information
118  TypeName("fvModel");
119 
120 
121  //- The keywords read by this class
122  static const wordHashSet keywords;
123 
124 
125  // Declare run-time constructor selection table
126 
128  (
129  autoPtr,
130  fvModel,
131  dictionary,
132  (
133  const word& name,
134  const word& modelType,
135  const fvMesh& mesh,
136  const dictionary& dict
137  ),
138  (name, modelType, mesh, dict)
139  );
140 
141 
142  // Static Member Functions
143 
144  //- Return the coefficients sub-dictionary for a given model type
145  inline static const dictionary& coeffs
146  (
147  const word& modelType,
148  const dictionary&
149  );
150 
151  //- Return the dimensions of the matrix of a source term
152  template<class AlphaRhoFieldType, class ... AlphaRhoFieldTypes>
153  static dimensionSet sourceDims
154  (
155  const dimensionSet& ds,
156  const AlphaRhoFieldType& alphaRhoField,
157  const AlphaRhoFieldTypes& ... alphaRhoFields
158  );
159 
160  //- Return the dimensions of the matrix of a source term (base
161  // condition for the above)
162  inline static const dimensionSet& sourceDims(const dimensionSet& ds);
163 
164  //- Return the name of the field associated with a source term
165  template<class AlphaRhoFieldType, class ... AlphaRhoFieldTypes>
166  static const word& fieldName
167  (
168  const AlphaRhoFieldType& alphaRhoField,
169  const AlphaRhoFieldTypes& ... alphaRhoFields
170  );
171 
172  //- Return the name of the field associated with a source term (base
173  // condition for the above)
174  template<class AlphaRhoFieldType>
175  static const word& fieldName(const AlphaRhoFieldType& alphaRhoField);
176 
177  //- Return the name of the field associated with a source term. Special
178  // overload for volume sources with no associated field.
179  static const word& fieldName();
180 
181 
182  // Constructors
183 
184  //- Construct from components
185  fvModel
186  (
187  const word& name,
188  const word& modelType,
189  const fvMesh& mesh,
190  const dictionary& dict
191  );
192 
193  //- Return clone
194  autoPtr<fvModel> clone() const
195  {
197  return autoPtr<fvModel>(nullptr);
198  }
199 
200  //- Return pointer to new fvModel object created
201  // on the freestore from an Istream
202  class iNew
203  {
204  //- Reference to the mesh
205  const fvMesh& mesh_;
206 
207  const word& name_;
208 
209  public:
210 
211  iNew
212  (
213  const fvMesh& mesh,
214  const word& name
215  )
216  :
217  mesh_(mesh),
218  name_(name)
219  {}
220 
222  {
223  // const word name(is);
224  const dictionary dict(is);
225 
226  return autoPtr<fvModel>
227  (
228  fvModel::New(name_, mesh_, dict)
229  );
230  }
231  };
232 
233 
234  // Selectors
235 
236  //- Return a reference to the selected fvModel
237  static autoPtr<fvModel> New
238  (
239  const word& name,
240  const fvMesh& mesh,
241  const dictionary& dict
242  );
243 
244 
245  //- Destructor
246  virtual ~fvModel();
247 
248 
249  // Member Functions
250 
251  // Access
252 
253  //- Return const access to the source name
254  inline const word& name() const;
255 
256  //- Return name as the keyword
257  inline const word& keyword() const;
258 
259  //- Return const access to the mesh database
260  inline const fvMesh& mesh() const;
261 
262  //- Return the coefficients sub-dictionary
263  inline const dictionary& coeffs(const dictionary&) const;
264 
265 
266  // Checks
267 
268  //- Return the list of fields for which the fvModel adds source term
269  // to the transport equation
270  virtual wordList addSupFields() const;
271 
272  //- Return true if the fvModel adds a source term to the given
273  // field's transport equation
274  virtual bool addsSupToField(const word& fieldName) const;
275 
276  //- Return the maximum time-step for stable operation
277  virtual scalar maxDeltaT() const;
278 
279 
280  // Sources
281 
282  //- Add a source term to a field-less proxy equation
283  virtual void addSup(fvMatrix<scalar>& eqn) const;
284 
285  //- Add a source term to an equation
287 
288  //- Add a source term to a compressible equation
290 
291  //- Add a source term to a phase equation
293 
294  //- Return source for an equation
295  template<class Type>
296  tmp<fvMatrix<Type>> sourceProxy
297  (
298  const VolField<Type>& eqnField
299  ) const;
300 
301  //- Return source for an equation
302  template<class Type>
303  tmp<fvMatrix<Type>> source
304  (
305  const VolField<Type>& field
306  ) const;
307 
308  //- Return source for an equation
309  template<class Type>
310  tmp<fvMatrix<Type>> sourceProxy
311  (
312  const VolField<Type>& field,
313  const VolField<Type>& eqnField
314  ) const;
315 
316  //- Return source for a compressible equation
317  template<class Type>
318  tmp<fvMatrix<Type>> source
319  (
320  const volScalarField& rho,
321  const VolField<Type>& field
322  ) const;
323 
324  //- Return source for a compressible equation
325  template<class Type>
326  tmp<fvMatrix<Type>> sourceProxy
327  (
328  const volScalarField& rho,
329  const VolField<Type>& field,
330  const VolField<Type>& eqnField
331  ) const;
332 
333  //- Return source for a phase equation
334  template<class Type>
335  tmp<fvMatrix<Type>> source
336  (
337  const volScalarField& alpha,
338  const volScalarField& rho,
339  const VolField<Type>& field
340  ) const;
341 
342  //- Return source for a phase equation
343  template<class Type>
344  tmp<fvMatrix<Type>> sourceProxy
345  (
346  const volScalarField& alpha,
347  const volScalarField& rho,
348  const VolField<Type>& field,
349  const VolField<Type>& eqnField
350  ) const;
351 
352  //- Return source for a phase equation
353  template<class Type>
354  tmp<fvMatrix<Type>> source
355  (
356  const volScalarField& alpha,
357  const geometricOneField& rho,
358  const VolField<Type>& field
359  ) const;
360 
361  //- Return source for a phase equation
362  template<class Type>
363  tmp<fvMatrix<Type>> source
364  (
365  const geometricOneField& alpha,
366  const volScalarField& rho,
367  const VolField<Type>& field
368  ) const;
369 
370  //- Return source for a phase equation
371  template<class Type>
372  tmp<fvMatrix<Type>> source
373  (
374  const geometricOneField& alpha,
375  const geometricOneField& rho,
376  const VolField<Type>& field
377  ) const;
378 
379  //- Return source for an equation with a second time derivative
380  template<class Type>
381  tmp<fvMatrix<Type>> d2dt2
382  (
383  const VolField<Type>& field
384  ) const;
385 
386 
387  // Mesh changes
388 
389  //- Prepare for mesh update
390  virtual void preUpdateMesh();
391 
392  //- Update for mesh motion
393  virtual bool movePoints() = 0;
394 
395  //- Update topology using the given map
396  virtual void topoChange(const polyTopoChangeMap&) = 0;
397 
398  //- Update from another mesh using the given map
399  virtual void mapMesh(const polyMeshMap&) = 0;
400 
401  //- Redistribute or update using the given distribution map
402  virtual void distribute(const polyDistributionMap&) = 0;
403 
404 
405  //- Correct the fvModel
406  // e.g. solve equations, update model, for film, Lagrangian etc.
407  virtual void correct();
408 
409 
410  // IO
411 
412  //- Read source dictionary
413  virtual bool read(const dictionary& dict);
414 
415  //- Write fvModel data
416  virtual bool write(const bool write = true) const;
417 };
418 
419 
420 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
421 
422 } // End namespace Foam
423 
424 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
425 
426 #include "fvModelI.H"
427 
428 #ifdef NoRepository
429  #include "fvModelTemplates.C"
430 #endif
431 
432 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
433 
434 #endif
435 
436 // ************************************************************************* //
Generic GeometricField class.
A HashTable with keys but without contents.
Definition: HashSet.H:62
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
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:96
Return pointer to new fvModel object created.
Definition: fvModel.H:202
autoPtr< fvModel > operator()(Istream &is) const
Definition: fvModel.H:220
iNew(const fvMesh &mesh, const word &name)
Definition: fvModel.H:211
Finite volume model abstract base class.
Definition: fvModel.H:60
tmp< fvMatrix< Type > > sourceTerm(const VolField< Type > &eqnField, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhoFields) const
Return a source for an equation.
virtual void mapMesh(const polyMeshMap &)=0
Update from another mesh using the given map.
static autoPtr< fvModel > New(const word &name, const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected fvModel.
Definition: fvModel.C:95
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
Definition: fvModel.C:161
static const dictionary & coeffs(const word &modelType, const dictionary &)
Return the coefficients sub-dictionary for a given model type.
Definition: fvModelI.H:31
virtual bool movePoints()=0
Update for mesh motion.
TypeName("fvModel")
Runtime type information.
tmp< fvMatrix< Type > > d2dt2(const VolField< Type > &field) const
Return source for an equation with a second time derivative.
virtual void distribute(const polyDistributionMap &)=0
Redistribute or update using the given distribution map.
void addSupType(const VolField< Type > &field, fvMatrix< Type > &eqn) const
Add a source term to an equation.
Definition: fvModel.C:45
declareRunTimeSelectionTable(autoPtr, fvModel, dictionary,(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict),(name, modelType, mesh, dict))
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
virtual ~fvModel()
Destructor.
Definition: fvModel.C:155
virtual void correct()
Correct the fvModel.
Definition: fvModel.C:196
static dimensionSet sourceDims(const dimensionSet &ds, const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields)
Return the dimensions of the matrix of a source term.
virtual scalar maxDeltaT() const
Return the maximum time-step for stable operation.
Definition: fvModel.C:173
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:200
static const wordHashSet keywords
The keywords read by this class.
Definition: fvModel.H:121
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModel.C:192
tmp< fvMatrix< Type > > sourceProxy(const VolField< Type > &eqnField) const
Add a source term to an equation.
virtual bool write(const bool write=true) const
Write fvModel data.
Definition: fvModel.C:206
const word & keyword() const
Return name as the keyword.
Definition: fvModelI.H:63
virtual void addSup(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
Definition: fvModel.C:179
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:69
fvModel(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
Definition: fvModel.C:76
autoPtr< fvModel > clone() const
Return clone.
Definition: fvModel.H:193
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:57
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
Definition: fvModel.C:167
static const word & fieldName()
Return the name of the field associated with a source term. Special.
Definition: fvModelI.H:49
virtual void topoChange(const polyTopoChangeMap &)=0
Update topology using the given map.
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
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
Forward declarations of fvMatrix specialisations.
#define DEFINE_FV_MODEL_ADD_ALPHA_RHO_FIELD_SUP(Type, nullArg)
Definition: fvModelM.H:62
#define DEFINE_FV_MODEL_ADD_RHO_FIELD_SUP(Type, nullArg)
Definition: fvModelM.H:43
#define DEFINE_FV_MODEL_ADD_FIELD_SUP(Type, nullArg)
Definition: fvModelM.H:26
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
dictionary dict