oscillatingVelocityPointPatchVectorField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "pointPatchFields.H"
29 #include "Time.H"
30 #include "polyMesh.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
41 (
42  const pointPatch& p,
44 )
45 :
47  amplitude_(Zero),
48  omega_(0.0),
49  p0_(p.localPoints())
50 {}
51 
52 
55 (
56  const pointPatch& p,
58  const dictionary& dict
59 )
60 :
62  amplitude_(dict.lookup("amplitude")),
63  omega_(dict.lookup<scalar>("omega"))
64 {
65  if (!dict.found("value"))
66  {
67  updateCoeffs();
68  }
69 
70  if (dict.found("p0"))
71  {
72  p0_ = vectorField("p0", dict , p.size());
73  }
74  else
75  {
76  p0_ = p.localPoints();
77  }
78 }
79 
80 
83 (
85  const pointPatch& p,
87  const pointPatchFieldMapper& mapper
88 )
89 :
90  fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
91  amplitude_(ptf.amplitude_),
92  omega_(ptf.omega_),
93  p0_(mapper(ptf.p0_))
94 {}
95 
96 
99 (
102 )
103 :
105  amplitude_(ptf.amplitude_),
106  omega_(ptf.omega_),
107  p0_(ptf.p0_)
108 {}
109 
110 
111 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112 
114 (
115  const pointPatchFieldMapper& m
116 )
117 {
119 
120  m(p0_, p0_);
121 }
122 
123 
125 (
126  const pointPatchField<vector>& ptf,
127  const labelList& addr
128 )
129 {
131  refCast<const oscillatingVelocityPointPatchVectorField>(ptf);
132 
134 
135  p0_.rmap(oVptf.p0_, addr);
136 }
137 
138 
140 {
141  if (this->updated())
142  {
143  return;
144  }
145 
146  const polyMesh& mesh = this->internalField().mesh()();
147  const Time& t = mesh.time();
148  const pointPatch& p = this->patch();
149 
151  (
152  (p0_ + amplitude_*sin(omega_*t.value()) - p.localPoints())
153  /t.deltaTValue()
154  );
155 
157 }
158 
159 
161 {
163  writeEntry(os, "amplitude", amplitude_);
164  writeEntry(os, "omega", omega_);
165  writeEntry(os, "p0", p0_);
166  writeEntry(os, "value", *this);
167 }
168 
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
173 (
176 );
177 
178 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179 
180 } // End namespace Foam
181 
182 // ************************************************************************* //
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
oscillatingVelocityPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
Foam::pointPatchFieldMapper.
label size() const
Return size.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Macros for easy insertion into run-time selection tables.
dynamicFvMesh & mesh
Pre-declare SubField and related Field type.
Definition: Field.H:56
const pointPatch & patch() const
Return patch.
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
const Type & value() const
Return const reference to value.
virtual void rmap(const pointPatchField< vector > &, const labelList &)
Reverse map the given pointPatchField onto this pointPatchField.
static const zero Zero
Definition: zero.H:97
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:54
dimensionedScalar sin(const dimensionedScalar &ds)
const Time & time() const
Return time.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
virtual const vectorField & localPoints() const =0
Return mesh points.
const DimensionedField< vector, pointMesh > & internalField() const
Return dimensioned internal field reference.
virtual void write(Ostream &) const
Write.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:56
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
bool updated() const
Return true if the boundary condition has already been updated.
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.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
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 this pointPatchField.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844