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 | Copyright (C) 2011-2015 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 \*---------------------------------------------------------------------------*/
25 
26 #include "SRFModel.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33  namespace SRF
34  {
35  defineTypeNameAndDebug(SRFModel, 0);
36  defineRunTimeSelectionTable(SRFModel, dictionary);
37  }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const word& type,
46  const volVectorField& Urel
47 )
48 :
50  (
51  IOobject
52  (
53  "SRFProperties",
54  Urel.time().constant(),
55  Urel.db(),
56  IOobject::MUST_READ_IF_MODIFIED,
57  IOobject::NO_WRITE
58  )
59  ),
60  Urel_(Urel),
61  mesh_(Urel_.mesh()),
62  origin_("origin", dimLength, lookup("origin")),
63  axis_(lookup("axis")),
64  SRFModelCoeffs_(subDict(type + "Coeffs")),
65  omega_(dimensionedVector("omega", dimless/dimTime, vector::zero))
66 {
67  // Normalise the axis
68  axis_ /= mag(axis_);
69 }
70 
71 
72 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
73 
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
81 {
82  if (regIOobject::read())
83  {
84  // Re-read origin
85  lookup("origin") >> origin_;
86 
87  // Re-read axis
88  lookup("axis") >> axis_;
89  axis_ /= mag(axis_);
90 
91  // Re-read sub-model coeffs
92  SRFModelCoeffs_ = subDict(type() + "Coeffs");
93 
94  return true;
95  }
96  else
97  {
98  return false;
99  }
100 }
101 
102 
104 {
105  return origin_;
106 }
107 
108 
110 {
111  return axis_;
112 }
113 
114 
116 {
117  return omega_;
118 }
119 
120 
123 {
125  (
127  (
128  IOobject
129  (
130  "Fcoriolis",
131  mesh_.time().timeName(),
132  mesh_,
135  ),
136  2.0*omega_ ^ Urel_
137  )
138  );
139 }
140 
141 
144 {
146  (
148  (
149  IOobject
150  (
151  "Fcentrifugal",
152  mesh_.time().timeName(),
153  mesh_,
156  ),
157  omega_ ^ (omega_ ^ (mesh_.C() - origin_))
158  )
159  );
160 }
161 
162 
165 {
166  return Fcoriolis() + Fcentrifugal();
167 }
168 
169 
171 (
172  const vectorField& positions
173 ) const
174 {
175  tmp<vectorField> tfld =
176  omega_.value()
177  ^ (
178  (positions - origin_.value())
179  - axis_*(axis_ & (positions - origin_.value()))
180  );
181 
182  return tfld();
183 }
184 
185 
187 {
188  return tmp<volVectorField>
189  (
190  new volVectorField
191  (
192  IOobject
193  (
194  "Usrf",
195  mesh_.time().timeName(),
196  mesh_,
199  ),
200  omega_
201  ^ ((mesh_.C() - origin_) - axis_*(axis_ & (mesh_.C() - origin_)))
202  )
203  );
204 }
205 
206 
208 {
209  tmp<volVectorField> Usrf = U();
210 
211  tmp<volVectorField> tUabs
212  (
213  new volVectorField
214  (
215  IOobject
216  (
217  "Uabs",
218  mesh_.time().timeName(),
219  mesh_,
222  false
223  ),
224  Usrf
225  )
226  );
227 
228  // Add SRF contribution to internal field
229  tUabs().internalField() += Urel_.internalField();
230 
231  // Add Urel boundary contributions
232  const volVectorField::GeometricBoundaryField& 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  tUabs().boundaryField()[i] += Urel_.boundaryField()[i];
245  }
246  }
247  else
248  {
249  tUabs().boundaryField()[i] += Urel_.boundaryField()[i];
250  }
251  }
252 
253  return tUabs;
254 }
255 
256 
257 // ************************************************************************* //
SRF
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, simple.dict(), pRefCell, pRefValue);mesh.setFluxRequired(p.name());Info<< "Creating SRF model\n"<< endl;autoPtr< SRF::SRFModel > SRF(SRF::SRFModel::New(Urel))
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::SRF::SRFModel::~SRFModel
virtual ~SRFModel()
Destructor.
Definition: SRFModel.C:74
Foam::SRF::SRFModel::Uabs
tmp< volVectorField > Uabs() const
Return absolute velocity for complete mesh.
Definition: SRFModel.C:207
SRFVelocityFvPatchVectorField.H
Foam::regIOobject::read
virtual bool read()
Read object.
Definition: regIOobjectRead.C:171
Foam::SRF::SRFModel::SRFModel
SRFModel(const SRFModel &)
Disallow default bitwise copy construct.
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::SRF::SRFModel::velocity
vectorField velocity(const vectorField &positions) const
Return velocity vector from positions.
Definition: SRFModel.C:171
Foam::SRF::defineRunTimeSelectionTable
defineRunTimeSelectionTable(SRFModel, dictionary)
Foam::SRF::SRFModel::Fcoriolis
tmp< DimensionedField< vector, volMesh > > Fcoriolis() const
Return the coriolis force.
Definition: SRFModel.C:122
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
U
U
Definition: pEqn.H:46
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:48
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::SRF::SRFModel::axis
const vector & axis() const
Return the axis of rotation.
Definition: SRFModel.C:109
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::SRF::SRFModel::read
virtual bool read()
Read radiationProperties dictionary.
Definition: SRFModel.C:80
Foam::SRF::SRFModel::Su
tmp< DimensionedField< vector, volMesh > > Su() const
Source term component for momentum equation.
Definition: SRFModel.C:164
Urel
Urel
Definition: pEqn.H:45
Foam::SRF::SRFModel::omega
const dimensionedVector & omega() const
Return the angular velocity field [rad/s].
Definition: SRFModel.C:115
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
SRFModel.H
Foam::SRF::SRFModel::origin
const dimensionedVector & origin() const
Return the origin of rotation.
Definition: SRFModel.C:103
Foam::Vector< scalar >
Foam::SRFVelocityFvPatchVectorField
Velocity condition to be used in conjunction with the single rotating frame (SRF) model (see: SRFMode...
Definition: SRFVelocityFvPatchVectorField.H:130
Foam::SRF::SRFModel::Fcentrifugal
tmp< DimensionedField< vector, volMesh > > Fcentrifugal() const
Return the centrifugal force.
Definition: SRFModel.C:143
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::GeometricField::GeometricBoundaryField
Definition: GeometricField.H:105
Foam::SRF::defineTypeNameAndDebug
defineTypeNameAndDebug(rpm, 0)
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
Foam::SRF::SRFModel::U
tmp< volVectorField > U() const
Return velocity of SRF for complete mesh.
Definition: SRFModel.C:186
Foam::SRFVelocityFvPatchVectorField::relative
const Switch & relative() const
Return const access to the relative flag.
Definition: SRFVelocityFvPatchVectorField.H:216
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress