kShellIntegration.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "kShellIntegration.H"
29 #include "mathematicalConstants.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
34 (
35  const complexVectorField& Ek,
36  const Kmesh& K
37 )
38 {
39  // evaluate the radial component of the spectra as an average
40  // over the shells of thickness dk
41 
42  graph kShellMeanEk = kShellMean(Ek, K);
43  const scalarField& x = kShellMeanEk.x();
44  scalarField& y = *kShellMeanEk.begin()();
45 
46  // now multiply by 4pi k^2 (the volume of each shell) to get the
47  // spectra E(k). int E(k) dk is now the total energy in a box
48  // of side 2pi
49 
51 
52  // now scale this to get the energy in a box of side l0
53 
54  scalar l0(K.sizeOfBox()[0]*(scalar(K.nn()[0])/(scalar(K.nn()[0])-1.0)));
55  scalar factor = pow((l0/(2.0*constant::mathematical::pi)),3.0);
56 
57  y *= factor;
58 
59  // and divide by the number of points in the box, to give the
60  // energy density.
61 
62  y /= scalar(K.size());
63 
64  return kShellMeanEk;
65 }
66 
67 
68 // kShellMean : average over the points in a k-shell to evaluate the
69 // radial part of the energy spectrum.
70 
72 (
73  const complexVectorField& Ek,
74  const Kmesh& K
75 )
76 {
77  const label tnp = Ek.size();
78  const label NoSubintervals = label
79  (
80  pow(scalar(tnp), 1.0/vector::dim)*pow(1.0/vector::dim, 0.5) - 0.5
81  );
82 
83  scalarField k1D(NoSubintervals);
84  scalarField Ek1D(NoSubintervals);
85  scalarField EWeight(NoSubintervals);
86 
87  scalar kmax = K.max()*pow(1.0/vector::dim,0.5);
88  scalar delta_k = kmax/(NoSubintervals);
89 
90  forAll(Ek1D, a)
91  {
92  k1D[a] = (a + 1)*delta_k;
93  Ek1D[a] = 0.0;
94  EWeight[a] = 0;
95  }
96 
97  forAll(K, l)
98  {
99  scalar kmag = mag(K[l]);
100 
101  for (label a=0; a<NoSubintervals; a++)
102  {
103  if
104  (
105  kmag <= ((a + 1)*delta_k + delta_k/2.0)
106  && kmag > ((a + 1)*delta_k - delta_k/2.0)
107  )
108  {
109  scalar dist = delta_k/2.0 - mag((a + 1)*delta_k - kmag);
110 
111  Ek1D[a] += dist*
112  magSqr
113  (
114  vector
115  (
116  mag(Ek[l].x()),
117  mag(Ek[l].y()),
118  mag(Ek[l].z())
119  )
120  );
121 
122  EWeight[a] += dist;
123  }
124  }
125  }
126 
127  for (label a=0; a<NoSubintervals; a++)
128  {
129  if (EWeight[a] > 0)
130  {
131  Ek1D[a] /= EWeight[a];
132  }
133  }
134 
135  return graph("E(k)", "k", "E(k)", k1D, Ek1D);
136 }
137 
138 
139 // ************************************************************************* //
mathematicalConstants.H
Foam::graph
Class to create, store and output qgraph files.
Definition: graph.H:61
Ek
autoPtr< Function1< scalar > > Ek(Function1< scalar >::New("Ek", dict, &runTime))
Foam::graph::x
const scalarField & x() const
Definition: graph.H:168
Foam::Kmesh
Calculate the wavenumber vector field corresponding to the space vector field of a finite volume mesh...
Definition: Kmesh.H:51
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:58
Foam::kShellMean
graph kShellMean(const complexVectorField &Ek, const Kmesh &K)
Definition: kShellIntegration.C:72
Foam::Field< complexVector >
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
kShellIntegration.H
Integrate a multi-dimensional complexVectorField in k-shells to create the 1D.
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::Vector< scalar >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::kShellIntegration
graph kShellIntegration(const complexVectorField &Ek, const Kmesh &K)
Definition: kShellIntegration.C:34
y
scalar y
Definition: LISASMDCalcMethod1.H:14