37 template<
class CloudType>
41 const Field<point>&
points,
42 const Field<scalar>& area
50 Info<<
"Creating output file" <<
endl;
53 if (Pstream::master())
56 mkDir(this->writeTimeDir());
61 new OFstream(this->writeTimeDir()/(
type() +
".dat"))
65 <<
"# Source : " <<
type() <<
nl
66 <<
"# Bins : " << faces.size() <<
nl
67 <<
"# Total area : " <<
sum(area) <<
nl;
70 <<
"# Geometry :" <<
nl
73 <<
tab <<
"(Centre_x Centre_y Centre_z)"
89 <<
"# Output format:" <<
nl;
94 word binId =
"bin_" + id;
100 <<
tab <<
"mass[" <<
id <<
"]"
101 <<
tab <<
"massFlowRate[" <<
id <<
"]"
109 template<
class CloudType>
112 const List<Field<point>>& polygons
120 label np = polygons[polyI].size();
124 <<
"polygons must consist of at least 3 points"
131 label pointOffset = 0;
133 faces_.setSize(polygons.size());
134 area_.setSize(polygons.size());
137 const Field<point>& polyPoints = polygons[facei];
139 UIndirectList<point>(points_,
f) = polyPoints;
140 area_[facei] =
f.mag(points_);
141 faces_[facei].transfer(
f);
143 pointOffset += polyPoints.size();
148 template<
class CloudType>
151 mode_ = mtConcentricCircle;
153 vector origin(this->coeffDict().lookup(
"origin"));
155 this->coeffDict().lookup(
"radius") >> radius_;
156 nSector_ = this->coeffDict().template lookup<label>(
"nSector");
163 refDir = this->coeffDict().lookup(
"refDir");
164 refDir -= normal_[0]*(normal_[0] & refDir);
165 refDir /=
mag(refDir);
175 scalar dThetaSector =
degToRad(360.0)/scalar(nS);
176 label intervalPerSector =
max(1, ceil(dThetaSector/dTheta));
177 dTheta = dThetaSector/scalar(intervalPerSector);
179 label nPointPerSector = intervalPerSector + 1;
181 label nPointPerRadius = nS*(nPointPerSector - 1);
182 label nPoint = radius_.size()*nPointPerRadius;
183 label nFace = radius_.size()*nS;
188 points_.setSize(nPoint);
189 faces_.setSize(nFace);
190 area_.setSize(nFace);
193 coordinateSystems::cylindrical(
"coordSys", origin, normal_[0], refDir);
202 label pointOffset = radI*nPointPerRadius + 1;
204 for (
label i = 0; i < nPointPerRadius; i++)
206 label pI = i + pointOffset;
207 point pCyl(radius_[radI], i*dTheta, 0.0);
208 points_[pI] = coordSys_.globalPosition(pCyl);
213 DynamicList<label> facePts(2*nPointPerSector);
218 for (
label secI = 0; secI < nS; secI++)
225 for (
label ptI = 0; ptI < nPointPerSector; ptI++)
227 label i = ptI + secI*(nPointPerSector - 1);
228 label id = ptIDs.fcIndex(i - 1) + 1;
232 label facei = secI + radI*nS;
234 faces_[facei] = face(facePts);
235 area_[facei] = faces_[facei].mag(points_);
240 for (
label secI = 0; secI < nS; secI++)
246 for (
label ptI = 0; ptI < nPointPerSector; ptI++)
248 label i = ptI + secI*(nPointPerSector - 1);
252 for (
label ptI = nPointPerSector-1; ptI >= 0; ptI--)
254 label i = ptI + secI*(nPointPerSector - 1);
255 label id =
offset + nPointPerRadius + ptIDs.fcIndex(i - 1);
259 label facei = secI + radI*nS;
261 faces_[facei] = face(facePts);
262 area_[facei] = faces_[facei].mag(points_);
269 template<
class CloudType>
278 const label facePoint0 = faces_[facei][0];
280 const point& pf = points_[facePoint0];
282 const scalar d1 = normal_[facei] & (p1 - pf);
283 const scalar d2 = normal_[facei] & (p2 - pf);
292 const point pIntersect = p1 + (d1/(d1 - d2))*(p2 - p1);
300 const face&
f = faces_[facei];
301 const vector a =
f.area(points_);
303 for (
label i = 0; i <
f.size(); ++ i)
305 const label j =
f.fcIndex(i);
306 const triPointRef t(pIntersect, points_[
f[i]], points_[
f[j]]);
307 if ((a & t.area()) < 0)
317 hitFaceIndices_.append(facei);
323 template<
class CloudType>
332 const scalar d1 = normal_[0] & (p1 - coordSys_.origin());
333 const scalar d2 = normal_[0] & (p2 - coordSys_.origin());
342 const point pCyl = coordSys_.localPosition(p1 + (d1/(d1 - d2))*(p2 - p1));
346 if (r < radius_.last())
349 while (r > radius_[radI])
373 hitFaceIndices_.append(secI);
380 template<
class CloudType>
385 scalar timeNew = time.
value();
386 scalar timeElapsed = timeNew - timeOld_;
388 totalTime_ += timeElapsed;
390 const scalar
alpha = (totalTime_ - timeElapsed)/totalTime_;
391 const scalar beta = timeElapsed/totalTime_;
395 massFlowRate_[facei] =
396 alpha*massFlowRate_[facei] + beta*mass_[facei]/timeElapsed;
397 massTotal_[facei] += mass_[facei];
400 const label proci = Pstream::myProcNo();
405 this->getModelProperty(
"massTotal", faceMassTotal);
408 this->getModelProperty(
"massFlowRate", faceMassFlowRate);
411 scalar sumTotalMass = 0.0;
412 scalar sumAverageMFR = 0.0;
416 allProcMass[proci] = massTotal_[facei];
417 Pstream::gatherList(allProcMass);
418 faceMassTotal[facei] +=
sum(allProcMass);
420 scalarList allProcMassFlowRate(Pstream::nProcs());
421 allProcMassFlowRate[proci] = massFlowRate_[facei];
422 Pstream::gatherList(allProcMassFlowRate);
423 faceMassFlowRate[facei] +=
sum(allProcMassFlowRate);
425 sumTotalMass += faceMassTotal[facei];
426 sumAverageMFR += faceMassFlowRate[facei];
428 if (outputFilePtr_.valid())
433 <<
tab << faceMassTotal[facei]
434 <<
tab << faceMassFlowRate[facei]
439 Info<<
" sum(total mass) = " << sumTotalMass <<
nl
440 <<
" sum(average mass flow rate) = " << sumAverageMFR <<
nl
444 if (surfaceFormat_ !=
"none")
446 if (Pstream::master())
455 this->writeTimeDir(),
472 this->setModelProperty(
"massTotal", dummy);
473 this->setModelProperty(
"massFlowRate", dummy);
480 this->setModelProperty(
"massTotal", faceMassTotal);
481 this->setModelProperty(
"massFlowRate", faceMassFlowRate);
487 massTotal_[facei] = 0.0;
488 massFlowRate_[facei] = 0.0;
495 template<
class CloudType>
500 const word& modelName
505 parcelType_(this->coeffDict().lookupOrDefault(
"parcelType", -1)),
506 removeCollected_(this->coeffDict().lookup(
"removeCollected")),
513 negateParcelsOppositeNormal_
515 readBool(this->coeffDict().lookup(
"negateParcelsOppositeNormal"))
517 surfaceFormat_(this->coeffDict().lookup(
"surfaceFormat")),
518 resetOnWrite_(this->coeffDict().lookup(
"resetOnWrite")),
523 log_(this->coeffDict().lookup(
"log")),
525 timeOld_(owner.mesh().time().value()),
528 normal_ /=
mag(normal_);
531 if (
mode ==
"polygon")
535 initPolygons(polygons);
540 else if (
mode ==
"polygonWithNormal")
552 polygons[polyI] = polygonAndNormal[polyI].
first();
553 normal_[polyI] = polygonAndNormal[polyI].second();
554 normal_[polyI] /=
mag(normal_[polyI]) + rootVSmall;
557 initPolygons(polygons);
559 else if (
mode ==
"concentricCircle")
564 initConcentricCircles();
569 <<
"Unknown mode " << mode <<
". Available options are "
570 <<
"polygon, polygonWithNormal and concentricCircle"
578 makeLogFile(faces_, points_, area_);
582 template<
class CloudType>
590 parcelType_(pc.parcelType_),
591 removeCollected_(pc.removeCollected_),
594 nSector_(pc.nSector_),
596 coordSys_(pc.coordSys_),
598 negateParcelsOppositeNormal_(pc.negateParcelsOppositeNormal_),
599 surfaceFormat_(pc.surfaceFormat_),
600 resetOnWrite_(pc.resetOnWrite_),
601 totalTime_(pc.totalTime_),
603 massTotal_(pc.massTotal_),
604 massFlowRate_(pc.massFlowRate_),
614 template<
class CloudType>
621 template<
class CloudType>
626 const point& position0,
630 if ((parcelType_ != -1) && (parcelType_ !=
p.typeId()))
635 hitFaceIndices_.clear();
644 p.position(this->owner().mesh())
648 case mtConcentricCircle:
650 collectParcelConcentricCircles
653 p.position(this->owner().mesh())
661 forAll(hitFaceIndices_, i)
663 label facei = hitFaceIndices_[i];
664 scalar m =
p.nParticle()*
p.mass();
666 if (negateParcelsOppositeNormal_)
674 Unormal = Uhat & normal_[facei];
677 case mtConcentricCircle:
679 Unormal = Uhat & normal_[0];
686 Uhat /=
mag(Uhat) + rootVSmall;
699 mass_[facei + 1] += m;
700 mass_[facei + 2] += m;
701 mass_[facei + 3] += m;
704 if (removeCollected_)
706 keepParticle =
false;
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
#define forAll(list, i)
Loop across all elements in list.
Templated cloud function object base class.
Templated base class for dsmc cloud.
void size(const label)
Override size to be inconsistent with allocated storage.
void setSize(const label)
Reset size of List.
Function object to collect the parcel mass- and mass flow rate over a set of polygons....
virtual void postMove(parcelType &p, const scalar dt, const point &position0, bool &keepParticle)
Post-move hook.
virtual ~ParticleCollector()
Destructor.
ParticleCollector(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
void write()
Write post-processing info.
Templated 3D SphericalTensor derived from VectorSpace adding construction from 1 component,...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
T & first()
Return the first element of the list.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
A coordinate rotation specified using global axis.
A list of keyword definitions, which are a keyword followed by any number of values (e....
const Type & value() const
Return const reference to value.
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
const Time & time() const
Return the top-level database.
const polyMesh & mesh() const
Return reference to polyMesh.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
A class for handling words, derived from string.
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const scalar twoPi(2 *pi)
errorManipArg< error, int > exit(error &err, const int errNo=1)
mode_t mode(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file mode.
scalar degToRad(const scalar deg)
Convert degrees to radians.
dimensionedScalar sign(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
word name(const bool)
Return a word representation of a bool.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
static const Identity< scalar > I
vector point
Point is a vector.
Vector< scalar > vector
A scalar version of the templated Vector.
dimensionSet normalised(const dimensionSet &)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Field< vector > vectorField
Specialisation of Field<T> for vector.
triangle< point, const point & > triPointRef
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
dimensionSet perpendicular(const dimensionSet &)
void offset(label &lst, const label o)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.