Loading...
Searching...
No Matches
SRFFreestreamVelocityFvPatchVectorField.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 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
31#include "volFields.H"
33#include "SRFModel.H"
35
36
37// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38
41(
42 const fvPatch& p,
44)
46 inletOutletFvPatchVectorField(p, iF),
47 relative_(false),
48 UInf_(Zero)
49{}
50
51
54(
56 const fvPatch& p,
58 const fvPatchFieldMapper& mapper
59)
61 inletOutletFvPatchVectorField(ptf, p, iF, mapper),
62 relative_(ptf.relative_),
63 UInf_(ptf.UInf_)
64{}
65
66
69(
70 const fvPatch& p,
72 const dictionary& dict
73)
74:
75 inletOutletFvPatchVectorField(p, iF),
76 relative_(dict.getOrDefault("relative", false)),
77 UInf_(dict.get<vector>("UInf"))
78{
79 this->phiName_ = dict.getOrDefault<word>("phi", "phi");
80 this->readValueEntry(dict, IOobjectOption::MUST_READ);
81}
82
83
86(
88)
90 inletOutletFvPatchVectorField(srfvpvf),
91 relative_(srfvpvf.relative_),
92 UInf_(srfvpvf.UInf_)
93{}
94
95
98(
101)
102:
103 inletOutletFvPatchVectorField(srfvpvf, iF),
104 relative_(srfvpvf.relative_),
105 UInf_(srfvpvf.UInf_)
106{}
107
108
109// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110
112{
113 if (updated())
114 {
115 return;
116 }
117
118 // Get reference to the SRF model
119 const SRF::SRFModel& srf =
120 db().lookupObject<SRF::SRFModel>("SRFProperties");
121
122
123 word ddtScheme
124 (
125 this->internalField().mesh()
126 .ddtScheme(this->internalField().name())
127 );
128
129 if (ddtScheme == fv::steadyStateDdtScheme<scalar>::typeName)
130 {
131 // If not relative to the SRF include the effect of the SRF
132 if (!relative_)
133 {
134 refValue() = UInf_ - srf.velocity(patch().Cf());
135 }
136 // If already relative to the SRF simply supply the inlet value
137 // as a fixed value
138 else
139 {
140 refValue() = UInf_;
141 }
142 }
143 else
144 {
145 scalar time = this->db().time().value();
146 scalar theta = time*mag(srf.omega().value());
147
148 refValue() =
149 cos(theta)*UInf_ + sin(theta)*(srf.axis() ^ UInf_)
150 - srf.velocity(patch().Cf());
151 }
152
153 // Set the inlet-outlet choice based on the direction of the freestream
154 valueFraction() = neg(refValue() & patch().Sf());
155
157}
158
159
161{
163 os.writeEntry("relative", relative_);
164 os.writeEntry("UInf", UInf_);
165 os.writeEntry("phi", this->phiName_);
168
169
170// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171
172namespace Foam
173{
175 (
177 SRFFreestreamVelocityFvPatchVectorField
178 );
179}
180
181
182// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
@ MUST_READ
Reading required.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Freestream velocity condition to be used in conjunction with the single rotating frame (SRF) model (s...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
SRFFreestreamVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Top level model for single rotating frame.
Definition SRFModel.H:64
vectorField velocity(const vectorField &positions) const
Return velocity vector from positions.
Definition SRFModel.C:148
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
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const Type & value() const noexcept
Return const reference to value.
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
SteadyState implicit/explicit ddt which returns 0.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
auto & name
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
dimensionedScalar sin(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar neg(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Vector< scalar > vector
Definition vector.H:57
fvPatchField< vector > fvPatchVectorField
dimensionedScalar cos(const dimensionedScalar &ds)
dictionary dict
Foam::surfaceFields.