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 -------------------------------------------------------------------------------
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);
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  (
76  ComplexField(U, vector::zero),
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 = newFileAtTime(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 :
128  fvMeshFunctionObject(name, runTime, dict),
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  boundBox cellBb(mesh_.cellBb(0));
154 
155  const vector L(meshBb.span());
156  const vector nCellXYZ(cmptDivide(L, cellBb.span()));
157 
158  N_ = Vector<int>
159  (
160  round(nCellXYZ.x()),
161  round(nCellXYZ.z()),
162  round(nCellXYZ.z())
163  );
164 
165  // Check that the mesh is a structured box
166  vector cellDx(cellBb.span());
167  vector expectedMax(N_.x()*cellDx.x(), N_.y()*cellDx.y(), N_.z()*cellDx.z());
168  vector relativeSize(cmptDivide(L, expectedMax));
169  for (direction i = 0; i < 3; ++i)
170  {
171  if (mag(relativeSize[i] - 1) > 1e-3)
172  {
174  << name() << " function object is only appropriate for "
175  << "isotropic structured IJK meshes. Mesh extents: " << L
176  << ", computed IJK mesh extents: " << expectedMax
177  << exit(FatalError);
178  }
179  }
180  Log << "Mesh extents (deltax,deltay,deltaz): " << L << endl;
181  Log << "Number of cells (Nx,Ny,Nz): " << N_ << endl;
182 
183  // Map into i-j-k co-ordinates
184  const vectorField& C = mesh_.C();
185  c0_ = returnReduce(min(C), minOp<vector>());
186  deltaC_ = returnReduce(max(C), maxOp<vector>());
187  deltaC_ -= c0_;
188 
189  forAll(C, celli)
190  {
191  label i = round((C[celli].x() - c0_.x())/(deltaC_.x())*(N_.x() - 1));
192  label j = round((C[celli].y() - c0_.y())/(deltaC_.y())*(N_.y() - 1));
193  label k = round((C[celli].z() - c0_.z())/(deltaC_.z())*(N_.z() - 1));
194 
195  cellAddr_[celli] = k + j*N_.y() + i*N_.y()*N_.z();
196  }
199 
200  return true;
201 }
202 
205 {
206  return true;
207 }
208 
209 
211 {
212  const auto& U = mesh_.lookupObject<volVectorField>(UName_);
213  const vectorField& Uc = U.primitiveField();
214  const vectorField& C = mesh_.C();
215 
216  if (Pstream::parRun())
217  {
219 
220  {
221  UOPstream toMaster(Pstream::masterNo(), pBufs);
222  toMaster << Uc << C << cellAddr_;
223  }
224 
225  pBufs.finishedGathers();
226 
227  if (Pstream::master())
228  {
229  const label nGlobalCells(cmptProduct(N_));
230  vectorField Uijk(nGlobalCells);
231  vectorField Cijk(nGlobalCells);
232 
233  for (const int proci : Pstream::allProcs())
234  {
235  UIPstream fromProc(proci, pBufs);
236  vectorField Up;
237  vectorField Cp;
238  labelList cellAddrp;
239  fromProc >> Up >> Cp >> cellAddrp;
240  UIndirectList<vector>(Uijk, cellAddrp) = Up;
241  UIndirectList<vector>(Cijk, cellAddrp) = Cp;
242  }
243 
244  calcAndWriteSpectrum(Uijk, Cijk, c0_, deltaC_, N_, kappaNorm_);
245  }
246 
247  }
248  else
249  {
250  vectorField Uijk(Uc, cellAddr_);
251  vectorField Cijk(C, cellAddr_);
252 
253  calcAndWriteSpectrum(Uijk, Cijk, c0_, deltaC_, N_, kappaNorm_);
254  }
255 
256  return true;
257 }
258 
259 
260 // ************************************************************************* //
virtual bool read(const dictionary &)
Read the field min/max data.
dictionary dict
Graphite solid properties.
Definition: C.H:46
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
defineTypeNameAndDebug(ObukhovLength, 0)
uint8_t direction
Definition: direction.H:46
Cmpt cmptProduct(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:597
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:339
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const vector L(dict.get< vector >("L"))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition: UPstream.H:1176
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
virtual bool execute()
Execute, currently does nothing.
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:49
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
const Cmpt & y() const noexcept
Access to the vector y component.
Definition: Vector.H:140
label k
Boltzmann constant.
Abstract base-class for Time/database function objects.
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
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.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Input inter-processor communications stream using MPI send/recv etc. - operating on external buffer...
Definition: UIPstream.H:287
scalar y
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:313
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
virtual bool write()
Write the energySpectrum.
C()
Construct null.
Definition: C.C:36
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
constexpr scalar twoPi(2 *M_PI)
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static constexpr int masterNo() noexcept
Relative rank for the master process - is always 0.
Definition: UPstream.H:1059
complexField ComplexField(const UList< scalar > &realValues, const UList< scalar > &imagValues)
Create complex field by zipping two lists of real/imag values.
Definition: complexField.C:150
Calculates the energy spectrum for a structured IJK mesh.
autoPtr< surfaceVectorField > Uf
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
Output inter-processor communications stream using MPI send/recv etc. - operating on external buffer...
Definition: UOPstream.H:395
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const Cmpt & x() const noexcept
Access to the vector x component.
Definition: Vector.H:135
const volScalarField & Cp
Definition: EEqn.H:7
OBJstream os(runTime.globalPath()/outputName)
volScalarField & C
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
static tmp< complexField > forwardTransform(const tmp< complexField > &field, const UList< int > &nn)
Definition: fft.C:236
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:241
Buffers for inter-processor communications streams (UOPstream, UIPstream).
virtual void writeFileHeader(Ostream &os)
Output file header information.
const Vector< label > N(dict.get< Vector< label >>("N"))
const Cmpt & z() const noexcept
Access to the vector z component.
Definition: Vector.H:145
U
Definition: pEqn.H:72
vector point
Point is a vector.
Definition: point.H:37
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:56
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
#define Log
Definition: PDRblock.C:28
"nonBlocking" : (MPI_Isend, MPI_Irecv)
virtual bool read(const dictionary &dict)
Read optional controls.
Field< vector > vectorField
Specialisation of Field<T> for vector.
energySpectrum(const energySpectrum &)=delete
No copy construct.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127