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 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 mapPolyMesh;
51 
52 /*---------------------------------------------------------------------------*\
53  Class fvModel Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 class fvModel
57 {
58  // Private Member Data
59 
60  //- Model name
61  const word name_;
62 
63  //- Model type
64  const word modelType_;
65 
66  //- Reference to the mesh database
67  const fvMesh& mesh_;
68 
69  //- Top level source dictionary
70  dictionary dict_;
71 
72  //- Dictionary containing source coefficients
73  dictionary coeffs_;
74 
75 
76 protected:
77 
78  // Protected Member Functions
79 
80  //- Add a source term to an equation
81  template<class Type>
82  void addSupType
83  (
84  fvMatrix<Type>& eqn,
85  const word& fieldName
86  ) const;
87 
88  //- Add a source term to a compressible equation
89  template<class Type>
90  void addSupType
91  (
92  const volScalarField& rho,
93  fvMatrix<Type>& eqn,
94  const word& fieldName
95  ) const;
96 
97  //- Add a source term to a phase equation
98  template<class Type>
99  void addSupType
100  (
101  const volScalarField& alpha,
102  const volScalarField& rho,
103  fvMatrix<Type>& eqn,
104  const word& fieldName
105  ) const;
106 
107 
108  //- Return source for equation with specified name and dimensions
109  template<class Type, class ... AlphaRhoFieldTypes>
111  (
113  const word& fieldName,
114  const dimensionSet& ds,
115  const AlphaRhoFieldTypes& ... alphaRhos
116  ) const;
117 
118 
119 public:
120 
121  //- Runtime type information
122  TypeName("fvModel");
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 dictionary& dict,
136  const fvMesh& mesh
137  ),
138  (name, modelType, dict, mesh)
139  );
140 
141 
142  // Static Member Functions
143 
144  //- Return the dimensions of the matrix of a source term
145  template
146  <
147  class Type,
148  class AlphaRhoFieldType,
149  class ... AlphaRhoFieldTypes
150  >
151  inline static dimensionSet sourceDims
152  (
154  const dimensionSet& ds,
155  const AlphaRhoFieldType& alphaRho,
156  const AlphaRhoFieldTypes& ... alphaRhos
157  );
158 
159  //- Return the dimensions of the matrix of a source term (base
160  // condition for the above)
161  template<class Type>
162  inline static dimensionSet sourceDims
163  (
165  const dimensionSet& ds
166  );
167 
168 
169  // Constructors
170 
171  //- Construct from components
172  fvModel
173  (
174  const word& name,
175  const word& modelType,
176  const dictionary& dict,
177  const fvMesh& mesh
178  );
179 
180  //- Return clone
181  autoPtr<fvModel> clone() const
182  {
184  return autoPtr<fvModel>(nullptr);
185  }
186 
187  //- Return pointer to new fvModel object created
188  // on the freestore from an Istream
189  class iNew
190  {
191  //- Reference to the mesh
192  const fvMesh& mesh_;
193 
194  const word& name_;
195 
196  public:
197 
199  (
200  const fvMesh& mesh,
201  const word& name
202  )
203  :
204  mesh_(mesh),
205  name_(name)
206  {}
209  {
210  // const word name(is);
211  const dictionary dict(is);
212 
213  return autoPtr<fvModel>
214  (
215  fvModel::New(name_, dict, mesh_)
216  );
217  }
218  };
219 
220 
221  // Selectors
222 
223  //- Return a reference to the selected fvModel
224  static autoPtr<fvModel> New
225  (
226  const word& name,
227  const dictionary& dict,
228  const fvMesh& mesh
229  );
230 
231 
232  //- Destructor
233  virtual ~fvModel();
234 
235 
236  // Member Functions
237 
238  // Access
239 
240  //- Return const access to the source name
241  inline const word& name() const;
242 
243  //- Return const access to the mesh database
244  inline const fvMesh& mesh() const;
245 
246  //- Return dictionary
247  inline const dictionary& coeffs() const;
248 
249 
250  // Checks
251 
252  //- Return the list of fields for which the fvModel adds source term
253  // to the transport equation
254  virtual wordList addSupFields() const;
255 
256  //- Return true if the fvModel adds a source term to the given
257  // field's transport equation
258  virtual bool addsSupToField(const word& fieldName) const;
259 
260 
261  // Sources
262 
263  //- Add a source term to an equation
265 
266  //- Add a source term to a compressible equation
268 
269  //- Add a source term to a phase equation
271 
272  //- Return source for an equation
273  template<class Type>
275  (
277  ) const;
278 
279  //- Return source for an equation with a specified name
280  template<class Type>
282  (
284  const word& fieldName
285  ) const;
286 
287  //- Return source for a compressible equation
288  template<class Type>
290  (
291  const volScalarField& rho,
293  ) const;
294 
295  //- Return source for a compressible equation with a specified name
296  template<class Type>
298  (
299  const volScalarField& rho,
301  const word& fieldName
302  ) const;
303 
304  //- Return source for a phase equation
305  template<class Type>
307  (
308  const volScalarField& alpha,
309  const volScalarField& rho,
311  ) const;
312 
313  //- Return source for a phase equation with a specified name
314  template<class Type>
316  (
317  const volScalarField& alpha,
318  const volScalarField& rho,
320  const word& fieldName
321  ) const;
322 
323 
324  // Mesh changes
325 
326  //- Prepare for mesh update
327  virtual void preUpdateMesh();
328 
329  //- Update for mesh changes
330  virtual void updateMesh(const mapPolyMesh&);
331 
332  //- Update for mesh motion
333  virtual bool movePoints();
334 
335 
336  //- Correct the fvModel
337  // e.g. solve equations, update model, for film, Lagrangian etc.
338  virtual void correct();
339 
340 
341  // IO
342 
343  //- Read source dictionary
344  virtual bool read(const dictionary& dict);
345 };
346 
347 
348 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
349 
350 } // End namespace Foam
351 
352 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
353 
354 #include "fvModelI.H"
355 
356 #ifdef NoRepository
357  #include "fvModelTemplates.C"
358 #endif
359 
360 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
361 
362 #endif
363 
364 // ************************************************************************* //
iNew(const fvMesh &mesh, const word &name)
Definition: fvModel.H:198
dictionary dict
virtual void updateMesh(const mapPolyMesh &)
Update for mesh changes.
Definition: fvModel.C:170
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:158
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
Finite volume model abstract base class.
Definition: fvModel.H:55
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModel.C:166
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
TypeName("fvModel")
Runtime type information.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Return pointer to new fvModel object created.
Definition: fvModel.H:188
Dimension set for the base types.
Definition: dimensionSet.H:120
autoPtr< fvModel > operator()(Istream &is) const
Definition: fvModel.H:207
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
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
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
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:180
Forward declarations of fvMatrix specialisations.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual void correct()
Correct the fvModel.
Definition: fvModel.C:180
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
virtual bool movePoints()
Update for mesh motion.
Definition: fvModel.C:174
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:370
Namespace for OpenFOAM.