SRFModel.C
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-2017 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
26\*---------------------------------------------------------------------------*/
27
28#include "SRFModel.H"
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33namespace Foam
34{
35 namespace SRF
36 {
39 }
40}
41
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
46(
47 const word& type,
48 const volVectorField& Urel
49)
50:
52 (
54 (
55 "SRFProperties",
56 Urel.time().constant(),
57 Urel.db(),
58 IOobject::MUST_READ_IF_MODIFIED,
59 IOobject::NO_WRITE
60 )
61 ),
62 Urel_(Urel),
63 mesh_(Urel_.mesh()),
64 origin_("origin", dimLength, get<vector>("origin")),
65 axis_(normalised(get<vector>("axis"))),
66 SRFModelCoeffs_(optionalSubDict(type + "Coeffs")),
67 omega_(dimensionedVector("omega", dimless/dimTime, Zero))
68{}
69
70
71// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
72
74{}
75
76
77// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78
80{
82 {
83 // Re-read origin
84 readEntry("origin", origin_);
85
86 // Re-read axis
87 readEntry("axis", axis_);
88 axis_.normalise();
89
90 // Re-read sub-model coeffs
91 SRFModelCoeffs_ = optionalSubDict(type() + "Coeffs");
92
93 return true;
94 }
95
96 return false;
97}
98
99
101{
102 return origin_;
103}
104
105
107{
108 return axis_;
109}
110
111
113{
114 return omega_;
115}
116
117
120{
122 (
124 (
126 (
127 "Fcoriolis",
128 mesh_.time().timeName(),
129 mesh_,
132 ),
133 2.0*omega_ ^ Urel_
134 )
135 );
136}
137
138
141{
143 (
145 (
147 (
148 "Fcentrifugal",
149 mesh_.time().timeName(),
150 mesh_,
153 ),
154 omega_ ^ (omega_ ^ (mesh_.C() - origin_))
155 )
156 );
157}
158
159
162{
163 return Fcoriolis() + Fcentrifugal();
164}
165
166
168(
169 const vectorField& positions
170) const
171{
172 tmp<vectorField> tfld =
173 omega_.value()
174 ^ (
175 (positions - origin_.value())
176 - axis_*(axis_ & (positions - origin_.value()))
177 );
178
179 return tfld();
180}
181
182
184{
186 (
188 (
190 (
191 "Usrf",
192 mesh_.time().timeName(),
193 mesh_,
196 ),
197 omega_
198 ^ ((mesh_.C() - origin_) - axis_*(axis_ & (mesh_.C() - origin_)))
199 )
200 );
201}
202
203
205{
206 tmp<volVectorField> Usrf = U();
207
209 (
211 (
213 (
214 "Uabs",
215 mesh_.time().timeName(),
216 mesh_,
219 false
220 ),
221 Usrf
222 )
223 );
224
225 volVectorField& Uabs = tUabs.ref();
226
227 // Add SRF contribution to internal field
228 Uabs.primitiveFieldRef() += Urel_.primitiveField();
229
230 // Add Urel boundary contributions
232 const volVectorField::Boundary& bvf = Urel_.boundaryField();
233
234 forAll(bvf, i)
235 {
236 if (isA<SRFVelocityFvPatchVectorField>(bvf[i]))
237 {
238 // Only include relative contributions from
239 // SRFVelocityFvPatchVectorField's
240 const SRFVelocityFvPatchVectorField& UrelPatch =
241 refCast<const SRFVelocityFvPatchVectorField>(bvf[i]);
242 if (UrelPatch.relative())
243 {
244 Uabsbf[i] += Urel_.boundaryField()[i];
245 }
246 }
247 else
248 {
249 Uabsbf[i] += Urel_.boundaryField()[i];
250 }
251 }
252
253 return tUabs;
254}
255
256
257// ************************************************************************* //
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Velocity condition to be used in conjunction with the single rotating frame (SRF) model (see: SRFMode...
bool relative() const
Is supplied inlet value relative to the SRF?
Top level model for single rotating frame.
Definition: SRFModel.H:67
tmp< volVectorField::Internal > Fcentrifugal() const
Return the centrifugal force.
Definition: SRFModel.C:140
tmp< volVectorField > Uabs() const
Return absolute velocity for complete mesh.
Definition: SRFModel.C:204
tmp< volVectorField::Internal > Su() const
Source term component for momentum equation.
Definition: SRFModel.C:161
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
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
virtual ~SRFModel()
Destructor.
Definition: SRFModel.C:73
virtual bool read()
Read radiationProperties dictionary.
Definition: SRFModel.C:79
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
constant condensation/saturation model.
virtual bool read()
Read object.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
dynamicFvMesh & mesh
Urel
Definition: pEqn.H:56
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:680
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
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
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
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 forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333