Loading...
Searching...
No Matches
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"
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33namespace Foam
34{
35 namespace SRF
36 {
39 }
40}
41
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
45Foam::SRF::SRFModel::SRFModel
46(
47 const word& type,
48 const volVectorField& Urel
49)
50:
52 (
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"))),
69{}
70
71
72// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
75{}
76
77
78// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79
81{
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
104}
105
108{
109 return axis_;
110}
111
112
124 "Fcoriolis",
126 (
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 auto oldConsistency = FieldBase::localBoundaryConsistency(0);
173 tmp<volVectorField> relPos(mesh_.C() - origin_);
174
175 auto tU = volVectorField::New
176 (
177 "Usrf",
179 (
180 omega_ ^ (relPos() - axis_*(axis_ & relPos()))
181 )
182 );
184
185 return tU;
186}
187
188
190{
191 tmp<volVectorField> Usrf = U();
192
193 auto tUabs = volVectorField::New
194 (
195 "Uabs",
197 Usrf
198 );
199 auto& Uabs = tUabs.ref();
200
201 // Add SRF contribution to internal field
202 Uabs.primitiveFieldRef() += Urel_.primitiveField();
203
204 // Add Urel boundary contributions
205 volVectorField::Boundary& Uabsbf = Uabs.boundaryFieldRef();
206 const volVectorField::Boundary& bvf = Urel_.boundaryField();
207
208 forAll(bvf, i)
209 {
211 {
212 // Only include relative contributions from
213 // SRFVelocityFvPatchVectorField's
214 const SRFVelocityFvPatchVectorField& UrelPatch =
216 if (UrelPatch.relative())
217 {
218 Uabsbf[i] += Urel_.boundaryField()[i];
219 }
220 }
221 else
222 {
223 Uabsbf[i] += Urel_.boundaryField()[i];
224 }
225 }
226
227 return tUabs;
228}
229
230
231// ************************************************************************* //
static int localBoundaryConsistency() noexcept
Get flag for local boundary consistency checks.
Definition Field.H:133
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())
GeometricBoundaryField< vector, fvPatchField, volMesh > Boundary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
IOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
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:64
tmp< volVectorField::Internal > Fcentrifugal() const
Return the centrifugal force.
Definition SRFModel.C:127
tmp< volVectorField > Uabs() const
Return absolute velocity for complete mesh.
Definition SRFModel.C:182
dimensionedVector omega_
Angular velocity of the frame (rad/s).
Definition SRFModel.H:98
const fvMesh & mesh_
Reference to the mesh.
Definition SRFModel.H:78
tmp< volVectorField::Internal > Su() const
Source term component for momentum equation.
Definition SRFModel.C:141
vector axis_
Axis of rotation, a direction vector which passes through the origin.
Definition SRFModel.H:88
dictionary SRFModelCoeffs_
SRF model coefficients dictionary.
Definition SRFModel.H:93
const dimensionedVector & origin() const
Return the origin of rotation.
Definition SRFModel.C:94
vectorField velocity(const vectorField &positions) const
Return velocity vector from positions.
Definition SRFModel.C:148
tmp< volVectorField::Internal > Fcoriolis() const
Return the coriolis force.
Definition SRFModel.C:113
tmp< volVectorField > U() const
Return velocity of SRF for complete mesh.
Definition SRFModel.C:163
dimensionedVector origin_
Origin of the axis.
Definition SRFModel.H:83
const vector & axis() const
Return the axis of rotation.
Definition SRFModel.C:100
const dimensionedVector & omega() const
Return the angular velocity field [rad/s].
Definition SRFModel.C:106
const volVectorField & Urel_
Reference to the relative velocity field.
Definition SRFModel.H:73
virtual ~SRFModel()
Destructor.
Definition SRFModel.C:67
virtual bool read()
Read radiationProperties dictionary.
Definition SRFModel.C:73
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
Definition dictionary.C:560
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
dictionary()
Default construct, a top-level empty dictionary.
Definition dictionary.C:68
virtual bool read()
Read object.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
dynamicFvMesh & mesh
Urel
Definition pEqn.H:56
Namespace for single rotating frame (SRF) models.
Definition rpm.C:31
Different types of constants.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Vector< scalar > vector
Definition vector.H:57
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299