oscillatingVelocityPointPatchVectorField.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  amplitude_(Zero),
50  omega_(0.0),
51  p0_(p.localPoints())
52 {}
53 
54 
57 (
58  const pointPatch& p,
60  const dictionary& dict
61 )
62 :
64  amplitude_(dict.lookup("amplitude")),
65  omega_(dict.get<scalar>("omega"))
66 {
67  if (!dict.found("value"))
68  {
69  updateCoeffs();
70  }
71 
72  if (dict.found("p0"))
73  {
74  p0_ = vectorField("p0", dict , p.size());
75  }
76  else
77  {
78  p0_ = p.localPoints();
79  }
80 }
81 
82 
85 (
86  const oscillatingVelocityPointPatchVectorField& ptf,
87  const pointPatch& p,
88  const DimensionedField<vector, pointMesh>& iF,
89  const pointPatchFieldMapper& mapper
90 )
91 :
92  fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
93  amplitude_(ptf.amplitude_),
94  omega_(ptf.omega_),
95  p0_(ptf.p0_, mapper)
96 {}
97 
98 
101 (
104 )
105 :
107  amplitude_(ptf.amplitude_),
108  omega_(ptf.omega_),
109  p0_(ptf.p0_)
110 {}
111 
112 
113 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 
116 (
117  const pointPatchFieldMapper& m
118 )
119 {
121 
122  p0_.autoMap(m);
123 }
124 
125 
127 (
128  const pointPatchField<vector>& ptf,
129  const labelList& addr
130 )
131 {
133  refCast<const oscillatingVelocityPointPatchVectorField>(ptf);
136 
137  p0_.rmap(oVptf.p0_, addr);
138 }
139 
140 
142 {
143  if (this->updated())
144  {
145  return;
146  }
147 
148  const polyMesh& mesh = this->internalField().mesh()();
149  const Time& t = mesh.time();
150  const pointPatch& p = this->patch();
151 
153  (
154  (p0_ + amplitude_*sin(omega_*t.value()) - p.localPoints())
155  /t.deltaTValue()
156  );
157 
159 }
160 
161 
163 {
165  os.writeEntry("amplitude", amplitude_);
166  os.writeEntry("omega", omega_);
167  p0_.writeEntry("p0", os);
168  this->writeValueEntry(os);
169 }
170 
171 
172 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173 
175 (
177  oscillatingVelocityPointPatchVectorField
178 );
179 
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 
182 } // End namespace Foam
183 
184 // ************************************************************************* //
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.
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
oscillatingVelocityPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
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
virtual void rmap(const pointPatchField< vector > &, const labelList &)
Reverse map the given pointPatchField onto this pointPatchField.
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...
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)
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField & p
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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