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-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 #include "OneConstant.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
37 (
38  const pointPatch& p,
40 )
41 :
43  amplitude_(Zero),
44  omega_(0.0),
45  waveNumber_(Zero)
46 {}
47 
48 
51 (
52  const pointPatch& p,
54  const dictionary& dict
55 )
56 :
58  amplitude_(dict.lookup("amplitude")),
59  omega_(dict.lookup<scalar>("omega")),
60  waveNumber_(dict.lookupOrDefault<vector>("waveNumber", Zero)),
61  startRamp_
62  (
63  dict.found("startRamp")
64  ? Function1<scalar>::New("startRamp", dict)
65  : autoPtr<Function1<scalar>>
66  (
67  new Function1s::OneConstant<scalar>("startRamp")
68  )
69  ),
70  endRamp_
71  (
72  dict.found("endRamp")
73  ? Function1<scalar>::New("endRamp", dict)
74  : autoPtr<Function1<scalar>>
75  (
76  new Function1s::OneConstant<scalar>("endRamp")
77  )
78  ),
79  timeRamp_
80  (
81  dict.found("timeRamp")
82  ? Function1<scalar>::New("timeRamp", dict)
83  : autoPtr<Function1<scalar>>
84  (
85  new Function1s::OneConstant<scalar>("timeRamp")
86  )
87  )
88 {
89  if (!dict.found("value"))
90  {
91  updateCoeffs();
92  }
93 }
94 
95 
98 (
100  const pointPatch& p,
102  const pointPatchFieldMapper& mapper
103 )
104 :
105  fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
106  amplitude_(ptf.amplitude_),
107  omega_(ptf.omega_),
108  waveNumber_(ptf.waveNumber_)
109 {}
110 
111 
114 (
117 )
118 :
120  amplitude_(ptf.amplitude_),
121  omega_(ptf.omega_),
122  waveNumber_(ptf.waveNumber_)
123 {}
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
129 {
130  if (this->updated())
131  {
132  return;
133  }
134 
135  const polyMesh& mesh = this->internalField().mesh()();
136  const Time& t = mesh.time();
137 
138  const scalarField points(waveNumber_ & patch().localPoints());
139 
140  const scalar timeRamp = timeRamp_->value(t.value());
141 
142  const scalarField startRamp(startRamp_->value(points));
143 
144  const scalarField endRamp
145  (
146  endRamp_->value(points[points.size() - 1] - points)
147  );
148 
150  (
151  timeRamp*startRamp*endRamp*amplitude_*cos(omega_*t.value() - points)
152  );
153 
155 }
156 
157 
159 {
161  writeEntry(os, "amplitude", amplitude_);
162  writeEntry(os, "omega", omega_);
163  writeEntry(os, "waveNumber", waveNumber_);
164 
165  if (!isType<Function1s::OneConstant<scalar>>(startRamp_()))
166  {
167  writeEntry(os, startRamp_());
168  }
169 
170  if (!isType<Function1s::OneConstant<scalar>>(endRamp_()))
171  {
172  writeEntry(os, endRamp_());
173  }
174 
175  if (!isType<Function1s::OneConstant<scalar>>(timeRamp_()))
176  {
177  writeEntry(os, timeRamp_());
178  }
179 
180  writeEntry(os, "value", *this);
181 }
182 
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 namespace Foam
187 {
189  (
192  );
193 }
194 
195 // ************************************************************************* //
Run-time selectable general function of one variable.
Definition: Function1.H:52
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
waveDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Foam::pointPatchFieldMapper.
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.
Templated function that returns the corresponding 1 (one).
Definition: OneConstant.H:56
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
const pointField & points
virtual Type value(const scalar x) const =0
Return value as a function of scalar x.
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
const pointPatch & patch() const
Return patch.
bool isType(const Type &t)
Check the typeid.
Definition: typeInfo.H:126
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:54
const Time & time() const
Return time.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
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)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
bool found
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