41namespace functionObjects
68 const scalar kappaNorm
83 Log <<
"Computing wave numbers" <<
endl;
88 kappa[kappai] = kappai*kappaNorm;
91 Log <<
"Computing energy spectrum" <<
endl;
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)));
102 E[kappai] += Ec[celli];
107 Log <<
"Writing spectrum" <<
endl;
114 os << kappa[kappai] <<
tab << E[kappai] <<
nl;
130 cellAddr_(mesh_.nCells()),
148 dict.readIfPresent(
"U", UName_);
153 const cell& c = mesh_.cells()[0];
157 const face&
f = mesh_.faces()[c[facei]];
158 cellBb.
add(mesh_.points(),
f);
173 vector expectedMax(N_.x()*cellDx.
x(), N_.y()*cellDx.
y(), N_.z()*cellDx.
z());
177 if (
mag(relativeSize[i] - 1) > 1
e-3)
180 <<
name() <<
" function object is only appropriate for "
181 <<
"isotropic structured IJK meshes. Mesh extents: " <<
L
182 <<
", computed IJK mesh extents: " << expectedMax
186 Log <<
"Mesh extents (deltax,deltay,deltaz): " <<
L <<
endl;
187 Log <<
"Number of cells (Nx,Ny,Nz): " << N_ <<
endl;
193 deltaC_ = cMax - c0_;
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));
201 cellAddr_[celli] =
k + j*N_.y() + i*N_.y()*N_.z();
228 toMaster << Uc <<
C << cellAddr_;
245 fromProc >> Up >>
Cp >> cellAddrp;
250 calcAndWriteSpectrum(Uijk, Cijk, c0_, deltaC_, N_, kappaNorm_);
259 calcAndWriteSpectrum(Uijk, Cijk, c0_, deltaC_, N_, kappaNorm_);
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.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Output to file stream, using an OSstream.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
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.
A List with indirect addressing. Like IndirectList but does not store addressing.
@ nonBlocking
"nonBlocking"
static constexpr int masterNo() noexcept
Process index of the master (always 0)
static bool & parRun() noexcept
Test if this a parallel run.
const Cmpt & z() const
Access to the vector z component.
const Cmpt & y() const
Access to the vector y component.
const Cmpt & x() const
Access to the vector x component.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
T & ref()
Return reference to the managed object without nullptr checking.
A bounding box defined in terms of min/max extrema points.
const point & min() const
Minimum describing the bounding box.
const point & max() const
Maximum describing the bounding box.
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
void add(const boundBox &bb)
Extend to include the second box.
A cell is defined as a list of faces with extra functionality.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
A face is a list of labels corresponding to mesh vertices.
static tmp< complexField > forwardTransform(const tmp< complexField > &field, const UList< int > &nn)
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.
Computes the magnitude of an input field.
Base class for writing single files from the function objects.
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
splitCell * master() const
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
const volScalarField & Cp
autoPtr< surfaceVectorField > Uf
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
constexpr scalar twoPi(2 *M_PI)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
Cmpt cmptProduct(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
complexField ReComplexField(const UList< scalar > &re)
Create complex field from a list of real (using imag == 0)
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
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)
constexpr char nl
The newline '\n' character (0x0a)
constexpr char tab
The tab '\t' character(0x09)
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(rndGen, 1e-3))
#define forAll(list, i)
Loop across all elements in list.
const vector L(dict.get< vector >("L"))
const Vector< label > N(dict.get< Vector< label > >("N"))