angularOscillatingVelocityPointPatchVectorField.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 -------------------------------------------------------------------------------
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 
29 #include "pointPatchFields.H"
31 #include "Time.H"
32 #include "polyMesh.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace 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.lookup("axis")),
68  origin_(dict.lookup("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  }
82  else
83  {
84  p0_ = p.localPoints();
85  }
86 }
87 
88 
91 (
92  const angularOscillatingVelocityPointPatchVectorField& ptf,
93  const pointPatch& p,
94  const DimensionedField<vector, pointMesh>& iF,
95  const pointPatchFieldMapper& mapper
96 )
97 :
98  fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
99  axis_(ptf.axis_),
100  origin_(ptf.origin_),
101  angle0_(ptf.angle0_),
102  amplitude_(ptf.amplitude_),
103  omega_(ptf.omega_),
104  p0_(ptf.p0_)
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 )
131 {
133 
134  p0_.autoMap(m);
135 }
136 
137 
139 (
140  const pointPatchField<vector>& ptf,
141  const labelList& addr
142 )
143 {
145  refCast<const angularOscillatingVelocityPointPatchVectorField>(ptf);
148 
149  p0_.rmap(aOVptf.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  const pointPatch& p = this->patch();
163 
164  scalar angle = angle0_ + amplitude_*sin(omega_*t.value());
165  vector axisHat = axis_/mag(axis_);
166  vectorField p0Rel(p0_ - origin_);
167 
169  (
170  (
171  p0_
172  + p0Rel*(cos(angle) - 1)
173  + (axisHat ^ p0Rel*sin(angle))
174  + (axisHat & p0Rel)*(1 - cos(angle))*axisHat
175  - p.localPoints()
176  )/t.deltaTValue()
177  );
178 
180 }
181 
182 
184 (
185  Ostream& os
186 ) const
187 {
189  os.writeEntry("axis", axis_);
190  os.writeEntry("origin", origin_);
191  os.writeEntry("angle0", angle0_);
192  os.writeEntry("amplitude", amplitude_);
193  os.writeEntry("omega", omega_);
194  p0_.writeEntry("p0", os);
195  this->writeValueEntry(os);
196 }
197 
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
202 (
204  angularOscillatingVelocityPointPatchVectorField
205 );
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
209 } // End namespace Foam
210 
211 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const pointPatch & patch() const noexcept
Return the patch.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Foam::pointPatchFieldMapper.
virtual void write(Ostream &os) const
Write.
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:720
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Lookup type of boundary radiation properties.
Definition: lookup.H:57
Macros for easy insertion into run-time selection tables.
friend Ostream & operator(Ostream &, const Field< vector > &)
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
angularOscillatingVelocityPointPatchVectorField(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.
Vector< scalar > vector
Definition: vector.H:57
bool updated() const noexcept
True if the boundary condition has already been updated.
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
dimensionedScalar sin(const dimensionedScalar &ds)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
virtual void rmap(const pointPatchField< vector > &, const labelList &)
Reverse map the given pointPatchField onto this pointPatchField.
OBJstream os(runTime.globalPath()/outputName)
pointPatchField< vector > pointPatchVectorField
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition: Field.C:466
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:529
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:61
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
makePointPatchTypeField(pointPatchVectorField, solidBodyMotionDisplacementPointPatchVectorField)
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField & p
virtual void rmap(const pointPatchField< Type > &, const labelList &)
Reverse map the given PointPatchField onto.
const DimensionedField< vector, pointMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127