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 "volFields.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  fvMatrix<Type>& eqn,
88  const word& fieldName
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  fvMatrix<Type>& eqn,
97  const word& fieldName
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  fvMatrix<Type>& eqn,
107  const word& fieldName
108  ) const;
109 
110 
111  //- Return source for equation with specified name and dimensions
112  template<class Type, class ... AlphaRhoFieldTypes>
114  (
115  const VolField<Type>& field,
116  const word& fieldName,
117  const dimensionSet& ds,
118  const AlphaRhoFieldTypes& ... alphaRhos
119  ) const;
120 
121 
122 public:
123 
124  //- Runtime type information
125  TypeName("fvModel");
126 
127 
128  // Declare run-time constructor selection table
129 
131  (
132  autoPtr,
133  fvModel,
134  dictionary,
135  (
136  const word& name,
137  const word& modelType,
138  const fvMesh& mesh,
139  const dictionary& dict
140  ),
141  (name, modelType, mesh, dict)
142  );
143 
144 
145  // Static Member Functions
146 
147  //- Return the dimensions of the matrix of a source term
148  template
149  <
150  class Type,
151  class AlphaRhoFieldType,
152  class ... AlphaRhoFieldTypes
153  >
154  inline static dimensionSet sourceDims
155  (
156  const VolField<Type>& field,
157  const dimensionSet& ds,
158  const AlphaRhoFieldType& alphaRho,
159  const AlphaRhoFieldTypes& ... alphaRhos
160  );
161 
162  //- Return the dimensions of the matrix of a source term (base
163  // condition for the above)
164  template<class Type>
165  inline static dimensionSet sourceDims
166  (
167  const VolField<Type>& field,
168  const dimensionSet& ds
169  );
170 
171 
172  // Constructors
173 
174  //- Construct from components
175  fvModel
176  (
177  const word& name,
178  const word& modelType,
179  const fvMesh& mesh,
180  const dictionary& dict
181  );
182 
183  //- Return clone
184  autoPtr<fvModel> clone() const
185  {
187  return autoPtr<fvModel>(nullptr);
188  }
189 
190  //- Return pointer to new fvModel object created
191  // on the freestore from an Istream
192  class iNew
193  {
194  //- Reference to the mesh
195  const fvMesh& mesh_;
196 
197  const word& name_;
198 
199  public:
200 
201  iNew
202  (
203  const fvMesh& mesh,
204  const word& name
205  )
206  :
207  mesh_(mesh),
208  name_(name)
209  {}
210 
212  {
213  // const word name(is);
214  const dictionary dict(is);
215 
216  return autoPtr<fvModel>
217  (
218  fvModel::New(name_, mesh_, dict)
219  );
220  }
221  };
222 
223 
224  // Selectors
225 
226  //- Return a reference to the selected fvModel
227  static autoPtr<fvModel> New
228  (
229  const word& name,
230  const fvMesh& mesh,
231  const dictionary& dict
232  );
233 
234 
235  //- Destructor
236  virtual ~fvModel();
237 
238 
239  // Member Functions
240 
241  // Access
242 
243  //- Return const access to the source name
244  inline const word& name() const;
245 
246  //- Return const access to the mesh database
247  inline const fvMesh& mesh() const;
248 
249  //- Return dictionary
250  inline const dictionary& coeffs() const;
251 
252 
253  // Checks
254 
255  //- Return the list of fields for which the fvModel adds source term
256  // to the transport equation
257  virtual wordList addSupFields() const;
258 
259  //- Return true if the fvModel adds a source term to the given
260  // field's transport equation
261  virtual bool addsSupToField(const word& fieldName) const;
262 
263  //- Return the maximum time-step for stable operation
264  virtual scalar maxDeltaT() const;
265 
266 
267  // Sources
268 
269  //- Add a source term to an equation
271 
272  //- Add a source term to a compressible equation
274 
275  //- Add a source term to a phase equation
277 
278  //- Return source for an equation
279  template<class Type>
281  (
282  const VolField<Type>& field
283  ) const;
284 
285  //- Return source for an equation with a specified name
286  template<class Type>
288  (
289  const VolField<Type>& field,
290  const word& fieldName
291  ) const;
292 
293  //- Return source for a compressible equation
294  template<class Type>
296  (
297  const volScalarField& rho,
298  const VolField<Type>& field
299  ) const;
300 
301  //- Return source for a compressible equation with a specified name
302  template<class Type>
304  (
305  const volScalarField& rho,
306  const VolField<Type>& field,
307  const word& fieldName
308  ) const;
309 
310  //- Return source for a phase equation
311  template<class Type>
313  (
314  const volScalarField& alpha,
315  const volScalarField& rho,
316  const VolField<Type>& field
317  ) const;
318 
319  //- Return source for a phase equation with a specified name
320  template<class Type>
322  (
323  const volScalarField& alpha,
324  const volScalarField& rho,
325  const VolField<Type>& field,
326  const word& fieldName
327  ) const;
328 
329  //- Return source for a phase equation
330  template<class Type>
332  (
333  const volScalarField& alpha,
334  const geometricOneField& rho,
335  const VolField<Type>& field
336  ) const;
337 
338  //- Return source for a phase equation
339  template<class Type>
341  (
342  const geometricOneField& alpha,
343  const volScalarField& rho,
344  const VolField<Type>& field
345  ) const;
346 
347  //- Return source for a phase equation
348  template<class Type>
350  (
351  const geometricOneField& alpha,
352  const geometricOneField& rho,
353  const VolField<Type>& field
354  ) const;
355 
356  //- Return source for an equation with a second time derivative
357  template<class Type>
359  (
360  const VolField<Type>& field
361  ) const;
362 
363  //- Return source for an equation with a second time derivative
364  template<class Type>
366  (
367  const VolField<Type>& field,
368  const word& fieldName
369  ) const;
370 
371 
372  // Mesh changes
373 
374  //- Prepare for mesh update
375  virtual void preUpdateMesh();
376 
377  //- Update for mesh motion
378  virtual bool movePoints() = 0;
379 
380  //- Update topology using the given map
381  virtual void topoChange(const polyTopoChangeMap&) = 0;
382 
383  //- Update from another mesh using the given map
384  virtual void mapMesh(const polyMeshMap&) = 0;
385 
386  //- Redistribute or update using the given distribution map
387  virtual void distribute(const polyDistributionMap&) = 0;
388 
389 
390  //- Correct the fvModel
391  // e.g. solve equations, update model, for film, Lagrangian etc.
392  virtual void correct();
393 
394 
395  // IO
396 
397  //- Read source dictionary
398  virtual bool read(const dictionary& dict);
399 };
400 
401 
402 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
403 
404 } // End namespace Foam
405 
406 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
407 
408 #include "fvModelI.H"
409 
410 #ifdef NoRepository
411  #include "fvModelTemplates.C"
412 #endif
413 
414 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
415 
416 #endif
417 
418 // ************************************************************************* //
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:160
Dimension set for the base types.
Definition: dimensionSet.H:122
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:101
Return pointer to new fvModel object created.
Definition: fvModel.H:192
autoPtr< fvModel > operator()(Istream &is) const
Definition: fvModel.H:210
iNew(const fvMesh &mesh, const word &name)
Definition: fvModel.H:201
Finite volume model abstract base class.
Definition: fvModel.H:59
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.
declareRunTimeSelectionTable(autoPtr, fvModel, dictionary,(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict),(name, modelType, mesh, dict))
virtual ~fvModel()
Destructor.
Definition: fvModel.C:154
virtual void correct()
Correct the fvModel.
Definition: fvModel.C:199
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:40
void addSupType(fvMatrix< Type > &eqn, const word &fieldName) const
Add a source term to an equation.
Definition: fvModel.C:42
tmp< fvMatrix< Type > > source(const VolField< Type > &field, const word &fieldName, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhos) const
Return source for equation with specified name and dimensions.
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:187
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModel.C:195
static dimensionSet sourceDims(const VolField< Type > &field, const dimensionSet &ds, const AlphaRhoFieldType &alphaRho, const AlphaRhoFieldTypes &... alphaRhos)
Return the dimensions of the matrix of a source term.
FOR_ALL_FIELD_TYPES(DEFINE_FV_MODEL_ADD_SUP)
Add a source term to an equation.
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
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:183
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:28
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
Definition: fvModel.C:166
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:353
Forward declarations of fvMatrix specialisations.
#define DEFINE_FV_MODEL_ADD_SUP(Type, nullArg)
Definition: fvModelM.H:26
#define DEFINE_FV_MODEL_ADD_RHO_SUP(Type, nullArg)
Definition: fvModelM.H:43
#define DEFINE_FV_MODEL_ADD_ALPHA_RHO_SUP(Type, nullArg)
Definition: fvModelM.H:62
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
dictionary dict