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