SRFModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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(),
58  )
59  ),
60  Urel_(Urel),
61  mesh_(Urel_.mesh()),
62  origin_("origin", dimLength, lookup("origin")),
63  axis_(lookup("axis")),
64  SRFModelCoeffs_(optionalSubDict(type + "Coeffs")),
65  omega_(dimensionedVector("omega", dimless/dimTime, 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_ = optionalSubDict(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  (
126  "Fcoriolis",
127  2.0*omega_ ^ Urel_
128  );
129 }
130 
131 
134 {
136  (
137  "Fcentrifugal",
138  omega_ ^ (omega_ ^ (mesh_.C() - origin_))
139  );
140 }
141 
142 
145 {
146  return Fcoriolis() + Fcentrifugal();
147 }
148 
149 
151 (
152  const vectorField& positions
153 ) const
154 {
155  tmp<vectorField> tfld =
156  omega_.value()
157  ^ (
158  (positions - origin_.value())
159  - axis_*(axis_ & (positions - origin_.value()))
160  );
161 
162  return tfld();
163 }
164 
165 
167 {
168  return volVectorField::New
169  (
170  "Usrf",
171  omega_ ^ ((mesh_.C() - origin_) - axis_*(axis_ & (mesh_.C() - origin_)))
172  );
173 }
174 
175 
177 {
178  tmp<volVectorField> Usrf = U();
179 
180  tmp<volVectorField> tUabs
181  (
182  volVectorField::New("Uabs", Usrf)
183  );
184 
185  volVectorField& Uabs = tUabs.ref();
186 
187  // Add SRF contribution to internal field
188  Uabs.primitiveFieldRef() += Urel_.primitiveField();
189 
190  // Add Urel boundary contributions
192  const volVectorField::Boundary& bvf = Urel_.boundaryField();
193 
194  forAll(bvf, i)
195  {
196  if (isA<SRFVelocityFvPatchVectorField>(bvf[i]))
197  {
198  // Only include relative contributions from
199  // SRFVelocityFvPatchVectorField's
200  const SRFVelocityFvPatchVectorField& UrelPatch =
201  refCast<const SRFVelocityFvPatchVectorField>(bvf[i]);
202  if (UrelPatch.relative())
203  {
204  Uabsbf[i] += Urel_.boundaryField()[i];
205  }
206  }
207  else
208  {
209  Uabsbf[i] += Urel_.boundaryField()[i];
210  }
211  }
212 
213  return tUabs;
214 }
215 
216 
217 // ************************************************************************* //
defineTypeNameAndDebug(rpm, 0)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
const Switch & relative() const
Return const access to the relative flag.
U
Definition: pEqn.H:72
virtual bool read()
Read object.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< vector >> &)
Return a temporary field constructed from name,.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const dimensionSet dimless
defineRunTimeSelectionTable(SRFModel, dictionary)
const vector & axis() const
Return the axis of rotation.
Definition: SRFModel.C:109
const dimensionSet dimLength
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
const dimensionSet dimTime
virtual bool read()
Read SRFProperties dictionary.
Definition: SRFModel.C:80
stressControl lookup("compactNormalStress") >> compactNormalStress
A class for handling words, derived from string.
Definition: word.H:59
const word & constant() const
Return constant name.
Definition: TimePaths.H:123
const dimensionedVector & omega() const
Return the angular velocity field [rad/s].
Definition: SRFModel.C:115
static const zero Zero
Definition: zero.H:97
const dimensionedVector & origin() const
Return the origin of rotation.
Definition: SRFModel.C:103
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
tmp< volVectorField::Internal > Fcoriolis() const
Return the coriolis force.
Definition: SRFModel.C:122
tmp< volVectorField::Internal > Su() const
Source term component for momentum equation.
Definition: SRFModel.C:144
tmp< volVectorField > U() const
Return velocity of SRF for complete mesh.
Definition: SRFModel.C:166
tmp< volVectorField::Internal > Fcentrifugal() const
Return the centrifugal force.
Definition: SRFModel.C:133
SRFModel(const word &type, const volVectorField &Urel)
Construct from components.
Definition: SRFModel.C:44
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const Time & time() const
Return time.
Definition: IOobject.C:318
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());pressureReference pressureReference(p, pimple.dict());mesh.schemes().setFluxRequired(p.name());Info<< "Creating SRF model\"<< endl;autoPtr< SRF::SRFModel > SRF(SRF::SRFModel::New(Urel))
tmp< volVectorField > Uabs() const
Return absolute velocity for complete mesh.
Definition: SRFModel.C:176
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual ~SRFModel()
Destructor.
Definition: SRFModel.C:74
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
vectorField velocity(const vectorField &positions) const
Return velocity vector from positions.
Definition: SRFModel.C:151
Namespace for OpenFOAM.
Velocity condition to be used in conjunction with the single rotating frame (SRF) model (see: SRFMode...