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-2018 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  "Water Wave Mechanics For Engineers And Scientists"
33  R G Dean and R A Dalrymple
34  Pages 314-317
35  \endverbatim
36 
37 SourceFiles
38  solitary.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef solitary_H
43 #define solitary_H
44 
45 #include "waveModel.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 namespace waveModels
52 {
53 
54 /*---------------------------------------------------------------------------*\
55  Class solitary Declaration
56 \*---------------------------------------------------------------------------*/
57 
58 class solitary
59 :
60  public waveModel
61 {
62  //- Private data
63 
64  //- Offset [m]
65  const scalar offset_;
66 
67  //- Depth [m]
68  const scalar depth_;
69 
70 
71  // Private Member Functions
72 
73  //- The wavenumber [1/m]
74  scalar k(const scalar t) const;
75 
76  //- The dimensionless amplitude [1]
77  scalar alpha(const scalar t) const;
78 
79  //- The wave celerity [m/s]
80  scalar celerity(const scalar t) const;
81 
82  //- The evolution parameter [1]
83  // This is analogous to the oscillation angle of a periodic wave
84  tmp<scalarField> parameter
85  (
86  const scalar t,
87  const scalar u,
88  const scalarField& x
89  ) const;
90 
91  //- The dimensionless elevation [1]
93  (
94  const scalar t,
95  const scalar u,
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, mean velocity and local
142  // coordinates. Local x is aligned with the mean velocity.
144  (
145  const scalar t,
146  const scalar u,
147  const scalarField& x
148  ) const;
149 
150  //- Get the wave velocity at a given time, mean velocity and local
151  // coordinates. Local x is aligned with the mean velocity, and z with
152  // negative gravity.
154  (
155  const scalar t,
156  const scalar u,
157  const vector2DField& xz
158  ) const;
159 
160  //- Get the wave pressure at a given time, mean velocity and local
161  // coordinates. Local x is aligned with the mean velocity, and z with
162  // negative gravity.
163  virtual tmp<scalarField> pressure
164  (
165  const scalar t,
166  const scalar u,
167  const vector2DField& xz
168  ) const;
169 
170  //- Write
171  virtual void write(Ostream& os) const;
172 };
173 
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 } // End namespace waveModels
178 } // End namespace Foam
179 
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 
182 #endif
183 
184 // ************************************************************************* //
solitary(const solitary &wave)
Construct a copy.
Definition: solitary.C:88
virtual tmp< scalarField > elevation(const scalar t, const scalar u, const scalarField &x) const
Get the wave elevation at a given time, mean velocity and local.
Definition: solitary.C:117
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:137
virtual tmp< vector2DField > velocity(const scalar t, const scalar u, const vector2DField &xz) const
Get the wave velocity at a given time, mean velocity and local.
Definition: solitary.C:128
scalar offset() const
Get the offset.
Definition: solitary.H:129
virtual void write(Ostream &os) const
Write.
Definition: solitary.C:168
Generic base class for waves. Derived classes must implement field functions which return the elevati...
Definition: waveModel.H:55
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
virtual tmp< scalarField > pressure(const scalar t, const scalar u, const vector2DField &xz) const
Get the wave pressure at a given time, mean velocity and local.
Definition: solitary.C:157
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:110
TypeName("solitary")
Runtime type information.
Namespace for OpenFOAM.
Solitary wave model.
Definition: solitary.H:57