waveDisplacementPointPatchVectorField.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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const pointPatch& p,
39 )
40 :
42  amplitude_(Zero),
43  omega_(0.0),
44  waveNumber_(Zero)
45 {}
46 
47 
50 (
51  const pointPatch& p,
53  const dictionary& dict
54 )
55 :
57  amplitude_(dict.lookup("amplitude")),
58  omega_(readScalar(dict.lookup("omega"))),
59  waveNumber_(dict.lookupOrDefault<vector>("waveNumber", Zero))
60 {
61  if (!dict.found("value"))
62  {
63  updateCoeffs();
64  }
65 }
66 
67 
70 (
72  const pointPatch& p,
74  const pointPatchFieldMapper& mapper
75 )
76 :
77  fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
78  amplitude_(ptf.amplitude_),
79  omega_(ptf.omega_),
80  waveNumber_(ptf.waveNumber_)
81 {}
82 
83 
86 (
89 )
90 :
92  amplitude_(ptf.amplitude_),
93  omega_(ptf.omega_),
94  waveNumber_(ptf.waveNumber_)
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
101 {
102  if (this->updated())
103  {
104  return;
105  }
106 
107  const polyMesh& mesh = this->internalField().mesh()();
108  const Time& t = mesh.time();
109 
110  const scalarField points( waveNumber_ & patch().localPoints());
111 
113  (
114  amplitude_*cos(omega_*t.value() - points)
115  );
116 
118 }
119 
120 
122 {
124  os.writeKeyword("amplitude")
125  << amplitude_ << token::END_STATEMENT << nl;
126  os.writeKeyword("omega")
127  << omega_ << token::END_STATEMENT << nl;
128  os.writeKeyword("waveNumber")
129  << waveNumber_ << token::END_STATEMENT << nl;
130  writeEntry("value", os);
131 }
132 
133 
134 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
135 
136 namespace Foam
137 {
139  (
142  );
143 }
144 
145 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
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
waveDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
Foam::pointPatchFieldMapper.
Foam::waveDisplacementPointPatchVectorField.
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.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
const pointField & points
const pointPatch & patch() const
Return patch.
const Type & value() const
Return const reference to value.
static const zero Zero
Definition: zero.H:97
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:265
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.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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