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-2021 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 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 Foam::vector Foam::streamlinesParticle::interpolateFields
34 (
35  const trackingData& td,
36  const point& position,
37  const label celli,
38  const label facei
39 )
40 {
41  if (celli == -1)
42  {
44  << "Cell:" << celli << abort(FatalError);
45  }
46 
47  sampledScalars_.setSize(td.vsInterp_.size());
48  forAll(td.vsInterp_, scalarI)
49  {
50  sampledScalars_[scalarI].append
51  (
52  td.vsInterp_[scalarI].interpolate
53  (
54  position,
55  celli,
56  facei
57  )
58  );
59  }
60 
61  sampledVectors_.setSize(td.vvInterp_.size());
62  forAll(td.vvInterp_, vectorI)
63  {
64  sampledVectors_[vectorI].append
65  (
66  td.vvInterp_[vectorI].interpolate
67  (
68  position,
69  celli,
70  facei
71  )
72  );
73  }
74 
75  const DynamicList<vector>& U = sampledVectors_[td.UIndex_];
76 
77  return U.last();
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
82 
84 (
85  const polyMesh& mesh,
86  const vector& position,
87  const label celli,
88  const label lifeTime
89 )
90 :
91  particle(mesh, position, celli),
92  lifeTime_(lifeTime),
93  age_(0)
94 {}
95 
96 
98 (
99  const polyMesh& mesh,
100  Istream& is,
101  bool readFields
102 )
103 :
104  particle(mesh, is, readFields)
105 {
106  if (readFields)
107  {
108  List<scalarList> sampledScalars;
109  List<vectorList> sampledVectors;
110 
111  is >> lifeTime_ >> age_ >> sampledPositions_ >> sampledTimes_
112  >> sampledScalars >> sampledVectors;
113 
114  sampledScalars_.setSize(sampledScalars.size());
115  forAll(sampledScalars, i)
116  {
117  sampledScalars_[i].transfer(sampledScalars[i]);
118  }
119  sampledVectors_.setSize(sampledVectors.size());
120  forAll(sampledVectors, i)
121  {
122  sampledVectors_[i].transfer(sampledVectors[i]);
123  }
124  }
125 
126  // Check state of Istream
127  is.check
128  (
129  "streamlinesParticle::streamlinesParticle"
130  "(const Cloud<streamlinesParticle>&, Istream&, bool)"
131  );
132 }
133 
134 
136 (
137  const streamlinesParticle& p
138 )
139 :
140  particle(p),
141  lifeTime_(p.lifeTime_),
142  age_(p.age_),
143  sampledPositions_(p.sampledPositions_),
144  sampledTimes_(p.sampledTimes_),
145  sampledScalars_(p.sampledScalars_),
146  sampledVectors_(p.sampledVectors_)
147 {}
148 
149 
150 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 
153 (
154  streamlinesCloud& cloud,
155  trackingData& td,
156  const scalar
157 )
158 {
159  td.switchProcessor = false;
160  td.keepParticle = true;
161 
162  const scalar maxDt = mesh().bounds().mag();
163 
164  while (td.keepParticle && !td.switchProcessor && lifeTime_ > 0)
165  {
166  scalar dt = maxDt;
167 
168  // Cross cell in steps:
169  // - at subiter 0 calculate dt to cross cell in nSubCycle steps
170  // - at the last subiter do all of the remaining track
171  for (label subIter = 0; subIter < max(1, td.nSubCycle_); subIter++)
172  {
173  --lifeTime_;
174 
175  // Store current position and sampled velocity.
176  sampledPositions_.append(position());
177  sampledTimes_.append(age_);
178  vector U = interpolateFields(td, position(), cell(), face());
179 
180  if (!td.trackForward_)
181  {
182  U = -U;
183  }
184 
185  scalar magU = mag(U);
186 
187  if (magU < small)
188  {
189  // Stagnant particle. Might as well stop
190  lifeTime_ = 0;
191  break;
192  }
193 
194  U /= magU;
195 
196  if (td.trackLength_ < great)
197  {
198  // No sub-cycling. Track a set length on each step.
199  dt = td.trackLength_;
200  }
201  else if (subIter == 0)
202  {
203  // Sub-cycling. Cross the cell in nSubCycle steps.
204  particle copy(*this);
205  copy.trackToFace(maxDt*U, 1);
206  dt *= (copy.stepFraction() - stepFraction())/td.nSubCycle_;
207  }
208  else if (subIter == td.nSubCycle_ - 1)
209  {
210  // Sub-cycling. Track the whole cell on the last step.
211  dt = maxDt;
212  }
213 
214  age_ +=
215  (td.trackForward_ ? +1 : -1)
216  *dt
217  *(1 - trackToAndHitFace(dt*U, 0, cloud, td));
218 
219  if
220  (
221  onFace()
222  || !td.keepParticle
223  || td.switchProcessor
224  || lifeTime_ == 0
225  )
226  {
227  break;
228  }
229  }
230  }
231 
232  if (!td.keepParticle || lifeTime_ == 0)
233  {
234  if (lifeTime_ == 0)
235  {
236  // Failure exit. Particle stagnated or it's life ran out.
237  if (debug)
238  {
239  Pout<< "streamlinesParticle: Removing stagnant particle:"
240  << position() << " sampled positions:"
241  << sampledPositions_.size() << endl;
242  }
243  td.keepParticle = false;
244  }
245  else
246  {
247  // Normal exit. Store last position and fields
248  sampledPositions_.append(position());
249  sampledTimes_.append(age_);
250  interpolateFields(td, position(), cell(), face());
251 
252  if (debug)
253  {
254  Pout<< "streamlinesParticle: Removing particle:" << position()
255  << " sampled positions:" << sampledPositions_.size()
256  << endl;
257  }
258  }
259 
260  // Transfer particle data into trackingData.
261  {
262  td.allPositions_.append(vectorList());
263  vectorList& top = td.allPositions_.last();
264  top.transfer(sampledPositions_);
265  }
266  {
267  td.allTimes_.append(scalarList());
268  scalarList& top = td.allTimes_.last();
269  top.transfer(sampledTimes_);
270  }
271  forAll(sampledScalars_, i)
272  {
273  td.allScalars_[i].append(scalarList());
274  scalarList& top = td.allScalars_[i].last();
275  top.transfer(sampledScalars_[i]);
276  }
277  forAll(sampledVectors_, i)
278  {
279  td.allVectors_[i].append(vectorList());
280  vectorList& top = td.allVectors_[i].last();
281  top.transfer(sampledVectors_[i]);
282  }
283  }
284 
285  return td.keepParticle;
286 }
287 
288 
290 {
291  // Disable generic patch interaction
292  return false;
293 }
294 
295 
297 (
299  trackingData& td
300 )
301 {
302  // Remove particle
303  td.keepParticle = false;
304 }
305 
306 
308 (
310  trackingData& td
311 )
312 {
313  // Remove particle
314  td.keepParticle = false;
315 }
316 
317 
319 (
321  trackingData& td
322 )
323 {
324  // Remove particle
325  td.keepParticle = false;
326 }
327 
328 
330 (
332  trackingData& td
333 )
334 {
335  // Remove particle
336  td.keepParticle = false;
337 }
338 
339 
341 (
342  const vector&,
343  const scalar,
345  trackingData& td
346 )
347 {
348  // Remove particle
349  td.keepParticle = false;
350 }
351 
352 
354 (
355  const vector&,
356  const scalar,
358  trackingData& td
359 )
360 {
361  // Remove particle
362  td.keepParticle = false;
363 }
364 
365 
367 (
368  const vector&,
369  const scalar,
371  trackingData& td
372 )
373 {
374  // Remove particle
375  td.keepParticle = false;
376 }
377 
378 
380 (
382  trackingData& td
383 )
384 {
385  // Switch particle
386  td.switchProcessor = true;
387 }
388 
389 
391 (
393  trackingData& td
394 )
395 {
396  // Remove particle
397  td.keepParticle = false;
398 }
399 
400 
402 {
403 // if (!c.size())
404 // {
405 // return;
406 // }
407  bool valid = c.size();
408 
410 
411  IOField<label> lifeTime
412  (
413  c.fieldIOobject("lifeTime", IOobject::MUST_READ),
414  valid
415  );
416  c.checkFieldIOobject(c, lifeTime);
417 
418  vectorFieldIOField sampledPositions
419  (
420  c.fieldIOobject("sampledPositions", IOobject::MUST_READ),
421  valid
422  );
423  c.checkFieldIOobject(c, sampledPositions);
424 
425  scalarFieldIOField sampledTimes
426  (
427  c.fieldIOobject("sampledTimes", IOobject::MUST_READ),
428  valid
429  );
430  c.checkFieldIOobject(c, sampledTimes);
431 
432  label i = 0;
434  {
435  iter().lifeTime_ = lifeTime[i];
436  iter().sampledPositions_.transfer(sampledPositions[i]);
437  iter().sampledTimes_.transfer(sampledTimes[i]);
438  i++;
439  }
440 }
441 
442 
444 {
446 
447  label np = c.size();
448 
449  IOField<label> lifeTime
450  (
451  c.fieldIOobject("lifeTime", IOobject::NO_READ),
452  np
453  );
454  vectorFieldIOField sampledPositions
455  (
456  c.fieldIOobject("sampledPositions", IOobject::NO_READ),
457  np
458  );
459  scalarFieldIOField sampledTimes
460  (
461  c.fieldIOobject("sampledTimes", IOobject::NO_READ),
462  np
463  );
464 
465  label i = 0;
467  {
468  lifeTime[i] = iter().lifeTime_;
469  sampledPositions[i] = iter().sampledPositions_;
470  sampledTimes[i] = iter().sampledTimes_;
471  i++;
472  }
473 
474  lifeTime.write(np > 0);
475  sampledPositions.write(np > 0);
476  sampledTimes.write(np > 0);
477 }
478 
479 
480 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
481 
483 {
484  os << static_cast<const particle&>(p)
485  << token::SPACE << p.lifeTime_
486  << token::SPACE << p.age_
487  << token::SPACE << p.sampledPositions_
488  << token::SPACE << p.sampledTimes_
489  << token::SPACE << p.sampledScalars_
490  << token::SPACE << p.sampledVectors_;
491 
492  // Check state of Ostream
493  os.check("Ostream& operator<<(Ostream&, const streamlinesParticle&)");
494 
495  return os;
496 }
497 
498 
499 // ************************************************************************* //
void hitCyclicACMIPatch(const vector &, const scalar, streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
void hitProcessorPatch(streamlinesCloud &, trackingData &)
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
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
bool move(streamlinesCloud &, trackingData &, const scalar)
Track all particles to their end point.
void hitCyclicAMIPatch(const vector &, const scalar, streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
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:645
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
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.
label size() const
Return the number of particles in the cloud.
Definition: Cloud.H:148
const dimensionedScalar c
Speed of light in a vacuum.
List< vector > vectorList
A List of vectors.
Definition: vectorList.H:50
Base particle class.
Definition: particle.H:83
static void writeFields(const Cloud< streamlinesParticle > &)
Write.
label cell() const
Return current cell particle is in.
Definition: particleI.H:137
DynamicList< vectorList > & allPositions_
List< DynamicList< scalarList > > & allScalars_
label face() const
Return current face particle is on otherwise -1.
Definition: particleI.H:155
bool switchProcessor
Flag to switch processor.
Definition: particle.H:109
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
Definition: particle.C:513
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
A Cloud of streamlines particles.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
bool hitPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a patch.
void hitWedgePatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a wedge.
Particle class that samples fields as it passes through. Used in streamlines calculation.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
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
bool onFace() const
Is the particle on a face?
Definition: particleI.H:254
void hitCyclicPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a cyclic.
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:440
void hitSymmetryPatch(streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
U
Definition: pEqn.H:72
vector point
Point is a vector.
Definition: point.H:41
List< DynamicList< vectorList > > & allVectors_
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:112
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
void hitCyclicRepeatAMIPatch(const vector &, const scalar, streamlinesCloud &, trackingData &)
Overridable function to handle the particle hitting a.
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:187
dimensioned< scalar > mag(const dimensioned< Type > &)
static void readFields(Cloud< streamlinesParticle > &)
Read.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
virtual bool write(const bool write=true) const
Write using setting from DB.
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)
Construct from components.
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
Definition: CloudIO.C:167
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
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342