36 template<
class CloudType>
40 const Field<point>& points,
41 const Field<scalar>& area
49 Info<<
"Creating output file" <<
endl;
52 if (Pstream::master())
55 mkDir(this->outputTimeDir());
60 new OFstream(this->outputTimeDir()/(
type() +
".dat"))
64 <<
"# Source : " <<
type() <<
nl 65 <<
"# Bins : " << faces.size() <<
nl 66 <<
"# Total area : " <<
sum(area) <<
nl;
69 <<
"# Geometry :" << nl
72 <<
tab <<
"(Centre_x Centre_y Centre_z)" 81 <<
tab << faces[i].centre(points)
88 <<
"# Output format:" <<
nl;
93 word binId =
"bin_" + id;
99 <<
tab <<
"mass[" <<
id <<
"]" 100 <<
tab <<
"massFlowRate[" <<
id <<
"]" 108 template<
class CloudType>
111 const List<Field<point> >& polygons
119 label np = polygons[polyI].size();
124 "Foam::ParticleCollector<CloudType>::initPolygons()",
127 <<
"polygons must consist of at least 3 points" 134 label pointOffset = 0;
135 points_.setSize(nPoints);
136 faces_.setSize(polygons.size());
137 faceTris_.setSize(polygons.size());
138 area_.setSize(polygons.size());
141 const Field<point>& polyPoints = polygons[faceI];
142 face
f(
identity(polyPoints.size()) + pointOffset);
143 UIndirectList<point>(points_,
f) = polyPoints;
144 area_[faceI] =
f.mag(points_);
146 DynamicList<face> tris;
147 f.triangles(points_, tris);
148 faceTris_[faceI].transfer(tris);
150 faces_[faceI].transfer(
f);
152 pointOffset += polyPoints.size();
157 template<
class CloudType>
160 mode_ = mtConcentricCircle;
164 radius_ = this->coeffDict().lookup(
"radius");
172 refDir = this->coeffDict().lookup(
"refDir");
173 refDir -= normal_[0]*(normal_[0] & refDir);
174 refDir /=
mag(refDir);
181 vector tangent = vector::zero;
182 scalar magTangent = 0.0;
185 while (magTangent < SMALL)
187 vector v = rnd.vector01();
189 tangent = v - (v & normal_[0])*normal_[0];
190 magTangent =
mag(tangent);
193 refDir = tangent/magTangent;
197 scalar dThetaSector = 360.0/scalar(nS);
198 label intervalPerSector =
max(1, ceil(dThetaSector/dTheta));
199 dTheta = dThetaSector/scalar(intervalPerSector);
201 label nPointPerSector = intervalPerSector + 1;
203 label nPointPerRadius = nS*(nPointPerSector - 1);
204 label nPoint = radius_.size()*nPointPerRadius;
205 label nFace = radius_.size()*nS;
210 points_.setSize(nPoint);
211 faces_.setSize(nFace);
212 area_.setSize(nFace);
214 coordSys_ = cylindricalCS(
"coordSys", origin, normal_[0], refDir,
false);
216 List<label> ptIDs(
identity(nPointPerRadius));
223 label pointOffset = radI*nPointPerRadius + 1;
225 for (
label i = 0; i < nPointPerRadius; i++)
227 label pI = i + pointOffset;
229 points_[pI] = coordSys_.globalPosition(pCyl);
234 DynamicList<label> facePts(2*nPointPerSector);
239 for (
label secI = 0; secI < nS; secI++)
246 for (
label ptI = 0; ptI < nPointPerSector; ptI++)
248 label i = ptI + secI*(nPointPerSector - 1);
249 label id = ptIDs.fcIndex(i - 1) + 1;
253 label faceI = secI + radI*nS;
255 faces_[faceI] = face(facePts);
256 area_[faceI] = faces_[faceI].mag(points_);
261 for (
label secI = 0; secI < nS; secI++)
265 label offset = (radI - 1)*nPointPerRadius + 1;
267 for (
label ptI = 0; ptI < nPointPerSector; ptI++)
269 label i = ptI + secI*(nPointPerSector - 1);
270 label id = offset + ptIDs.fcIndex(i - 1);
273 for (
label ptI = nPointPerSector-1; ptI >= 0; ptI--)
275 label i = ptI + secI*(nPointPerSector - 1);
276 label id = offset + nPointPerRadius + ptIDs.fcIndex(i - 1);
280 label faceI = secI + radI*nS;
282 faces_[faceI] = face(facePts);
283 area_[faceI] = faces_[faceI].mag(points_);
290 template<
class CloudType>
297 label dummyNearType = -1;
298 label dummyNearLabel = -1;
302 const label facePoint0 = faces_[faceI][0];
304 const point& pf = points_[facePoint0];
306 const scalar d1 = normal_[faceI] & (p1 - pf);
307 const scalar d2 = normal_[faceI] & (p2 - pf);
316 const point pIntersect = p1 + (d1/(d1 - d2))*(p2 - p1);
318 const List<face>& tris = faceTris_[faceI];
323 const face& tri = tris[triI];
331 if (t.classify(pIntersect, dummyNearType, dummyNearLabel))
333 hitFaceIDs_.append(faceI);
340 template<
class CloudType>
349 const scalar d1 = normal_[0] & (p1 - coordSys_.origin());
350 const scalar d2 = normal_[0] & (p2 - coordSys_.origin());
359 const point pCyl = coordSys_.localPosition(p1 + (d1/(d1 - d2))*(p2 - p1));
363 if (r < radius_.last())
366 while (r > radius_[radI])
388 hitFaceIDs_.append(secI);
392 template<
class CloudType>
397 scalar timeNew = time.
value();
398 scalar timeElapsed = timeNew - timeOld_;
400 totalTime_ += timeElapsed;
402 const scalar
alpha = (totalTime_ - timeElapsed)/totalTime_;
403 const scalar
beta = timeElapsed/totalTime_;
407 massFlowRate_[faceI] =
408 alpha*massFlowRate_[faceI] + beta*mass_[faceI]/timeElapsed;
409 massTotal_[faceI] += mass_[faceI];
412 const label procI = Pstream::myProcNo();
417 this->getModelProperty(
"massTotal", faceMassTotal);
420 this->getModelProperty(
"massFlowRate", faceMassFlowRate);
423 scalar sumTotalMass = 0.0;
424 scalar sumAverageMFR = 0.0;
428 allProcMass[procI] = massTotal_[faceI];
429 Pstream::gatherList(allProcMass);
430 faceMassTotal[faceI] +=
sum(allProcMass);
432 scalarList allProcMassFlowRate(Pstream::nProcs());
433 allProcMassFlowRate[procI] = massFlowRate_[faceI];
434 Pstream::gatherList(allProcMassFlowRate);
435 faceMassFlowRate[faceI] +=
sum(allProcMassFlowRate);
437 sumTotalMass += faceMassTotal[faceI];
438 sumAverageMFR += faceMassFlowRate[faceI];
440 if (outputFilePtr_.valid())
445 <<
tab << faceMassTotal[faceI]
446 <<
tab << faceMassFlowRate[faceI]
451 Info<<
" sum(total mass) = " << sumTotalMass << nl
452 <<
" sum(average mass flow rate) = " << sumAverageMFR << nl
456 if (surfaceFormat_ !=
"none")
458 if (Pstream::master())
464 this->outputTimeDir(),
475 this->outputTimeDir(),
490 this->setModelProperty(
"massTotal", dummy);
491 this->setModelProperty(
"massFlowRate", dummy);
498 this->setModelProperty(
"massTotal", faceMassTotal);
499 this->setModelProperty(
"massFlowRate", faceMassFlowRate);
505 massTotal_[faceI] = 0.0;
506 massFlowRate_[faceI] = 0.0;
513 template<
class CloudType>
518 const word& modelName
523 parcelType_(this->coeffDict().lookupOrDefault(
"parcelType", -1)),
524 removeCollected_(this->coeffDict().
lookup(
"removeCollected")),
532 negateParcelsOppositeNormal_
536 surfaceFormat_(this->coeffDict().
lookup(
"surfaceFormat")),
537 resetOnWrite_(this->coeffDict().
lookup(
"resetOnWrite")),
542 log_(this->coeffDict().
lookup(
"log")),
544 timeOld_(owner.mesh().time().value()),
547 normal_ /=
mag(normal_);
550 if (mode ==
"polygon")
554 initPolygons(polygons);
559 else if (mode ==
"polygonWithNormal")
563 this->coeffDict().
lookup(
"polygons")
571 polygons[polyI] = polygonAndNormal[polyI].
first();
572 normal_[polyI] = polygonAndNormal[polyI].second();
573 normal_[polyI] /=
mag(normal_[polyI]) + ROOTVSMALL;
576 initPolygons(polygons);
578 else if (mode ==
"concentricCircle")
583 initConcentricCircles();
589 "Foam::ParticleCollector<CloudType>::ParticleCollector" 597 <<
"Unknown mode " << mode <<
". Available options are " 598 <<
"polygon, polygonWithNormal and concentricCircle" 602 mass_.setSize(faces_.size(), 0.0);
603 massTotal_.setSize(faces_.size(), 0.0);
604 massFlowRate_.setSize(faces_.size(), 0.0);
606 makeLogFile(faces_, points_, area_);
610 template<
class CloudType>
618 parcelType_(pc.parcelType_),
619 removeCollected_(pc.removeCollected_),
622 faceTris_(pc.faceTris_),
623 nSector_(pc.nSector_),
625 coordSys_(pc.coordSys_),
627 negateParcelsOppositeNormal_(pc.negateParcelsOppositeNormal_),
628 surfaceFormat_(pc.surfaceFormat_),
629 resetOnWrite_(pc.resetOnWrite_),
630 totalTime_(pc.totalTime_),
632 massTotal_(pc.massTotal_),
633 massFlowRate_(pc.massFlowRate_),
643 template<
class CloudType>
650 template<
class CloudType>
656 const point& position0,
660 if ((parcelType_ != -1) && (parcelType_ != p.typeId()))
666 const point position1 = position0 + 1.0001*(p.position() - position0);
674 collectParcelPolygon(position0, position1);
677 case mtConcentricCircle:
679 collectParcelConcentricCircles(position0, position1);
690 label faceI = hitFaceIDs_[i];
691 scalar m = p.nParticle()*p.mass();
693 if (negateParcelsOppositeNormal_)
696 Uhat /=
mag(Uhat) + ROOTVSMALL;
697 if ((Uhat & normal_[faceI]) < 0)
708 mass_[faceI + 1] += m;
709 mass_[faceI + 2] += m;
710 mass_[faceI + 3] += m;
713 if (removeCollected_)
715 keepParticle =
false;
Mesh data needed to do the Finite Volume discretisation.
vector point
Point is a vector.
mode_t mode(const fileName &)
Return the file mode.
dimensioned< scalar > mag(const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
virtual void postMove(parcelType &p, const label cellI, const scalar dt, const point &position0, bool &keepParticle)
Post-move hook.
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
A class for handling words, derived from string.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
void size(const label)
Override size to be inconsistent with allocated storage.
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensionedScalar sign(const dimensionedScalar &ds)
ParticleCollector(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Function object to collect the parcel mass- and mass flow rate over a set of polygons. The polygons can either be specified by sets of user- supplied points, or in a concentric circles arrangement. If a parcel is 'collected', it can be flagged to be removed from the domain using the removeCollected entry.
T & first()
Return the first element of the list.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
label readLabel(Istream &is)
triangle< point, const point & > triPointRef
const Time & time() const
Return the top-level database.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
void setSize(const label)
Reset size of List.
Ostream & endl(Ostream &os)
Add newline and flush stream.
stressControl lookup("compactNormalStress") >> compactNormalStress
volVectorField vectorField(fieldObject, mesh)
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Unit conversion functions.
Base class for graphics format writing. Entry points are.
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
void write()
Write post-processing info.
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
Vector< scalar > vector
A scalar version of the templated Vector.
Templated base class for dsmc cloud.
virtual ~ParticleCollector()
Destructor.
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Templated cloud function object base class.
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
const scalar twoPi(2 *pi)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
const Type & value() const
Return const reference to value.