wallBoundedStreamLineParticle.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-2013 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 
27 #include "vectorFieldIOField.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 Foam::vector Foam::wallBoundedStreamLineParticle::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  {
41  FatalErrorIn("wallBoundedStreamLineParticle::interpolateFields(..)")
42  << "Cell:" << cellI << abort(FatalError);
43  }
44 
45  const tetIndices ti = currentTetIndices();
46 
47  const vector U = td.vvInterp_[td.UIndex_].interpolate
48  (
49  position,
50  ti, //cellI,
51  faceI
52  );
53 
54  // Check if at different position
55  if
56  (
57  !sampledPositions_.size()
58  || magSqr(sampledPositions_.last()-position) > Foam::sqr(SMALL)
59  )
60  {
61  // Store the location
62  sampledPositions_.append(position);
63 
64  // Store the scalar fields
65  sampledScalars_.setSize(td.vsInterp_.size());
66  forAll(td.vsInterp_, scalarI)
67  {
68  sampledScalars_[scalarI].append
69  (
70  td.vsInterp_[scalarI].interpolate
71  (
72  position,
73  ti, //cellI,
74  faceI
75  )
76  );
77  }
78 
79  // Store the vector fields
80  sampledVectors_.setSize(td.vvInterp_.size());
81  forAll(td.vvInterp_, vectorI)
82  {
83  vector positionU;
84  if (vectorI == td.UIndex_)
85  {
86  positionU = U;
87  }
88  else
89  {
90  positionU = td.vvInterp_[vectorI].interpolate
91  (
92  position,
93  ti, //cellI,
94  faceI
95  );
96  }
97  sampledVectors_[vectorI].append(positionU);
98  }
99  }
100 
101  return U;
102 }
103 
104 
105 Foam::vector Foam::wallBoundedStreamLineParticle::sample
106 (
107  trackingData& td
108 )
109 {
110  vector U = interpolateFields(td, position(), cell(), tetFace());
111 
112  if (!td.trackForward_)
113  {
114  U = -U;
115  }
116 
117  scalar magU = mag(U);
118 
119  if (magU < SMALL)
120  {
121  // Stagnant particle. Might as well stop
122  lifeTime_ = 0;
123  }
124 
125  U /= magU;
126 
127  return U;
128 }
129 
130 
131 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
132 
134 (
135  const polyMesh& mesh,
136  const vector& position,
137  const label cellI,
138  const label tetFaceI,
139  const label tetPtI,
140  const label meshEdgeStart,
141  const label diagEdge,
142  const label lifeTime
143 )
144 :
146  (
147  mesh,
148  position,
149  cellI,
150  tetFaceI,
151  tetPtI,
152  meshEdgeStart,
153  diagEdge
154  ),
155  lifeTime_(lifeTime)
156 {}
157 
158 
160 (
161  const polyMesh& mesh,
162  Istream& is,
163  bool readFields
164 )
165 :
166  wallBoundedParticle(mesh, is, readFields)
167 {
168  if (readFields)
169  {
170  List<scalarList> sampledScalars;
171  List<vectorList> sampledVectors;
172 
173  is >> lifeTime_
174  >> sampledPositions_ >> sampledScalars >> sampledVectors;
175 
176  sampledScalars_.setSize(sampledScalars.size());
177  forAll(sampledScalars, i)
178  {
179  sampledScalars_[i].transfer(sampledScalars[i]);
180  }
181  sampledVectors_.setSize(sampledVectors.size());
182  forAll(sampledVectors, i)
183  {
184  sampledVectors_[i].transfer(sampledVectors[i]);
185  }
186  }
187 
188  // Check state of Istream
189  is.check
190  (
191  "wallBoundedStreamLineParticle::wallBoundedStreamLineParticle"
192  "(const Cloud<wallBoundedStreamLineParticle>&, Istream&, bool)"
193  );
194 }
195 
196 
198 (
200 )
201 :
203  lifeTime_(p.lifeTime_),
204  sampledPositions_(p.sampledPositions_),
205  sampledScalars_(p.sampledScalars_),
206  sampledVectors_(p.sampledVectors_)
207 {}
208 
209 
210 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
211 
213 (
214  trackingData& td,
215  const scalar trackTime
216 )
217 {
218  wallBoundedStreamLineParticle& p = static_cast
219  <
221  >(*this);
222 
223 
224  // Check position is inside tet
225  //checkInside();
226 
227  td.switchProcessor = false;
228  td.keepParticle = true;
229 
230  scalar tEnd = (1.0 - stepFraction())*trackTime;
231  scalar maxDt = mesh_.bounds().mag();
232 
233  while
234  (
235  td.keepParticle
236  && !td.switchProcessor
237  && lifeTime_ > 0
238  )
239  {
240  // set the lagrangian time-step
241  scalar dt = maxDt;
242 
243  --lifeTime_;
244 
245  // Get sampled velocity and fields. Store if position changed.
246  vector U = sample(td);
247 
248  // !user parameter!
249  if (dt < SMALL)
250  {
251  // Force removal
252  lifeTime_ = 0;
253  break;
254  }
255 
256 
257  if (td.trackLength_ < GREAT)
258  {
259  dt = td.trackLength_;
260  }
261 
262 
263  scalar fraction = trackToEdge(td, position() + dt*U);
264  dt *= fraction;
265 
266  tEnd -= dt;
267  stepFraction() = 1.0 - tEnd/trackTime;
268 
269 
270  if (tEnd <= ROOTVSMALL)
271  {
272  // Force removal
273  lifeTime_ = 0;
274  }
275 
276  if
277  (
278  !td.keepParticle
279  || td.switchProcessor
280  || lifeTime_ == 0
281  )
282  {
283  break;
284  }
285  }
286 
287 
288  if (!td.keepParticle || lifeTime_ == 0)
289  {
290  if (lifeTime_ == 0)
291  {
292  if (debug)
293  {
294  Pout<< "wallBoundedStreamLineParticle :"
295  << " Removing stagnant particle:"
296  << p.position()
297  << " sampled positions:" << sampledPositions_.size()
298  << endl;
299  }
300  td.keepParticle = false;
301  }
302  else
303  {
304  // Normal exit. Store last position and fields
305  sample(td);
306 
307  if (debug)
308  {
309  Pout<< "wallBoundedStreamLineParticle : Removing particle:"
310  << p.position()
311  << " sampled positions:" << sampledPositions_.size()
312  << endl;
313  }
314  }
315 
316  // Transfer particle data into trackingData.
317  {
318  //td.allPositions_.append(sampledPositions_);
319  td.allPositions_.append(vectorList());
320  vectorList& top = td.allPositions_.last();
321  top.transfer(sampledPositions_);
322  }
323 
324  forAll(sampledScalars_, i)
325  {
326  //td.allScalars_[i].append(sampledScalars_[i]);
327  td.allScalars_[i].append(scalarList());
328  scalarList& top = td.allScalars_[i].last();
329  top.transfer(sampledScalars_[i]);
330  }
331  forAll(sampledVectors_, i)
332  {
333  //td.allVectors_[i].append(sampledVectors_[i]);
334  td.allVectors_[i].append(vectorList());
335  vectorList& top = td.allVectors_[i].last();
336  top.transfer(sampledVectors_[i]);
337  }
338  }
339 
340  return td.keepParticle;
341 }
342 
343 
345 (
347 )
348 {
349  if (!c.size())
350  {
351  return;
352  }
353 
355 
356  IOField<label> lifeTime
357  (
358  c.fieldIOobject("lifeTime", IOobject::MUST_READ)
359  );
360  c.checkFieldIOobject(c, lifeTime);
361 
362  vectorFieldIOField sampledPositions
363  (
364  c.fieldIOobject("sampledPositions", IOobject::MUST_READ)
365  );
366  c.checkFieldIOobject(c, sampledPositions);
367 
368  label i = 0;
370  {
371  iter().lifeTime_ = lifeTime[i];
372  iter().sampledPositions_.transfer(sampledPositions[i]);
373  i++;
374  }
375 }
376 
377 
379 (
381 )
382 {
384 
385  label np = c.size();
386 
387  IOField<label> lifeTime
388  (
389  c.fieldIOobject("lifeTime", IOobject::NO_READ),
390  np
391  );
392  vectorFieldIOField sampledPositions
393  (
394  c.fieldIOobject("sampledPositions", IOobject::NO_READ),
395  np
396  );
397 
398  label i = 0;
400  {
401  lifeTime[i] = iter().lifeTime_;
402  sampledPositions[i] = iter().sampledPositions_;
403  i++;
404  }
405 
406  lifeTime.write();
407  sampledPositions.write();
408 }
409 
410 
411 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
412 
413 Foam::Ostream& Foam::operator<<
414 (
415  Ostream& os,
417 )
418 {
419  os << static_cast<const wallBoundedParticle&>(p)
420  << token::SPACE << p.lifeTime_
421  << token::SPACE << p.sampledPositions_
422  << token::SPACE << p.sampledScalars_
423  << token::SPACE << p.sampledVectors_;
424 
425  // Check state of Ostream
426  os.check
427  (
428  "Ostream& operator<<(Ostream&, const wallBoundedStreamLineParticle&)"
429  );
430 
431  return os;
432 }
433 
434 
435 // ************************************************************************* //
label & cell()
Return current cell particle is in.
Definition: particleI.H:621
vector point
Point is a vector.
Definition: point.H:41
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:390
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
dimensioned< scalar > magSqr(const dimensioned< Type > &)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
#define forAllIter(Container, container, iter)
Definition: UList.H:440
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 size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
Particle class that samples fields as it passes through. Used in streamline calculation.
wallBoundedStreamLineParticle(const polyMesh &c, const vector &position, const label cellI, const label tetFaceI, const label tetPtI, const label meshEdgeStart, const label diagEdge, const label lifeTime)
Construct from components.
Class used to pass tracking data to the trackToEdge function.
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
virtual bool write() const
Write using setting from DB.
static void writeFields(const CloudType &)
Write.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static void readFields(Cloud< wallBoundedStreamLineParticle > &)
Read.
volScalarField & p
Definition: createFields.H:51
static void writeFields(const Cloud< wallBoundedStreamLineParticle > &)
Write.
const vector & position() const
Return current particle position.
Definition: particleI.H:603
#define forAll(list, i)
Definition: UList.H:421
scalar trackToEdge(TrackData &td, const vector &endPosition)
Equivalent of trackToFace.
const polyMesh & mesh_
Reference to the polyMesh database.
Definition: particle.H:137
label & tetFace()
Return current tet face particle is in.
Definition: particleI.H:633
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
errorManip< error > abort(error &err)
Definition: errorManip.H:131
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
static void readFields(CloudType &)
Read.
void checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
Definition: CloudIO.C:215
scalar & stepFraction()
Return the fraction of time-step completed.
Definition: particleI.H:834
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< vector > vectorList
A List of vectors.
Definition: vectorList.H:50
label size() const
Definition: Cloud.H:175
error FatalError
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:427
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
wallBoundedParticle(const polyMesh &c, const vector &position, const label cellI, const label tetFaceI, const label tetPtI, const label meshEdgeStart, const label diagEdge)
Construct from components.
bool move(trackingData &, const scalar trackTime)
Track all particles to their end point.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tetIndices currentTetIndices() const
Return the indices of the current tet that the.
Definition: particleI.H:651
U
Definition: pEqn.H:82
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
Definition: CloudIO.C:195
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50