wallBoundedStreamLine.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "wallBoundedStreamLine.H"
31 #include "sampledSet.H"
32 #include "faceSet.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace functionObjects
40 {
41  defineTypeNameAndDebug(wallBoundedStreamLine, 0);
43  (
44  functionObject,
45  wallBoundedStreamLine,
46  dictionary
47  );
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
56 (
57  const bitSet& isWallPatch,
58  const point& seedPt,
59  const label celli
60 ) const
61 {
62  const cell& cFaces = mesh_.cells()[celli];
63 
64  label minFacei = -1;
65  label minTetPti = -1;
66  scalar minDistSqr = sqr(GREAT);
67  point nearestPt(GREAT, GREAT, GREAT);
68 
69  for (label facei : cFaces)
70  {
71  if (isWallPatch[facei])
72  {
73  const face& f = mesh_.faces()[facei];
74  const label fp0 = mesh_.tetBasePtIs()[facei];
75  const point& basePoint = mesh_.points()[f[fp0]];
76 
77  label fp = f.fcIndex(fp0);
78  for (label i = 2; i < f.size(); i++)
79  {
80  const point& thisPoint = mesh_.points()[f[fp]];
81  const label nextFp = f.fcIndex(fp);
82  const point& nextPoint = mesh_.points()[f[nextFp]];
83 
84  const triPointRef tri(basePoint, thisPoint, nextPoint);
85 
86  //const scalar d2 = magSqr(tri.centre() - seedPt);
87  const pointHit nearInfo(tri.nearestPoint(seedPt));
88  const scalar d2 = nearInfo.distance();
89  if (d2 < minDistSqr)
90  {
91  nearestPt = nearInfo.point();
92  minDistSqr = d2;
93  minFacei = facei;
94  minTetPti = i-1;
95  }
96  fp = nextFp;
97  }
98  }
99  }
100 
101  // Return tet and nearest point on wall triangle
103  (
104  tetIndices
105  (
106  celli,
107  minFacei,
108  minTetPti
109  ),
110  nearestPt
111  );
112 }
113 
114 
116 (
117  const triPointRef& tri,
118  const point& pt
119 ) const
120 {
121  //pointHit nearPt(tri.nearestPoint(pt));
122  //if (nearPt.distance() > ROOTSMALL)
123  //{
124  // FatalErrorInFunction << "tri:" << tri
125  // << " seed:" << pt << exit(FatalError);
126  //}
127 
128  return (1.0 - ROOTSMALL)*pt + ROOTSMALL*tri.centre();
129 }
130 
131 
133 {
134  // Determine the 'wall' patches
135  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
136  // These are the faces that need to be followed
137 
139  bitSet isWallPatch(mesh_.nFaces(), boundaryPatch().addressing());
140 
141 
142  // Find nearest wall particle for the seedPoints
143  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
144 
146  (
147  mesh_,
148  Foam::zero{},
149  cloudName_
150  );
151 
152  {
153  // Get the seed points
154  // ~~~~~~~~~~~~~~~~~~~
155 
156  const sampledSet& seedPoints = sampledSetPoints();
157 
158  forAll(seedPoints, seedi)
159  {
160  const label celli = seedPoints.cells()[seedi];
161 
162  if (celli != -1)
163  {
164  const point& seedPt = seedPoints[seedi];
165  Tuple2<tetIndices, point> nearestId
166  (
167  findNearestTet(isWallPatch, seedPt, celli)
168  );
169  const tetIndices& ids = nearestId.first();
170 
171  if (ids.face() != -1 && isWallPatch[ids.face()])
172  {
173  //Pout<< "Seeding particle :" << nl
174  // << " seedPt:" << seedPt << nl
175  // << " face :" << ids.face() << nl
176  // << " at :" << mesh_.faceCentres()[ids.face()]
177  // << nl
178  // << " cell :" << mesh_.cellCentres()[ids.cell()]
179  // << nl << endl;
180 
181  particles.addParticle
182  (
184  (
185  mesh_,
186  // Perturb seed point to be inside triangle
187  pushIn(ids.faceTri(mesh_), nearestId.second()),
188  ids.cell(),
189  ids.face(), // tetFace
190  ids.tetPt(),
191  -1, // not on a mesh edge
192  -1, // not on a diagonal edge
193  (trackDir_ == trackDirType::FORWARD), // forward?
194  lifeTime_ // lifetime
195  )
196  );
197 
198  if (trackDir_ == trackDirType::BIDIRECTIONAL)
199  {
200  // Additional particle for other half of track
201  particles.addParticle
202  (
204  (
205  mesh_,
206  // Perturb seed point to be inside triangle
207  pushIn(ids.faceTri(mesh_), nearestId.second()),
208  ids.cell(),
209  ids.face(), // tetFace
210  ids.tetPt(),
211  -1, // not on a mesh edge
212  -1, // not on a diagonal edge
213  true, // forward
214  lifeTime_ // lifetime
215  )
216  );
217  }
218  }
219  else
220  {
221  Pout<< type() << " : ignoring seed " << seedPt
222  << " since not in wall cell." << endl;
223  }
224  }
225  }
226  }
227 
228  const label nSeeds = returnReduce(particles.size(), sumOp<label>());
229 
230  Log << type() << " : seeded " << nSeeds << " particles." << endl;
231 
232 
233  // Field interpolators
234  PtrList<interpolation<scalar>> vsInterp;
235  PtrList<interpolation<vector>> vvInterp;
236 
237  refPtr<interpolation<vector>> UInterp
238  (
239  initInterpolations(nSeeds, vsInterp, vvInterp)
240  );
241 
242  // Additional particle info
243  wallBoundedStreamLineParticle::trackingData td
244  (
245  particles,
246  vsInterp,
247  vvInterp,
248  UInterp.cref(), // velocity interpolator (possibly within vvInterp)
249  trackLength_, // fixed track length
250  isWallPatch, // which faces are to follow
251 
252  allTracks_,
253  allScalars_,
254  allVectors_
255  );
256 
257 
258  // Set very large dt. Note: cannot use GREAT since 1/GREAT is SMALL
259  // which is a trigger value for the tracking...
260  const scalar trackTime = Foam::sqrt(GREAT);
261 
262  // Track
263  particles.move(particles, td, trackTime);
264 }
265 
266 
267 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
268 
270 (
271  const word& name,
272  const Time& runTime,
273  const dictionary& dict
274 )
275 :
277 {
278  read(dict_);
279 }
280 
281 
283 (
284  const word& name,
285  const Time& runTime,
286  const dictionary& dict,
287  const wordList& fieldNames
288 )
289 :
290  streamLineBase(name, runTime, dict, fieldNames)
291 {
292  read(dict_);
293 }
294 
295 
296 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
297 
299 {
301  {
302  // Make sure that the mesh is trackable
303  if (debug)
304  {
305  // 1. Positive volume decomposition tets
306  faceSet faces(mesh_, "lowQualityTetFaces", mesh_.nFaces()/100+1);
307 
308  if
309  (
311  (
312  mesh_,
314  true,
315  &faces
316  )
317  )
318  {
319  label nFaces = returnReduce(faces.size(), sumOp<label>());
320 
322  << "Found " << nFaces
323  <<" faces with low quality or negative volume "
324  << "decomposition tets. Writing to faceSet " << faces.name()
325  << endl;
326  }
327 
328  // 2. All edges on a cell having two faces
329  EdgeMap<label> numFacesPerEdge;
330  for (const cell& cFaces : mesh_.cells())
331  {
332  numFacesPerEdge.clear();
333 
334  for (const label facei : cFaces)
335  {
336  const face& f = mesh_.faces()[facei];
337  forAll(f, fp)
338  {
339  const edge e(f[fp], f.nextLabel(fp));
340 
341  ++(numFacesPerEdge(e, 0));
342  }
343  }
344 
345  forAllConstIters(numFacesPerEdge, iter)
346  {
347  if (iter() != 2)
348  {
350  << "problem cell:" << cFaces
351  << abort(FatalError);
352  }
353  }
354  }
355  }
356  }
357 
358  return true;
359 }
360 
361 
362 // ************************************************************************* //
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
A triangle primitive used to calculate face normals and swept volumes. Uses referred points...
Definition: triangle.H:217
wallBoundedStreamLine(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
A list of face labels.
Definition: faceSet.H:47
Like polyPatch but without reference to mesh. Used in boundaryMesh to hold data on patches...
Definition: boundaryPatch.H:51
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: stringOps.H:54
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:899
triPointRef faceTri(const polyMesh &mesh) const
The triangle geometry for the face for this tet. The normal of the tri points out of the cell...
Definition: tetIndicesI.H:182
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual void track()
Do the actual tracking to fill the track data.
engineTime & runTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual bool read(const dictionary &)
Read settings.
const cellList & cells() const
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
A Cloud of wall-bounded streamLine particles.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static bool checkFaceTets(const polyMesh &mesh, scalar tol=minTetQuality, const bool report=false, labelHashSet *setPtr=nullptr)
Check face-decomposition tet volume.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
pointHit nearestPoint(const point &p) const
Return nearest point to p on triangle.
Definition: triangleI.H:861
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
static Point centre(const Point &p0, const Point &p1, const Point &p2)
The centre (centroid) of three points.
Definition: triangleI.H:115
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
point pushIn(const triPointRef &tri, const point &pt) const
Push a point a tiny bit towards the centre of the triangle it is in to avoid tracking problems...
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:79
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition: edge.H:59
virtual bool read(const dictionary &)
Read the field average data.
label cell() const noexcept
Return the cell index.
Definition: tetIndices.H:146
A class for handling words, derived from Foam::string.
Definition: word.H:63
void clear()
Remove all entries from table.
Definition: HashTable.C:725
label tetPt() const noexcept
Return the characterising tet point index.
Definition: tetIndices.H:166
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:78
Describes the interaction of a object and a (templated) point. It carries the info of a successful hi...
Definition: pointHit.H:43
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Tuple2< tetIndices, point > findNearestTet(const bitSet &isWallPatch, const point &seedPt, const label celli) const
Find wall tet on cell.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
labelList f(nPoints)
label face() const noexcept
Return the face index.
Definition: tetIndices.H:156
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
#define WarningInFunction
Report a warning using Foam::Warning.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:53
scalar distance() const noexcept
Return distance to hit.
Definition: pointHit.H:169
#define Log
Definition: PDRblock.C:28
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
const labelList & cells() const noexcept
Definition: sampledSet.H:388
Particle class that samples fields as it passes through. Used in streamline calculation.
const fvMesh & mesh_
Reference to the fvMesh.
static const scalar minTetQuality
Minimum tetrahedron quality.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28