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-2023 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  }
123  #undef EndTrackType
124  }
125 }
126 
127 
128 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
129 
131 (
132  const polyMesh& mesh,
133  const vector& position,
134  const label celli,
135  label& nLocateBoundaryHits,
136  const label lifeTime,
137  const label trackIndex
138 )
139 :
140  particle(mesh, position, celli, nLocateBoundaryHits),
141  lifeTime_(lifeTime),
142  trackIndex_(trackIndex),
143  trackPartIndex_(0),
144  age_(0),
145  transform_(transformer::I)
146 {}
147 
148 
150 :
151  particle(is, readFields)
152 {
153  if (readFields)
154  {
155  is >> lifeTime_ >> trackIndex_ >> trackPartIndex_ >> age_
156  >> transform_ >> sampledPositions_ >> sampledAges_;
157 
158  #define ReadSampledTypes(Type, nullArg) \
159  List<List<Type>> sampled##Type##s; \
160  is >> sampled##Type##s; \
161  sampled##Type##s_.setSize(sampled##Type##s.size()); \
162  forAll(sampled##Type##s, i) \
163  { \
164  sampled##Type##s_[i].transfer(sampled##Type##s[i]); \
165  }
166  FOR_ALL_FIELD_TYPES(ReadSampledTypes);
167  #undef ReadSampledTypes
168  }
169 
170  // Check state of Istream
171  is.check
172  (
173  "streamlinesParticle::streamlinesParticle"
174  "(const Cloud<streamlinesParticle>&, Istream&, bool)"
175  );
176 }
177 
178 
180 (
181  const streamlinesParticle& p
182 )
183 :
184  particle(p),
185  lifeTime_(p.lifeTime_),
186  trackIndex_(p.trackIndex_),
187  trackPartIndex_(p.trackPartIndex_),
188  age_(p.age_),
189  transform_(p.transform_),
190  sampledPositions_(p.sampledPositions_),
191  sampledAges_(p.sampledAges_)
192  #define SampledTypesInit(Type, nullArg) \
193  , sampled##Type##s_(p.sampled##Type##s_)
195  #undef SampledTypesInit
196 {}
197 
198 
199 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
200 
202 (
204  trackingData& td
205 )
206 {
207  td.keepParticle = true;
208  td.sendToProc = -1;
209 
210  const scalar maxDt = td.mesh.bounds().mag();
211 
212  while (td.keepParticle && td.sendToProc == -1 && lifeTime_ > 0)
213  {
214  scalar dt = maxDt;
215 
216  // Cross cell in steps:
217  // - at subiter 0 calculate dt to cross cell in nSubCycle steps
218  // - at the last subiter do all of the remaining track
219  for (label subIter = 0; subIter < max(1, td.nSubCycle_); subIter++)
220  {
221  --lifeTime_;
222 
223  // Store current position and sampled velocity.
224  sampledPositions_.append
225  (
226  td.trackOutside_
227  ? transform_.invTransformPosition(position(td.mesh))
228  : position(td.mesh)
229  );
230  sampledAges_.append(age_);
231  vector U = interpolateFields(td, position(td.mesh), cell(), face());
232 
233  if (!td.trackForward_)
234  {
235  U = -U;
236  }
237 
238  scalar magU = mag(U);
239 
240  if (magU < small)
241  {
242  // Stagnant particle. Might as well stop
243  lifeTime_ = 0;
244  break;
245  }
246 
247  U /= magU;
248 
249  if (td.trackLength_ < great)
250  {
251  // No sub-cycling. Track a set length on each step.
252  dt = td.trackLength_;
253  }
254  else if (subIter == 0)
255  {
256  // Sub-cycling. Cross the cell in nSubCycle steps.
257  particle copy(*this);
258  copy.trackToFace(td.mesh, maxDt*U, 1);
259  dt *= (copy.stepFraction() - stepFraction())/td.nSubCycle_;
260  }
261  else if (subIter == td.nSubCycle_ - 1)
262  {
263  // Sub-cycling. Track the whole cell on the last step.
264  dt = maxDt;
265  }
266 
267  age_ +=
268  (td.trackForward_ ? +1 : -1)
269  *dt/magU
270  *(1 - trackToAndHitFace(dt*U, 0, cloud, td));
271 
272  if (!td.keepParticle || td.sendToProc != -1 || lifeTime_ == 0)
273  {
274  break;
275  }
276  }
277  }
278 
279  if (!td.keepParticle || lifeTime_ == 0)
280  {
281  if (lifeTime_ == 0)
282  {
283  // Failure exit. Particle stagnated or its life ran out.
284  if (debug)
285  {
286  Pout<< "streamlinesParticle: Removing stagnant particle:"
287  << position(td.mesh) << " sampled positions:"
288  << sampledPositions_.size() << endl;
289  }
290  td.keepParticle = false;
291  }
292  else
293  {
294  // Normal exit. Store last position and fields
295  sampledPositions_.append
296  (
297  td.trackOutside_
298  ? transform_.invTransformPosition(position(td.mesh))
299  : position(td.mesh)
300  );
301  sampledAges_.append(age_);
302  interpolateFields(td, position(td.mesh), cell(), face());
303 
304  if (debug)
305  {
306  Pout<< "streamlinesParticle: Removing particle:"
307  << position(td.mesh) << " sampled positions:"
308  << sampledPositions_.size() << endl;
309  }
310  }
311 
312  // End this track
313  endTrack(td);
314  }
315 
316  return td.keepParticle;
317 }
318 
319 
321 (
323  trackingData& td
324 )
325 {
326  // Remove particle
327  td.keepParticle = false;
328 }
329 
330 
332 (
334  trackingData& td
335 )
336 {
337  // Remove particle
338  td.keepParticle = false;
339 }
340 
341 
343 (
345  trackingData& td
346 )
347 {
348  // Remove particle
349  td.keepParticle = false;
350 }
351 
352 
354 (
356  trackingData& td
357 )
358 {
359  const cyclicPolyPatch& cpp =
360  static_cast<const cyclicPolyPatch&>
361  (
362  td.mesh.boundaryMesh()[patch(td.mesh)]
363  );
364 
365  // End this track
366  if (!td.trackOutside_ && cpp.transform().transformsPosition())
367  {
368  endTrack(td);
369  }
370 
371  // Move across the cyclic
373 }
374 
375 
377 (
378  const vector& displacement,
379  const scalar fraction,
380  const label patchi,
382  trackingData& td
383 )
384 {
385  // Move across the cyclic
386  const bool result =
388  (
389  displacement,
390  fraction,
391  patchi,
392  cloud,
393  td
394  );
395 
396  // End this track
397  if (result) endTrack(td);
398 
399  return result;
400 }
401 
403 (
405  trackingData& td
406 )
407 {
408  const processorPolyPatch& ppp =
409  static_cast<const processorPolyPatch&>
410  (
411  td.mesh.boundaryMesh()[patch(td.mesh)]
412  );
413 
414  // End this track
415  if (!td.trackOutside_ && ppp.transform().transformsPosition())
416  {
417  endTrack(td);
418  }
419 
421 }
422 
423 
425 (
427  trackingData& td
428 )
429 {
430  // Remove particle
431  td.keepParticle = false;
432 }
433 
434 
436 (
437  const transformer& transform
438 )
439 {
440  transform_ = transform & transform_;
441 }
442 
443 
445 {
446  bool valid = c.size();
447 
449 
450  IOField<label> lifeTime
451  (
452  c.fieldIOobject("lifeTime", IOobject::MUST_READ),
453  valid
454  );
455  c.checkFieldIOobject(c, lifeTime);
456 
457  IOField<label> trackIndex
458  (
459  c.fieldIOobject("trackIndex", IOobject::MUST_READ),
460  valid
461  );
462  c.checkFieldIOobject(c, trackIndex);
463 
464  IOField<label> trackPartIndex
465  (
466  c.fieldIOobject("trackPartIndex", IOobject::MUST_READ),
467  valid
468  );
469  c.checkFieldIOobject(c, trackPartIndex);
470 
471  IOField<scalar> age
472  (
473  c.fieldIOobject("age", IOobject::MUST_READ),
474  valid
475  );
476  c.checkFieldIOobject(c, age);
477 
479  (
480  c.fieldIOobject("transform", IOobject::MUST_READ),
481  valid
482  );
483  //c.checkFieldIOobject(c, transform);
484 
485  vectorFieldIOField sampledPositions
486  (
487  c.fieldIOobject("sampledPositions", IOobject::MUST_READ),
488  valid
489  );
490  c.checkFieldIOobject(c, sampledPositions);
491 
492  scalarFieldIOField sampledAges
493  (
494  c.fieldIOobject("sampledAges", IOobject::MUST_READ),
495  valid
496  );
497  c.checkFieldIOobject(c, sampledAges);
498 
499  label i = 0;
501  {
502  iter().lifeTime_ = lifeTime[i];
503  iter().trackIndex_ = trackIndex[i];
504  iter().trackPartIndex_ = trackPartIndex[i];
505  iter().age_ = age[i];
506  iter().transform_ = transform[i];
507  iter().sampledPositions_.transfer(sampledPositions[i]);
508  iter().sampledAges_.transfer(sampledAges[i]);
509  i++;
510  }
511 }
512 
513 
515 {
517 
518  label np = c.size();
519 
520  IOField<label> lifeTime
521  (
522  c.fieldIOobject("lifeTime", IOobject::NO_READ),
523  np
524  );
525 
526  IOList<label> trackIndex
527  (
528  c.fieldIOobject("trackIndex", IOobject::NO_READ),
529  np
530  );
531 
532  IOList<label> trackPartIndex
533  (
534  c.fieldIOobject("trackPartIndex", IOobject::NO_READ),
535  np
536  );
537 
538  IOField<scalar> age
539  (
540  c.fieldIOobject("age", IOobject::NO_READ),
541  np
542  );
543 
545  (
546  c.fieldIOobject("transform", IOobject::NO_READ),
547  np
548  );
549 
550  vectorFieldIOField sampledPositions
551  (
552  c.fieldIOobject("sampledPositions", IOobject::NO_READ),
553  np
554  );
555 
556  scalarFieldIOField sampledAges
557  (
558  c.fieldIOobject("sampledAges", IOobject::NO_READ),
559  np
560  );
561 
562  label i = 0;
564  {
565  lifeTime[i] = iter().lifeTime_;
566  trackIndex[i] = iter().trackIndex_;
567  trackPartIndex[i] = iter().trackPartIndex_;
568  age[i] = iter().age_;
569  transform[i] = iter().transform_;
570  sampledPositions[i] = iter().sampledPositions_;
571  sampledAges[i] = iter().sampledAges_;
572  i++;
573  }
574 
575  lifeTime.write(np > 0);
576  trackIndex.write(np > 0);
577  trackPartIndex.write(np > 0);
578  age.write(np > 0);
579  transform.write(np > 0);
580  sampledPositions.write(np > 0);
581  sampledAges.write(np > 0);
582 }
583 
584 
585 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
586 
588 {
589 
590  os << static_cast<const particle&>(p)
591  << token::SPACE << p.lifeTime_
592  << token::SPACE << p.trackIndex_
593  << token::SPACE << p.trackPartIndex_
594  << token::SPACE << p.age_
595  << token::SPACE << p.transform_
596  << token::SPACE << p.sampledPositions_
597  << token::SPACE << p.sampledAges_
598  #define WriteSampledTypes(Type, nullArg) \
599  << token::SPACE << p.sampled##Type##s_
601  #undef WriteSampledTypes
602 
603  // Check state of Ostream
604  os.check("Ostream& operator<<(Ostream&, const streamlinesParticle&)");
605 
606  return os;
607 }
608 
609 
610 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
A primitive field of type <Type> with automated input and output.
Definition: IOField.H:53
A List of objects of type <Type> with automated input and output.
Definition: IOList.H:124
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:168
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:96
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
Cyclic plane patch.
virtual const transformer & transform() const
Return transformation between the coupled patches.
Ostream & write(Ostream &os) const
Write.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
label sendToProc
Processor to send the particle to. -1 indicates that this.
Definition: particle.H:113
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:109
const polyMesh & mesh
Reference to the mesh.
Definition: particle.H:106
Base particle class.
Definition: particle.H:83
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
void hitCyclicPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
void hitProcessorPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
scalar trackToFace(const polyMesh &mesh, const vector &displacement, const scalar fraction)
As particle::track, but stops when a face is hit.
Definition: particle.C:565
scalar stepFraction() const
Return the fraction of time-step completed.
Definition: particleI.H:164
vector position(const polyMesh &mesh) const
Return current particle position.
Definition: particleI.H:255
bool hitNonConformalCyclicPatch(const vector &displacement, const scalar fraction, const label patchi, TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting an.
static void writeFields(const TrackCloudType &c)
Write the fields associated with the owner cloud.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:410
Neighbour processor patch.
virtual const transformer & transform() const
Return null transform between processor patches.
virtual bool write(const bool write=true) const
Write using setting from DB.
A Cloud of streamlines particles.
Particle class that samples fields as it passes through. Used in streamlines calculation.
bool move(streamlinesCloud &, trackingData &)
Track all particles to their end point.
static void writeFields(const Cloud< streamlinesParticle > &)
Write.
void hitProcessorPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
bool hitNonConformalCyclicPatch(const vector &displacement, const scalar fraction, const label patchi, streamlinesCloud &cloud, trackingData &td)
Overridable function to handle the particle hitting an.
void hitSymmetryPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
void hitCyclicPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a cyclic.
void hitWallPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
void hitSymmetryPlanePatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
static void readFields(Cloud< streamlinesParticle > &)
Read.
streamlinesParticle(const polyMesh &mesh, const vector &position, const label celli, label &nLocateBoundaryHits, const label lifeTime, const label trackIndex)
Construct from components.
virtual void transformProperties(const transformer &)
Transform the physical properties of the particle.
void hitWedgePatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a wedge.
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space.
Definition: transformer.H:84
bool transformsPosition() const
Return true if the transformer transforms a point.
Definition: transformerI.H:146
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
U
Definition: pEqn.H:72
bool valid(const PtrList< ModelType > &l)
const dimensionedScalar c
Speed of light in a vacuum.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
errorManip< error > abort(error &err)
Definition: errorManip.H:131
static const Identity< scalar > I
Definition: Identity.H:93
vector point
Point is a vector.
Definition: point.H:41
void reverse(UList< T > &, const label n)
Definition: UListI.H:334
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:504
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
error FatalError
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
volScalarField & p
#define InterpolateType(Type, nullArg)
#define WriteSampledTypes(Type, nullArg)
#define EndTrackType(Type, nullArg)
#define SampledTypesInit(Type, nullArg)
#define ReadSampledTypes(Type, nullArg)