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::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 {
81  if (regIOobject::read())
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 
107 {
108  return axis_;
109 }
110 
111 
113 {
114  return omega_;
115 }
116 
117 
120 {
122  (
124  (
125  IOobject
126  (
127  "Fcoriolis",
128  mesh_.time().timeName(),
129  mesh_,
132  ),
133  2.0*omega_ ^ Urel_
134  )
135  );
136 }
137 
138 
141 {
143  (
145  (
146  IOobject
147  (
148  "Fcentrifugal",
149  mesh_.time().timeName(),
150  mesh_,
153  ),
154  omega_ ^ (omega_ ^ (mesh_.C().internalField() - 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 {
185  const int oldLocal = volVectorField::Boundary::localConsistency;
187  tmp<volVectorField> relPos(mesh_.C() - origin_);
188 
190  (
191  new volVectorField
192  (
193  IOobject
194  (
195  "Usrf",
196  mesh_.time().timeName(),
197  mesh_,
200  ),
201  omega_ ^ (relPos() - axis_*(axis_ & relPos()))
202  )
203  );
205 
206  return tU;
207 }
208 
209 
211 {
212  tmp<volVectorField> Usrf = U();
213 
214  tmp<volVectorField> tUabs
215  (
216  new volVectorField
217  (
218  IOobject
219  (
220  "Uabs",
221  mesh_.time().timeName(),
222  mesh_,
226  ),
227  Usrf
228  )
229  );
230 
231  volVectorField& Uabs = tUabs.ref();
232 
233  // Add SRF contribution to internal field
234  Uabs.primitiveFieldRef() += Urel_.primitiveField();
235 
236  // Add Urel boundary contributions
237  volVectorField::Boundary& Uabsbf = Uabs.boundaryFieldRef();
238  const volVectorField::Boundary& bvf = Urel_.boundaryField();
239 
240  forAll(bvf, i)
241  {
242  if (isA<SRFVelocityFvPatchVectorField>(bvf[i]))
243  {
244  // Only include relative contributions from
245  // SRFVelocityFvPatchVectorField's
246  const SRFVelocityFvPatchVectorField& UrelPatch =
247  refCast<const SRFVelocityFvPatchVectorField>(bvf[i]);
248  if (UrelPatch.relative())
249  {
250  Uabsbf[i] += Urel_.boundaryField()[i];
251  }
252  }
253  else
254  {
255  Uabsbf[i] += Urel_.boundaryField()[i];
256  }
257  }
258 
259  return tUabs;
260 }
261 
262 
263 // ************************************************************************* //
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?
Ignore writing from objectRegistry::writeObject()
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.
DimensionedField< vector, volMesh > Internal
The internal field type from which this GeometricField is derived.
const vector & axis() const
Return the axis of rotation.
Definition: SRFModel.C:99
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:72
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:105
const dimensionedVector & origin() const
Return the origin of rotation.
Definition: SRFModel.C:93
Urel
Definition: pEqn.H:56
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:112
tmp< volVectorField::Internal > Su() const
Source term component for momentum equation.
Definition: SRFModel.C:154
tmp< volVectorField > U() const
Return velocity of SRF for complete mesh.
Definition: SRFModel.C:176
tmp< volVectorField::Internal > Fcentrifugal() const
Return the centrifugal force.
Definition: SRFModel.C:133
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:203
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual ~SRFModel()
Destructor.
Definition: SRFModel.C:66
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
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:172
Do not request registration (bool: false)
vectorField velocity(const vectorField &positions) const
Return velocity vector from positions.
Definition: SRFModel.C:161
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...