solitary.H
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) 2017-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 Class
25  Foam::waveModels::solitary
26 
27 Description
28  Solitary wave model.
29 
30  Reference:
31  \verbatim
32  Dean, R. G., & Dalrymple, R. A. (1991).
33  Water wave mechanics for engineers and scientists (Vol. 2).
34  World Scientific Publishing Company.
35  \endverbatim
36 
37  See pages 314-317.
38 
39 Usage
40  \table
41  Property | Description | Required? | Default
42  depth | The water depth [m] | no | great
43  amplitude | Peak amplitude [m] | yes |
44  offset | The positional offset [m] | yes |
45  \endtable
46 
47  Example specification in constant/waveProperties:
48  \verbatim
49  waves
50  (
51  solitary
52  {
53  length 40;
54  amplitude 0.5;
55  offset 0;
56  }
57  );
58  \endverbatim
59 
60 SourceFiles
61  solitary.C
62 
63 \*---------------------------------------------------------------------------*/
64 
65 #ifndef solitary_H
66 #define solitary_H
67 
68 #include "waveModel.H"
69 
70 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72 namespace Foam
73 {
74 namespace waveModels
75 {
76 
77 /*---------------------------------------------------------------------------*\
78  Class solitary Declaration
79 \*---------------------------------------------------------------------------*/
80 
81 class solitary
82 :
83  public waveModel
84 {
85  // Private Data
86 
87  //- Depth [m]
88  const scalar depth_;
89 
90  //- Amplitude [m]
91  const scalar amplitude_;
92 
93  //- Offset [m]
94  const scalar offset_;
95 
96 
97  // Private Member Functions
98 
99  //- The wavenumber [1/m]
100  scalar k() const;
101 
102  //- The dimensionless amplitude [1]
103  scalar alpha() const;
104 
105  //- The evolution parameter [1]
106  // This is analogous to the oscillation angle of a periodic wave
107  tmp<scalarField> parameter
108  (
109  const scalar t,
110  const scalarField& x
111  ) const;
112 
113  //- The dimensionless elevation [1]
115  (
116  const scalar t,
117  const scalarField& x
118  ) const;
119 
120 
121 public:
122 
123  //- Runtime type information
124  TypeName("solitary");
125 
126 
127  // Constructors
128 
129  //- Construct a copy
130  solitary(const solitary& wave);
131 
132  //- Construct from a dictionary and gravity
133  solitary(const dictionary& dict, const scalar g);
134 
135  //- Construct a clone
136  virtual autoPtr<waveModel> clone() const
137  {
138  return autoPtr<waveModel>(new solitary(*this));
139  }
140 
141 
142  //- Destructor
143  virtual ~solitary();
144 
145 
146  // Member Functions
147 
148  // Access
149 
150  //- Get the depth
151  scalar depth() const
152  {
153  return depth_;
154  }
155 
156  //- Get the amplitude
157  scalar amplitude() const
158  {
159  return amplitude_;
160  }
161 
162  //- Get the offset
163  scalar offset() const
164  {
165  return offset_;
166  }
167 
168  //- The wave celerity [m/s]
169  scalar celerity() const;
170 
171  //- Get the wave elevation at a given time and local coordinates. Local
172  // x is aligned with the direction of propagation.
174  (
175  const scalar t,
176  const scalarField& x
177  ) const;
178 
179  //- Get the wave velocity at a given time and local coordinates. Local
180  // x is aligned with the direction of propagation, and z with negative
181  // gravity.
183  (
184  const scalar t,
185  const vector2DField& xz
186  ) const;
187 
188  //- Write
189  virtual void write(Ostream& os) const;
190 };
191 
192 
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 
195 } // End namespace waveModels
196 } // End namespace Foam
197 
198 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
199 
200 #endif
201 
202 // ************************************************************************* //
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
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
A class for managing temporary objects.
Definition: tmp.H:55
scalar g() const
Get the value of gravity.
Definition: waveModel.H:117
Solitary wave model.
Definition: solitary.H:103
solitary(const solitary &wave)
Construct a copy.
Definition: solitary.C:86
scalar celerity() const
The wave celerity [m/s].
Definition: solitary.C:56
virtual tmp< scalarField > elevation(const scalar t, const scalarField &x) const
Get the wave elevation at a given time and local coordinates. Local.
Definition: solitary.C:117
TypeName("solitary")
Runtime type information.
virtual autoPtr< waveModel > clone() const
Construct a clone.
Definition: solitary.H:155
virtual void write(Ostream &os) const
Write.
Definition: solitary.C:154
virtual tmp< vector2DField > velocity(const scalar t, const vector2DField &xz) const
Get the wave velocity at a given time and local coordinates. Local.
Definition: solitary.C:127
virtual ~solitary()
Destructor.
Definition: solitary.C:110
scalar amplitude() const
Get the amplitude.
Definition: solitary.H:176
scalar offset() const
Get the offset.
Definition: solitary.H:182
scalar depth() const
Get the depth.
Definition: solitary.H:170
label wave(const fvMesh &mesh, const List< labelPair > &changedPatchAndFaces, const label nCorrections, GeometricField< scalar, PatchField, GeoMesh > &distance, TrackingData &td, GeometricField< DataType, PatchField, GeoMesh > &... data)
Wave distance (and maybe additional) data from faces. If nCorrections is.
Namespace for OpenFOAM.
dictionary dict