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 | Copyright (C) 2011-2017 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"
27 #include "vectorFieldIOField.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 Foam::vector Foam::streamLineParticle::interpolateFields
32 (
33  const trackingData& td,
34  const point& position,
35  const label celli,
36  const label facei
37 )
38 {
39  if (celli == -1)
40  {
42  << "Cell:" << celli << abort(FatalError);
43  }
44 
45  sampledScalars_.setSize(td.vsInterp_.size());
46  forAll(td.vsInterp_, scalarI)
47  {
48  sampledScalars_[scalarI].append
49  (
50  td.vsInterp_[scalarI].interpolate
51  (
52  position,
53  celli,
54  facei
55  )
56  );
57  }
58 
59  sampledVectors_.setSize(td.vvInterp_.size());
60  forAll(td.vvInterp_, vectorI)
61  {
62  sampledVectors_[vectorI].append
63  (
64  td.vvInterp_[vectorI].interpolate
65  (
66  position,
67  celli,
68  facei
69  )
70  );
71  }
72 
73  const DynamicList<vector>& U = sampledVectors_[td.UIndex_];
74 
75  return U.last();
76 }
77 
78 
79 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
80 
82 (
83  const polyMesh& mesh,
84  const vector& position,
85  const label celli,
86  const label lifeTime
87 )
88 :
89  particle(mesh, position, celli),
90  lifeTime_(lifeTime)
91 {}
92 
93 
95 (
96  const polyMesh& mesh,
97  Istream& is,
98  bool readFields
99 )
100 :
101  particle(mesh, is, readFields)
102 {
103  if (readFields)
104  {
105  List<scalarList> sampledScalars;
106  List<vectorList> sampledVectors;
107 
108  is >> lifeTime_ >> sampledPositions_ >> sampledScalars
109  >> sampledVectors;
110 
111  sampledScalars_.setSize(sampledScalars.size());
112  forAll(sampledScalars, i)
113  {
114  sampledScalars_[i].transfer(sampledScalars[i]);
115  }
116  sampledVectors_.setSize(sampledVectors.size());
117  forAll(sampledVectors, i)
118  {
119  sampledVectors_[i].transfer(sampledVectors[i]);
120  }
121  }
122 
123  // Check state of Istream
124  is.check
125  (
126  "streamLineParticle::streamLineParticle"
127  "(const Cloud<streamLineParticle>&, Istream&, bool)"
128  );
129 }
130 
131 
133 (
134  const streamLineParticle& p
135 )
136 :
137  particle(p),
138  lifeTime_(p.lifeTime_),
139  sampledPositions_(p.sampledPositions_),
140  sampledScalars_(p.sampledScalars_)
141 {}
142 
143 
144 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
145 
147 (
148  trackingData& td,
149  const scalar
150 )
151 {
152  td.switchProcessor = false;
153  td.keepParticle = true;
154 
155  const scalar maxDt = mesh().bounds().mag();
156 
157  while (td.keepParticle && !td.switchProcessor && lifeTime_ > 0)
158  {
159  scalar dt = maxDt;
160 
161  // Cross cell in steps:
162  // - at subiter 0 calculate dt to cross cell in nSubCycle steps
163  // - at the last subiter do all of the remaining track
164  for (label subIter = 0; subIter < max(1, td.nSubCycle_); subIter++)
165  {
166  --lifeTime_;
167 
168  // Store current position and sampled velocity.
169  sampledPositions_.append(position());
170  vector U = interpolateFields(td, position(), cell(), face());
171 
172  if (!td.trackForward_)
173  {
174  U = -U;
175  }
176 
177  scalar magU = mag(U);
178 
179  if (magU < SMALL)
180  {
181  // Stagnant particle. Might as well stop
182  lifeTime_ = 0;
183  break;
184  }
185 
186  U /= magU;
187 
188  if (td.trackLength_ < GREAT)
189  {
190  // No sub-cycling. Track a set length on each step.
191  dt = td.trackLength_;
192  }
193  else if (subIter == 0)
194  {
195  // Sub-cycling. Cross the cell in nSubCycle steps.
196  particle copy(*this);
197  copy.trackToFace(maxDt*U, 1);
198  dt *= (copy.stepFraction() - stepFraction())/td.nSubCycle_;
199  }
200  else if (subIter == td.nSubCycle_ - 1)
201  {
202  // Sub-cycling. Track the whole cell on the last step.
203  dt = maxDt;
204  }
205 
206  trackToFace(dt*U, 0, td);
207 
208  if
209  (
210  onFace()
211  || !td.keepParticle
212  || td.switchProcessor
213  || lifeTime_ == 0
214  )
215  {
216  break;
217  }
218  }
219  }
220 
221  if (!td.keepParticle || lifeTime_ == 0)
222  {
223  if (lifeTime_ == 0)
224  {
225  // Failure exit. Particle stagnated or it's life ran out.
226  if (debug)
227  {
228  Pout<< "streamLineParticle: Removing stagnant particle:"
229  << position() << " sampled positions:"
230  << sampledPositions_.size() << endl;
231  }
232  td.keepParticle = false;
233  }
234  else
235  {
236  // Normal exit. Store last position and fields
237  sampledPositions_.append(position());
238  interpolateFields(td, position(), cell(), face());
239 
240  if (debug)
241  {
242  Pout<< "streamLineParticle: Removing particle:" << position()
243  << " sampled positions:" << sampledPositions_.size()
244  << endl;
245  }
246  }
247 
248  // Transfer particle data into trackingData.
249  td.allPositions_.append(vectorList());
250  vectorList& top = td.allPositions_.last();
251  top.transfer(sampledPositions_);
252 
253  forAll(sampledScalars_, i)
254  {
255  td.allScalars_[i].append(scalarList());
256  scalarList& top = td.allScalars_[i].last();
257  top.transfer(sampledScalars_[i]);
258  }
259  forAll(sampledVectors_, i)
260  {
261  td.allVectors_[i].append(vectorList());
262  vectorList& top = td.allVectors_[i].last();
263  top.transfer(sampledVectors_[i]);
264  }
265  }
266 
267  return td.keepParticle;
268 }
269 
270 
272 (
273  const polyPatch&,
274  trackingData& td,
275  const label patchi,
276  const scalar trackFraction,
277  const tetIndices& tetIs
278 )
279 {
280  // Disable generic patch interaction
281  return false;
282 }
283 
284 
286 (
287  const wedgePolyPatch& pp,
288  trackingData& td
289 )
290 {
291  // Remove particle
292  td.keepParticle = false;
293 }
294 
295 
297 (
298  const symmetryPlanePolyPatch& pp,
299  trackingData& td
300 )
301 {
302  // Remove particle
303  td.keepParticle = false;
304 }
305 
306 
308 (
309  const symmetryPolyPatch& pp,
310  trackingData& td
311 )
312 {
313  // Remove particle
314  td.keepParticle = false;
315 }
316 
317 
319 (
320  const cyclicPolyPatch& pp,
321  trackingData& td
322 )
323 {
324  // Remove particle
325  td.keepParticle = false;
326 }
327 
328 
330 (
331  const processorPolyPatch&,
332  trackingData& td
333 )
334 {
335  // Switch particle
336  td.switchProcessor = true;
337 }
338 
339 
341 (
342  const wallPolyPatch& wpp,
343  trackingData& td,
344  const tetIndices&
345 )
346 {
347  // Remove particle
348  td.keepParticle = false;
349 }
350 
351 
353 (
354  const polyPatch& wpp,
355  trackingData& td
356 )
357 {
358  // Remove particle
359  td.keepParticle = false;
360 }
361 
362 
364 {
365 // if (!c.size())
366 // {
367 // return;
368 // }
369  bool valid = c.size();
370 
372 
373  IOField<label> lifeTime
374  (
375  c.fieldIOobject("lifeTime", IOobject::MUST_READ),
376  valid
377  );
378  c.checkFieldIOobject(c, lifeTime);
379 
380  vectorFieldIOField sampledPositions
381  (
382  c.fieldIOobject("sampledPositions", IOobject::MUST_READ),
383  valid
384  );
385  c.checkFieldIOobject(c, sampledPositions);
386 
387  label i = 0;
389  {
390  iter().lifeTime_ = lifeTime[i];
391  iter().sampledPositions_.transfer(sampledPositions[i]);
392  i++;
393  }
394 }
395 
396 
398 {
400 
401  label np = c.size();
402 
403  IOField<label> lifeTime
404  (
405  c.fieldIOobject("lifeTime", IOobject::NO_READ),
406  np
407  );
408  vectorFieldIOField sampledPositions
409  (
410  c.fieldIOobject("sampledPositions", IOobject::NO_READ),
411  np
412  );
413 
414  label i = 0;
416  {
417  lifeTime[i] = iter().lifeTime_;
418  sampledPositions[i] = iter().sampledPositions_;
419  i++;
420  }
421 
422  lifeTime.write(np > 0);
423  sampledPositions.write(np > 0);
424 }
425 
426 
427 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
428 
430 {
431  os << static_cast<const particle&>(p)
432  << token::SPACE << p.lifeTime_
433  << token::SPACE << p.sampledPositions_
434  << token::SPACE << p.sampledScalars_
435  << token::SPACE << p.sampledVectors_;
436 
437  // Check state of Ostream
438  os.check("Ostream& operator<<(Ostream&, const streamLineParticle&)");
439 
440  return os;
441 }
442 
443 
444 // ************************************************************************* //
Symmetry patch for non-planar or multi-plane patches.
static void writeFields(const Cloud< streamLineParticle > &)
Write.
bool move(trackingData &, const scalar)
Track all particles to their end point.
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.H:45
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
List< DynamicList< vectorList > > & allVectors_
U
Definition: pEqn.H:83
void hitSymmetryPlanePatch(const symmetryPlanePolyPatch &, trackingData &td)
Overridable function to handle the particle hitting a.
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
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:453
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 also stops on internal faces.
Definition: particle.C:622
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
void hitCyclicPatch(const cyclicPolyPatch &, trackingData &td)
Overridable function to handle the particle hitting a cyclic.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
label size() const
Definition: Cloud.H:162
static void readFields(CloudType &c)
Read the fields associated with the owner cloud.
List< DynamicList< scalarList > > & allScalars_
List< vector > vectorList
A List of vectors.
Definition: vectorList.H:50
Base particle class.
Definition: particle.H:81
bool hitPatch(const polyPatch &, trackingData &td, const label patchi, const scalar trackFraction, const tetIndices &tetIs)
Overridable function to handle the particle hitting a patch.
Neighbour processor patch.
label cell() const
Return current cell particle is in.
Definition: particleI.H:57
void hitSymmetryPatch(const symmetryPolyPatch &, trackingData &td)
Overridable function to handle the particle hitting a.
label face() const
Return current face particle is on otherwise -1.
Definition: particleI.H:75
static void writeFields(const CloudType &c)
Write the fields associated with the owner cloud.
DynamicList< vectorList > & allPositions_
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:524
Particle class that samples fields as it passes through. Used in streamline calculation.
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:48
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
Wedge front and back plane patch.
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:125
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.
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:180
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:81
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:430
vector point
Point is a vector.
Definition: point.H:41
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
const dimensionedScalar c
Speed of light in a vacuum.
Ostream & operator<<(Ostream &, const ensightPart &)
void checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
Definition: CloudIO.C:186
dimensioned< scalar > mag(const dimensioned< Type > &)
void hitProcessorPatch(const processorPolyPatch &, trackingData &td)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
bool switchProcessor
Flag to switch processor.
Definition: particle.H:122
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
virtual bool write(const bool valid=true) const
Write using setting from DB.
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
Definition: CloudIO.C:166
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
void hitWedgePatch(const wedgePolyPatch &, trackingData &td)
Overridable function to handle the particle hitting a wedge.
vector position() const
Return current particle position.
Definition: particleI.H:204
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void hitWallPatch(const wallPolyPatch &, trackingData &td, const tetIndices &)
Overridable function to handle the particle hitting a wallPatch.