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 
31 (
32  sizeof(wallBoundedParticle) - sizeof(particle)
33 );
34 
35 
36 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37 
39 {
40  if ((meshEdgeStart_ != -1) == (diagEdge_ != -1))
41  {
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
73  particle::crossEdgeConnectedFace(cell(), tetFace(), tetPt(), meshEdge);
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  << "Problem :"
98  << " particle:"
99  << info()
100  << "face:" << tetFace()
101  << " verts:" << f
102  << " meshEdge:" << meshEdge
103  << abort(FatalError);
104  }
105  }
106 
107  diagEdge_ = -1;
108 
109  // Check that still on same mesh edge
110  const edge eNew(f[meshEdgeStart_], f.nextLabel(meshEdgeStart_));
111  if (eNew != meshEdge)
112  {
114  << "Problem" << abort(FatalError);
115  }
116 }
117 
118 
120 {
121  if (diagEdge_ == -1)
122  {
124  << "Particle:"
125  << info()
126  << "not on a diagonal edge" << abort(FatalError);
127  }
128  if (meshEdgeStart_ != -1)
129  {
131  << "Particle:"
132  << info()
133  << "meshEdgeStart_:" << meshEdgeStart_ << abort(FatalError);
134  }
135 
136  const Foam::face& f = mesh_.faces()[tetFace()];
137 
138  // tetPtI starts from 1, goes up to f.size()-2
139 
140  if (tetPt() == diagEdge_)
141  {
142  tetPt() = f.rcIndex(tetPt());
143  }
144  else
145  {
146  label nextTetPt = f.fcIndex(tetPt());
147  if (diagEdge_ == nextTetPt)
148  {
149  tetPt() = nextTetPt;
150  }
151  else
152  {
154  << "Particle:"
155  << info()
156  << "tetPt:" << tetPt()
157  << " diagEdge:" << diagEdge_ << abort(FatalError);
158  }
159  }
160 
161  meshEdgeStart_ = -1;
162 }
163 
164 
166 (
167  const vector& endPosition,
168  label& minEdgeI
169 )
170 {
171  // Track p from position to endPosition
172  const triFace tri(currentTetIndices().faceTriIs(mesh_));
173  vector n = tri.normal(mesh_.points());
174  n /= mag(n)+VSMALL;
175 
176  // Check which edge intersects the trajectory.
177  // Project trajectory onto triangle
178  minEdgeI = -1;
179  scalar minS = 1; // end position
180 
181  edge currentE(-1, -1);
182  if (meshEdgeStart_ != -1 || diagEdge_ != -1)
183  {
184  currentE = currentEdge();
185  }
186 
187  // Determine path along line position+s*d to see where intersections are.
188  forAll(tri, i)
189  {
190  label j = tri.fcIndex(i);
191 
192  const point& pt0 = mesh_.points()[tri[i]];
193  const point& pt1 = mesh_.points()[tri[j]];
194 
195  if (edge(tri[i], tri[j]) == currentE)
196  {
197  // Do not check particle is on
198  continue;
199  }
200 
201  // Outwards pointing normal
202  vector edgeNormal = (pt1-pt0)^n;
203 
204  edgeNormal /= mag(edgeNormal)+VSMALL;
205 
206  // Determine whether position and end point on either side of edge.
207  scalar sEnd = (endPosition-pt0)&edgeNormal;
208  if (sEnd >= 0)
209  {
210  // endPos is outside triangle. position() should always be
211  // inside.
212  scalar sStart = (position()-pt0)&edgeNormal;
213  if (mag(sEnd-sStart) > VSMALL)
214  {
215  scalar s = sStart/(sStart-sEnd);
216 
217  if (s >= 0 && s < minS)
218  {
219  minS = s;
220  minEdgeI = i;
221  }
222  }
223  }
224  }
225 
226  if (minEdgeI != -1)
227  {
228  position() += minS*(endPosition-position());
229  }
230  else
231  {
232  // Did not hit any edge so tracked to the end position. Set position
233  // without any calculation to avoid truncation errors.
234  position() = endPosition;
235  minS = 1.0;
236  }
237 
238  // Project position() back onto plane of triangle
239  const point& triPt = mesh_.points()[tri[0]];
240  position() -= ((position()-triPt)&n)*n;
241 
242  return minS;
243 }
244 
245 
247 (
248  const point& endPosition
249 ) const
250 {
251  const triFace triVerts(currentTetIndices().faceTriIs(mesh_));
252  const edge currentE = currentEdge();
253 
254  if
255  (
256  currentE[0] == currentE[1]
257  || findIndex(triVerts, currentE[0]) == -1
258  || findIndex(triVerts, currentE[1]) == -1
259  )
260  {
262  << "Edge " << currentE << " not on triangle " << triVerts
263  << info()
264  << abort(FatalError);
265  }
266 
267 
268  const vector dir = endPosition-position();
269 
270  // Get normal of currentE
271  vector n = triVerts.normal(mesh_.points());
272  n /= mag(n);
273 
274  forAll(triVerts, i)
275  {
276  label j = triVerts.fcIndex(i);
277  const point& pt0 = mesh_.points()[triVerts[i]];
278  const point& pt1 = mesh_.points()[triVerts[j]];
279 
280  if (edge(triVerts[i], triVerts[j]) == currentE)
281  {
282  vector edgeNormal = (pt1-pt0)^n;
283  return (dir&edgeNormal) < 0;
284  }
285  }
286 
288  << "Problem" << abort(FatalError);
289 
290  return false;
291 }
292 
293 
294 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
295 
297 (
298  const polyMesh& mesh,
299  const vector& position,
300  const label cellI,
301  const label tetFaceI,
302  const label tetPtI,
303  const label meshEdgeStart,
304  const label diagEdge
305 )
306 :
307  particle(mesh, position, cellI, tetFaceI, tetPtI),
308  meshEdgeStart_(meshEdgeStart),
309  diagEdge_(diagEdge)
310 {}
311 
312 
314 (
315  const polyMesh& mesh,
316  Istream& is,
317  bool readFields
318 )
319 :
320  particle(mesh, is, readFields)
321 {
322  if (readFields)
323  {
324  if (is.format() == IOstream::ASCII)
325  {
326  is >> meshEdgeStart_ >> diagEdge_;
327  }
328  else
329  {
330  is.read(reinterpret_cast<char*>(&meshEdgeStart_), sizeofFields_);
331  }
332  }
333 
334  // Check state of Istream
335  is.check
336  (
337  "wallBoundedParticle::wallBoundedParticle"
338  "(const Cloud<wallBoundedParticle>&, Istream&, bool)"
339  );
340 }
341 
342 
344 (
345  const wallBoundedParticle& p
346 )
347 :
348  particle(p),
349  meshEdgeStart_(p.meshEdgeStart_),
350  diagEdge_(p.diagEdge_)
351 {}
352 
353 
354 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
355 
356 Foam::Ostream& Foam::operator<<
357 (
358  Ostream& os,
359  const wallBoundedParticle& p
360 )
361 {
362  if (os.format() == IOstream::ASCII)
363  {
364  os << static_cast<const particle&>(p)
365  << token::SPACE << p.meshEdgeStart_
366  << token::SPACE << p.diagEdge_;
367  }
368  else
369  {
370  os << static_cast<const particle&>(p);
371  os.write
372  (
373  reinterpret_cast<const char*>(&p.meshEdgeStart_),
375  );
376  }
377 
378  return os;
379 }
380 
381 
382 Foam::Ostream& Foam::operator<<
383 (
384  Ostream& os,
385  const InfoProxy<wallBoundedParticle>& ip
386 )
387 {
388  const wallBoundedParticle& p = ip.t_;
389 
390  tetPointRef tpr(p.currentTet());
391 
392  os << " " << static_cast<const particle&>(p) << nl
393  << " tet:" << nl;
394  os << " ";
395  meshTools::writeOBJ(os, tpr.a());
396  os << " ";
397  meshTools::writeOBJ(os, tpr.b());
398  os << " ";
399  meshTools::writeOBJ(os, tpr.c());
400  os << " ";
401  meshTools::writeOBJ(os, tpr.d());
402  os << " l 1 2" << nl
403  << " l 1 3" << nl
404  << " l 1 4" << nl
405  << " l 2 3" << nl
406  << " l 2 4" << nl
407  << " l 3 4" << nl;
408  os << " ";
409  meshTools::writeOBJ(os, p.position());
410 
411  return os;
412 }
413 
414 
415 // ************************************************************************* //
Foam::triFace::normal
vector normal(const pointField &) const
Vector normal; magnitude is equal to area of face.
Definition: triFaceI.H:181
Foam::IOstream::format
streamFormat format() const
Return current stream format.
Definition: IOstream.H:377
Foam::wallBoundedParticle::sizeofFields_
static const std::size_t sizeofFields_
Size in bytes of the fields.
Definition: wallBoundedParticle.H:59
p
p
Definition: pEqn.H:62
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::particle::crossEdgeConnectedFace
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:458
Foam::findIndex
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
Foam::polyMesh::tetBasePtIs
const labelList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:818
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::wallBoundedParticle::wallBoundedParticle
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.
Definition: wallBoundedParticle.C:297
Foam::wallBoundedParticle::meshEdgeStart_
label meshEdgeStart_
Particle is on mesh edge:
Definition: wallBoundedParticle.H:102
Foam::wallBoundedParticle::crossDiagonalEdge
void crossDiagonalEdge()
Cross diagonal edge into different triangle on same face,cell.
Definition: wallBoundedParticle.C:119
Foam::wallBoundedParticle::info
InfoProxy< wallBoundedParticle > info() const
Return info proxy.
Definition: wallBoundedParticle.H:299
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::readFields
This function object reads fields from the time directories and adds them to the mesh database for fu...
Definition: readFields.H:104
Foam::wallBoundedParticle::trackFaceTri
scalar trackFaceTri(const vector &endPosition, label &minEdgeI)
Track through single triangle.
Definition: wallBoundedParticle.C:166
Foam::FixedList::fcIndex
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: FixedListI.H:115
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::wallBoundedParticle::crossEdgeConnectedFace
void crossEdgeConnectedFace(const edge &meshEdge)
Cross mesh edge into different face on same cell.
Definition: wallBoundedParticle.C:68
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::IOstream::ASCII
@ ASCII
Definition: IOstream.H:88
wallBoundedParticle.H
Foam::particle::tetFace
label & tetFace()
Return current tet face particle is in.
Definition: particleI.H:610
Foam::label
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
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::wallBoundedParticle
Particle class that tracks on triangles of boundary faces. Use trackToEdge similar to trackToFace on ...
Definition: wallBoundedParticle.H:52
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::tetPointRef
tetrahedron< point, const point & > tetPointRef
Definition: tetrahedron.H:78
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
s
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))
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::particle
Base particle class.
Definition: particle.H:78
Foam::wallBoundedParticle::currentEdge
edge currentEdge() const
Construct current edge.
Definition: wallBoundedParticle.C:38
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Foam::wallBoundedParticle::isTriAlongTrack
bool isTriAlongTrack(const point &endPosition) const
Is current triangle in the track direction.
Definition: wallBoundedParticle.C:247
Foam::particle::mesh_
const polyMesh & mesh_
Reference to the polyMesh database.
Definition: particle.H:137
Foam::wallBoundedParticle::diagEdge_
label diagEdge_
Particle is on diagonal edge:
Definition: wallBoundedParticle.H:109
Foam::token::SPACE
@ SPACE
Definition: token.H:95
Foam::Istream::read
virtual Istream & read(token &)=0
Return next token from stream.