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-2023 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  const dictionary& dict
41 )
42 :
44  amplitude_(dict.lookup("amplitude")),
45  omega_(dict.lookup<scalar>("omega")),
46  waveNumber_(dict.lookupOrDefault<vector>("waveNumber", Zero)),
47  startRamp_
48  (
49  dict.found("startRamp")
50  ? Function1<scalar>::New("startRamp", dict)
51  : autoPtr<Function1<scalar>>
52  (
53  new Function1s::OneConstant<scalar>("startRamp")
54  )
55  ),
56  endRamp_
57  (
58  dict.found("endRamp")
59  ? Function1<scalar>::New("endRamp", dict)
60  : autoPtr<Function1<scalar>>
61  (
62  new Function1s::OneConstant<scalar>("endRamp")
63  )
64  ),
65  timeRamp_
66  (
67  dict.found("timeRamp")
68  ? Function1<scalar>::New("timeRamp", dict)
69  : autoPtr<Function1<scalar>>
70  (
71  new Function1s::OneConstant<scalar>("timeRamp")
72  )
73  )
74 {
75  if (!dict.found("value"))
76  {
77  updateCoeffs();
78  }
79 }
80 
81 
84 (
86  const pointPatch& p,
88  const pointPatchFieldMapper& mapper
89 )
90 :
91  fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
92  amplitude_(ptf.amplitude_),
93  omega_(ptf.omega_),
94  waveNumber_(ptf.waveNumber_)
95 {}
96 
97 
100 (
103 )
104 :
106  amplitude_(ptf.amplitude_),
107  omega_(ptf.omega_),
108  waveNumber_(ptf.waveNumber_)
109 {}
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
115 {
116  if (this->updated())
117  {
118  return;
119  }
120 
121  const polyMesh& mesh = this->internalField().mesh()();
122  const Time& t = mesh.time();
123 
124  const scalarField points(waveNumber_ & patch().localPoints());
125 
126  const scalar timeRamp = timeRamp_->value(t.value());
127 
128  const scalarField startRamp(startRamp_->value(points));
129 
130  const scalarField endRamp
131  (
132  endRamp_->value(points[points.size() - 1] - points)
133  );
134 
136  (
137  timeRamp*startRamp*endRamp*amplitude_*cos(omega_*t.value() - points)
138  );
139 
141 }
142 
143 
145 {
147  writeEntry(os, "amplitude", amplitude_);
148  writeEntry(os, "omega", omega_);
149  writeEntry(os, "waveNumber", waveNumber_);
150 
151  if (!isType<Function1s::OneConstant<scalar>>(startRamp_()))
152  {
153  writeEntry(os, startRamp_());
154  }
155 
156  if (!isType<Function1s::OneConstant<scalar>>(endRamp_()))
157  {
158  writeEntry(os, endRamp_());
159  }
160 
161  if (!isType<Function1s::OneConstant<scalar>>(timeRamp_()))
162  {
163  writeEntry(os, timeRamp_());
164  }
165 
166  writeEntry(os, "value", *this);
167 }
168 
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
172 namespace Foam
173 {
175  (
178  );
179 }
180 
181 // ************************************************************************* //
bool found
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Run-time selectable general function of one variable.
Definition: Function1.H:64
Templated function that returns the corresponding 1 (one).
Definition: OneConstant.H:59
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const Type & value() const
Return const reference to value.
A FixedValue boundary condition for pointField.
const Time & time() const
Return time.
Foam::pointPatchFieldMapper.
Abstract base class for point-mesh patch fields.
virtual void write(Ostream &) const
Write.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:57
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
waveDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
const pointField & points
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
makePointPatchTypeField(pointPatchVectorField, solidBodyMotionDisplacementPointPatchVectorField)
bool isType(const Type &t)
Check the typeid.
Definition: typeInfo.H:131
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar cos(const dimensionedScalar &ds)
dictionary dict
volScalarField & p