waveSpectrum.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) 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::waveSpectrum
26 
27 Description
28  Base class for wave spectra
29 
30 SourceFiles
31  waveSpectrum.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef waveSpectrum_H
36 #define waveSpectrum_H
37 
38 #include "dictionary.H"
39 #include "scalarField.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 /*---------------------------------------------------------------------------*\
47  Class waveSpectrum Declaration
48 \*---------------------------------------------------------------------------*/
49 
50 class waveSpectrum
51 {
52  // Private Data
53 
54  //- Gravitational acceleration [m/s^2]
55  const scalar g_;
56 
57  //- Number of intervals to use for integration
58  const label n_;
59 
60 
61  // Private Member Functions
62 
63  //- Refine or coarsen the given field to the given size, using linear
64  // interpolation to construct intermediate values. Guaranteed
65  // symmetric, so refined(x.size(), refined(n, x)) will return x (if n
66  // is larger than x.size()).
67  static tmp<scalarField> refined(const label n, const scalarField& x);
68 
69  //- Tmp variant of the above
70  static tmp<scalarField> refined
71  (
72  const label n,
73  const tmp<scalarField>& x
74  );
75 
76 
77 protected:
78 
79  // Protected Member Functions
80 
81  //- Return the frequency below which a given fraction of the spectrum's
82  // total energy falls, given also the maximum frequency with
83  // appreciable energy
84  scalar fFraction(const scalar fraction, const scalar f1) const;
85 
86 
87 public:
88 
89  //- Runtime type information
90  TypeName("waveSpectrum");
91 
92 
93  // Declare runtime construction
95  (
96  autoPtr,
98  dictionary,
99  (const dictionary& dict, const scalar g),
100  (dict, g)
101  );
102 
103 
104  // Constructors
105 
106  //- Construct a copy
107  waveSpectrum(const waveSpectrum& spectrum);
108 
109  //- Construct from a dictionary and gravity
110  waveSpectrum(const dictionary& dict, const scalar g);
111 
112  //- Construct a clone
113  virtual autoPtr<waveSpectrum> clone() const = 0;
114 
115 
116  // Selectors
117 
118  //- Select given a dictionary and gravity
120  (
121  const dictionary& dict,
122  const scalar g
123  );
124 
125 
126  //- Destructor
127  virtual ~waveSpectrum();
128 
129 
130  // Member Functions
131 
132  //- Access the gravitation acceleration [m/s^2]
133  inline scalar g() const
134  {
135  return g_;
136  }
137 
138  //- Evaluate the wave spectral density at the given frequencies [m^2/Hz]
139  virtual tmp<scalarField> S(const scalarField& f) const = 0;
140 
141  //- Evaluate the integral of the wave spectral density at the given
142  // frequencies [m^2]
143  virtual tmp<scalarField> integralS(const scalarField& f) const;
144 
145  //- Evaluate the integral of the wave spectral density multiplied by
146  // the frequency at the given frequencies [Hz m^2]
147  virtual tmp<scalarField> integralFS(const scalarField& f) const;
148 
149  //- Return the frequency below which a given fraction of the spectrum's
150  // total energy falls
151  virtual scalar fFraction(const scalar fraction) const = 0;
152 
153  //- Write
154  virtual void write(Ostream& os) const;
155 
156 
157  // Member Operators
158 
159  //- Disallow default bitwise assignment
160  void operator=(const waveSpectrum&) = delete;
161 };
162 
163 
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165 
166 } // End namespace Foam
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 #endif
171 
172 // ************************************************************************* //
label n
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
Base class for wave spectra.
Definition: waveSpectrum.H:50
virtual autoPtr< waveSpectrum > clone() const =0
Construct a clone.
scalar g() const
Access the gravitation acceleration [m/s^2].
Definition: waveSpectrum.H:132
declareRunTimeSelectionTable(autoPtr, waveSpectrum, dictionary,(const dictionary &dict, const scalar g),(dict, g))
TypeName("waveSpectrum")
Runtime type information.
waveSpectrum(const waveSpectrum &spectrum)
Construct a copy.
Definition: waveSpectrum.C:127
virtual tmp< scalarField > S(const scalarField &f) const =0
Evaluate the wave spectral density at the given frequencies [m^2/Hz].
virtual ~waveSpectrum()
Destructor.
Definition: waveSpectrum.C:171
void operator=(const waveSpectrum &)=delete
Disallow default bitwise assignment.
static autoPtr< waveSpectrum > New(const dictionary &dict, const scalar g)
Select given a dictionary and gravity.
Definition: waveSpectrum.C:144
scalar fFraction(const scalar fraction, const scalar f1) const
Return the frequency below which a given fraction of the spectrum's.
Definition: waveSpectrum.C:105
virtual void write(Ostream &os) const
Write.
Definition: waveSpectrum.C:207
virtual tmp< scalarField > integralFS(const scalarField &f) const
Evaluate the integral of the wave spectral density multiplied by.
Definition: waveSpectrum.C:193
virtual tmp< scalarField > integralS(const scalarField &f) const
Evaluate the integral of the wave spectral density at the given.
Definition: waveSpectrum.C:178
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
labelList f(nPoints)
dictionary dict