waveSpectrum.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) 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 \*---------------------------------------------------------------------------*/
25 
26 #include "waveSpectrum.H"
27 #include "unintegrable.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
35 }
36 
37 
38 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39 
40 Foam::tmp<Foam::scalarField> Foam::waveSpectrum::refined
41 (
42  const label n,
43  const scalarField& x
44 )
45 {
46  if (n == x.size())
47  {
48  return x;
49  }
50 
51  tmp<scalarField> txx(new scalarField(n, -vGreat));
52  scalarField& xx = txx.ref();
53 
54  if (n > x.size())
55  {
56  for (label i = 0; i < x.size() - 1; ++ i)
57  {
58  const label j0 = i*(n - 1)/(x.size() - 1);
59  const label j1 = (i + 1)*(n - 1)/(x.size() - 1);
60 
61  xx[j0] = x[i];
62 
63  for (label j = j0 + 1; j < j1; ++ j)
64  {
65  const scalar f = scalar(j - j0)/(j1 - j0);
66 
67  xx[j] = (1 - f)*x[i] + f*x[i + 1];
68  }
69  }
70 
71  xx.last() = x.last();
72  }
73  else
74  {
75  forAll(xx, j)
76  {
77  const label i = j*(x.size() - 1)/(n - 1);
78 
79  xx[j] = x[i];
80  }
81  }
82 
83  return txx;
84 }
85 
86 
87 Foam::tmp<Foam::scalarField> Foam::waveSpectrum::refined
88 (
89  const label n,
90  const tmp<scalarField>& x
91 )
92 {
93  if (n == x().size())
94  {
95  return x;
96  }
97 
98  return refined(n, x());
99 }
100 
101 
102 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
103 
105 (
106  const scalar fraction,
107  const scalar f1
108 ) const
109 {
110  // Integrate from zero to the maximum frequency
111  const scalarField ff(refined(n_, scalarField(scalarList({small*f1, f1}))));
112  const scalarField integralSS(integralS(ff));
113 
114  // Sample the integrated values for the frequency at the given fraction
115  return
117  (
118  ff,
119  integralSS/integralSS.last(),
120  fraction
121  );
122 }
123 
124 
125 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
126 
128 :
129  g_(spectrum.g_),
130  n_(spectrum.n_)
131 {}
132 
133 
135 :
136  g_(g),
137  n_((1<<dict.lookupOrDefault<label>("level", 16)) + 1)
138 {}
139 
140 
141 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
142 
144 Foam::waveSpectrum::New(const dictionary& dict, const scalar g)
145 {
146  const word type = dict.lookup<word>("spectrum");
147 
148  if (debug)
149  {
150  Info<< "Selecting " << waveSpectrum::typeName << " " << type << endl;
151  }
152 
153  dictionaryConstructorTable::iterator cstrIter =
154  dictionaryConstructorTablePtr_->find(type);
155 
156  if (cstrIter == dictionaryConstructorTablePtr_->end())
157  {
159  << "Unknown " << waveSpectrum::typeName << " " << type << nl << nl
160  << "Valid spectrum types are:" << nl
161  << dictionaryConstructorTablePtr_->sortedToc()
162  << exit(FatalError);
163  }
164 
165  return cstrIter()(dict.optionalSubDict(type + "Coeffs"), g);
166 }
167 
168 
169 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
170 
172 {}
173 
174 
175 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
176 
178 (
179  const scalarField& f
180 ) const
181 {
182  const scalarField ff(refined(n_, f));
183  return
184  refined
185  (
186  f.size(),
188  );
189 }
190 
191 
193 (
194  const scalarField& f
195 ) const
196 {
197  const scalarField ff(refined(n_, f));
198  return
199  refined
200  (
201  f.size(),
203  );
204 }
205 
206 
208 {}
209 
210 
211 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
T & last()
Return the last element of the list.
Definition: UListI.H:128
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:921
scalar sample() const
Sample the distribution.
Definition: unintegrable.C:351
static tmp< scalarField > integrate(const scalarField &x, const scalarField &y)
Integrate the values y with respect to the coordinates x.
Definition: unintegrable.C:34
static tmp< scalarField > integrateX(const scalarField &x, const scalarField &y)
Integrate the values x*y with respect to the coordinates x.
Definition: unintegrable.C:54
A class for managing temporary objects.
Definition: tmp.H:55
Base class for wave spectra.
Definition: waveSpectrum.H:50
waveSpectrum(const waveSpectrum &spectrum)
Construct a copy.
Definition: waveSpectrum.C:127
virtual ~waveSpectrum()
Destructor.
Definition: waveSpectrum.C:171
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
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
volScalarField scalarField(fieldObject, mesh)
const FieldField< fvPatchField, Type > & ff(const FieldField< fvPatchField, Type > &bf)
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
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
dimensionedScalar j1(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
messageStream Info
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
defineTypeNameAndDebug(combustionModel, 0)
error FatalError
static const char nl
Definition: Ostream.H:266
dimensionedScalar j0(const dimensionedScalar &ds)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
labelList f(nPoints)
dictionary dict