Go to the documentation of this file.
34 template<
class CloudType>
48 Info<<
"Creating output file." <<
endl;
51 if (Pstream::master())
54 mkDir(this->outputTimeDir());
62 this->outputTimeDir()/(
type() +
'_' + zoneName +
".dat")
67 <<
"# Source : " <<
type() <<
nl
68 <<
"# Face zone : " << zoneName <<
nl
69 <<
"# Faces : " << nFaces <<
nl
70 <<
"# Area : " << totArea <<
nl
71 <<
"# Time" <<
tab <<
"mass" <<
tab <<
"massFlowRate" <<
endl;
77 template<
class CloudType>
83 scalar timeNew = time.
value();
84 scalar timeElapsed = timeNew - timeOld_;
86 totalTime_ += timeElapsed;
88 const scalar
alpha = (totalTime_ - timeElapsed)/totalTime_;
89 const scalar
beta = timeElapsed/totalTime_;
91 forAll(faceZoneIDs_, zoneI)
93 massFlowRate_[zoneI] =
94 alpha*massFlowRate_[zoneI] +
beta*mass_[zoneI]/timeElapsed;
95 massTotal_[zoneI] += mass_[zoneI];
98 const label procI = Pstream::myProcNo();
104 forAll(faceZoneIDs_, zoneI)
106 const word& zoneName = fzm[faceZoneIDs_[zoneI]].name();
109 allProcMass[procI] = massTotal_[zoneI];
110 Pstream::gatherList(allProcMass);
111 zoneMassTotal[zoneI] =
112 ListListOps::combine<scalarList>
116 const scalar sumMassTotal =
sum(zoneMassTotal[zoneI]);
119 allProcMassFlowRate[procI] = massFlowRate_[zoneI];
120 Pstream::gatherList(allProcMassFlowRate);
121 zoneMassFlowRate[zoneI] =
122 ListListOps::combine<scalarList>
126 const scalar sumMassFlowRate =
sum(zoneMassFlowRate[zoneI]);
128 Info<<
" " << zoneName
129 <<
": total mass = " << sumMassTotal
130 <<
"; average mass flow rate = " << sumMassFlowRate
133 if (outputFilePtr_.set(zoneI))
135 OFstream& os = outputFilePtr_[zoneI];
136 os << time.
timeName() << token::TAB << sumMassTotal << token::TAB
137 << sumMassFlowRate<<
endl;
144 if (surfaceFormat_ !=
"none")
146 forAll(faceZoneIDs_, zoneI)
148 const faceZone& fZone = fzm[faceZoneIDs_[zoneI]];
153 mesh.globalData().mergePoints
155 fZone().meshPoints(),
156 fZone().meshPointMap(),
158 uniqueMeshPointLabels
163 allProcPoints[procI] = uniquePoints;
164 Pstream::gatherList(allProcPoints);
166 faceList faces(fZone().localFaces());
172 allProcFaces[procI] = faces;
173 Pstream::gatherList(allProcFaces);
175 if (Pstream::master())
179 ListListOps::combine<pointField>
187 ListListOps::combine<faceList>
198 this->coeffDict().subOrEmptyDict(
"formatOptions").
199 subOrEmptyDict(surfaceFormat_)
205 this->outputTimeDir(),
210 zoneMassTotal[zoneI],
216 this->outputTimeDir(),
221 zoneMassFlowRate[zoneI],
231 forAll(faceZoneIDs_, zoneI)
233 massFlowRate_[zoneI] = 0.0;
250 template<
class CloudType>
255 const word& modelName
260 surfaceFormat_(this->coeffDict().
lookup(
"surfaceFormat")),
261 resetOnWrite_(this->coeffDict().
lookup(
"resetOnWrite")),
266 log_(this->coeffDict().
lookup(
"log")),
271 mass_.setSize(faceZoneNames.
size());
272 massTotal_.setSize(faceZoneNames.
size());
273 massFlowRate_.setSize(faceZoneNames.
size());
275 outputFilePtr_.setSize(faceZoneNames.
size());
283 const word& zoneName = faceZoneNames[i];
289 mass_[i].setSize(fz.
size(), 0.0);
290 massTotal_[i].setSize(fz.
size(), 0.0);
291 massFlowRate_[i].setSize(fz.
size(), 0.0);
294 Info<<
" " << zoneName <<
" faces: " << nFaces <<
nl;
296 scalar totArea = 0.0;
302 totArea += magSf[fz[j]];
313 || refCast<const coupledPolyPatch>(pp).owner()
323 makeLogFile(zoneName, i, nFaces, totArea);
327 faceZoneIDs_.transfer(zoneIDs);
333 template<
class CloudType>
355 template<
class CloudType>
362 template<
class CloudType>
373 || this->owner().
solution().
transient()
380 const faceZone& fz = fzm[faceZoneIDs_[i]];
394 mass_[i][
faceId] +=
p.mass()*
p.nParticle();
void makeLogFile(const word &zoneName, const label zoneI, const label nFaces, const scalar totArea)
Helper function to create log files.
List< scalarField > mass_
Mass storage.
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Selector class for relaxation factors, solver type and solution.
A class for handling words, derived from string.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
#define forAll(list, i)
Loop across all elements in list.
Switch log_
Flag to indicate whether data should be written to file.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual void write(const coordSet &, const wordList &, const List< const Field< Type > * > &, Ostream &) const =0
General entry point for writing.
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
List< scalarField > massTotal_
Mass total storage.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type & value() const
Return const reference to value.
Records particle face quantities on used-specified face zone.
const faceZoneMesh & faceZones() const
Return face zone mesh.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
const labelList & patchID() const
Per boundary face label the patch index.
virtual ~FacePostProcessing()
Destructor.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
A subset of mesh faces organised as a primitive patch.
const fvMesh & mesh() const
Return refernce to the mesh.
A patch is a list of labels that address the faces in the global face list.
labelList faceZoneIDs_
Face zone IDs.
Templated base class for dsmc cloud.
A list of keyword definitions, which are a keyword followed by any number of values (e....
label nInternalFaces() const
Base class for graphics format writing. Entry points are.
Mesh data needed to do the Finite Volume discretisation.
const word & name() const
Return name.
label findZoneID(const word &zoneName) const
Find zone index given a name.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
FacePostProcessing(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
const word surfaceFormat_
Surface output format.
Switch resetOnWrite_
Flag to indicate whether data should be reset/cleared on writing.
label whichFace(const label l) const
Return label of face in patch from global face label.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
virtual void postFace(const parcelType &p, const label faceI, bool &keepParticle)
Post-face hook.
Templated cloud function object base class.
List< scalarField > massFlowRate_
Mass flow rate storage.
void write()
Write post-processing info.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
scalar totalTime_
Total time.
const Time & time() const
Return the top-level database.
CloudType::parcelType parcelType
Convenience typedef for parcel type.
void size(const label)
Override size to be inconsistent with allocated storage.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Generic GeometricField class.
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
stressControl lookup("compactNormalStress") >> compactNormalStress