waveSpectra::GodaJONSWAP: New irregular wave spectrum

This is an alternative, approximate parameterisation of the JONSWAP
spectrum, in which the significant wave height and period are specified
instead of the wind speed and fetch.

Example specification, in constant/waveProperties:

    waves
    (
        irregular
        {
            spectrum GodaJONSWAP;
            Hs      2;  // <- significant wave height [m]
            Ts      6;  // <- significant wave period [s]
            n       12; // <- number of samples []
            angle   0;
        }
    );
This commit is contained in:
Will Bainbridge
2024-02-27 09:58:41 +00:00
parent 8ba48acb7f
commit cbec00456f
3 changed files with 272 additions and 0 deletions

View File

@ -8,6 +8,7 @@ waveModels/solitary/solitary.C
waveModels/irregular/waveSpectra/waveSpectrum/waveSpectrum.C
waveModels/irregular/waveSpectra/PiersonMoskowitz/PiersonMoskowitz.C
waveModels/irregular/waveSpectra/JONSWAP/JONSWAP.C
waveModels/irregular/waveSpectra/GodaJONSWAP/GodaJONSWAP.C
waveModels/irregular/irregular.C
waveSuperpositions/waveSuperposition/waveSuperposition.C

View File

@ -0,0 +1,116 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2024 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "GodaJONSWAP.H"
#include "mathematicalConstants.H"
#include "addToRunTimeSelectionTable.H"
using namespace Foam::constant::mathematical;
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace waveSpectra
{
defineTypeNameAndDebug(GodaJONSWAP, 0);
addToRunTimeSelectionTable(waveSpectrum, GodaJONSWAP, dictionary);
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::waveSpectra::GodaJONSWAP::GodaJONSWAP
(
const GodaJONSWAP& spectrum
)
:
waveSpectrum(spectrum),
Hs_(spectrum.Hs_),
Ts_(spectrum.Ts_),
gamma_(spectrum.gamma_)
{}
Foam::waveSpectra::GodaJONSWAP::GodaJONSWAP
(
const dictionary& dict,
const scalar g
)
:
waveSpectrum(dict, g),
Hs_(dict.lookup<scalar>("Hs")),
Ts_(dict.lookup<scalar>("Ts")),
gamma_(dict.lookupOrDefault<scalar>("gamma", 3.3))
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::waveSpectra::GodaJONSWAP::~GodaJONSWAP()
{}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
Foam::tmp<Foam::scalarField> Foam::waveSpectra::GodaJONSWAP::S
(
const scalarField& f
) const
{
const scalar betaJ =
0.0624
/(0.230 + 0.0336*gamma_ - 0.185/(1.9 + gamma_))
*(1.094 - 0.01915*log(gamma_));
const scalar Tp = Ts_/(1 - 0.132*pow(gamma_ + 0.2, -0.559));
const scalarField sigma(0.07 + 0.02*pos(f - 1/Tp));
const scalarField r(exp(- sqr(Tp*f - 1)/(2*sqr(sigma))));
return betaJ*sqr(Hs_)/f/pow4(Tp*f)*exp(- 1.25/pow4(Tp*f))*pow(gamma_, r);
}
Foam::scalar Foam::waveSpectra::GodaJONSWAP::fFraction
(
const scalar fraction
) const
{
const scalar Tp = Ts_/(1 - 0.132*pow(gamma_ + 0.2, -0.559));
const scalar f1 = 1/Tp/pow025(rootSmall/1.25);
return waveSpectrum::fFraction(fraction, f1);
}
void Foam::waveSpectra::GodaJONSWAP::write(Ostream& os) const
{
waveSpectrum::write(os);
writeEntry(os, "Hs", Hs_);
writeEntry(os, "Ts", Ts_);
writeEntryIfDifferent(os, "gamma", scalar(3.3), gamma_);
}
// ************************************************************************* //

View File

@ -0,0 +1,155 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2024 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::GodaJONSWAP
Description
GodaJONSWAP wave spectrum. This is an alternative, approximate
parameterisation of the JONSWAP spectrum, in which the significant wave
height and period are specified instead of the wind speed and fetch.
References:
\verbatim
Goda, Y. (1988).
Statistical variability of sea state parameters as a function of wave
spectrum.
Coastal Engineering in Japan, 31(1), 39-52.
\endverbatim
\verbatim
Goda, Y. (2010).
Random seas and design of maritime structures.
World Scientific Publishing Company.
\endverbatim
See page 29 of the second reference for a convenient formulation.
Usage
\table
Property | Description | Required? | Default
Hs | The significant wave height [m] | yes |
Tp | The significant wave period [s] | yes |
gamma | Peaked-ness parameter | no | 3.3
\endtable
Example specification:
\verbatim
spectrum GodaJONSWAP;
GodaJONSWAPCoeffs
{
Hs 2;
Ts 6;
}
\endverbatim
See also
Foam::waveSpectra::JONSWAP
SourceFiles
GodaJONSWAP.C
\*---------------------------------------------------------------------------*/
#ifndef GodaJONSWAP_H
#define GodaJONSWAP_H
#include "waveSpectrum.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace waveSpectra
{
/*---------------------------------------------------------------------------*\
Class GodaJONSWAP Declaration
\*---------------------------------------------------------------------------*/
class GodaJONSWAP
:
public waveSpectrum
{
// Private Data
//- Significant Wave Height [m]
const scalar Hs_;
//- Significant wave period [s]
const scalar Ts_;
//- Peaked-ness parameter
const scalar gamma_;
public:
//- Runtime type information
TypeName("GodaJONSWAP");
// Constructors
//- Construct a copy
GodaJONSWAP(const GodaJONSWAP& spectrum);
//- Construct from a dictionary and gravity
GodaJONSWAP(const dictionary& dict, const scalar g);
//- Construct a clone
virtual autoPtr<waveSpectrum> clone() const
{
return autoPtr<waveSpectrum>(new GodaJONSWAP(*this));
}
//- Destructor
virtual ~GodaJONSWAP();
// Member Functions
//- Evaluate the wave spectral density at the given frequencies [m^2/Hz]
virtual tmp<scalarField> S(const scalarField& f) const;
//- Return the frequency below which a given fraction of the spectrum's
// total energy falls []
virtual scalar fFraction(const scalar fraction) const;
//- Write
virtual void write(Ostream& os) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace waveSpectra
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //