Loading...
Searching...
No Matches
streamLineParticle.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "streamLineParticle.H"
31#include "vectorFieldIOField.H"
32
33// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34
35Foam::vector Foam::streamLineParticle::interpolateFields
36(
37 const trackingData& td,
38 const barycentric& tetCoords,
39 const tetIndices& tetIs
40)
41{
42 if (tetIs.cell() < 0)
43 {
45 << "Invalid cell (-1)" << abort(FatalError);
46 }
47
48 const point position
49 (
50 tetIs.barycentricToPoint(this->mesh(), tetCoords)
51 );
52
53 bool foundU = false;
54 vector U(Zero);
55
56 // If current position is different
57 if
58 (
59 sampledPositions_.empty()
60 || sampledPositions_.back().distSqr(position) > Foam::sqr(SMALL)
61 )
62 {
63 // Store new location
64 sampledPositions_.push_back(position);
65
66 // Scalar fields
67 sampledScalars_.resize(td.vsInterp_.size());
68 forAll(td.vsInterp_, i)
69 {
70 sampledScalars_[i].push_back
71 (
72 td.vsInterp_[i].interpolate(tetCoords, tetIs, tetIs.face())
73 );
74 }
75
76 // Vector fields
77 sampledVectors_.resize(td.vvInterp_.size());
78 forAll(td.vvInterp_, i)
79 {
80 sampledVectors_[i].push_back
81 (
82 td.vvInterp_[i].interpolate(tetCoords, tetIs, tetIs.face())
83 );
84
85 if (td.vvInterp_.get(i) == &(td.UInterp_))
86 {
87 foundU = true;
88 U = sampledVectors_[i].back();
89 }
90 }
91 }
92
93 if (!foundU)
94 {
95 U = td.UInterp_.interpolate(tetCoords, tetIs, tetIs.face());
96 }
98 return U;
99}
100
101
102// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
103
105(
106 const polyMesh& mesh,
107 const vector& position,
108 const label celli,
109 const bool trackForward,
110 const label lifeTime
111)
113 particle(mesh, position, celli),
114 trackForward_(trackForward),
115 lifeTime_(lifeTime)
116{}
117
118
120(
121 const polyMesh& mesh,
122 Istream& is,
123 bool readFields,
124 bool newFormat
125)
126:
127 particle(mesh, is, readFields, newFormat)
128{
129 if (readFields)
130 {
131 List<scalarList> sampledScalars;
132 List<vectorList> sampledVectors;
133
134 is >> trackForward_ >> lifeTime_
135 >> sampledPositions_ >> sampledScalars
136 >> sampledVectors;
137
138 sampledScalars_.resize(sampledScalars.size());
139 forAll(sampledScalars, i)
140 {
141 sampledScalars_[i].transfer(sampledScalars[i]);
142 }
143 sampledVectors_.resize(sampledVectors.size());
144 forAll(sampledVectors, i)
145 {
146 sampledVectors_[i].transfer(sampledVectors[i]);
148 }
149
151}
152
153
155(
156 const streamLineParticle& p
157)
158:
159 particle(p),
160 trackForward_(p.trackForward_),
161 lifeTime_(p.lifeTime_),
162 sampledPositions_(p.sampledPositions_),
163 sampledScalars_(p.sampledScalars_),
164 sampledVectors_(p.sampledVectors_)
165{}
166
167
168// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
169
171(
173 trackingData& td,
174 const scalar
175)
176{
177 td.switchProcessor = false;
178 td.keepParticle = true;
179
180 const scalar maxDt = mesh().bounds().mag();
181
182 while (td.keepParticle && !td.switchProcessor && lifeTime_ > 0)
183 {
184 // Set the lagrangian time-step
185 scalar dt = maxDt;
186
187 // Cross cell in steps:
188 // - at subiter 0 calculate dt to cross cell in nSubCycle steps
189 // - at the last subiter do all of the remaining track
190 for (label subIter = 0; subIter < max(1, td.nSubCycle_); subIter++)
191 {
192 --lifeTime_;
193
194 // Store current position, return sampled velocity
195 vector U =
196 interpolateFields(td, coordinates(), currentTetIndices());
197
198 const scalar magU = mag(U);
199
200 if (magU < SMALL)
201 {
202 // Stagnant particle. Might as well stop
203 lifeTime_ = 0;
204 break;
205 }
206
207 if (!trackForward_)
208 {
209 U = -U;
210 }
211
212 U /= magU;
213
214 if (td.trackLength_ < GREAT)
215 {
216 // No sub-cycling. Track a set length on each step.
217 dt = td.trackLength_;
218 }
219 else if (subIter == 0)
220 {
221 // Sub-cycling. Cross the cell in nSubCycle steps.
222 particle copy(*this);
223 copy.trackToFace(maxDt*U, 1);
224 dt *= (copy.stepFraction() - stepFraction())/td.nSubCycle_;
225 }
226 else if (subIter == td.nSubCycle_ - 1)
227 {
228 // Sub-cycling. Track the whole cell on the last step.
229 dt = maxDt;
230 }
231
232 trackToAndHitFace(dt*U, 0, cloud, td);
233
234 if
235 (
236 onFace()
237 || !td.keepParticle
238 || td.switchProcessor
239 || lifeTime_ == 0
240 )
241 {
242 break;
243 }
244 }
245 }
246
247 if (!td.keepParticle || lifeTime_ == 0)
248 {
249 if (lifeTime_ == 0)
250 {
251 // Failure exit. Particle stagnated or it's life ran out.
252 if (debug)
253 {
254 Pout<< "streamLineParticle: Removing stagnant particle:"
255 << position() << " sampled positions:"
256 << sampledPositions_.size() << endl;
257 }
258 td.keepParticle = false;
259 }
260 else
261 {
262 // Normal exit. Store last position and fields
263 (void)interpolateFields(td, coordinates(), currentTetIndices());
264
265 if (debug)
266 {
267 Pout<< "streamLineParticle: Removing particle:" << position()
268 << " sampled positions:" << sampledPositions_.size()
269 << endl;
270 }
271 }
272
273 // Transfer particle data into trackingData.
274 {
275 td.allPositions_.emplace_back().transfer(sampledPositions_);
276 }
277
278 forAll(sampledScalars_, i)
279 {
280 td.allScalars_[i].emplace_back().transfer(sampledScalars_[i]);
281 }
282 forAll(sampledVectors_, i)
283 {
284 td.allVectors_[i].emplace_back().transfer(sampledVectors_[i]);
286 }
287
288 return td.keepParticle;
289}
290
291
293{
294 // Disable generic patch interaction
295 return false;
296}
297
298
300(
302 trackingData& td
304{
305 // Remove particle
306 td.keepParticle = false;
307}
308
309
311(
313 trackingData& td
315{
316 // Remove particle
317 td.keepParticle = false;
318}
319
320
322(
324 trackingData& td
326{
327 // Remove particle
328 td.keepParticle = false;
329}
330
331
333(
335 trackingData& td
337{
338 // Remove particle
339 td.keepParticle = false;
340}
341
342
344(
346 trackingData& td,
347 const vector&
349{
350 // Remove particle
351 td.keepParticle = false;
352}
353
354
356(
358 trackingData& td,
359 const vector&
361{
362 // Remove particle
363 td.keepParticle = false;
364}
365
366
368(
370 trackingData& td
372{
373 // Switch particle
374 td.switchProcessor = true;
375}
376
377
379(
381 trackingData& td
383{
384 // Remove particle
385 td.keepParticle = false;
386}
387
388
390{
391 const bool readOnProc = c.size();
392
394
395 IOField<label> lifeTime
396 (
397 c.newIOobject("lifeTime", IOobject::MUST_READ),
398 readOnProc
399 );
400 c.checkFieldIOobject(c, lifeTime);
401
402 vectorFieldIOField sampledPositions
403 (
404 c.newIOobject("sampledPositions", IOobject::MUST_READ),
405 readOnProc
406 );
407 c.checkFieldIOobject(c, sampledPositions);
408
409 label i = 0;
410 for (streamLineParticle& p : c)
411 {
412 p.lifeTime_ = lifeTime[i];
413 p.sampledPositions_.transfer(sampledPositions[i]);
414 ++i;
415 }
416}
417
418
419void Foam::streamLineParticle::writeFields(const Cloud<streamLineParticle>& c)
420{
422
423 const label np = c.size();
424 const bool writeOnProc = c.size();
425
426 IOField<label> lifeTime
427 (
428 c.newIOobject("lifeTime", IOobject::NO_READ),
429 np
430 );
431 vectorFieldIOField sampledPositions
432 (
433 c.newIOobject("sampledPositions", IOobject::NO_READ),
434 np
435 );
436
437 label i = 0;
438 for (const streamLineParticle& p : c)
439 {
440 lifeTime[i] = p.lifeTime_;
441 sampledPositions[i] = p.sampledPositions_;
442 ++i;
443 }
444
445 lifeTime.write(writeOnProc);
446 sampledPositions.write(writeOnProc);
447}
448
449
450// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
451
452Foam::Ostream& Foam::operator<<(Ostream& os, const streamLineParticle& p)
453{
455 << token::SPACE << p.trackForward_
456 << token::SPACE << p.lifeTime_
457 << token::SPACE << p.sampledPositions_
458 << token::SPACE << p.sampledScalars_
459 << token::SPACE << p.sampledVectors_;
460
461 os.check(FUNCTION_NAME);
462 return os;
463}
464
465
466// ************************************************************************* //
Base cloud calls templated on particle type.
Definition Cloud.H:64
A primitive field of type <T> with automated input and output.
Definition IOField.H:53
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition IOstream.C:45
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
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 resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
scalar mag() const
The magnitude/length of the bounding box diagonal.
Definition boundBoxI.H:198
A cloud is a registry collection of lagrangian particles.
Definition cloud.H:56
Base particle class.
Definition particle.H:72
tetIndices currentTetIndices() const noexcept
Return indices of the current tet that the particle occupies.
Definition particleI.H:236
scalar trackToFace(const vector &displacement, const scalar fraction)
As particle::track, but also stops on internal faces.
Definition particle.C:639
vector position() const
Return current particle position.
Definition particleI.H:283
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
bool onFace() const noexcept
Is the particle on a face?
Definition particleI.H:259
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition particleI.H:110
scalar stepFraction() const noexcept
Return the fraction of time-step completed.
Definition particleI.H:170
void trackToAndHitFace(const vector &direction, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Convenience function. Combines trackToFace and hitFace.
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
Definition particle.C:507
static void writeFields(const TrackCloudType &c)
Write the fields associated with the owner cloud.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const boundBox & bounds() const noexcept
Return mesh bounding box.
Definition polyMesh.H:617
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
A Cloud of streamLine particles.
Particle class that samples fields as it passes through. Used in streamline calculation.
static void writeFields(const Cloud< streamLineParticle > &)
Write.
void hitCyclicACMIPatch(streamLineParticleCloud &, trackingData &, const vector &direction)
Overridable function to handle the particle hitting a.
void hitSymmetryPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a.
void hitSymmetryPlanePatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a.
void hitCyclicPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a cyclic.
void hitCyclicAMIPatch(streamLineParticleCloud &, trackingData &, const vector &direction)
Overridable function to handle the particle hitting a.
static void readFields(Cloud< streamLineParticle > &)
Read.
streamLineParticle(const polyMesh &mesh, const vector &position, const label celli, const bool trackForward, const label lifeTime)
Construct from components.
void hitWedgePatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a wedge.
void hitProcessorPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a processorPatch.
bool move(streamLineParticleCloud &, trackingData &, const scalar)
Track all particles to their end point.
bool hitPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a patch.
void hitWallPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
@ SPACE
Space [isspace].
Definition token.H:144
U
Definition pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
#define FUNCTION_NAME
Namespace for handling debugging switches.
Definition debug.C:45
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
IOField< vectorField > vectorFieldIOField
IO for a Field of vectorField.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Vector< scalar > vector
Definition vector.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299