streamlinesParticle.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "streamlinesParticle.H"
27 #include "streamlinesCloud.H"
28 #include "vectorFieldIOField.H"
29 #include "scalarFieldIOField.H"
30 #include "transformerIOList.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 Foam::vector Foam::streamlinesParticle::interpolateFields
35 (
36  const trackingData& td,
37  const point& position,
38  const label celli,
39  const label facei
40 )
41 {
42  if (celli == -1)
43  {
45  << "Cell:" << celli << abort(FatalError);
46  }
47 
48  bool interpolatedU = false;
49  vector U = vector::uniform(NaN);
50 
51  forAll(td.scalarInterp_, fieldi)
52  {
53  #define InterpolateType(Type, nullArg) \
54  if (td.Type##Interp_.set(fieldi)) \
55  { \
56  const Type s = \
57  td.Type##Interp_[fieldi].interpolate \
58  ( \
59  position, \
60  celli, \
61  facei \
62  ); \
63  \
64  sampled##Type##s_.setSize(td.Type##Interp_.size()); \
65  sampled##Type##s_[fieldi].append \
66  ( \
67  td.trackOutside_ ? transform_.invTransform(s) : s \
68  ); \
69  }
70  FOR_ALL_FIELD_TYPES(InterpolateType);
71  #undef InterpolateType
72 
73  if
74  (
75  td.vectorInterp_.set(fieldi)
76  && &td.vectorInterp_[fieldi] == &td.UInterp_
77  )
78  {
79  interpolatedU = true;
80  U = sampledvectors_[fieldi].last();
81  }
82  }
83 
84  // Interpolate the velocity if it has not already been done
85  if (!interpolatedU)
86  {
87  U = td.UInterp_.interpolate(position, celli, facei);
88  }
89 
90  return U;
91 }
92 
93 
94 void Foam::streamlinesParticle::endTrack(trackingData& td)
95 {
96  const label n = sampledPositions_.size();
97 
98  const label trackPartIndex =
99  td.trackForward_ ? trackPartIndex_ : -1 - trackPartIndex_;
100 
101  if (!td.trackForward_) reverse(sampledPositions_);
102  td.allPositions_.append(sampledPositions_);
103  sampledPositions_.clear();
104 
105  td.allTracks_.append(List<label>(n, trackIndex_));
106  td.allTrackParts_.append(List<label>(n, trackPartIndex));
107  trackPartIndex_ ++;
108 
109  if (!td.trackForward_) reverse(sampledAges_);
110  td.allAges_.append(sampledAges_);
111  sampledAges_.clear();
112 
113  forAll(td.scalarInterp_, fieldi)
114  {
115  #define EndTrackType(Type, nullArg) \
116  if (td.Type##Interp_.set(fieldi)) \
117  { \
118  if (!td.trackForward_) reverse(sampled##Type##s_[fieldi]); \
119  td.all##Type##s_[fieldi].append(sampled##Type##s_[fieldi]); \
120  sampled##Type##s_[fieldi].clear(); \
121  }
122  FOR_ALL_FIELD_TYPES(EndTrackType);
123  #undef EndTrackType
124  }
125 }
126 
127 
128 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
129 
131 (
132  const polyMesh& mesh,
133  const vector& position,
134  const label celli,
135  const label lifeTime,
136  const label trackIndex
137 )
138 :
139  particle(mesh, position, celli),
140  lifeTime_(lifeTime),
141  trackIndex_(trackIndex),
142  trackPartIndex_(0),
143  age_(0),
144  transform_(transformer::I)
145 {}
146 
147 
149 (
150  const polyMesh& mesh,
151  Istream& is,
152  bool readFields
153 )
154 :
155  particle(mesh, is, readFields)
156 {
157  if (readFields)
158  {
159  is >> lifeTime_ >> trackIndex_ >> trackPartIndex_ >> age_
160  >> transform_ >> sampledPositions_ >> sampledAges_;
161 
162  #define ReadSampledTypes(Type, nullArg) \
163  List<List<Type>> sampled##Type##s; \
164  is >> sampled##Type##s; \
165  sampled##Type##s_.setSize(sampled##Type##s.size()); \
166  forAll(sampled##Type##s, i) \
167  { \
168  sampled##Type##s_[i].transfer(sampled##Type##s[i]); \
169  }
170  FOR_ALL_FIELD_TYPES(ReadSampledTypes);
171  #undef ReadSampledTypes
172  }
173 
174  // Check state of Istream
175  is.check
176  (
177  "streamlinesParticle::streamlinesParticle"
178  "(const Cloud<streamlinesParticle>&, Istream&, bool)"
179  );
180 }
181 
182 
184 (
185  const streamlinesParticle& p
186 )
187 :
188  particle(p),
189  lifeTime_(p.lifeTime_),
190  trackIndex_(p.trackIndex_),
191  trackPartIndex_(p.trackPartIndex_),
192  age_(p.age_),
193  transform_(p.transform_),
194  sampledPositions_(p.sampledPositions_),
195  sampledAges_(p.sampledAges_)
196  #define SampledTypesInit(Type, nullArg) \
197  , sampled##Type##s_(p.sampled##Type##s_)
198  FOR_ALL_FIELD_TYPES(SampledTypesInit)
199  #undef SampledTypesInit
200 {}
201 
202 
203 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
204 
206 (
207  streamlinesCloud& cloud,
208  trackingData& td,
209  const scalar
210 )
211 {
212  td.keepParticle = true;
213  td.sendToProc = -1;
214 
215  const scalar maxDt = mesh().bounds().mag();
216 
217  while (td.keepParticle && td.sendToProc == -1 && lifeTime_ > 0)
218  {
219  scalar dt = maxDt;
220 
221  // Cross cell in steps:
222  // - at subiter 0 calculate dt to cross cell in nSubCycle steps
223  // - at the last subiter do all of the remaining track
224  for (label subIter = 0; subIter < max(1, td.nSubCycle_); subIter++)
225  {
226  --lifeTime_;
227 
228  // Store current position and sampled velocity.
229  sampledPositions_.append
230  (
231  td.trackOutside_
232  ? transform_.invTransformPosition(position())
233  : position()
234  );
235  sampledAges_.append(age_);
236  vector U = interpolateFields(td, position(), cell(), face());
237 
238  if (!td.trackForward_)
239  {
240  U = -U;
241  }
242 
243  scalar magU = mag(U);
244 
245  if (magU < small)
246  {
247  // Stagnant particle. Might as well stop
248  lifeTime_ = 0;
249  break;
250  }
251 
252  U /= magU;
253 
254  if (td.trackLength_ < great)
255  {
256  // No sub-cycling. Track a set length on each step.
257  dt = td.trackLength_;
258  }
259  else if (subIter == 0)
260  {
261  // Sub-cycling. Cross the cell in nSubCycle steps.
262  particle copy(*this);
263  copy.trackToFace(maxDt*U, 1);
264  dt *= (copy.stepFraction() - stepFraction())/td.nSubCycle_;
265  }
266  else if (subIter == td.nSubCycle_ - 1)
267  {
268  // Sub-cycling. Track the whole cell on the last step.
269  dt = maxDt;
270  }
271 
272  age_ +=
273  (td.trackForward_ ? +1 : -1)
274  *dt/magU
275  *(1 - trackToAndHitFace(dt*U, 0, cloud, td));
276 
277  if (!td.keepParticle || td.sendToProc != -1 || lifeTime_ == 0)
278  {
279  break;
280  }
281  }
282  }
283 
284  if (!td.keepParticle || lifeTime_ == 0)
285  {
286  if (lifeTime_ == 0)
287  {
288  // Failure exit. Particle stagnated or it's life ran out.
289  if (debug)
290  {
291  Pout<< "streamlinesParticle: Removing stagnant particle:"
292  << position() << " sampled positions:"
293  << sampledPositions_.size() << endl;
294  }
295  td.keepParticle = false;
296  }
297  else
298  {
299  // Normal exit. Store last position and fields
300  sampledPositions_.append
301  (
302  td.trackOutside_
303  ? transform_.invTransformPosition(position())
304  : position()
305  );
306  sampledAges_.append(age_);
307  interpolateFields(td, position(), cell(), face());
308 
309  if (debug)
310  {
311  Pout<< "streamlinesParticle: Removing particle:" << position()
312  << " sampled positions:" << sampledPositions_.size()
313  << endl;
314  }
315  }
316 
317  // End this track
318  endTrack(td);
319  }
320 
321  return td.keepParticle;
322 }
323 
324 
326 (
328  trackingData& td
329 )
330 {
331  // Remove particle
332  td.keepParticle = false;
333 }
334 
335 
337 (
339  trackingData& td
340 )
341 {
342  // Remove particle
343  td.keepParticle = false;
344 }
345 
346 
348 (
350  trackingData& td
351 )
352 {
353  // Remove particle
354  td.keepParticle = false;
355 }
356 
357 
359 (
360  streamlinesCloud& cloud,
361  trackingData& td
362 )
363 {
364  const cyclicPolyPatch& cpp =
365  static_cast<const cyclicPolyPatch&>(mesh().boundaryMesh()[patch()]);
366 
367  // End this track
368  if (!td.trackOutside_ && cpp.transform().transformsPosition())
369  {
370  endTrack(td);
371  }
372 
373  // Move across the cyclic
374  particle::hitCyclicPatch(cloud, td);
375 }
376 
377 
379 (
380  const vector& displacement,
381  const scalar fraction,
382  streamlinesCloud& cloud,
383  trackingData& td
384 )
385 {
386  // End this track
387  endTrack(td);
388 
389  // Move across the cyclic
390  particle::hitCyclicAMIPatch(displacement, fraction, cloud, td);
391 }
392 
393 
395 (
396  const vector& displacement,
397  const scalar fraction,
398  const label patchi,
399  streamlinesCloud& cloud,
400  trackingData& td
401 )
402 {
403  // Move across the cyclic
404  const bool result =
406  (
407  displacement,
408  fraction,
409  patchi,
410  cloud,
411  td
412  );
413 
414  // End this track
415  if (result) endTrack(td);
416 
417  return result;
418 }
419 
421 (
422  streamlinesCloud& cloud,
423  trackingData& td
424 )
425 {
426  const processorPolyPatch& ppp =
427  static_cast<const processorPolyPatch&>(mesh().boundaryMesh()[patch()]);
428 
429  // End this track
430  if (!td.trackOutside_ && ppp.transform().transformsPosition())
431  {
432  endTrack(td);
433  }
434 
435  particle::hitProcessorPatch(cloud, td);
436 }
437 
438 
440 (
442  trackingData& td
443 )
444 {
445  // Remove particle
446  td.keepParticle = false;
447 }
448 
449 
451 (
452  const transformer& transform
453 )
454 {
455  transform_ = transform & transform_;
456 }
457 
458 
460 {
461  bool valid = c.size();
462 
464 
465  IOField<label> lifeTime
466  (
467  c.fieldIOobject("lifeTime", IOobject::MUST_READ),
468  valid
469  );
470  c.checkFieldIOobject(c, lifeTime);
471 
472  IOField<label> trackIndex
473  (
474  c.fieldIOobject("trackIndex", IOobject::MUST_READ),
475  valid
476  );
477  c.checkFieldIOobject(c, trackIndex);
478 
479  IOField<label> trackPartIndex
480  (
481  c.fieldIOobject("trackPartIndex", IOobject::MUST_READ),
482  valid
483  );
484  c.checkFieldIOobject(c, trackPartIndex);
485 
486  IOField<scalar> age
487  (
489  valid
490  );
491  c.checkFieldIOobject(c, age);
492 
494  (
495  c.fieldIOobject("transform", IOobject::MUST_READ),
496  valid
497  );
498  //c.checkFieldIOobject(c, transform);
499 
500  vectorFieldIOField sampledPositions
501  (
502  c.fieldIOobject("sampledPositions", IOobject::MUST_READ),
503  valid
504  );
505  c.checkFieldIOobject(c, sampledPositions);
506 
507  scalarFieldIOField sampledAges
508  (
509  c.fieldIOobject("sampledAges", IOobject::MUST_READ),
510  valid
511  );
512  c.checkFieldIOobject(c, sampledAges);
513 
514  label i = 0;
516  {
517  iter().lifeTime_ = lifeTime[i];
518  iter().trackIndex_ = trackIndex[i];
519  iter().trackPartIndex_ = trackPartIndex[i];
520  iter().age_ = age[i];
521  iter().transform_ = transform[i];
522  iter().sampledPositions_.transfer(sampledPositions[i]);
523  iter().sampledAges_.transfer(sampledAges[i]);
524  i++;
525  }
526 }
527 
528 
530 {
532 
533  label np = c.size();
534 
535  IOField<label> lifeTime
536  (
537  c.fieldIOobject("lifeTime", IOobject::NO_READ),
538  np
539  );
540 
541  IOList<label> trackIndex
542  (
543  c.fieldIOobject("trackIndex", IOobject::NO_READ),
544  np
545  );
546 
547  IOList<label> trackPartIndex
548  (
549  c.fieldIOobject("trackPartIndex", IOobject::NO_READ),
550  np
551  );
552 
553  IOField<scalar> age
554  (
556  np
557  );
558 
560  (
561  c.fieldIOobject("transform", IOobject::NO_READ),
562  np
563  );
564 
565  vectorFieldIOField sampledPositions
566  (
567  c.fieldIOobject("sampledPositions", IOobject::NO_READ),
568  np
569  );
570 
571  scalarFieldIOField sampledAges
572  (
573  c.fieldIOobject("sampledAges", IOobject::NO_READ),
574  np
575  );
576 
577  label i = 0;
579  {
580  lifeTime[i] = iter().lifeTime_;
581  trackIndex[i] = iter().trackIndex_;
582  trackPartIndex[i] = iter().trackPartIndex_;
583  age[i] = iter().age_;
584  transform[i] = iter().transform_;
585  sampledPositions[i] = iter().sampledPositions_;
586  sampledAges[i] = iter().sampledAges_;
587  i++;
588  }
589 
590  lifeTime.write(np > 0);
591  trackIndex.write(np > 0);
592  trackPartIndex.write(np > 0);
593  age.write(np > 0);
594  transform.write(np > 0);
595  sampledPositions.write(np > 0);
596  sampledAges.write(np > 0);
597 }
598 
599 
600 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
601 
603 {
604 
605  os << static_cast<const particle&>(p)
606  << token::SPACE << p.lifeTime_
607  << token::SPACE << p.trackIndex_
608  << token::SPACE << p.trackPartIndex_
609  << token::SPACE << p.age_
610  << token::SPACE << p.transform_
611  << token::SPACE << p.sampledPositions_
612  << token::SPACE << p.sampledAges_
613  #define WriteSampledTypes(Type, nullArg) \
614  << token::SPACE << p.sampled##Type##s_
615  FOR_ALL_FIELD_TYPES(WriteSampledTypes);
616  #undef WriteSampledTypes
617 
618  // Check state of Ostream
619  os.check("Ostream& operator<<(Ostream&, const streamlinesParticle&)");
620 
621  return os;
622 }
623 
624 
625 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
void hitProcessorPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.H:125
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
bool move(streamlinesCloud &, trackingData &, const scalar)
Track all particles to their end point.
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space...
Definition: transformer.H:83
#define EndTrackType(Type, nullArg)
U
Definition: pEqn.H:72
void hitCyclicAMIPatch(const vector &, const scalar, streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
void hitWallPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
scalar trackToFace(const vector &displacement, const scalar fraction)
As particle::track, but stops when a face is hit.
Definition: particle.C:606
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
scalar trackToAndHitFace(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Convenience function. Combines trackToFace and hitFace.
A List of objects of type <T> with automated input and output.
Definition: IOList.H:50
label size() const
Return the number of particles in the cloud.
Definition: Cloud.H:190
virtual const transformer & transform() const
Return transformation between the coupled patches.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void clear()
Clear the addressed list, i.e. set the size to zero.
const dimensionedScalar c
Speed of light in a vacuum.
#define WriteSampledTypes(Type, nullArg)
Base particle class.
Definition: particle.H:81
Neighbour processor patch.
static void writeFields(const Cloud< streamlinesParticle > &)
Write.
label cell() const
Return current cell particle is in.
Definition: particleI.H:137
vector invTransformPosition(const vector &v) const
Inverse transform the given position.
Definition: transformerI.H:177
label face() const
Return current face particle is on otherwise -1.
Definition: particleI.H:155
void hitProcessorPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a processorPatch.
Cyclic plane patch.
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
Definition: particle.C:474
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
A Cloud of streamlines particles.
void hitWedgePatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a wedge.
void hitCyclicPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a cyclicPatch.
static const transformer I
Definition: transformer.H:125
Particle class that samples fields as it passes through. Used in streamlines calculation.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void reverse(UList< T > &, const label n)
Definition: UListI.H:334
virtual void transformProperties(const transformer &)
Transform the physical properties of the particle.
void hitCyclicPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a cyclic.
#define ReadSampledTypes(Type, nullArg)
bool transformsPosition() const
Return true if the transformer transforms a point.
Definition: transformerI.H:146
bool hitNonConformalCyclicPatch(const vector &displacement, const scalar fraction, const label patchi, streamlinesCloud &cloud, trackingData &td)
Overridable function to handle the particle hitting an.
scalar stepFraction() const
Return the fraction of time-step completed.
Definition: particleI.H:161
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:459
virtual const transformer & transform() const
Return null transform between processor patches.
void hitSymmetryPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
vector point
Point is a vector.
Definition: point.H:41
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:107
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
static void writeFields(const TrackCloudType &c)
Write the fields associated with the owner cloud.
Ostream & operator<<(Ostream &, const ensightPart &)
void checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
Definition: CloudIO.C:190
static Vector< scalar > uniform(const scalar &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:168
dimensioned< scalar > mag(const dimensioned< Type > &)
static void readFields(Cloud< streamlinesParticle > &)
Read.
void hitCyclicAMIPatch(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a cyclicAMIPatch.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
volScalarField & p
void hitSymmetryPlanePatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
streamlinesParticle(const polyMesh &c, const vector &position, const label celli, const label lifeTime, const label trackIndex)
Construct from components.
bool hitNonConformalCyclicPatch(const vector &displacement, const scalar fraction, const label patchi, TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting an.
label sendToProc
Processor to send the particle to. -1 indicates that this.
Definition: particle.H:111
#define SampledTypesInit(Type, nullArg)
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
Definition: CloudIO.C:170
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
vector position() const
Return current particle position.
Definition: particleI.H:278
DynamicField< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
label patch() const
Return the index of patch that the particle is on.
Definition: particleI.H:272
#define InterpolateType(Type, nullArg)
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:483