Loading...
Searching...
No Matches
angularOscillatingDisplacementPointPatchVectorField.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-2016 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
29#include "pointPatchFields.H"
31#include "Time.H"
32#include "polyMesh.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38
39// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40
43(
44 const pointPatch& p,
46)
47:
49 axis_(Zero),
50 origin_(Zero),
51 angle0_(0.0),
52 amplitude_(0.0),
53 omega_(0.0),
54 p0_(p.localPoints())
55{}
56
57
60(
61 const pointPatch& p,
63 const dictionary& dict
64)
65:
67 axis_(dict.get<vector>("axis")),
68 origin_(dict.get<vector>("origin")),
69 angle0_(dict.get<scalar>("angle0")),
70 amplitude_(dict.get<scalar>("amplitude")),
71 omega_(dict.get<scalar>("omega"))
72{
73 if (!dict.found("value"))
74 {
75 updateCoeffs();
76 }
77
78 if (dict.found("p0"))
79 {
80 p0_ = vectorField("p0", dict , p.size());
81 }
99 axis_(ptf.axis_),
100 origin_(ptf.origin_),
101 angle0_(ptf.angle0_),
102 amplitude_(ptf.amplitude_),
103 omega_(ptf.omega_),
104 p0_(ptf.p0_, mapper)
105{}
106
107
110(
113)
114:
116 axis_(ptf.axis_),
117 origin_(ptf.origin_),
118 angle0_(ptf.angle0_),
119 amplitude_(ptf.amplitude_),
120 omega_(ptf.omega_),
121 p0_(ptf.p0_)
122{}
123
124
125// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
126
128(
129 const pointPatchFieldMapper& m
130)
133
134 p0_.autoMap(m);
135}
136
137
139(
140 const pointPatchField<vector>& ptf,
141 const labelList& addr
142)
143{
148
149 p0_.rmap(aODptf.p0_, addr);
150}
151
152
154{
155 if (this->updated())
156 {
157 return;
158 }
159
160 const polyMesh& mesh = this->internalField().mesh()();
161 const Time& t = mesh.time();
162
163 scalar angle = angle0_ + amplitude_*sin(omega_*t.value());
164 vector axisHat = axis_/mag(axis_);
165 vectorField p0Rel(p0_ - origin_);
166
167 vectorField::operator=
168 (
169 p0Rel*(cos(angle) - 1)
170 + (axisHat ^ p0Rel*sin(angle))
171 + (axisHat & p0Rel)*(1 - cos(angle))*axisHat
172 );
173
175}
176
177
179(
180 Ostream& os
181) const
182{
184 os.writeEntry("axis", axis_);
185 os.writeEntry("origin", origin_);
186 os.writeEntry("angle0", angle0_);
187 os.writeEntry("amplitude", amplitude_);
188 os.writeEntry("omega", omega_);
189 p0_.writeEntry("p0", os);
190 this->writeValueEntry(os);
191}
192
193
194// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195
197(
199 angularOscillatingDisplacementPointPatchVectorField
200);
201
202// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203
204} // End namespace Foam
205
206// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
if(patchID !=-1)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const noexcept
Return const reference to mesh.
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition Field.C:465
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition Field.C:754
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition Field.C:528
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
bool get(const label i) const
Definition UList.H:868
angularOscillatingDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void rmap(const pointPatchField< vector > &, const labelList &)
Reverse map the given pointPatchField onto this pointPatchField.
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 FixedValue boundary condition for pointField.
fixedValuePointPatchField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
bool updated() const noexcept
True if the boundary condition has already been updated.
Foam::pointPatchFieldMapper.
Abstract base class for point-mesh patch fields.
virtual void write(Ostream &os) const
Write.
const DimensionedField< vector, pointMesh > & internalField() const noexcept
Basic pointPatch represents a set of points from the mesh.
Definition pointPatch.H:67
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedScalar sin(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
pointPatchField< vector > pointPatchVectorField
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Vector< scalar > vector
Definition vector.H:57
dimensionedScalar cos(const dimensionedScalar &ds)
#define makePointPatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete pointPatchField type and add to run-time tables Example, (pointPatchScalarField,...
dictionary dict