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-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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"
34#include "volFields.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
41namespace functionObjects
42{
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 {
227 UOPstream toMaster(Pstream::masterNo(), pBufs);
228 toMaster << Uc << C << cellAddr_;
229 }
230
231 pBufs.finishedGathers();
232
233 if (Pstream::master())
234 {
235 const label nGlobalCells(cmptProduct(N_));
236 vectorField Uijk(nGlobalCells);
237 vectorField Cijk(nGlobalCells);
238
239 for (const int proci : Pstream::allProcs())
240 {
241 UIPstream fromProc(proci, pBufs);
242 vectorField Up;
244 labelList cellAddrp;
245 fromProc >> Up >> Cp >> cellAddrp;
246 UIndirectList<vector>(Uijk, cellAddrp) = Up;
247 UIndirectList<vector>(Cijk, cellAddrp) = Cp;
248 }
249
250 calcAndWriteSpectrum(Uijk, Cijk, c0_, deltaC_, N_, kappaNorm_);
251 }
252
253 }
254 else
255 {
256 vectorField Uijk(Uc, cellAddr_);
257 vectorField Cijk(C, cellAddr_);
258
259 calcAndWriteSpectrum(Uijk, Cijk, c0_, deltaC_, N_, kappaNorm_);
260 }
261
262 return true;
263}
264
265
266// ************************************************************************* //
scalar y
label k
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Graphite solid properties.
Definition: C.H:53
C()
Construct null.
Definition: C.C:43
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
Output to file stream, using an OSstream.
Definition: OFstream.H:57
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Buffers for inter-processor communications streams (UOPstream, UIPstream).
UPstream::rangeType allProcs() const noexcept
Range of ranks indices associated with PstreamBuffers.
void finishedGathers(const bool wait=true)
Mark all sends to master as done.
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
A List with indirect addressing. Like IndirectList but does not store addressing.
Definition: IndirectList.H:79
@ nonBlocking
"nonBlocking"
static constexpr int masterNo() noexcept
Process index of the master (always 0)
Definition: UPstream.H:451
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
T & ref()
Return reference to the managed object without nullptr checking.
Definition: autoPtr.H:161
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition: boundBox.H:86
void add(const boundBox &bb)
Extend to include the second box.
Definition: boundBoxI.H:191
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
static tmp< complexField > forwardTransform(const tmp< complexField > &field, const UList< int > &nn)
Definition: fft.C:243
Abstract base-class for Time/database function objects.
Calculates the energy spectrum for a structured IJK mesh.
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.
virtual void writeFileHeader(Ostream &os)
Output file header information.
virtual bool execute()
Execute, currently does nothing.
virtual bool write()
Write the energySpectrum.
virtual bool read(const dictionary &)
Read the field min/max data.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Computes the magnitude of the square of an input field.
Definition: magSqr.H:153
Computes the magnitude of an input field.
Definition: mag.H:153
Base class for writing single files from the function objects.
Definition: writeFile.H:120
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:295
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:269
splitCell * master() const
Definition: splitCell.H:113
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
const volScalarField & Cp
Definition: EEqn.H:7
engineTime & runTime
autoPtr< surfaceVectorField > Uf
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
constexpr scalar twoPi(2 *M_PI)
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
Cmpt cmptProduct(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
complexField ReComplexField(const UList< scalar > &re)
Create complex field from a list of real (using imag == 0)
Definition: complexField.C:118
uint8_t direction
Definition: direction.H:56
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
constexpr char tab
The tab '\t' character(0x09)
Definition: Ostream.H:52
labelList f(nPoints)
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(rndGen, 1e-3))
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
const vector L(dict.get< vector >("L"))
const Vector< label > N(dict.get< Vector< label > >("N"))