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