oscillatingDisplacementPointPatchVectorField.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 {}
52 
53 
56 (
57  const pointPatch& p,
59  const dictionary& dict
60 )
61 :
63  amplitude_(dict.lookup("amplitude")),
64  omega_(dict.get<scalar>("omega"))
65 {
66  if (!dict.found("value"))
67  {
69  }
70 }
71 
72 
75 (
76  const oscillatingDisplacementPointPatchVectorField& ptf,
77  const pointPatch& p,
78  const DimensionedField<vector, pointMesh>& iF,
79  const pointPatchFieldMapper& mapper
80 )
81 :
82  fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
83  amplitude_(ptf.amplitude_),
84  omega_(ptf.omega_)
85 {}
86 
87 
90 (
93 )
94 :
96  amplitude_(ptf.amplitude_),
97  omega_(ptf.omega_)
98 {}
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
104 {
105  if (this->updated())
106  {
107  return;
108  }
109 
110  const polyMesh& mesh = this->internalField().mesh()();
111  const Time& t = mesh.time();
113  Field<vector>::operator=(amplitude_*sin(omega_*t.value()));
114 
116 }
117 
118 
120 {
122  os.writeEntry("amplitude", amplitude_);
123  os.writeEntry("omega", omega_);
124  this->writeValueEntry(os);
125 }
126 
127 
128 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
129 
131 (
133  oscillatingDisplacementPointPatchVectorField
134 );
135 
136 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137 
138 } // End namespace Foam
139 
140 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
dictionary dict
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void write(Ostream &os) const
Write.
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.
dynamicFvMesh & mesh
Vector< scalar > vector
Definition: vector.H:57
bool updated() const noexcept
True if the boundary condition has already been updated.
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)
void operator=(const Field< Type > &)
Copy assignment.
Definition: Field.C:747
pointPatchField< vector > pointPatchVectorField
oscillatingDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
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)
volScalarField & p
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