wallBoundedParticle.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-2015 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 "wallBoundedParticle.H"
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 const std::size_t Foam::wallBoundedParticle::sizeofFields_
31 (
32  sizeof(wallBoundedParticle) - sizeof(particle)
33 );
34 
35 
36 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37 
39 {
40  if ((meshEdgeStart_ != -1) == (diagEdge_ != -1))
41  {
42  FatalErrorIn("wallBoundedParticle::currentEdge() const")
43  << "Particle:"
44  << info()
45  << "cannot both be on a mesh edge and a face-diagonal edge."
46  << " meshEdgeStart_:" << meshEdgeStart_
47  << " diagEdge_:" << diagEdge_
48  << abort(FatalError);
49  }
50 
51  const Foam::face& f = mesh_.faces()[tetFace()];
52 
53  if (meshEdgeStart_ != -1)
54  {
55  return edge(f[meshEdgeStart_], f.nextLabel(meshEdgeStart_));
56  }
57  else
58  {
59  label faceBasePtI = mesh_.tetBasePtIs()[tetFace()];
60  label diagPtI = (faceBasePtI+diagEdge_)%f.size();
61 
62  return edge(f[faceBasePtI], f[diagPtI]);
63  }
64 }
65 
66 
68 (
69  const edge& meshEdge
70 )
71 {
72  // Update tetFace, tetPt
74 
75  // Update face to be same as tracking one
76  face() = tetFace();
77 
78  // And adapt meshEdgeStart_.
79  const Foam::face& f = mesh_.faces()[tetFace()];
80  label fp = findIndex(f, meshEdge[0]);
81 
82  if (f.nextLabel(fp) == meshEdge[1])
83  {
84  meshEdgeStart_ = fp;
85  }
86  else
87  {
88  label fpMin1 = f.rcIndex(fp);
89 
90  if (f[fpMin1] == meshEdge[1])
91  {
92  meshEdgeStart_ = fpMin1;
93  }
94  else
95  {
97  (
98  "wallBoundedParticle::crossEdgeConnectedFace"
99  "(const edge&)"
100  ) << "Problem :"
101  << " particle:"
102  << info()
103  << "face:" << tetFace()
104  << " verts:" << f
105  << " meshEdge:" << meshEdge
106  << abort(FatalError);
107  }
108  }
109 
110  diagEdge_ = -1;
111 
112  // Check that still on same mesh edge
113  const edge eNew(f[meshEdgeStart_], f.nextLabel(meshEdgeStart_));
114  if (eNew != meshEdge)
115  {
117  (
118  "wallBoundedParticle::crossEdgeConnectedFace"
119  "(const edge&)"
120  ) << "Problem" << abort(FatalError);
121  }
122 }
123 
124 
126 {
127  if (diagEdge_ == -1)
128  {
129  FatalErrorIn("wallBoundedParticle::crossDiagonalEdge()")
130  << "Particle:"
131  << info()
132  << "not on a diagonal edge" << abort(FatalError);
133  }
134  if (meshEdgeStart_ != -1)
135  {
136  FatalErrorIn("wallBoundedParticle::crossDiagonalEdge()")
137  << "Particle:"
138  << info()
139  << "meshEdgeStart_:" << meshEdgeStart_ << abort(FatalError);
140  }
141 
142  const Foam::face& f = mesh_.faces()[tetFace()];
143 
144  // tetPtI starts from 1, goes up to f.size()-2
145 
146  if (tetPt() == diagEdge_)
147  {
148  tetPt() = f.rcIndex(tetPt());
149  }
150  else
151  {
152  label nextTetPt = f.fcIndex(tetPt());
153  if (diagEdge_ == nextTetPt)
154  {
155  tetPt() = nextTetPt;
156  }
157  else
158  {
159  FatalErrorIn("wallBoundedParticle::crossDiagonalEdge()")
160  << "Particle:"
161  << info()
162  << "tetPt:" << tetPt()
163  << " diagEdge:" << diagEdge_ << abort(FatalError);
164  }
165  }
166 
167  meshEdgeStart_ = -1;
168 }
169 
170 
172 (
173  const vector& endPosition,
174  label& minEdgeI
175 )
176 {
177  // Track p from position to endPosition
178  const triFace tri(currentTetIndices().faceTriIs(mesh_));
179  vector n = tri.normal(mesh_.points());
180  n /= mag(n)+VSMALL;
181 
182  // Check which edge intersects the trajectory.
183  // Project trajectory onto triangle
184  minEdgeI = -1;
185  scalar minS = 1; // end position
186 
187  edge currentE(-1, -1);
188  if (meshEdgeStart_ != -1 || diagEdge_ != -1)
189  {
190  currentE = currentEdge();
191  }
192 
193  // Determine path along line position+s*d to see where intersections
194  // are.
195 
196  forAll(tri, i)
197  {
198  label j = tri.fcIndex(i);
199 
200  const point& pt0 = mesh_.points()[tri[i]];
201  const point& pt1 = mesh_.points()[tri[j]];
202 
203  if (edge(tri[i], tri[j]) == currentE)
204  {
205  // Do not check particle is on
206  continue;
207  }
208 
209  // Outwards pointing normal
210  vector edgeNormal = (pt1-pt0)^n;
211 
212  edgeNormal /= mag(edgeNormal)+VSMALL;
213 
214  // Determine whether position and end point on either side of edge.
215  scalar sEnd = (endPosition-pt0)&edgeNormal;
216  if (sEnd >= 0)
217  {
218  // endPos is outside triangle. position() should always be
219  // inside.
220  scalar sStart = (position()-pt0)&edgeNormal;
221  if (mag(sEnd-sStart) > VSMALL)
222  {
223  scalar s = sStart/(sStart-sEnd);
224 
225  if (s >= 0 && s < minS)
226  {
227  minS = s;
228  minEdgeI = i;
229  }
230  }
231  }
232  }
233 
234  if (minEdgeI != -1)
235  {
236  position() += minS*(endPosition-position());
237  }
238  else
239  {
240  // Did not hit any edge so tracked to the end position. Set position
241  // without any calculation to avoid truncation errors.
242  position() = endPosition;
243  minS = 1.0;
244  }
245 
246  // Project position() back onto plane of triangle
247  const point& triPt = mesh_.points()[tri[0]];
248  position() -= ((position()-triPt)&n)*n;
249 
250  return minS;
251 }
252 
253 
255 (
256  const point& endPosition
257 ) const
258 {
259  const triFace triVerts(currentTetIndices().faceTriIs(mesh_));
260  const edge currentE = currentEdge();
261 
262  //if (debug)
263  {
264  if
265  (
266  currentE[0] == currentE[1]
267  || findIndex(triVerts, currentE[0]) == -1
268  || findIndex(triVerts, currentE[1]) == -1
269  )
270  {
272  (
273  "wallBoundedParticle::isTriAlongTrack"
274  "(const point&)"
275  ) << "Edge " << currentE << " not on triangle " << triVerts
276  << info()
277  << abort(FatalError);
278  }
279  }
280 
281 
282  const vector dir = endPosition-position();
283 
284  // Get normal of currentE
285  vector n = triVerts.normal(mesh_.points());
286  n /= mag(n);
287 
288  forAll(triVerts, i)
289  {
290  label j = triVerts.fcIndex(i);
291  const point& pt0 = mesh_.points()[triVerts[i]];
292  const point& pt1 = mesh_.points()[triVerts[j]];
293 
294  if (edge(triVerts[i], triVerts[j]) == currentE)
295  {
296  vector edgeNormal = (pt1-pt0)^n;
297  return (dir&edgeNormal) < 0;
298  }
299  }
300 
302  (
303  "wallBoundedParticle::isTriAlongTrack"
304  "(const point&)"
305  ) << "Problem" << abort(FatalError);
306 
307  return false;
308 }
309 
310 
311 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
312 
314 (
315  const polyMesh& mesh,
316  const vector& position,
317  const label cellI,
318  const label tetFaceI,
319  const label tetPtI,
320  const label meshEdgeStart,
321  const label diagEdge
322 )
323 :
324  particle(mesh, position, cellI, tetFaceI, tetPtI),
325  meshEdgeStart_(meshEdgeStart),
326  diagEdge_(diagEdge)
327 {}
328 
329 
331 (
332  const polyMesh& mesh,
333  Istream& is,
334  bool readFields
335 )
336 :
337  particle(mesh, is, readFields)
338 {
339  if (readFields)
340  {
341  if (is.format() == IOstream::ASCII)
342  {
343  is >> meshEdgeStart_ >> diagEdge_;
344  }
345  else
346  {
347  is.read
348  (
349  reinterpret_cast<char*>(&meshEdgeStart_),
350  sizeof(meshEdgeStart_)
351  + sizeof(diagEdge_)
352  );
353  }
354  }
355 
356  // Check state of Istream
357  is.check
358  (
359  "wallBoundedParticle::wallBoundedParticle"
360  "(const Cloud<wallBoundedParticle>&, Istream&, bool)"
361  );
362 }
363 
364 
366 (
367  const wallBoundedParticle& p
368 )
369 :
370  particle(p),
373 {}
374 
375 
376 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
377 
378 Foam::Ostream& Foam::operator<<
379 (
380  Ostream& os,
381  const wallBoundedParticle& p
382 )
383 {
384  if (os.format() == IOstream::ASCII)
385  {
386  os << static_cast<const particle&>(p)
388  << token::SPACE << p.diagEdge_;
389  }
390  else
391  {
392  os << static_cast<const particle&>(p);
393  os.write
394  (
395  reinterpret_cast<const char*>(&p.meshEdgeStart_),
396  wallBoundedParticle::sizeofFields_
397  );
398  }
399 
400  return os;
401 }
402 
403 
404 Foam::Ostream& Foam::operator<<
405 (
406  Ostream& os,
408 )
409 {
410  const wallBoundedParticle& p = ip.t_;
411 
412  tetPointRef tpr(p.currentTet());
413 
414  os << " " << static_cast<const particle&>(p) << nl
415  << " tet:" << nl;
416  os << " ";
417  meshTools::writeOBJ(os, tpr.a());
418  os << " ";
419  meshTools::writeOBJ(os, tpr.b());
420  os << " ";
421  meshTools::writeOBJ(os, tpr.c());
422  os << " ";
423  meshTools::writeOBJ(os, tpr.d());
424  os << " l 1 2" << nl
425  << " l 1 3" << nl
426  << " l 1 4" << nl
427  << " l 2 3" << nl
428  << " l 2 4" << nl
429  << " l 3 4" << nl;
430  os << " ";
431  meshTools::writeOBJ(os, p.position());
432 
433  return os;
434 }
435 
436 
437 // ************************************************************************* //
label & cell()
Return current cell particle is in.
Definition: particleI.H:621
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject( name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE ))
label & face()
Return current face particle is on otherwise -1.
Definition: particleI.H:681
dimensioned< scalar > mag(const dimensioned< Type > &)
labelList f(nPoints)
const labelList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:864
vector normal(const pointField &) const
Vector normal; magnitude is equal to area of face.
Definition: triFaceI.H:181
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
virtual Istream & read(token &)=0
Return next token from stream.
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
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:117
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
void crossDiagonalEdge()
Cross diagonal edge into different triangle on same face,cell.
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
label n
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
static const char nl
Definition: Ostream.H:260
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1035
volScalarField & p
Definition: createFields.H:51
A helper class for outputting values to Ostream.
Definition: InfoProxy.H:45
const vector & position() const
Return current particle position.
Definition: particleI.H:603
tetPointRef currentTet() const
Return the geometry of the current tet that the.
Definition: particleI.H:657
scalar trackFaceTri(const vector &endPosition, label &minEdgeI)
Track through single triangle.
#define forAll(list, i)
Definition: UList.H:421
streamFormat format() const
Return current stream format.
Definition: IOstream.H:377
label diagEdge_
Particle is on diagonal edge:
Base particle class.
Definition: particle.H:78
const polyMesh & mesh_
Reference to the polyMesh database.
Definition: particle.H:137
A tetrahedron primitive.
Definition: tetrahedron.H:62
label & tetFace()
Return current tet face particle is in.
Definition: particleI.H:633
errorManip< error > abort(error &err)
Definition: errorManip.H:131
InfoProxy< wallBoundedParticle > info() const
Return info proxy.
label & tetPt()
Return current tet face particle is in.
Definition: particleI.H:645
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
edge currentEdge() const
Construct current edge.
void crossEdgeConnectedFace(const edge &meshEdge)
Check if inside current tet.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
error FatalError
bool isTriAlongTrack(const point &endPosition) const
Is current triangle in the track direction.
particle(const polyMesh &mesh, const vector &position, const label cellI, const label tetFaceI, const label tetPtI)
Construct from components.
Definition: particle.C:48
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: FixedListI.H:115
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
void crossEdgeConnectedFace(const label &cellI, label &tetFaceI, label &tetPtI, const edge &e)
Cross the from the given face across the given edge of the.
Definition: particleI.H:465
label rcIndex(const label i) const
Return the reverse circular index, i.e. the previous index.
Definition: UListI.H:65
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.
Particle class that tracks on triangles of boundary faces. Use trackToEdge similar to trackToFace on ...
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1060
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:209
tetIndices currentTetIndices() const
Return the indices of the current tet that the.
Definition: particleI.H:651
label meshEdgeStart_
Particle is on mesh edge: