oscillatingDisplacementPointPatchVectorField.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-2018 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 {}
50 
51 
54 (
55  const pointPatch& p,
57  const dictionary& dict
58 )
59 :
61  amplitude_(dict.lookup("amplitude")),
62  omega_(readScalar(dict.lookup("omega")))
63 {
64  if (!dict.found("value"))
65  {
66  updateCoeffs();
67  }
68 }
69 
70 
73 (
75  const pointPatch& p,
77  const pointPatchFieldMapper& mapper
78 )
79 :
80  fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
81  amplitude_(ptf.amplitude_),
82  omega_(ptf.omega_)
83 {}
84 
85 
88 (
91 )
92 :
94  amplitude_(ptf.amplitude_),
95  omega_(ptf.omega_)
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
102 {
103  if (this->updated())
104  {
105  return;
106  }
107 
108  const polyMesh& mesh = this->internalField().mesh()();
109  const Time& t = mesh.time();
110 
111  Field<vector>::operator=(amplitude_*sin(omega_*t.value()));
112 
114 }
115 
116 
118 {
120  os.writeKeyword("amplitude")
121  << amplitude_ << token::END_STATEMENT << nl;
122  os.writeKeyword("omega")
123  << omega_ << token::END_STATEMENT << nl;
124  writeEntry("value", os);
125 }
126 
127 
128 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
129 
131 (
134 );
135 
136 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137 
138 } // End namespace Foam
139 
140 // ************************************************************************* //
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Foam::pointPatchFieldMapper.
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:726
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
const Type & value() const
Return const reference to value.
static const zero Zero
Definition: zero.H:97
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
dimensionedScalar sin(const dimensionedScalar &ds)
static const char nl
Definition: Ostream.H:265
void operator=(const Field< Type > &)
Definition: Field.C:764
const Time & time() const
Return time.
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
const DimensionedField< vector, pointMesh > & internalField() const
Return dimensioned internal field reference.
virtual void write(Ostream &) const
Write.
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: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)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576