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-2024 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<vector>("amplitude", dimLength)),
45  omega_(dict.lookup<scalar>("omega", unitRadians/dimTime)),
46  waveNumber_
47  (
48  dict.lookupOrDefault<vector>("waveNumber", dimless/dimLength, Zero)
49  ),
50  startRamp_
51  (
52  dict.found("startRamp")
53  ? Function1<scalar>::New("startRamp", dimless, dimless, dict)
54  : autoPtr<Function1<scalar>>
55  (
56  new Function1s::OneConstant<scalar>("startRamp")
57  )
58  ),
59  endRamp_
60  (
61  dict.found("endRamp")
62  ? Function1<scalar>::New("endRamp", dimless, dimless, dict)
63  : autoPtr<Function1<scalar>>
64  (
65  new Function1s::OneConstant<scalar>("endRamp")
66  )
67  ),
68  timeRamp_
69  (
70  dict.found("timeRamp")
71  ? Function1<scalar>::New("timeRamp", dimTime, dimless, dict)
72  : autoPtr<Function1<scalar>>
73  (
74  new Function1s::OneConstant<scalar>("timeRamp")
75  )
76  )
77 {
78  if (!dict.found("value"))
79  {
80  updateCoeffs();
81  }
82 }
83 
84 
87 (
89  const pointPatch& p,
91  const fieldMapper& mapper
92 )
93 :
94  fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
95  amplitude_(ptf.amplitude_),
96  omega_(ptf.omega_),
97  waveNumber_(ptf.waveNumber_)
98 {}
99 
100 
103 (
106 )
107 :
109  amplitude_(ptf.amplitude_),
110  omega_(ptf.omega_),
111  waveNumber_(ptf.waveNumber_)
112 {}
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
118 {
119  if (this->updated())
120  {
121  return;
122  }
123 
124  const polyMesh& mesh = this->internalField().mesh()();
125  const Time& t = mesh.time();
126 
127  const scalarField points(waveNumber_ & patch().localPoints());
128 
129  const scalar timeRamp = timeRamp_->value(t.value());
130 
131  const scalarField startRamp(startRamp_->value(points));
132 
133  const scalarField endRamp
134  (
135  endRamp_->value(points[points.size() - 1] - points)
136  );
137 
139  (
140  timeRamp*startRamp*endRamp*amplitude_*cos(omega_*t.value() - points)
141  );
142 
144 }
145 
146 
148 {
150  writeEntry(os, "amplitude", amplitude_);
151  writeEntry(os, "omega", omega_);
152  writeEntry(os, "waveNumber", waveNumber_);
153 
154  if (!isType<Function1s::OneConstant<scalar>>(startRamp_()))
155  {
156  writeEntry(os, startRamp_());
157  }
158 
159  if (!isType<Function1s::OneConstant<scalar>>(endRamp_()))
160  {
161  writeEntry(os, endRamp_());
162  }
163 
164  if (!isType<Function1s::OneConstant<scalar>>(timeRamp_()))
165  {
166  writeEntry(os, timeRamp_());
167  }
168 
169  writeEntry(os, "value", *this);
170 }
171 
172 
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174 
175 namespace Foam
176 {
178  (
181  );
182 }
183 
184 // ************************************************************************* //
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:125
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:162
const Type & value() const
Return const reference to value.
Abstract base class for field mapping.
Definition: fieldMapper.H:48
A FixedValue boundary condition for pointField.
const Time & time() const
Return time.
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, angularOscillatingDisplacementPointPatchVectorField)
const dimensionSet dimless
bool isType(const Type &t)
Check the typeid.
Definition: typeInfo.H:158
const dimensionSet dimLength
const dimensionSet dimTime
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const unitConversion unitRadians
dimensionedScalar cos(const dimensionedScalar &ds)
dictionary dict
volScalarField & p