rpm.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 OpenFOAM Foundation
9 Copyright (C) 2018 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::SRF::rpm
29
30Description
31 Basic SRF model whereby angular velocity is specified in terms of
32 a (global) axis and revolutions-per-minute [rpm]
33
34SourceFiles
35 rpm.C
36
37\*---------------------------------------------------------------------------*/
38
39#ifndef SRF_rpm_H
40#define SRF_rpm_H
41
42#include "SRFModel.H"
43
44// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45
46namespace Foam
47{
48namespace SRF
49{
50
51/*---------------------------------------------------------------------------*\
52 Class rpm Declaration
53\*---------------------------------------------------------------------------*/
55class rpm
56:
57 public SRFModel
58{
59
60 // Private data
61
62 //- Revolutions per minute
63 scalar rpm_;
64
65
66 // Private Member Functions
67
68 //- No copy construct
69 rpm(const rpm&) = delete;
70
71 //- No copy assignment
72 void operator=(const rpm&) = delete;
73
74
75public:
76
77 //- Runtime type information
78 TypeName("rpm");
79
80
81 // Constructors
82
83 //- Construct from components
84 rpm(const volVectorField& U);
85
86
87 //- Destructor
88 ~rpm() = default;
89
90
91 // Member functions
92
93 //- Read coefficients
94 bool read();
95};
96
97
98// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
99
100} // End namespace SRF
101} // End namespace Foam
102
103// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
104
105#endif
106
107// ************************************************************************* //
Top level model for single rotating frame.
Definition: SRFModel.H:67
tmp< volVectorField > U() const
Return velocity of SRF for complete mesh.
Definition: SRFModel.C:183
Basic SRF model whereby angular velocity is specified in terms of a (global) axis and revolutions-per...
Definition: rpm.H:57
TypeName("rpm")
Runtime type information.
~rpm()=default
Destructor.
bool read()
Read coefficients.
Definition: rpm.C:64
Namespace for OpenFOAM.
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