36Foam::label Foam::turbulentDFSEMInletFvPatchVectorField::seedIterMax_ = 1000;
40void Foam::turbulentDFSEMInletFvPatchVectorField::writeEddyOBJ()
const
44 OFstream
os(db().time().
path()/
"eddyBox.obj");
46 const polyPatch&
pp = this->
patch().patch();
47 const labelList& boundaryPoints =
pp.boundaryPoints();
50 const vector offset(patchNormal_*maxSigmaX_);
53 point p = localPoints[boundaryPoints[i]];
55 os <<
"v " <<
p.
x() <<
" " <<
p.y() <<
" " <<
p.z() <<
nl;
60 point p = localPoints[boundaryPoints[i]];
62 os <<
"v " <<
p.
x() <<
" " <<
p.y() <<
" " <<
p.z() <<
nl;
67 const Time& time = db().time();
70 time.path()/
"eddies_" +
Foam::name(time.timeIndex()) +
".obj"
73 label pointOffset = 0;
76 const eddy&
e = eddies_[eddyI];
77 pointOffset +=
e.writeSurfaceOBJ(pointOffset, patchNormal_,
os);
83void Foam::turbulentDFSEMInletFvPatchVectorField::writeLumleyCoeffs()
const
91 const scalar t = db().time().timeOutputValue();
107 const scalar xi =
cbrt(0.5*iii);
108 const scalar eta =
sqrt(-ii/3.0);
115void Foam::turbulentDFSEMInletFvPatchVectorField::initialisePatch()
128 <<
"Patch " <<
patch().name() <<
" is not planar"
132 patchNormal_ /=
mag(patchNormal_) + ROOTVSMALL;
141 for (
const face&
f : patch)
143 nTris +=
f.nTriangles();
156 for (
const auto& t : tris)
158 dynTriFace.emplace_back(t[0], t[1], t[2], facei);
163 triFace_.transfer(dynTriFace);
172 triCumulativeMagSf_.resize_nocopy(triFace_.size()+1);
174 auto iter = triCumulativeMagSf_.begin();
177 scalar patchArea = 0;
181 for (
const auto& t : triFace_)
183 patchArea += t.mag(
points);
187 sumTriMagSf_.resize_nocopy(numProc+1);
192 slice[myProci] = patchArea;
197 for (label i = 1; i < sumTriMagSf_.size(); ++i)
199 sumTriMagSf_[i] += sumTriMagSf_[i-1];
204 patchArea_ = sumTriMagSf_.back();
208 patchBounds_.inflate(0.1);
218void Foam::turbulentDFSEMInletFvPatchVectorField::initialiseEddyBox()
222 const scalarField L(L_->value(db().time().timeOutputValue())/Lref_);
230 scalar&
s = sigmax_[faceI];
235 s =
min(
mag(
L[faceI]), kappa_*delta_);
236 s =
max(
s, nCellPerEddy_*cellDx[faceI]);
240 maxSigmaX_ =
gMax(sigmax_);
243 v0_ = 2*
gSum(magSf)*maxSigmaX_;
247 Info<<
"Patch: " <<
patch().patch().name() <<
" eddy box:" <<
nl
248 <<
" volume : " << v0_ <<
nl
249 <<
" maxSigmaX : " << maxSigmaX_ <<
nl
267 const scalar areaFraction =
268 rndGen_.globalPosition<scalar>(0, patchArea_);
274 if (areaFraction >= sumTriMagSf_[i])
285 const scalar offset = sumTriMagSf_[proci];
288 if (areaFraction > triCumulativeMagSf_[i] + offset)
300 const scalar maxAreaLimit = triCumulativeMagSf_.back();
301 const scalar areaFraction = rndGen_.position<scalar>(0, maxAreaLimit);
305 if (areaFraction > triCumulativeMagSf_[i])
320 triFace_[triI].tri(
points).randomPoint(rndGen_),
321 triFace_[triI].index()
330void Foam::turbulentDFSEMInletFvPatchVectorField::initialiseEddies()
332 const scalar t = db().time().timeOutputValue();
338 scalar sumVolEddy = 0;
339 scalar sumVolEddyAllProc = 0;
341 while (sumVolEddyAllProc/v0_ < d_)
346 while (
search && iter++ < seedIterMax_)
350 const label patchFaceI =
pos.index();
353 if (patchFaceI != -1)
359 rndGen_.position<scalar>(-maxSigmaX_, maxSigmaX_),
366 if (
e.patchFaceI() != -1)
369 sumVolEddy +=
e.volume();
381 eddies_.transfer(eddies);
383 nEddy_ = eddies_.size();
387 Pout<<
"Patch:" <<
patch().patch().name();
394 Pout<<
" seeded:" << nEddy_ <<
" eddies" <<
endl;
401 Info<<
"Turbulent DFSEM patch: " <<
patch().name()
402 <<
" seeded " << nEddy_ <<
" eddies with total volume "
409 <<
"Patch: " <<
patch().patch().name()
410 <<
" on field " << internalField().name()
411 <<
": No eddies seeded - please check your set-up"
417void Foam::turbulentDFSEMInletFvPatchVectorField::convectEddies
423 const scalar t = db().time().timeOutputValue();
432 eddy&
e = eddies_[eddyI];
433 e.move(deltaT*(UBulk & patchNormal_));
435 const scalar position0 =
e.x();
438 if (position0 > maxSigmaX_)
443 while (
search && iter++ < seedIterMax_)
447 const label patchFaceI =
pos.index();
453 position0 - floor(position0/maxSigmaX_)*maxSigmaX_,
459 if (
e.patchFaceI() != -1)
475 Info<<
"Patch: " <<
patch().patch().name()
476 <<
" recycled " << nRecycled <<
" eddies"
483Foam::vector Foam::turbulentDFSEMInletFvPatchVectorField::uPrimeEddy
486 const point& patchFaceCf
493 const eddy&
e = eddies[
k];
494 uPrime +=
e.uPrime(patchFaceCf, patchNormal_);
501void Foam::turbulentDFSEMInletFvPatchVectorField::calcOverlappingProcEddies
518 const eddy&
e = eddies_[i];
521 const point x(
e.position(patchNormal_));
531 if (ebb.overlaps(patchBBs[proci]))
533 sendMap[proci].push_back(i);
554 for (
const int proci : pBufs.allProcs())
560 is >> overlappingEddies[proci];
592 triCumulativeMagSf_(),
600 sigmax_(size(),
Zero),
619 U_(ptf.U_.clone(patch().patch())),
620 R_(ptf.R_.clone(patch().patch())),
621 L_(ptf.L_.clone(patch().patch())),
629 nCellPerEddy_(ptf.nCellPerEddy_),
631 patchArea_(ptf.patchArea_),
632 triFace_(ptf.triFace_),
633 triCumulativeMagSf_(ptf.triCumulativeMagSf_),
634 sumTriMagSf_(ptf.sumTriMagSf_),
635 patchNormal_(ptf.patchNormal_),
636 patchBounds_(ptf.patchBounds_),
638 eddies_(ptf.eddies_),
640 rndGen_(ptf.rndGen_),
641 sigmax_(ptf.sigmax_, mapper),
642 maxSigmaX_(ptf.maxSigmaX_),
645 singleProc_(ptf.singleProc_),
646 writeEddies_(ptf.writeEddies_)
669 nCellPerEddy_(
dict.getOrDefault<label>(
"nCellPerEddy", 5)),
673 triCumulativeMagSf_(),
681 sigmax_(size(),
Zero),
686 writeEddies_(
dict.getOrDefault(
"writeEddies", false))
690 const scalar t = db().time().timeOutputValue();
705 U_(ptf.U_.clone(patch().patch())),
706 R_(ptf.R_.clone(patch().patch())),
707 L_(ptf.L_.clone(patch().patch())),
715 nCellPerEddy_(ptf.nCellPerEddy_),
717 patchArea_(ptf.patchArea_),
718 triFace_(ptf.triFace_),
719 triCumulativeMagSf_(ptf.triCumulativeMagSf_),
720 sumTriMagSf_(ptf.sumTriMagSf_),
721 patchNormal_(ptf.patchNormal_),
722 patchBounds_(ptf.patchBounds_),
724 eddies_(ptf.eddies_),
726 rndGen_(ptf.rndGen_),
727 sigmax_(ptf.sigmax_),
728 maxSigmaX_(ptf.maxSigmaX_),
731 singleProc_(ptf.singleProc_),
732 writeEddies_(ptf.writeEddies_)
753 constexpr label maxDiffs = 5;
761 if (maxDiffs < nDiffs)
763 Info<<
"More than " << maxDiffs <<
" times"
764 <<
" Reynolds-stress realizability checks failed."
765 <<
" Skipping further comparisons." <<
endl;
774 <<
"Reynolds stress " << r <<
" at index " << i
775 <<
" does not obey the constraint: Rxx >= 0"
780 if ((r.xx()*r.yy() -
sqr(r.xy())) < 0)
783 <<
"Reynolds stress " << r <<
" at index " << i
784 <<
" does not obey the constraint: Rxx*Ryy - sqr(Rxy) >= 0"
792 <<
"Reynolds stress " << r <<
" at index " << i
793 <<
" does not obey the constraint: det(R) >= 0"
814 <<
"Reynolds stresses contain at least one "
815 <<
"nonpositive element. min(R) = " <<
min(
R)
823 const fvPatchFieldMapper& m
853 const auto& dfsemptf =
858 U_->rmap(dfsemptf.U_(), addr);
862 R_->rmap(dfsemptf.R_(), addr);
866 L_->rmap(dfsemptf.L_(), addr);
869 sigmax_.rmap(dfsemptf.sigmax_, addr);
880 if (curTimeIndex_ == -1)
890 if (curTimeIndex_ != db().time().
timeIndex())
893 U_->value(db().time().timeOutputValue())/Uref_;
902 const scalar deltaT = db().time().deltaTValue();
903 convectEddies(UBulk, deltaT);
922 U[faceI] +=
c*uPrimeEddy(eddies_, Cf[faceI]);
930 U[faceI] +=
c*uPrimeEddy(eddies_, Cf[faceI]);
935 calcOverlappingProcEddies(overlappingEddies);
937 for (
const List<eddy>& eddies : overlappingEddies)
943 U[faceI] +=
c*uPrimeEddy(eddies, Cf[faceI]);
951 gSum((UBulk & patchNormal_)*
patch().magSf())
956 curTimeIndex_ = db().time().timeIndex();
967 Info<<
"Magnitude of bulk velocity: " << UBulk <<
endl;
969 Info<<
"Number of eddies: "
973 Info<<
"Patch:" <<
patch().patch().name()
974 <<
" min/max(U):" <<
limits.min() <<
", " <<
limits.max()
977 if (db().time().writeTime())
984 fixedValueFvPatchVectorField::updateCoeffs();
1023 turbulentDFSEMInletFvPatchVectorField
#define R(A, B, C, D, E, F, K, M)
Macros for easy insertion into run-time selection tables.
volVectorField UMean(UMeanHeader, mesh)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
SubList< scalar > subList
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void finishedSends(const bool wait=true)
Mark the send phase as being finished.
Inter-processor communications stream.
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual communication.
const Cmpt & yy() const noexcept
const Cmpt & xx() const noexcept
const Cmpt & xy() const noexcept
Input inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Output inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
static int & msgType() noexcept
Message tag of standard messages.
static void reduceAnd(bool &value, const int communicator=worldComm)
Logical (and) reduction (MPI_AllReduce).
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
const Cmpt & x() const noexcept
Access to the vector x component.
A bounding box defined in terms of min/max extrema points.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Class to describe eddies for the turbulentDFSEMInletFvPatchVectorField boundary condition.
static int debug
Flag to activate debug statements.
Class to handle errors and exceptions in a simple, consistent stream-based manner.
A face is a list of labels corresponding to mesh vertices.
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
A FieldMapper for finite-volume patch fields.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
A patch is a list of labels that address the faces in the global face list.
A class for managing temporary objects.
The turbulentDFSEMInlet is a synthesised-eddy based velocity inlet boundary condition to generate syn...
static void checkStresses(const symmTensorField &R)
Check if input Reynolds stresses are valid.
virtual void rmap(const fvPatchVectorField &ptf, const labelList &addr)
Reverse map the given fvPatchField onto this fvPatchField.
turbulentDFSEMInletFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void write(Ostream &) const
Write.
virtual tmp< fvPatchField< vector > > clone() const
Return a clone.
virtual void autoMap(const fvPatchFieldMapper &m)
Map (and resize as needed) from self given a mapping object.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for handling debugging switches.
const std::string patch
OpenFOAM patch number as a std::string.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar det(const dimensionedSphericalTensor &dt)
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
dimensionedScalar pos(const dimensionedScalar &ds)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Type gSum(const FieldField< Field, Type > &f)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Type gWeightedAverage(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted average of a field, using the mag() of the weights.
List< label > labelList
A List of labels.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Cmpt invariantII(const SymmTensor2D< Cmpt > &st)
Return the 2nd invariant of a SymmTensor2D.
messageStream Info
Information stream (stdout output on master, null elsewhere).
MinMax< scalar > scalarMinMax
A scalar min/max range.
Cmpt invariantIII(const SymmTensor< Cmpt > &st)
Return the 3rd invariant of a SymmTensor.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Type gMax(const FieldField< Field, Type > &f)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
fvPatchField< vector > fvPatchVectorField
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
fileName search(const word &file, const fileName &directory)
Recursively search the given directory for the file.
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.
#define forAllReverse(list, i)
Reverse loop across all elements in list.
const vector L(dict.get< vector >("L"))