SRFModel.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) 2011-2020 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 Namespace
25  Foam::SRF
26 
27 Description
28  Namespace for single rotating frame (SRF) models
29 
30 Class
31  Foam::SRF::SRFModel
32 
33 Description
34  Top level model for single rotating frame
35  - Steady state only - no time derivatives included
36 
37 SourceFiles
38  SRFModel.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef SRFModel_H
43 #define SRFModel_H
44 
45 #include "IOdictionary.H"
46 #include "autoPtr.H"
47 #include "runTimeSelectionTables.H"
48 #include "fvMesh.H"
49 #include "volFields.H"
50 #include "vectorField.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 namespace SRF
57 {
58 
59 /*---------------------------------------------------------------------------*\
60  Class SRFModel Declaration
61 \*---------------------------------------------------------------------------*/
62 
63 class SRFModel
64 :
65  public IOdictionary
66 {
67 
68 protected:
69 
70  // Protected data
71 
72  //- Reference to the relative velocity field
73  const volVectorField& Urel_;
74 
75  //- Reference to the mesh
76  const fvMesh& mesh_;
77 
78  //- Origin of the axis
80 
81  //- Axis of rotation, a direction vector which passes through the origin
82  vector axis_;
83 
84  //- SRF model coefficients dictionary
86 
87  //- Angular velocity of the frame (rad/s)
89 
90 
91 public:
92 
93  //- Runtime type information
94  TypeName("SRFModel");
95 
96 
97  // Declare runtime constructor selection table
98 
100  (
101  autoPtr,
102  SRFModel,
103  dictionary,
104  (
105  const volVectorField& Urel
106  ),
107  (Urel)
108  );
109 
110 
111  // Constructors
112 
113  //- Construct from components
114  SRFModel
115  (
116  const word& type,
117  const volVectorField& Urel
118  );
119 
120  //- Disallow default bitwise copy construction
121  SRFModel(const SRFModel&) = delete;
122 
123 
124  // Selectors
125 
126  //- Return a reference to the selected SRF model
127  static autoPtr<SRFModel> New
128  (
129  const volVectorField& Urel
130  );
131 
132 
133  //- Destructor
134  virtual ~SRFModel();
135 
136 
137  // Member Functions
138 
139  // Edit
140 
141  //- Read SRFProperties dictionary
142  virtual bool read();
143 
144 
145  // Access
146 
147  //- Return the origin of rotation
148  const dimensionedVector& origin() const;
149 
150  //- Return the axis of rotation
151  const vector& axis() const;
152 
153  //- Return the angular velocity field [rad/s]
154  const dimensionedVector& omega() const;
155 
156  //- Return the coriolis force
158 
159  //- Return the centrifugal force
161 
162  //- Source term component for momentum equation
164 
165  //- Return velocity vector from positions
166  vectorField velocity(const vectorField& positions) const;
167 
168  //- Return velocity of SRF for complete mesh
169  tmp<volVectorField> U() const;
170 
171  //- Return absolute velocity for complete mesh
172  tmp<volVectorField> Uabs() const;
173 
174 
175  // Member Operators
176 
177  //- Disallow default bitwise assignment
178  void operator=(const SRFModel&) = delete;
179 };
180 
181 
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
183 
184 } // End namespace SRF
185 } // End namespace Foam
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 #endif
190 
191 // ************************************************************************* //
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
void operator=(const SRFModel &)=delete
Disallow default bitwise assignment.
const vector & axis() const
Return the axis of rotation.
Definition: SRFModel.C:109
const fvMesh & mesh_
Reference to the mesh.
Definition: SRFModel.H:75
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
dimensionedVector omega_
Angular velocity of the frame (rad/s)
Definition: SRFModel.H:87
static autoPtr< SRFModel > New(const volVectorField &Urel)
Return a reference to the selected SRF model.
Definition: SRFModelNew.C:31
virtual bool read()
Read SRFProperties dictionary.
Definition: SRFModel.C:80
Info<< "Reading field p\"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field Urel\"<< endl;volVectorField Urel(IOobject("Urel", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading/calculating face flux field phi\"<< endl;surfaceScalarField phi(IOobject("phi", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), linearInterpolate(Urel) &mesh.Sf());pressureReference pressureReference(p, pimple.dict());mesh.setFluxRequired(p.name());Info<< "Creating SRF model\"<< endl;autoPtr< SRF::SRFModel > SRF(SRF::SRFModel::New(Urel))
dictionary SRFModelCoeffs_
SRF model coefficients dictionary.
Definition: SRFModel.H:84
A class for handling words, derived from string.
Definition: word.H:59
vector axis_
Axis of rotation, a direction vector which passes through the origin.
Definition: SRFModel.H:81
const dimensionedVector & omega() const
Return the angular velocity field [rad/s].
Definition: SRFModel.C:115
declareRunTimeSelectionTable(autoPtr, SRFModel, dictionary,(const volVectorField &Urel),(Urel))
const dimensionedVector & origin() const
Return the origin of rotation.
Definition: SRFModel.C:103
Urel
Definition: pEqn.H:60
dimensionedVector origin_
Origin of the axis.
Definition: SRFModel.H:78
tmp< volVectorField::Internal > Fcoriolis() const
Return the coriolis force.
Definition: SRFModel.C:122
const volVectorField & Urel_
Reference to the relative velocity field.
Definition: SRFModel.H:72
tmp< volVectorField::Internal > Su() const
Source term component for momentum equation.
Definition: SRFModel.C:144
tmp< volVectorField > U() const
Return velocity of SRF for complete mesh.
Definition: SRFModel.C:166
tmp< volVectorField::Internal > Fcentrifugal() const
Return the centrifugal force.
Definition: SRFModel.C:133
SRFModel(const word &type, const volVectorField &Urel)
Construct from components.
Definition: SRFModel.C:44
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
tmp< volVectorField > Uabs() const
Return absolute velocity for complete mesh.
Definition: SRFModel.C:176
TypeName("SRFModel")
Runtime type information.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual ~SRFModel()
Destructor.
Definition: SRFModel.C:74
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: PtrList.H:53
Top level model for single rotating frame.
Definition: SRFModel.H:62
vectorField velocity(const vectorField &positions) const
Return velocity vector from positions.
Definition: SRFModel.C:151
Namespace for OpenFOAM.