solitary.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) 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 \*---------------------------------------------------------------------------*/
25 
26 #include "solitary.H"
27 #include "mathematicalConstants.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace waveModels
35 {
36  defineTypeNameAndDebug(solitary, 0);
37  addToRunTimeSelectionTable(waveModel, solitary, objectRegistry);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 Foam::scalar Foam::waveModels::solitary::k(const scalar t) const
45 {
46  return sqrt(0.75*amplitude(t)/pow3(depth()));
47 }
48 
49 
50 Foam::scalar Foam::waveModels::solitary::alpha(const scalar t) const
51 {
52  return amplitude(t)/depth();
53 }
54 
55 
56 Foam::scalar Foam::waveModels::solitary::celerity(const scalar t) const
57 {
58  return sqrt(g()*depth()/(1 - alpha(t)));
59 }
60 
61 
62 Foam::tmp<Foam::scalarField> Foam::waveModels::solitary::parameter
63 (
64  const scalar t,
65  const scalarField& x
66 ) const
67 {
68  return k(t)*(x - offset_ - celerity(t)*t);
69 }
70 
71 
72 Foam::tmp<Foam::scalarField> Foam::waveModels::solitary::Pi
73 (
74  const scalar t,
75  const scalarField& x
76 ) const
77 {
78  const scalar clip = log(great);
79 
80  return 1/sqr(cosh(max(- clip, min(clip, parameter(t, x)))));
81 }
82 
83 
84 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
85 
87 :
88  waveModel(wave),
89  offset_(wave.offset_),
90  depth_(wave.depth_)
91 {}
92 
93 
95 (
96  const objectRegistry& db,
97  const dictionary& dict
98 )
99 :
100  waveModel(db, dict),
101  offset_(dict.lookup<scalar>("offset")),
102  depth_(dict.lookup<scalar>("depth"))
103 {}
104 
105 
106 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
107 
109 {}
110 
111 
112 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
113 
115 (
116  const scalar t,
117  const scalarField& x
118 ) const
119 {
120  return amplitude(t)*Pi(t, x);
121 }
122 
123 
125 (
126  const scalar t,
127  const vector2DField& xz
128 ) const
129 {
130  const scalar A = alpha(t);
131  const scalarField Z(max(scalar(0), 1 + xz.component(1)/depth()));
132  const scalarField P(Pi(t, xz.component(0)));
133 
134  return
135  celerity(t)
136  *zip
137  (
138  A/4
139  *(
140  (4 + 2*A - 6*A*sqr(Z))*P
141  + (- 7*A + 9*A*sqr(Z))*sqr(P)
142  ),
143  A*Z*depth()*k(t)*tanh(parameter(t, xz.component(0)))
144  *(
145  (2 + A - A*sqr(Z))*P
146  + (- 7*A + 3*A*sqr(Z))*sqr(P)
147  )
148  );
149 }
150 
151 
153 (
154  const scalar t,
155  const vector2DField& xz
156 ) const
157 {
159  return tmp<scalarField>(nullptr);
160 }
161 
162 
164 {
165  waveModel::write(os);
166 
167  writeEntry(os, "offset", offset_);
168  writeEntry(os, "depth", depth_);
169 }
170 
171 
172 // ************************************************************************* //
solitary(const solitary &wave)
Construct a copy.
Definition: solitary.C:86
tmp< vector2DField > zip(const tmp< scalarField > &x, const tmp< scalarField > &y)
Definition: vector2DField.C:31
dimensionedScalar tanh(const dimensionedScalar &ds)
scalar depth() const
Get the depth.
Definition: solitary.H:135
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
defineTypeNameAndDebug(Airy, 0)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
dimensionedScalar sqrt(const dimensionedScalar &ds)
label k
Boltzmann constant.
Macros for easy insertion into run-time selection tables.
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:56
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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
addToRunTimeSelectionTable(waveModel, Airy, objectRegistry)
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:441
waveModel(const waveModel &wave)
Construct a copy.
Definition: waveModel.C:41
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar pow3(const dimensionedScalar &ds)
scalar amplitude(const scalar t) const
Get the amplitude.
Definition: waveModel.H:123
dimensionedScalar cosh(const dimensionedScalar &ds)
A class for managing temporary objects.
Definition: PtrList.H:53
Registry of regIOobjects.
const dimensionedVector & g
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:370
virtual ~solitary()
Destructor.
Definition: solitary.C:108
virtual void write(Ostream &os) const
Write.
Definition: waveModel.C:73
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844
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