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 -------------------------------------------------------------------------------
10 License
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 
33 namespace Foam
34 {
35  namespace SRF
36  {
37  defineTypeNameAndDebug(SRFModel, 0);
38  defineRunTimeSelectionTable(SRFModel, dictionary);
39  }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 Foam::SRF::SRFModel::SRFModel
46 (
47  const word& type,
48  const volVectorField& Urel
49 )
50 :
52  (
53  IOobject
54  (
55  "SRFProperties",
56  Urel.time().constant(),
57  Urel.db(),
58  IOobject::READ_MODIFIED,
59  IOobject::NO_WRITE,
60  IOobject::REGISTER
61  )
62  ),
63  Urel_(Urel),
64  mesh_(Urel_.mesh()),
65  origin_("origin", dimLength, get<vector>("origin")),
66  axis_(normalised(get<vector>("axis"))),
67  SRFModelCoeffs_(optionalSubDict(type + "Coeffs")),
68  omega_(dimensionedVector("omega", dimless/dimTime, Zero))
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
73 
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
81 {
82  if (regIOobject::read())
83  {
84  // Re-read origin
85  readEntry("origin", origin_);
86 
87  // Re-read axis
88  readEntry("axis", axis_);
89  axis_.normalise();
90 
91  // Re-read sub-model coeffs
92  SRFModelCoeffs_ = optionalSubDict(type() + "Coeffs");
93 
94  return true;
95  }
96 
97  return false;
98 }
99 
102 {
103  return origin_;
104 }
105 
108 {
109  return axis_;
110 }
111 
112 
114 {
115  return omega_;
116 }
117 
118 
121 {
123  (
124  "Fcoriolis",
126  (
127  2.0*omega_ ^ Urel_.internalField()
128  )
129  );
130 }
131 
132 
135 {
137  (
138  "Fcentrifugal",
140  (
141  omega_ ^ (omega_ ^ (mesh_.C().internalField() - origin_))
142  )
143  );
144 }
145 
146 
149 {
150  return Fcoriolis() + Fcentrifugal();
151 }
152 
153 
155 (
156  const vectorField& positions
157 ) const
158 {
159  return vectorField
160  (
161  omega_.value()
162  ^ (
163  (positions - origin_.value())
164  - axis_*(axis_ & (positions - origin_.value()))
165  )
166  );
167 }
168 
169 
171 {
172  const int oldLocal = volVectorField::Boundary::localConsistency;
174  tmp<volVectorField> relPos(mesh_.C() - origin_);
175 
176  auto tU = volVectorField::New
177  (
178  "Usrf",
180  (
181  omega_ ^ (relPos() - axis_*(axis_ & relPos()))
182  )
183  );
185 
186  return tU;
187 }
188 
189 
191 {
192  tmp<volVectorField> Usrf = U();
193 
194  auto tUabs = volVectorField::New
195  (
196  "Uabs",
198  Usrf
199  );
200  auto& Uabs = tUabs.ref();
201 
202  // Add SRF contribution to internal field
203  Uabs.primitiveFieldRef() += Urel_.primitiveField();
204 
205  // Add Urel boundary contributions
206  volVectorField::Boundary& Uabsbf = Uabs.boundaryFieldRef();
207  const volVectorField::Boundary& bvf = Urel_.boundaryField();
208 
209  forAll(bvf, i)
210  {
211  if (isA<SRFVelocityFvPatchVectorField>(bvf[i]))
212  {
213  // Only include relative contributions from
214  // SRFVelocityFvPatchVectorField's
215  const SRFVelocityFvPatchVectorField& UrelPatch =
216  refCast<const SRFVelocityFvPatchVectorField>(bvf[i]);
217  if (UrelPatch.relative())
218  {
219  Uabsbf[i] += Urel_.boundaryField()[i];
220  }
221  }
222  else
223  {
224  Uabsbf[i] += Urel_.boundaryField()[i];
225  }
226  }
227 
228  return tUabs;
229 }
230 
231 
232 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
defineTypeNameAndDebug(rpm, 0)
type
Types of root.
Definition: Roots.H:52
virtual bool read()
Read object.
bool relative() const
Is supplied inlet value relative to the SRF?
defineRunTimeSelectionTable(SRFModel, dictionary)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
const dimensionSet dimless
Dimensionless.
const vector & axis() const
Return the axis of rotation.
Definition: SRFModel.C:100
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
virtual bool read()
Read radiationProperties dictionary.
Definition: SRFModel.C:73
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
const dimensionedVector & omega() const
Return the angular velocity field [rad/s].
Definition: SRFModel.C:106
const dimensionedVector & origin() const
Return the origin of rotation.
Definition: SRFModel.C:94
Urel
Definition: pEqn.H:56
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< vector >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
tmp< volVectorField::Internal > Fcoriolis() const
Return the coriolis force.
Definition: SRFModel.C:113
tmp< volVectorField::Internal > Su() const
Source term component for momentum equation.
Definition: SRFModel.C:141
tmp< volVectorField > U() const
Return velocity of SRF for complete mesh.
Definition: SRFModel.C:163
tmp< volVectorField::Internal > Fcentrifugal() const
Return the centrifugal force.
Definition: SRFModel.C:127
U
Definition: pEqn.H:72
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());label pRefCell=0;scalar pRefValue=0.0;setRefCell(p, pimple.dict(), pRefCell, pRefValue);mesh.setFluxRequired(p.name());Info<< "Creating SRF model\"<< endl;autoPtr< SRF::SRFModel > SRF(SRF::SRFModel::New(Urel))
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
tmp< volVectorField > Uabs() const
Return absolute velocity for complete mesh.
Definition: SRFModel.C:183
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const Field< Type > &iField)
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions, copy of internal field...
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual ~SRFModel()
Destructor.
Definition: SRFModel.C:67
Field< vector > vectorField
Specialisation of Field<T> for vector.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
Do not request registration (bool: false)
vectorField velocity(const vectorField &positions) const
Return velocity vector from positions.
Definition: SRFModel.C:148
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
Velocity condition to be used in conjunction with the single rotating frame (SRF) model (see: SRFMode...