SRFModel.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Namespace
27 Foam::SRF
28
29Description
30 Namespace for single rotating frame (SRF) models
31
32Class
33 Foam::SRF::SRFModel
34
35Description
36 Top level model for single rotating frame
37 - Steady state only - no time derivatives included
38
39SourceFiles
40 SRFModel.C
41
42\*---------------------------------------------------------------------------*/
43
44#ifndef SRFModel_H
45#define SRFModel_H
46
47#include "IOdictionary.H"
48#include "autoPtr.H"
50#include "fvMesh.H"
51#include "volFields.H"
52#include "vectorField.H"
53
54// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55
56namespace Foam
57{
58namespace SRF
59{
60
61/*---------------------------------------------------------------------------*\
62 Class SRFModel Declaration
63\*---------------------------------------------------------------------------*/
65class SRFModel
66:
67 public IOdictionary
68{
69
70protected:
71
72 // Protected data
73
74 //- Reference to the relative velocity field
75 const volVectorField& Urel_;
76
77 //- Reference to the mesh
78 const fvMesh& mesh_;
79
80 //- Origin of the axis
82
83 //- Axis of rotation, a direction vector which passes through the origin
85
86 //- SRF model coefficients dictionary
88
89 //- Angular velocity of the frame (rad/s)
91
92
93private:
94
95 // Private Member Functions
96
97 //- No copy construct
98 SRFModel(const SRFModel&) = delete;
99
100 //- No copy assignment
101 void operator=(const SRFModel&) = delete;
102
103
104public:
105
106 //- Runtime type information
107 TypeName("SRFModel");
108
109
110 // Declare runtime constructor selection table
113 (
114 autoPtr,
115 SRFModel,
117 (
118 const volVectorField& Urel
119 ),
120 (Urel)
121 );
122
123
124 // Constructors
125
126 //- Construct from components
128 (
129 const word& type,
130 const volVectorField& Urel
131 );
132
133
134 // Selectors
135
136 //- Return a reference to the selected SRF model
138 (
139 const volVectorField& Urel
140 );
141
142
143 //- Destructor
144 virtual ~SRFModel();
145
146
147 // Member Functions
148
149 // Edit
150
151 //- Read radiationProperties dictionary
152 virtual bool read();
153
154
155 // Access
156
157 //- Return the origin of rotation
158 const dimensionedVector& origin() const;
159
160 //- Return the axis of rotation
161 const vector& axis() const;
162
163 //- Return the angular velocity field [rad/s]
164 const dimensionedVector& omega() const;
165
166 //- Return the coriolis force
168
169 //- Return the centrifugal force
171
172 //- Source term component for momentum equation
174
175 //- Return velocity vector from positions
176 vectorField velocity(const vectorField& positions) const;
177
178 //- Return velocity of SRF for complete mesh
179 tmp<volVectorField> U() const;
180
181 //- Return absolute velocity for complete mesh
183};
184
185
186// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
187
188} // End namespace SRF
189} // End namespace Foam
190
191// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192
193#endif
194
195// ************************************************************************* //
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Top level model for single rotating frame.
Definition: SRFModel.H:67
tmp< volVectorField::Internal > Fcentrifugal() const
Return the centrifugal force.
Definition: SRFModel.C:140
TypeName("SRFModel")
Runtime type information.
tmp< volVectorField > Uabs() const
Return absolute velocity for complete mesh.
Definition: SRFModel.C:204
dimensionedVector omega_
Angular velocity of the frame (rad/s)
Definition: SRFModel.H:89
const fvMesh & mesh_
Reference to the mesh.
Definition: SRFModel.H:77
declareRunTimeSelectionTable(autoPtr, SRFModel, dictionary,(const volVectorField &Urel),(Urel))
tmp< volVectorField::Internal > Su() const
Source term component for momentum equation.
Definition: SRFModel.C:161
vector axis_
Axis of rotation, a direction vector which passes through the origin.
Definition: SRFModel.H:83
dictionary SRFModelCoeffs_
SRF model coefficients dictionary.
Definition: SRFModel.H:86
const dimensionedVector & origin() const
Return the origin of rotation.
Definition: SRFModel.C:100
vectorField velocity(const vectorField &positions) const
Return velocity vector from positions.
Definition: SRFModel.C:168
tmp< volVectorField::Internal > Fcoriolis() const
Return the coriolis force.
Definition: SRFModel.C:119
tmp< volVectorField > U() const
Return velocity of SRF for complete mesh.
Definition: SRFModel.C:183
dimensionedVector origin_
Origin of the axis.
Definition: SRFModel.H:80
const vector & axis() const
Return the axis of rotation.
Definition: SRFModel.C:106
const dimensionedVector & omega() const
Return the angular velocity field [rad/s].
Definition: SRFModel.C:112
const volVectorField & Urel_
Reference to the relative velocity field.
Definition: SRFModel.H:74
virtual ~SRFModel()
Destructor.
Definition: SRFModel.C:73
static autoPtr< SRFModel > New(const volVectorField &Urel)
Return a reference to the selected SRF model.
Definition: SRFModelNew.C:34
virtual bool read()
Read radiationProperties dictionary.
Definition: SRFModel.C:79
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
Urel
Definition: pEqn.H:56
Namespace for OpenFOAM.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes)
Info<< "Reading field p\n"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field Urel\n"<< endl;volVectorField Urel(IOobject("Urel", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading/calculating face flux field phi\n"<< endl;surfaceScalarField phi(IOobject("phi", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), linearInterpolate(Urel) &mesh.Sf());label pRefCell=0;scalar pRefValue=0.0;setRefCell(p, pimple.dict(), pRefCell, pRefValue);mesh.setFluxRequired(p.name());Info<< "Creating SRF model\n"<< endl;autoPtr< SRF::SRFModel > SRF(SRF::SRFModel::New(Urel))
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73