37template<
class CloudType>
38void Foam::FacePostProcessing<CloudType>::makeLogFile
51 if (Pstream::master())
54 mkDir(this->writeTimeDir());
62 this->writeTimeDir()/(
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;
79template<
class CloudType>
82 const fvMesh&
mesh = this->owner().mesh();
83 const Time& time =
mesh.time();
84 const faceZoneMesh& fzm =
mesh.faceZones();
85 scalar timeNew = time.value();
86 scalar timeElapsed = timeNew - timeOld_;
88 totalTime_ += timeElapsed;
90 const scalar
alpha = (totalTime_ - timeElapsed)/totalTime_;
91 const scalar
beta = timeElapsed/totalTime_;
93 forAll(faceZoneIDs_, zoneI)
95 massFlowRate_[zoneI] =
96 alpha*massFlowRate_[zoneI] +
beta*mass_[zoneI]/timeElapsed;
97 massTotal_[zoneI] += mass_[zoneI];
100 const label proci = Pstream::myProcNo();
104 List<scalarField> zoneMassTotal(mass_.size());
105 List<scalarField> zoneMassFlowRate(massFlowRate_.size());
106 forAll(faceZoneIDs_, zoneI)
108 const word& zoneName = fzm[faceZoneIDs_[zoneI]].name();
111 allProcMass[proci] = massTotal_[zoneI];
112 Pstream::gatherList(allProcMass);
113 zoneMassTotal[zoneI] =
114 ListListOps::combine<scalarList>
116 allProcMass, accessOp<scalarList>()
118 const scalar sumMassTotal =
sum(zoneMassTotal[zoneI]);
121 allProcMassFlowRate[proci] = massFlowRate_[zoneI];
122 Pstream::gatherList(allProcMassFlowRate);
123 zoneMassFlowRate[zoneI] =
124 ListListOps::combine<scalarList>
126 allProcMassFlowRate, accessOp<scalarList>()
128 const scalar sumMassFlowRate =
sum(zoneMassFlowRate[zoneI]);
130 Log_<<
" " << zoneName
131 <<
": total mass = " << sumMassTotal
132 <<
"; average mass flow rate = " << sumMassFlowRate
135 if (outputFilePtr_.set(zoneI))
137 OFstream&
os = outputFilePtr_[zoneI];
138 os << time.timeName() << token::TAB << sumMassTotal << token::TAB
139 << sumMassFlowRate<<
endl;
146 if (surfaceFormat_ !=
"none")
148 forAll(faceZoneIDs_, zoneI)
150 const faceZone& fZone = fzm[faceZoneIDs_[zoneI]];
154 autoPtr<globalIndex> globalPointsPtr =
155 mesh.globalData().mergePoints
157 fZone().meshPoints(),
158 fZone().meshPointMap(),
160 uniqueMeshPointLabels
164 List<pointField> allProcPoints(Pstream::nProcs());
165 allProcPoints[proci] = uniquePoints;
166 Pstream::gatherList(allProcPoints);
168 faceList faces(fZone().localFaces());
173 List<faceList> allProcFaces(Pstream::nProcs());
174 allProcFaces[proci] = faces;
175 Pstream::gatherList(allProcFaces);
177 if (Pstream::master())
181 ListListOps::combine<pointField>
183 allProcPoints, accessOp<pointField>()
189 ListListOps::combine<faceList>
191 allProcFaces, accessOp<faceList>()
195 auto writer = surfaceWriter::New
198 surfaceWriter::formatOptions
214 (this->writeTimeDir() / fZone.name()),
219 writer->write(
"massTotal", zoneMassTotal[zoneI]);
220 writer->write(
"massFlowRate", zoneMassFlowRate[zoneI]);
228 forAll(faceZoneIDs_, zoneI)
230 massFlowRate_[zoneI] = 0.0;
247template<
class CloudType>
252 const word& modelName
257 surfaceFormat_(this->coeffDict().getWord(
"surfaceFormat")),
258 resetOnWrite_(this->coeffDict().getBool(
"resetOnWrite")),
259 logToFile_(this->coeffDict().getBool(
"log")),
265 timeOld_(owner.
mesh().time().value())
272 outputFilePtr_.setSize(faceZoneNames.
size());
281 const word& zoneName = faceZoneNames[i];
293 Info<<
" " << zoneName <<
" faces: " << nFaces <<
nl;
295 scalar totArea = 0.0;
296 for (
const label facei : fz)
300 totArea += magSf[facei];
304 const label patchi =
pbm.patchID(facei);
313 label localFacei =
pp.whichFace(facei);
320 makeLogFile(zoneName, i, nFaces, totArea);
330template<
class CloudType>
337 faceZoneIDs_(pff.faceZoneIDs_),
338 surfaceFormat_(pff.surfaceFormat_),
339 resetOnWrite_(pff.resetOnWrite_),
340 logToFile_(pff.logToFile_),
341 totalTime_(pff.totalTime_),
343 massTotal_(pff.massTotal_),
344 massFlowRate_(pff.massFlowRate_),
352template<
class CloudType>
356 const typename parcelType::trackingData&
td
362 || this->owner().
solution().transient()
369 const faceZone& fz = fzm[faceZoneIDs_[i]];
374 if (fz[j] ==
p.face())
383 mass_[i][
faceId] +=
p.mass()*
p.nParticle();
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edgesCentres")))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
Templated cloud function object base class.
fileName writeTimeDir() const
Return the output time path.
CloudFunctionObject(CloudType &owner)
Construct null from owner.
const CloudType & owner() const
Return const access to the owner cloud.
const fvMesh & mesh() const
Return reference to the mesh.
particle::trackingData trackingData
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Records particle face quantities on used-specified face zone.
void write()
Write post-processing info.
FacePostProcessing(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
virtual bool postFace(const parcelType &p, const typename parcelType::trackingData &td)
Post-face hook.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
const word & name() const noexcept
Return the object name.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void setSize(label n)
Alias for resize().
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
static void gatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate.
void setSize(const label n)
Same as resize().
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
void size(const label n)
Older name for setAddressableSize.
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 master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
const Type & value() const noexcept
Return const reference to value.
A subset of mesh faces organised as a primitive patch.
Mesh data needed to do the Finite Volume discretisation.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
A patch is a list of labels that address the faces in the global face list.
label nInternalFaces() const noexcept
Number of internal faces.
Lookup type of boundary radiation properties.
Selector class for relaxation factors, solver type and solution.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
const word & modelName() const
Return const access to the name of the sub-model.
static autoPtr< surfaceWriter > New(const word &writeType)
Select construct a surfaceWriter.
static dictionary formatOptions(const dictionary &dict, const word &formatName, const word &entryName="formatOptions")
Same as fileFormats::getFormatOptions.
A class for handling words, derived from Foam::string.
const word & name() const noexcept
The zone name.
OBJstream os(runTime.globalPath()/outputName)
const labelIOList & zoneIDs
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
#define Log_
Report write to Foam::Info if the class log switch is true.
#define DebugInfo
Report an information message using Foam::Info.
void inplaceRenumber(const labelUList &oldToNew, IntListType &lists)
Inplace renumber the values (not the indices) of a list of lists.
AccessType combine(const UList< T > &lists, AccessOp aop=accessOp< T >())
Combines sub-lists into a single list.
Namespace for handling debugging switches.
List< scalarList > scalarListList
List of scalarList.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
List< word > wordList
List of word.
DSMCCloud< dsmcParcel > CloudType
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values within a list.
List< label > labelList
A List of labels.
messageStream Info
Information stream (stdout output on master, null elsewhere).
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
List< face > faceList
List of faces.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
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.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
vectorField pointField
pointField is a vectorField.
constexpr char nl
The newline '\n' character (0x0a).
constexpr char tab
The tab '\t' character(0x09).
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Object access operator or list access operator (default is pass-through).