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-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 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 SourceFiles
40  solitary.C
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #ifndef solitary_H
45 #define solitary_H
46 
47 #include "waveModel.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 namespace waveModels
54 {
55 
56 /*---------------------------------------------------------------------------*\
57  Class solitary Declaration
58 \*---------------------------------------------------------------------------*/
59 
60 class solitary
61 :
62  public waveModel
63 {
64  // Private Data
65 
66  //- Offset [m]
67  const scalar offset_;
68 
69  //- Depth [m]
70  const scalar depth_;
71 
72 
73  // Private Member Functions
74 
75  //- The wavenumber [1/m]
76  scalar k(const scalar t) const;
77 
78  //- The dimensionless amplitude [1]
79  scalar alpha(const scalar t) const;
80 
81  //- The wave celerity [m/s]
82  scalar celerity(const scalar t) const;
83 
84  //- The evolution parameter [1]
85  // This is analogous to the oscillation angle of a periodic wave
86  tmp<scalarField> parameter
87  (
88  const scalar t,
89  const scalarField& x
90  ) const;
91 
92  //- The dimensionless elevation [1]
94  (
95  const scalar t,
96  const scalarField& x
97  ) const;
98 
99 
100 public:
101 
102  //- Runtime type information
103  TypeName("solitary");
104 
105 
106  // Constructors
107 
108  //- Construct a copy
109  solitary(const solitary& wave);
110 
111  //- Construct from a database and a dictionary
112  solitary(const objectRegistry& db, const dictionary& dict);
113 
114  //- Construct a clone
115  virtual autoPtr<waveModel> clone() const
116  {
117  return autoPtr<waveModel>(new solitary(*this));
118  }
119 
120 
121  //- Destructor
122  virtual ~solitary();
123 
124 
125  // Member Functions
126 
127  // Access
128 
129  //- Get the offset
130  scalar offset() const
131  {
132  return offset_;
133  }
134 
135  //- Get the depth
136  scalar depth() const
137  {
138  return depth_;
139  }
140 
141  //- Get the wave elevation at a given time and local coordinates. Local
142  // x is aligned with the direction of propagation.
144  (
145  const scalar t,
146  const scalarField& x
147  ) const;
148 
149  //- Get the wave velocity at a given time and local coordinates. Local
150  // x is aligned with the direction of propagation, and z with negative
151  // gravity.
153  (
154  const scalar t,
155  const vector2DField& xz
156  ) const;
157 
158  //- Get the wave pressure at a given time and local coordinates. Local
159  // x is aligned with the direction of propagation, and z with negative
160  // gravity.
161  virtual tmp<scalarField> pressure
162  (
163  const scalar t,
164  const vector2DField& xz
165  ) const;
166 
167  //- Write
168  virtual void write(Ostream& os) const;
169 };
170 
171 
172 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173 
174 } // End namespace waveModels
175 } // End namespace Foam
176 
177 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 
179 #endif
180 
181 // ************************************************************************* //
solitary(const solitary &wave)
Construct a copy.
Definition: solitary.C:86
dictionary dict
scalar depth() const
Get the depth.
Definition: solitary.H:135
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
scalar offset() const
Get the offset.
Definition: solitary.H:129
virtual void write(Ostream &os) const
Write.
Definition: solitary.C:163
Generic base class for waves. Derived classes must implement field functions which return the elevati...
Definition: waveModel.H:55
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:125
virtual tmp< scalarField > pressure(const scalar t, const vector2DField &xz) const
Get the wave pressure at a given time and local coordinates. Local.
Definition: solitary.C:153
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
Registry of regIOobjects.
virtual autoPtr< waveModel > clone() const
Construct a clone.
Definition: solitary.H:114
virtual ~solitary()
Destructor.
Definition: solitary.C:108
TypeName("solitary")
Runtime type information.
Namespace for OpenFOAM.
Solitary wave model.
Definition: solitary.H:59
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:115