energySpectrum.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) 2018-2020 OpenCFD Ltd.
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 "energySpectrum.H"
29 #include "fft.H"
30 #include "fvMesh.H"
31 #include "boundBox.H"
32 #include "OFstream.H"
33 #include "mathematicalConstants.H"
34 #include "volFields.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace functionObjects
42 {
43  defineTypeNameAndDebug(energySpectrum, 0);
44  addToRunTimeSelectionTable(functionObject, energySpectrum, dictionary);
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
50 
52 {
53  writeHeader(os, "Turbulence energy spectra");
54 
55  writeCommented(os, "kappa E(kappa)");
56 
57  os << endl;
58 }
59 
60 
62 (
63  const vectorField& U,
64  const vectorField& C,
65  const vector& c0,
66  const vector& deltaC,
67  const Vector<int>& N,
68  const scalar kappaNorm
69 )
70 {
71  Log << "Computing FFT" << endl;
73  (
75  (
77  List<int>({N.x(), N.y(), N.z()})
78  )
79  /scalar(cmptProduct(N))
80  );
81 
82 
83  Log << "Computing wave numbers" << endl;
84  const label nMax = cmptMax(N);
85  scalarField kappa(nMax);
86  forAll(kappa, kappai)
87  {
88  kappa[kappai] = kappai*kappaNorm;
89  }
90 
91  Log << "Computing energy spectrum" << endl;
92  scalarField E(nMax, Zero);
93  const scalarField Ec(0.5*magSqr(Uf));
94  forAll(C, celli)
95  {
96  point Cc(C[celli] - c0);
97  label i = round((Cc.x() - 0.5*deltaC.x())/(deltaC.x())*(N.x() - 1));
98  label j = round((Cc.y() - 0.5*deltaC.y())/(deltaC.y())*(N.y() - 1));
99  label k = round((Cc.z() - 0.5*deltaC.z())/(deltaC.z())*(N.z() - 1));
100  label kappai = round(Foam::sqrt(scalar(i*i + j*j + k*k)));
101 
102  E[kappai] += Ec[celli];
103  }
104 
105  E /= kappaNorm;
106 
107  Log << "Writing spectrum" << endl;
108  autoPtr<OFstream> osPtr = createFile(name(), time_.value());
109  OFstream& os = osPtr.ref();
110  writeFileHeader(os);
111 
112  forAll(kappa, kappai)
113  {
114  os << kappa[kappai] << tab << E[kappai] << nl;
115  }
116 }
117 
118 
119 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
120 
122 (
123  const word& name,
124  const Time& runTime,
125  const dictionary& dict
126 )
127 :
129  writeFile(mesh_, name),
130  cellAddr_(mesh_.nCells()),
131  UName_("U"),
132  N_(Zero),
133  c0_(Zero),
134  deltaC_(Zero),
135  kappaNorm_(0)
136 {
137  read(dict);
138 }
139 
140 
141 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
142 
144 {
147 
148  dict.readIfPresent("U", UName_);
149 
150  const boundBox meshBb(mesh_.bounds());
151 
152  // Assume all cells are the same size...
153  const cell& c = mesh_.cells()[0];
155  forAll(c, facei)
156  {
157  const face& f = mesh_.faces()[c[facei]];
158  cellBb.add(mesh_.points(), f);
159  }
160 
161  const vector L(meshBb.max() - meshBb.min());
162  const vector nCellXYZ(cmptDivide(L, cellBb.max() - cellBb.min()));
163 
164  N_ = Vector<int>
165  (
166  round(nCellXYZ.x()),
167  round(nCellXYZ.z()),
168  round(nCellXYZ.z())
169  );
170 
171  // Check that the mesh is a structured box
172  vector cellDx(cellBb.max() - cellBb.min());
173  vector expectedMax(N_.x()*cellDx.x(), N_.y()*cellDx.y(), N_.z()*cellDx.z());
174  vector relativeSize(cmptDivide(L, expectedMax));
175  for (direction i = 0; i < 3; ++i)
176  {
177  if (mag(relativeSize[i] - 1) > 1e-3)
178  {
180  << name() << " function object is only appropriate for "
181  << "isotropic structured IJK meshes. Mesh extents: " << L
182  << ", computed IJK mesh extents: " << expectedMax
183  << exit(FatalError);
184  }
185  }
186  Log << "Mesh extents (deltax,deltay,deltaz): " << L << endl;
187  Log << "Number of cells (Nx,Ny,Nz): " << N_ << endl;
188 
189  // Map into i-j-k co-ordinates
190  const vectorField& C = mesh_.C();
191  c0_ = returnReduce(min(C), minOp<vector>());
192  const vector cMax = returnReduce(max(C), maxOp<vector>());
193  deltaC_ = cMax - c0_;
194 
195  forAll(C, celli)
196  {
197  label i = round((C[celli].x() - c0_.x())/(deltaC_.x())*(N_.x() - 1));
198  label j = round((C[celli].y() - c0_.y())/(deltaC_.y())*(N_.y() - 1));
199  label k = round((C[celli].z() - c0_.z())/(deltaC_.z())*(N_.z() - 1));
200 
201  cellAddr_[celli] = k + j*N_.y() + i*N_.y()*N_.z();
202  }
203 
205 
206  return true;
207 }
208 
209 
211 {
212  return true;
213 }
214 
215 
217 {
218  const auto& U = mesh_.lookupObject<volVectorField>(UName_);
219  const vectorField& Uc = U.primitiveField();
220  const vectorField& C = mesh_.C();
221 
222  if (Pstream::parRun())
223  {
225 
226  UOPstream toProc(0, pBufs);
227 
228  toProc << Uc << C << cellAddr_;
229 
230  pBufs.finishedSends();
231 
232  if (Pstream::master())
233  {
234  const label nGlobalCells(cmptProduct(N_));
235  vectorField Uijk(nGlobalCells);
236  vectorField Cijk(nGlobalCells);
237 
238  for (const int proci : Pstream::allProcs())
239  {
240  UIPstream fromProc(proci, pBufs);
241  vectorField Up;
242  vectorField Cp;
243  labelList cellAddrp;
244  fromProc >> Up >> Cp >> cellAddrp;
245  UIndirectList<vector>(Uijk, cellAddrp) = Up;
246  UIndirectList<vector>(Cijk, cellAddrp) = Cp;
247  }
248 
249  calcAndWriteSpectrum(Uijk, Cijk, c0_, deltaC_, N_, kappaNorm_);
250  }
251 
252  }
253  else
254  {
255  vectorField Uijk(Uc, cellAddr_);
256  vectorField Cijk(C, cellAddr_);
257 
258  calcAndWriteSpectrum(Uijk, Cijk, c0_, deltaC_, N_, kappaNorm_);
259  }
260 
261  return true;
262 }
263 
264 
265 // ************************************************************************* //
volFields.H
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::maxOp
Definition: ops.H:223
L
const vector L(dict.get< vector >("L"))
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
mathematicalConstants.H
Log
#define Log
Definition: PDRblock.C:35
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::functionObjects::energySpectrum::calcAndWriteSpectrum
void calcAndWriteSpectrum(const vectorField &U, const vectorField &C, const vector &c0, const vector &deltaC, const Vector< int > &N, const scalar kappaNorm)
Calculate and write the spectrum.
Definition: energySpectrum.C:62
Foam::UOPstream
Output inter-processor communications stream operating on external buffer.
Definition: UOPstream.H:57
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::minOp
Definition: ops.H:224
Foam::boundBox::invertedBox
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition: boundBox.H:86
Foam::functionObjects::energySpectrum::read
virtual bool read(const dictionary &)
Read the field min/max data.
Definition: energySpectrum.C:143
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::PstreamBuffers
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Definition: PstreamBuffers.H:88
Uf
autoPtr< surfaceVectorField > Uf
Definition: createUfIfPresent.H:33
Foam::cmptProduct
Cmpt cmptProduct(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:596
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
meshBb
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(rndGen, 1e-3))
Foam::constant::electromagnetic::kappa
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
OFstream.H
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::functionObjects::writeFile::writeHeader
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:298
Foam::ReComplexField
complexField ReComplexField(const UList< scalar > &re)
Create complex field from a list of real (using imag == 0)
Definition: complexField.C:118
Foam::cmptMax
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:253
Foam::functionObjects::writeFile::read
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:213
Foam::functionObjects::energySpectrum::execute
virtual bool execute()
Execute, currently does nothing.
Definition: energySpectrum.C:210
Foam::functionObjects::energySpectrum::writeFileHeader
virtual void writeFileHeader(Ostream &os)
Output file header information.
Definition: energySpectrum.C:51
Foam::Field< vector >
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::PstreamBuffers::finishedSends
void finishedSends(const bool block=true)
Mark all sends as having been done.
Definition: PstreamBuffers.C:73
Foam::constant::mathematical::twoPi
constexpr scalar twoPi(2 *M_PI)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:173
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::C::C
C()
Construct null.
Definition: C.C:43
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::functionObjects::energySpectrum::energySpectrum
energySpectrum(const energySpectrum &)=delete
No copy construct.
Foam::functionObjects::energySpectrum::write
virtual bool write()
Write the energySpectrum.
Definition: energySpectrum.C:216
boundBox.H
Foam::cmptDivide
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
U
U
Definition: pEqn.H:72
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::UPstream::allProcs
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition: UPstream.H:508
Foam::functionObjects::writeFile::writeCommented
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:272
Foam::UPstream::commsTypes::nonBlocking
Foam::tab
constexpr char tab
Definition: Ostream.H:403
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::autoPtr::ref
T & ref()
Return reference to the managed object without nullptr checking.
Definition: autoPtr.H:161
f
labelList f(nPoints)
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::Vector< scalar >
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List< int >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::fft::forwardTransform
static tmp< complexField > forwardTransform(const tmp< complexField > &field, const UList< int > &nn)
Definition: fft.C:243
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::functionObjects::writeFile
Base class for writing single files from the function objects.
Definition: writeFile.H:119
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
N
const Vector< label > N(dict.get< Vector< label >>("N"))
Foam::UIPstream
Input inter-processor communications stream operating on external buffer.
Definition: UIPstream.H:56
Foam::C
Graphite solid properties.
Definition: C.H:50
fft.H
Foam::GeometricField< vector, fvPatchField, volMesh >
energySpectrum.H
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::boundBox::add
void add(const boundBox &bb)
Extend to include the second box.
Definition: boundBoxI.H:191