Loading...
Searching...
No Matches
sampledSet.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-2017 OpenFOAM Foundation
9 Copyright (C) 2018-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "sampledSet.H"
30#include "polyMesh.H"
31#include "primitiveMesh.H"
32#include "meshSearch.H"
33#include "particle.H"
34#include "globalIndex.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
42}
43
44
45// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46
48{
49 if
50 (
51 (cells_.size() != size())
52 || (faces_.size() != size())
53 || (segments_.size() != size())
54 || (distance_.size() != size())
55 )
56 {
58 << "Sizes not equal : "
59 << " points:" << size()
60 << " cells:" << cells_.size()
61 << " faces:" << faces_.size()
62 << " segments:" << segments_.size()
63 << " distance:" << distance_.size()
64 << abort(FatalError);
65 }
66}
67
69Foam::label Foam::sampledSet::getBoundaryCell(const label facei) const
70{
71 return mesh().faceOwner()[facei];
72}
73
74
75Foam::label Foam::sampledSet::getNeighbourCell(const label facei) const
76{
77 if (facei < mesh().nInternalFaces())
78 {
79 return mesh().faceNeighbour()[facei];
80 }
81
82 return mesh().faceOwner()[facei];
83}
84
85
87(
88 const point& p,
89 const label samplei
90) const
91{
92 // Collect the face owner and neighbour cells of the sample into an array
93 // for convenience
94 const label cells[4] =
95 {
96 mesh().faceOwner()[faces_[samplei]],
97 getNeighbourCell(faces_[samplei]),
98 mesh().faceOwner()[faces_[samplei+1]],
99 getNeighbourCell(faces_[samplei+1])
100 };
101
102 // Find the sampled cell by checking the owners and neighbours of the
103 // sampled faces
104 label cellm =
105 (cells[0] == cells[2] || cells[0] == cells[3]) ? cells[0]
106 : (cells[1] == cells[2] || cells[1] == cells[3]) ? cells[1]
107 : -1;
108
109 if (cellm != -1)
110 {
111 // If found the sampled cell check the point is in the cell
112 // otherwise ignore
113 if (!mesh().pointInCell(p, cellm, searchEngine_.decompMode()))
114 {
115 cellm = -1;
116
117 if (debug)
118 {
120 << "Could not find mid-point " << p
121 << " cell " << cellm << endl;
122 }
123 }
124 }
125 else
126 {
127 // If the sample does not pass through a single cell check if the point
128 // is in any of the owners or neighbours otherwise ignore
129 for (label i=0; i<4; ++i)
130 {
131 if (mesh().pointInCell(p, cells[i], searchEngine_.decompMode()))
132 {
133 return cells[i];
134 }
135 }
136
137 if (debug)
138 {
140 << "Could not find cell for mid-point" << nl
141 << " samplei: " << samplei
142 << " pts[samplei]: " << operator[](samplei)
143 << " face[samplei]: " << faces_[samplei]
144 << " pts[samplei+1]: " << operator[](samplei+1)
145 << " face[samplei+1]: " << faces_[samplei+1]
146 << " cellio: " << cells[0]
147 << " cellin: " << cells[1]
148 << " celljo: " << cells[2]
149 << " celljn: " << cells[3]
150 << endl;
152 }
153
154 return cellm;
155}
156
157
159(
160 const label facei,
161 const point& sample
162) const
163{
164 vector vec = sample - mesh().faceCentres()[facei];
165
166 scalar magVec = mag(vec);
167
168 if (magVec < VSMALL)
169 {
170 // Sample on face centre. Regard as inside
171 return -1;
172 }
173
174 vec /= magVec;
176 const vector n = normalised(mesh().faceAreas()[facei]);
177
178 return n & vec;
179}
180
181
183(
184 const label celli,
185 const point& sample,
186 const scalar smallDist
187) const
188{
189 const cell& myFaces = mesh().cells()[celli];
190
191 forAll(myFaces, myFacei)
192 {
193 const face& f = mesh().faces()[myFaces[myFacei]];
194
195 pointHit inter = f.nearestPoint(sample, mesh().points());
196
197 scalar dist = inter.point().dist(sample);
198
199 if (dist < smallDist)
200 {
201 return myFaces[myFacei];
202 }
203 }
204 return -1;
205}
206
207
209(
210 const point& facePt,
211 const label facei
212) const
213{
214 label celli = mesh().faceOwner()[facei];
215 const point& cC = mesh().cellCentres()[celli];
216
217 point newPosition = facePt;
218
219 // Taken from particle::initCellFacePt()
220 label tetFacei;
221 label tetPtI;
222 mesh().findTetFacePt(celli, facePt, tetFacei, tetPtI);
223
224 // This is the tolerance that was defined as a static constant of the
225 // particle class. It is no longer used by particle, following the switch to
226 // barycentric tracking. The only place in which the tolerance is now used
227 // is here. I'm not sure what the purpose of this code is, but it probably
228 // wants removing. It is doing tet-searches for a particle position. This
229 // should almost certainly be left to the particle class.
230 const scalar trackingCorrectionTol = 1e-5;
231
232 if (tetFacei == -1 || tetPtI == -1)
233 {
234 newPosition = facePt;
235
236 label trap(1.0/trackingCorrectionTol + 1);
237
238 label iterNo = 0;
239
240 do
241 {
242 newPosition += trackingCorrectionTol*(cC - facePt);
243
245 (
246 celli,
247 newPosition,
248 tetFacei,
249 tetPtI
250 );
251
252 ++iterNo;
253
254 } while (tetFacei < 0 && iterNo <= trap);
255 }
256
257 if (tetFacei == -1)
258 {
260 << "After pushing " << facePt << " to " << newPosition
261 << " it is still outside face " << facei
262 << " at " << mesh().faceCentres()[facei]
263 << " of cell " << celli
264 << " at " << cC << endl
265 << "Please change your starting point"
267 }
268
269 return newPosition;
270}
271
272
274(
275 const point& samplePt,
276 const point& bPoint,
277 const label bFacei,
278 const scalar smallDist,
279
280 point& trackPt,
281 label& trackCelli,
282 label& trackFacei
283) const
284{
285 bool isGoodSample = false;
286
287 if (bFacei == -1)
288 {
289 // No boundary intersection. Try and find cell samplePt is in
290 trackCelli = mesh().findCell(samplePt, searchEngine_.decompMode());
291
292 if
293 (
294 (trackCelli == -1)
295 || !mesh().pointInCell
296 (
297 samplePt,
298 trackCelli,
299 searchEngine_.decompMode()
300 )
301 )
302 {
303 // Line samplePt - end_ does not intersect domain at all.
304 // (or is along edge)
305
306 trackCelli = -1;
307 trackFacei = -1;
308
309 isGoodSample = false;
310 }
311 else
312 {
313 // Start is inside. Use it as tracking point
314
315 trackPt = samplePt;
316 trackFacei = -1;
317
318 isGoodSample = true;
319 }
320 }
321 else if (mag(samplePt - bPoint) < smallDist)
322 {
323 // samplePt close to bPoint. Snap to it
324 trackPt = pushIn(bPoint, bFacei);
325 trackFacei = bFacei;
326 trackCelli = getBoundaryCell(trackFacei);
327
328 isGoodSample = true;
329 }
330 else
331 {
332 scalar sign = calcSign(bFacei, samplePt);
333
334 if (sign < 0)
335 {
336 // samplePt inside or marginally outside.
337 trackPt = samplePt;
338 trackFacei = -1;
339 trackCelli = mesh().findCell(trackPt, searchEngine_.decompMode());
340
341 isGoodSample = true;
342 }
343 else
344 {
345 // samplePt outside. use bPoint
346 trackPt = pushIn(bPoint, bFacei);
347 trackFacei = bFacei;
348 trackCelli = getBoundaryCell(trackFacei);
349
350 isGoodSample = false;
351 }
352 }
353
355 << " samplePt:" << samplePt
356 << " bPoint:" << bPoint
357 << " bFacei:" << bFacei << nl
358 << " Calculated first tracking point :"
359 << " trackPt:" << trackPt
360 << " trackCelli:" << trackCelli
361 << " trackFacei:" << trackFacei
362 << " isGoodSample:" << isGoodSample
363 << endl;
364
365 return isGoodSample;
366}
367
368
370(
371 const List<point>& samplingPts,
372 const labelList& samplingCells,
373 const labelList& samplingFaces,
374 const labelList& samplingSegments,
375 const scalarList& samplingDistance
376)
377{
378 setPoints(samplingPts);
379 setDistance(samplingDistance, false); // check=false
380
381 segments_ = samplingSegments;
382 cells_ = samplingCells;
383 faces_ = samplingFaces;
384
386}
387
388
390(
391 List<point>&& samplingPts,
392 labelList&& samplingCells,
393 labelList&& samplingFaces,
394 labelList&& samplingSegments,
395 scalarList&& samplingDistance
396)
397{
398 setPoints(std::move(samplingPts));
399 setDistance(std::move(samplingDistance), false); // check=false
400
401 segments_ = std::move(samplingSegments);
402 cells_ = std::move(samplingCells);
403 faces_ = std::move(samplingFaces);
406}
407
408
409// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
410
412(
413 const word& name,
414 const polyMesh& mesh,
415 const meshSearch& searchEngine,
416 const coordSet::coordFormat axisType
417)
418:
419 coordSet(name, axisType),
420 mesh_(mesh),
421 searchEngine_(searchEngine),
422 segments_(),
423 cells_(),
424 faces_()
425{}
426
427
429(
430 const word& name,
431 const polyMesh& mesh,
432 const meshSearch& searchEngine,
433 const word& axis
434)
435:
436 coordSet(name, axis),
437 mesh_(mesh),
438 searchEngine_(searchEngine),
439 segments_(),
440 cells_(),
441 faces_()
442{}
443
444
446(
447 const word& name,
448 const polyMesh& mesh,
449 const meshSearch& searchEngine,
450 const dictionary& dict
451)
452:
453 coordSet(name, dict.get<word>("axis")),
454 mesh_(mesh),
455 searchEngine_(searchEngine),
456 segments_(),
457 cells_(),
458 faces_()
459{}
460
461
462// Foam::autoPtr<Foam::coordSet> Foam::sampledSet::gather
463// (
464// labelList& sortOrder,
465// labelList& allSegments
466// ) const
467// {
468// autoPtr<coordSet> result(coordSet::gatherSort(sortOrder));
469//
470// // Optional
471// if (notNull(allSegments))
472// {
473// globalIndex::gatherOp(segments(), allSegments);
474// allSegments = UIndirectList<label>(allSegments, sortOrder)();
475// }
477// return result;
478// }
479
480
481// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
482
484(
485 const word& name,
486 const polyMesh& mesh,
488 const dictionary& dict
489)
490{
491 const word sampleType(dict.get<word>("type"));
492
493 auto* ctorPtr = wordConstructorTable(sampleType);
494
495 if (!ctorPtr)
496 {
498 (
499 dict,
500 "sample",
501 sampleType,
502 *wordConstructorTablePtr_
503 ) << exit(FatalIOError);
504 }
505
507 (
508 ctorPtr
509 (
510 name,
511 mesh,
513 dict.optionalSubDict(sampleType + "Coeffs")
514 )
515 );
516}
517
518
520{
522
523 os << nl << "\t(celli)\t(facei)" << nl;
524
525 forAll(*this, samplei)
526 {
527 os << '\t' << cells_[samplei]
528 << '\t' << faces_[samplei]
529 << nl;
530 }
531
532 return os;
533}
534
535
536// ************************************************************************* //
label n
Minimal example by using system/controlDict.functions:
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
const point_type & point() const noexcept
Return the point, no checks.
Definition pointHit.H:161
bool get(const label i) const
Definition UList.H:868
label size() const noexcept
Definition UList.H:706
void size(const label n)
Definition UList.H:118
vector & operator[](const label i)
scalar dist(const Vector< Cmpt > &v2) const
The L2-norm distance from another vector. The mag() of the difference.
Definition VectorI.H:107
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
Holds list of sampling positions.
Definition coordSet.H:52
scalarList distance_
Cumulative distance for "distance" write specifier.
Definition coordSet.H:89
const word & axis() const
The sort axis name.
Definition coordSet.H:160
Ostream & write(Ostream &os) const
Write to stream.
Definition coordSet.C:187
coordSet(const word &name, const coordFormat axisType)
Default construct with name and axis type.
Definition coordSet.C:70
void setPoints(const List< point > &newPoints)
Copy assign new points.
Definition coordSet.H:203
const word & name() const noexcept
The coord-set name.
Definition coordSet.H:152
void setDistance(const scalarList &dist, const bool check=true)
Copy assign the cumulative distance.
Definition coordSet.H:219
coordFormat
Enumeration defining the output format for coordinates.
Definition coordSet.H:61
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition meshSearch.H:57
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition polyMesh.C:1496
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
void findTetFacePt(const label celli, const point &p, label &tetFacei, label &tetPti) const
Find the tetFacei and tetPti for point p in celli.
Definition polyMesh.C:1370
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition polyMesh.C:1107
const vectorField & faceCentres() const
const vectorField & cellCentres() const
const cellList & cells() const
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition sampledSet.H:82
labelList faces_
Face numbers (-1 if not known).
Definition sampledSet.H:113
const meshSearch & searchEngine() const noexcept
Definition sampledSet.H:378
void setSamples(const List< point > &samplingPts, const labelList &samplingCells, const labelList &samplingFaces, const labelList &samplingSegments, const scalarList &samplingDistance)
Set sample data. Copy list contents.
Definition sampledSet.C:363
label pointInCell(const point &p, const label samplei) const
Return the cell in which the point on the sample line.
Definition sampledSet.C:80
label getNeighbourCell(const label) const
Returns the neighbour cell or the owner if face in on the boundary.
Definition sampledSet.C:68
scalar calcSign(const label facei, const point &sample) const
Calculates inproduct of face normal and vector sample-face centre.
Definition sampledSet.C:152
bool getTrackingPoint(const point &samplePt, const point &bPoint, const label bFacei, const scalar smallDist, point &trackPt, label &trackCelli, label &trackFacei) const
Calculates start of tracking given samplePt and first boundary.
Definition sampledSet.C:267
labelList segments_
Segment numbers.
Definition sampledSet.H:103
labelList cells_
Cell numbers.
Definition sampledSet.H:108
Ostream & write(Ostream &) const
Output for debugging.
Definition sampledSet.C:512
sampledSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const coordSet::coordFormat axisType)
Construct from components.
Definition sampledSet.C:405
label getBoundaryCell(const label) const
Returns cell next to boundary face.
Definition sampledSet.C:62
void checkDimensions() const
Check for consistent sizing.
Definition sampledSet.C:40
point pushIn(const point &sample, const label facei) const
Moves sample in direction of -n to it is 'inside' of facei.
Definition sampledSet.C:202
const polyMesh & mesh() const noexcept
Definition sampledSet.H:373
const labelList & cells() const noexcept
Definition sampledSet.H:388
static autoPtr< sampledSet > New(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Return a reference to the selected sampledSet.
Definition sampledSet.C:477
label findNearFace(const label celli, const point &sample, const scalar smallDist) const
Returns face label (or -1) of face which is close to sample.
Definition sampledSet.C:176
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
const cellShapeList & cells
#define WarningInFunction
Report a warning using Foam::Warning.
#define DebugInFunction
Report an information message using Foam::Info.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
PointHit< point > pointHit
A PointHit with a 3D point.
Definition pointHit.H:51
dimensionedScalar sign(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
vector point
Point is a vector.
Definition point.H:37
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299