Loading...
Searching...
No Matches
RemoveParcels.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2020-2023 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "RemoveParcels.H"
29#include "fvMesh.H"
30#include "faceZone.H"
31#include "OFstream.H"
32#include "surfaceFields.H"
33
34// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35
36template<class CloudType>
37void Foam::RemoveParcels<CloudType>::makeLogFile
38(
39 const word& zoneName,
40 const label zoneI,
41 const label nFaces,
42 const scalar totArea
43)
44{
45 // Create the output file if not already created
46 if (logToFile_)
47 {
48 DebugInfo<< "Creating output file." << endl;
49
50 if (Pstream::master())
51 {
52 // Create directory if does not exist
53 mkDir(this->writeTimeDir());
54
55 // Open new file at start up
56 outputFilePtr_.set
57 (
58 zoneI,
59 new OFstream
60 (
61 this->writeTimeDir()/(type() + '_' + zoneName + ".dat")
62 )
63 );
64
65 outputFilePtr_[zoneI]
66 << "# Source : " << type() << nl
67 << "# Face zone : " << zoneName << nl
68 << "# Faces : " << nFaces << nl
69 << "# Area : " << totArea << nl
70 << "# Time" << tab << "nParcels" << tab << "mass" << endl;
71 }
72 }
73}
74
75
76// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
77
78template<class CloudType>
80(
81 const typename parcelType::trackingData& td
82)
83{
84 Log_<< this->modelName() << " output:" << nl;
85
86 const fvMesh& mesh = this->owner().mesh();
87 const faceZoneMesh& fzm = mesh.faceZones();
88
89 forAll(faceZoneIDs_, i)
90 {
91 const word& zoneName = fzm[faceZoneIDs_[i]].name();
92
93 scalar zoneMass = returnReduce(mass_[i], sumOp<scalar>());
94 label zoneNParcels = returnReduce(nParcels_[i], sumOp<label>());
95
96 Log_<< " faceZone " << zoneName
97 << ": removed " << zoneNParcels
98 << " parcels with mass " << zoneMass
99 << nl;
101
103}
104
105
106template<class CloudType>
108{
109 const fvMesh& mesh = this->owner().mesh();
110 const Time& time = mesh.time();
111
112
113 List<scalar> allZoneMass(faceZoneIDs_.size(), 0.0);
114 List<label> allZoneNParcels(faceZoneIDs_.size(), 0);
115
116 forAll(faceZoneIDs_, i)
117 {
118 allZoneMass[i] = returnReduce(mass_[i], sumOp<scalar>());
119 allZoneNParcels[i] = returnReduce(nParcels_[i], sumOp<label>());
120
121 if (outputFilePtr_.set(i))
122 {
123 OFstream& os = outputFilePtr_[i];
124 os << time.timeName() << token::TAB
125 << allZoneNParcels[i] << token::TAB
126 << allZoneMass[i] << endl;
127 }
128 }
129
130 Log_<< endl;
131
132 if (resetOnWrite_)
133 {
134 forAll(mass_, i)
135 {
136 mass_[i] = 0.0;
137 nParcels_[i] = 0.0;
138 }
139 }
140
141 this->setModelProperty("mass", allZoneMass);
142 this->setModelProperty("nParcels", allZoneNParcels);
143}
144
145
146// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
147
148template<class CloudType>
150(
151 const dictionary& dict,
152 CloudType& owner,
153 const word& modelName
154)
155:
156 CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
157 faceZoneIDs_(),
158 nParcels_(),
159 mass_(),
160 typeId_(this->coeffDict().template getOrDefault<label>("parcelType", -1)),
161 logToFile_(this->coeffDict().getBool("log")),
162 resetOnWrite_(this->coeffDict().getBool("resetOnWrite")),
163 resetOnStart_(this->coeffDict().getBool("resetOnStart")),
164 outputFilePtr_()
165{
166 const wordList faceZoneNames(this->coeffDict().lookup("faceZones"));
167
168 nParcels_.setSize(faceZoneNames.size(), 0);
169 mass_.setSize(faceZoneNames.size(), 0.0);
170
171 if (!resetOnStart_ && Pstream::master())
172 {
173 this->getModelProperty("mass", mass_);
174 this->getModelProperty("nParcels", nParcels_);
175 }
176
177 outputFilePtr_.setSize(faceZoneNames.size());
178
179 DynamicList<label> zoneIDs;
180 const faceZoneMesh& fzm = owner.mesh().faceZones();
181 const surfaceScalarField& magSf = owner.mesh().magSf();
182 const polyBoundaryMesh& pbm = owner.mesh().boundaryMesh();
183
184 forAll(faceZoneNames, i)
185 {
186 const word& zoneName = faceZoneNames[i];
187 label zonei = fzm.findZoneID(zoneName);
188
189 if (zonei != -1)
190 {
191 zoneIDs.append(zonei);
192 const faceZone& fz = fzm[zonei];
193
194 label nFaces = returnReduce(fz.size(), sumOp<label>());
195 Info<< " " << zoneName << " faces: " << nFaces << nl;
196
197 scalar totArea = 0.0;
198 for (const label facei : fz)
199 {
200 if (facei < owner.mesh().nInternalFaces())
201 {
202 totArea += magSf[facei];
203 }
204 else
205 {
206 const label patchi = pbm.patchID(facei);
207 const polyPatch& pp = pbm[patchi];
208
209 if
210 (
211 !magSf.boundaryField()[patchi].coupled()
213 )
214 {
215 label localFacei = pp.whichFace(facei);
216 totArea += magSf.boundaryField()[patchi][localFacei];
217 }
218 }
219 }
220 totArea = returnReduce(totArea, sumOp<scalar>());
221
222 makeLogFile(zoneName, i, nFaces, totArea);
223 }
225
226 faceZoneIDs_.transfer(zoneIDs);
227}
228
229
230template<class CloudType>
232(
233 const RemoveParcels<CloudType>& rpf
234)
235:
236 CloudFunctionObject<CloudType>(rpf),
237 faceZoneIDs_(rpf.faceZoneIDs_),
238 nParcels_(rpf.nParcels_),
239 mass_(rpf.mass_),
240 typeId_(rpf.typeId_),
241 logToFile_(rpf.logToFile_),
242 resetOnWrite_(rpf.resetOnWrite_),
243 resetOnStart_(rpf.resetOnStart_),
244 outputFilePtr_()
245{}
246
247
248// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
249
250template<class CloudType>
252(
253 const parcelType& p,
254 const typename parcelType::trackingData& td
255)
256{
257 bool keepParticle = true;
258
259 if ((typeId_ >= 0) && (p.typeId() != typeId_))
260 {
261 // Not processing this particle type
262 return keepParticle;
263 }
264
265 if
266 (
267 this->owner().solution().output()
268 || this->owner().solution().transient()
269 )
270 {
271 const faceZoneMesh& fzm = this->owner().mesh().faceZones();
272
273 forAll(faceZoneIDs_, i)
274 {
275 const faceZone& fz = fzm[faceZoneIDs_[i]];
276 if (fz.found(p.face()))
277 {
278 nParcels_[i]++;
279 mass_[i] += p.mass()*p.nParticle();
280 keepParticle = false;
281 break;
282 }
283 }
284 }
285
286 return keepParticle;
287}
288
289
290// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
Templated cloud function object base class.
CloudFunctionObject(CloudType &owner)
Construct null from owner.
virtual void postEvolve(const typename parcelType::trackingData &td)
Post-evolve hook.
const CloudType & owner() const
Return const access to the owner cloud.
const fvMesh & mesh() const
Return reference to the mesh.
Definition DSMCCloudI.H:37
particle::trackingData trackingData
Definition DSMCParcel.H:155
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
void setSize(label n)
Alias for resize().
Definition List.H:536
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
RemoveParcels(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
virtual void postEvolve(const typename parcelType::trackingData &td)
Post-evolve hook.
void write()
Write post-processing info.
virtual bool postFace(const parcelType &p, const typename parcelType::trackingData &td)
Post-face hook.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
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.
Definition Time.C:714
bool found(const T &val, label pos=0) const
Same as contains().
Definition UList.H:983
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition ZoneMesh.C:757
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
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.
Definition polyMesh.H:609
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition polyMesh.H:671
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
label nInternalFaces() const noexcept
Number of internal faces.
Lookup type of boundary radiation properties.
Definition lookup.H:60
Selector class for relaxation factors, solver type and solution.
Definition solution.H:95
bool getModelProperty(const word &entryName, Type &value) const
Retrieve generic property from the sub-model.
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.
void setModelProperty(const word &entryName, const Type &value)
Add generic property to the sub-model.
@ TAB
Tab [isspace].
Definition token.H:142
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
const labelIOList & zoneIDs
Definition correctPhi.H:59
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
return returnReduce(nRefine-oldNRefine, sumOp< label >())
#define Log_
Report write to Foam::Info if the class log switch is true.
#define DebugInfo
Report an information message using Foam::Info.
type
Types of root.
Definition Roots.H:53
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< word > wordList
List of word.
Definition fileName.H:60
DSMCCloud< dsmcParcel > CloudType
messageStream Info
Information stream (stdout output on master, null elsewhere).
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
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.
Definition Ostream.H:519
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.
mkDir(pdfPath)