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-2022 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 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 class fvMesh;
50 class polyTopoChangeMap;
51 class polyMeshMap;
52 class polyDistributionMap;
53 
54 /*---------------------------------------------------------------------------*\
55  Class fvModel Declaration
56 \*---------------------------------------------------------------------------*/
57 
58 class fvModel
59 {
60  // Private Member Data
61 
62  //- Model name
63  const word name_;
64 
65  //- Model type
66  const word modelType_;
67 
68  //- Reference to the mesh database
69  const fvMesh& mesh_;
70 
71  //- Top level source dictionary
72  dictionary dict_;
73 
74  //- Dictionary containing source coefficients
75  dictionary coeffs_;
76 
77 
78 protected:
79 
80  // Protected Member Functions
81 
82  //- Add a source term to an equation
83  template<class Type>
84  void addSupType
85  (
86  fvMatrix<Type>& eqn,
87  const word& fieldName
88  ) const;
89 
90  //- Add a source term to a compressible equation
91  template<class Type>
92  void addSupType
93  (
94  const volScalarField& rho,
95  fvMatrix<Type>& eqn,
96  const word& fieldName
97  ) const;
98 
99  //- Add a source term to a phase equation
100  template<class Type>
101  void addSupType
102  (
103  const volScalarField& alpha,
104  const volScalarField& rho,
105  fvMatrix<Type>& eqn,
106  const word& fieldName
107  ) const;
108 
109 
110  //- Return source for equation with specified name and dimensions
111  template<class Type, class ... AlphaRhoFieldTypes>
113  (
115  const word& fieldName,
116  const dimensionSet& ds,
117  const AlphaRhoFieldTypes& ... alphaRhos
118  ) const;
119 
120 
121 public:
122 
123  //- Runtime type information
124  TypeName("fvModel");
125 
126 
127  // Declare run-time constructor selection table
128 
130  (
131  autoPtr,
132  fvModel,
133  dictionary,
134  (
135  const word& name,
136  const word& modelType,
137  const dictionary& dict,
138  const fvMesh& mesh
139  ),
140  (name, modelType, dict, mesh)
141  );
142 
143 
144  // Static Member Functions
145 
146  //- Return the dimensions of the matrix of a source term
147  template
148  <
149  class Type,
150  class AlphaRhoFieldType,
151  class ... AlphaRhoFieldTypes
152  >
153  inline static dimensionSet sourceDims
154  (
156  const dimensionSet& ds,
157  const AlphaRhoFieldType& alphaRho,
158  const AlphaRhoFieldTypes& ... alphaRhos
159  );
160 
161  //- Return the dimensions of the matrix of a source term (base
162  // condition for the above)
163  template<class Type>
164  inline static dimensionSet sourceDims
165  (
167  const dimensionSet& ds
168  );
169 
170 
171  // Constructors
172 
173  //- Construct from components
174  fvModel
175  (
176  const word& name,
177  const word& modelType,
178  const dictionary& dict,
179  const fvMesh& mesh
180  );
181 
182  //- Return clone
183  autoPtr<fvModel> clone() const
184  {
186  return autoPtr<fvModel>(nullptr);
187  }
188 
189  //- Return pointer to new fvModel object created
190  // on the freestore from an Istream
191  class iNew
192  {
193  //- Reference to the mesh
194  const fvMesh& mesh_;
195 
196  const word& name_;
197 
198  public:
199 
201  (
202  const fvMesh& mesh,
203  const word& name
204  )
205  :
206  mesh_(mesh),
207  name_(name)
208  {}
211  {
212  // const word name(is);
213  const dictionary dict(is);
214 
215  return autoPtr<fvModel>
216  (
217  fvModel::New(name_, dict, mesh_)
218  );
219  }
220  };
221 
222 
223  // Selectors
224 
225  //- Return a reference to the selected fvModel
226  static autoPtr<fvModel> New
227  (
228  const word& name,
229  const dictionary& dict,
230  const fvMesh& mesh
231  );
232 
233 
234  //- Destructor
235  virtual ~fvModel();
236 
237 
238  // Member Functions
239 
240  // Access
241 
242  //- Return const access to the source name
243  inline const word& name() const;
244 
245  //- Return const access to the mesh database
246  inline const fvMesh& mesh() const;
247 
248  //- Return dictionary
249  inline const dictionary& coeffs() const;
250 
251 
252  // Checks
253 
254  //- Return the list of fields for which the fvModel adds source term
255  // to the transport equation
256  virtual wordList addSupFields() const;
257 
258  //- Return true if the fvModel adds a source term to the given
259  // field's transport equation
260  virtual bool addsSupToField(const word& fieldName) const;
261 
262  //- Return the maximum time-step for stable operation
263  virtual scalar maxDeltaT() const;
264 
265 
266  // Sources
267 
268  //- Add a source term to an equation
270 
271  //- Add a source term to a compressible equation
273 
274  //- Add a source term to a phase equation
276 
277  //- Return source for an equation
278  template<class Type>
280  (
282  ) const;
283 
284  //- Return source for an equation with a specified name
285  template<class Type>
287  (
289  const word& fieldName
290  ) const;
291 
292  //- Return source for a compressible equation
293  template<class Type>
295  (
296  const volScalarField& rho,
298  ) const;
299 
300  //- Return source for a compressible equation with a specified name
301  template<class Type>
303  (
304  const volScalarField& rho,
306  const word& fieldName
307  ) const;
308 
309  //- Return source for a phase equation
310  template<class Type>
312  (
313  const volScalarField& alpha,
314  const volScalarField& rho,
316  ) const;
317 
318  //- Return source for a phase equation with a specified name
319  template<class Type>
321  (
322  const volScalarField& alpha,
323  const volScalarField& rho,
325  const word& fieldName
326  ) const;
327 
328 
329  // Mesh changes
330 
331  //- Prepare for mesh update
332  virtual void preUpdateMesh();
333 
334  //- Update for mesh motion
335  virtual bool movePoints() = 0;
336 
337  //- Update topology using the given map
338  virtual void topoChange(const polyTopoChangeMap&) = 0;
339 
340  //- Update from another mesh using the given map
341  virtual void mapMesh(const polyMeshMap&) = 0;
342 
343  //- Redistribute or update using the given distribution map
344  virtual void distribute(const polyDistributionMap&) = 0;
345 
346 
347  //- Correct the fvModel
348  // e.g. solve equations, update model, for film, Lagrangian etc.
349  virtual void correct();
350 
351 
352  // IO
353 
354  //- Read source dictionary
355  virtual bool read(const dictionary& dict);
356 };
357 
358 
359 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
360 
361 } // End namespace Foam
362 
363 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
364 
365 #include "fvModelI.H"
366 
367 #ifdef NoRepository
368  #include "fvModelTemplates.C"
369 #endif
370 
371 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
372 
373 #endif
374 
375 // ************************************************************************* //
iNew(const fvMesh &mesh, const word &name)
Definition: fvModel.H:200
dictionary dict
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
tmp< fvMatrix< Type > > source(const GeometricField< Type, fvPatchField, volMesh > &field, const word &fieldName, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhos) const
Return source for equation with specified name and dimensions.
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:28
virtual ~fvModel()
Destructor.
Definition: fvModel.C:131
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define DEFINE_FV_MODEL_ADD_SUP(Type, nullArg)
Definition: fvModelM.H:26
static dimensionSet sourceDims(const GeometricField< Type, fvPatchField, volMesh > &field, const dimensionSet &ds, const AlphaRhoFieldType &alphaRho, const AlphaRhoFieldTypes &... alphaRhos)
Return the dimensions of the matrix of a source term.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
Definition: fvModel.C:137
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:164
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
#define DEFINE_FV_MODEL_ADD_ALPHA_RHO_SUP(Type, nullArg)
Definition: fvModelM.H:62
virtual void mapMesh(const polyMeshMap &)=0
Update from another mesh using the given map.
Finite volume model abstract base class.
Definition: fvModel.H:57
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModel.C:172
virtual void distribute(const polyDistributionMap &)=0
Redistribute or update using the given distribution map.
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
TypeName("fvModel")
Runtime type information.
virtual bool movePoints()=0
Update for mesh motion.
Return pointer to new fvModel object created.
Definition: fvModel.H:190
Dimension set for the base types.
Definition: dimensionSet.H:121
autoPtr< fvModel > operator()(Istream &is) const
Definition: fvModel.H:209
A class for handling words, derived from string.
Definition: word.H:59
#define DEFINE_FV_MODEL_ADD_RHO_SUP(Type, nullArg)
Definition: fvModelM.H:43
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
static autoPtr< fvModel > New(const word &name, const dictionary &dict, const fvMesh &mesh)
Return a reference to the selected fvModel.
Definition: fvModel.C:94
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
Definition: fvModel.C:143
FOR_ALL_FIELD_TYPES(DEFINE_FV_MODEL_ADD_SUP)
Add a source term to an equation.
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:40
virtual scalar maxDeltaT() const
Return the maximum time-step for stable operation.
Definition: fvModel.C:149
fvModel(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition: fvModel.C:73
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
virtual void topoChange(const polyTopoChangeMap &)=0
Update topology using the given map.
void addSupType(fvMatrix< Type > &eqn, const word &fieldName) const
Add a source term to an equation.
Definition: fvModel.C:42
declareRunTimeSelectionTable(autoPtr, fvModel, dictionary,(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh),(name, modelType, dict, mesh))
autoPtr< fvModel > clone() const
Return clone.
Definition: fvModel.H:182
Forward declarations of fvMatrix specialisations.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
virtual void correct()
Correct the fvModel.
Definition: fvModel.C:176
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
rDeltaTY field()
A class for managing temporary objects.
Definition: PtrList.H:53
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
Namespace for OpenFOAM.