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
68
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;
151 }
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;
175
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;
198
199 if (inter.hit())
200 {
201 dist = mag(inter.hitPoint() - sample);
202 }
203 else
204 {
205 dist = mag(inter.missPoint() - sample);
206 }
207
208 if (dist < smallDist)
209 {
210 return myFaces[myFacei];
211 }
212 }
213 return -1;
214}
215
216
218(
219 const point& facePt,
220 const label facei
221) const
222{
223 label celli = mesh().faceOwner()[facei];
224 const point& cC = mesh().cellCentres()[celli];
225
226 point newPosition = facePt;
227
228 // Taken from particle::initCellFacePt()
229 label tetFacei;
230 label tetPtI;
231 mesh().findTetFacePt(celli, facePt, tetFacei, tetPtI);
232
233 // This is the tolerance that was defined as a static constant of the
234 // particle class. It is no longer used by particle, following the switch to
235 // barycentric tracking. The only place in which the tolerance is now used
236 // is here. I'm not sure what the purpose of this code is, but it probably
237 // wants removing. It is doing tet-searches for a particle position. This
238 // should almost certainly be left to the particle class.
239 const scalar trackingCorrectionTol = 1e-5;
240
241 if (tetFacei == -1 || tetPtI == -1)
242 {
243 newPosition = facePt;
244
245 label trap(1.0/trackingCorrectionTol + 1);
246
247 label iterNo = 0;
248
249 do
250 {
251 newPosition += trackingCorrectionTol*(cC - facePt);
252
254 (
255 celli,
256 newPosition,
257 tetFacei,
258 tetPtI
259 );
260
261 ++iterNo;
262
263 } while (tetFacei < 0 && iterNo <= trap);
264 }
265
266 if (tetFacei == -1)
267 {
269 << "After pushing " << facePt << " to " << newPosition
270 << " it is still outside face " << facei
271 << " at " << mesh().faceCentres()[facei]
272 << " of cell " << celli
273 << " at " << cC << endl
274 << "Please change your starting point"
275 << abort(FatalError);
276 }
277
278 return newPosition;
279}
280
281
283(
284 const point& samplePt,
285 const point& bPoint,
286 const label bFacei,
287 const scalar smallDist,
288
289 point& trackPt,
290 label& trackCelli,
291 label& trackFacei
292) const
293{
294 bool isGoodSample = false;
295
296 if (bFacei == -1)
297 {
298 // No boundary intersection. Try and find cell samplePt is in
299 trackCelli = mesh().findCell(samplePt, searchEngine_.decompMode());
300
301 if
302 (
303 (trackCelli == -1)
304 || !mesh().pointInCell
305 (
306 samplePt,
307 trackCelli,
308 searchEngine_.decompMode()
309 )
310 )
311 {
312 // Line samplePt - end_ does not intersect domain at all.
313 // (or is along edge)
314
315 trackCelli = -1;
316 trackFacei = -1;
317
318 isGoodSample = false;
319 }
320 else
321 {
322 // Start is inside. Use it as tracking point
323
324 trackPt = samplePt;
325 trackFacei = -1;
326
327 isGoodSample = true;
328 }
329 }
330 else if (mag(samplePt - bPoint) < smallDist)
331 {
332 // samplePt close to bPoint. Snap to it
333 trackPt = pushIn(bPoint, bFacei);
334 trackFacei = bFacei;
335 trackCelli = getBoundaryCell(trackFacei);
336
337 isGoodSample = true;
338 }
339 else
340 {
341 scalar sign = calcSign(bFacei, samplePt);
342
343 if (sign < 0)
344 {
345 // samplePt inside or marginally outside.
346 trackPt = samplePt;
347 trackFacei = -1;
348 trackCelli = mesh().findCell(trackPt, searchEngine_.decompMode());
349
350 isGoodSample = true;
351 }
352 else
353 {
354 // samplePt outside. use bPoint
355 trackPt = pushIn(bPoint, bFacei);
356 trackFacei = bFacei;
357 trackCelli = getBoundaryCell(trackFacei);
358
359 isGoodSample = false;
360 }
361 }
362
364 << " samplePt:" << samplePt
365 << " bPoint:" << bPoint
366 << " bFacei:" << bFacei << nl
367 << " Calculated first tracking point :"
368 << " trackPt:" << trackPt
369 << " trackCelli:" << trackCelli
370 << " trackFacei:" << trackFacei
371 << " isGoodSample:" << isGoodSample
372 << endl;
373
374 return isGoodSample;
375}
376
377
379(
380 const List<point>& samplingPts,
381 const labelList& samplingCells,
382 const labelList& samplingFaces,
383 const labelList& samplingSegments,
384 const scalarList& samplingDistance
385)
386{
387 setPoints(samplingPts);
388 setDistance(samplingDistance, false); // check=false
389
390 segments_ = samplingSegments;
391 cells_ = samplingCells;
392 faces_ = samplingFaces;
393
394 checkDimensions();
395}
396
397
399(
400 List<point>&& samplingPts,
401 labelList&& samplingCells,
402 labelList&& samplingFaces,
403 labelList&& samplingSegments,
404 scalarList&& samplingDistance
405)
406{
407 setPoints(std::move(samplingPts));
408 setDistance(std::move(samplingDistance), false); // check=false
409
410 segments_ = std::move(samplingSegments);
411 cells_ = std::move(samplingCells);
412 faces_ = std::move(samplingFaces);
413
414 checkDimensions();
415}
416
417
418// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
419
421(
422 const word& name,
423 const polyMesh& mesh,
424 const meshSearch& searchEngine,
425 const coordSet::coordFormat axisType
426)
427:
428 coordSet(name, axisType),
429 mesh_(mesh),
430 searchEngine_(searchEngine),
431 segments_(),
432 cells_(),
433 faces_()
434{}
435
436
438(
439 const word& name,
440 const polyMesh& mesh,
441 const meshSearch& searchEngine,
442 const word& axis
443)
444:
445 coordSet(name, axis),
446 mesh_(mesh),
447 searchEngine_(searchEngine),
448 segments_(),
449 cells_(),
450 faces_()
451{}
452
453
455(
456 const word& name,
457 const polyMesh& mesh,
458 const meshSearch& searchEngine,
459 const dictionary& dict
460)
461:
462 coordSet(name, dict.get<word>("axis")),
463 mesh_(mesh),
464 searchEngine_(searchEngine),
465 segments_(),
466 cells_(),
467 faces_()
468{}
469
470
471// Foam::autoPtr<Foam::coordSet> Foam::sampledSet::gather
472// (
473// labelList& sortOrder,
474// labelList& allSegments
475// ) const
476// {
477// autoPtr<coordSet> result(coordSet::gatherSort(sortOrder));
478//
479// // Optional
480// if (notNull(allSegments))
481// {
482// globalIndex::gatherOp(segments(), allSegments);
483// allSegments = UIndirectList<label>(allSegments, sortOrder)();
484// }
485//
486// return result;
487// }
488
489
490// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
491
493(
494 const word& name,
495 const polyMesh& mesh,
496 const meshSearch& searchEngine,
497 const dictionary& dict
498)
499{
500 const word sampleType(dict.get<word>("type"));
501
502 auto* ctorPtr = wordConstructorTable(sampleType);
503
504 if (!ctorPtr)
505 {
507 (
508 dict,
509 "sample",
510 sampleType,
511 *wordConstructorTablePtr_
512 ) << exit(FatalIOError);
513 }
514
516 (
517 ctorPtr
518 (
519 name,
520 mesh,
521 searchEngine,
522 dict.optionalSubDict(sampleType + "Coeffs")
523 )
524 );
525}
526
527
529{
531
532 os << nl << "\t(celli)\t(facei)" << nl;
533
534 forAll(*this, samplei)
535 {
536 os << '\t' << cells_[samplei]
537 << '\t' << faces_[samplei]
538 << nl;
539 }
540
541 return os;
542}
543
544
545// ************************************************************************* //
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:77
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
Definition: PointHit.H:54
const point_type & missPoint() const
Return the miss point. Fatal if hit.
Definition: PointHit.H:158
bool hit() const noexcept
Is there a hit.
Definition: PointHit.H:121
const point_type & hitPoint() const
Return the hit point. Fatal if not hit.
Definition: PointHit.H:145
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
label size() const noexcept
The number of elements in the UList.
Definition: UListI.H:420
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
Holds list of sampling positions.
Definition: coordSet.H:56
scalarList distance_
Cumulative distance for "distance" write specifier.
Definition: coordSet.H:85
coordFormat
Enumeration defining the output format for coordinates.
Definition: coordSet.H:63
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
Definition: dictionary.C:577
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
virtual bool write()
Write the output fields.
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:61
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1522
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1121
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:1396
bool pointInCell(const point &p, label celli, const cellDecomposition=CELL_TETS) const
Test if point p is in the celli.
Definition: polyMesh.C:1412
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1127
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:86
labelList faces_
Face numbers (-1 if not known)
Definition: sampledSet.H:107
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:379
label pointInCell(const point &p, const label samplei) const
Return the cell in which the point on the sample line.
Definition: sampledSet.C:87
label getNeighbourCell(const label) const
Returns the neighbour cell or the owner if face in on the boundary.
Definition: sampledSet.C:75
scalar calcSign(const label facei, const point &sample) const
Calculates inproduct of face normal and vector sample-face centre.
Definition: sampledSet.C:159
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:283
labelList segments_
Segment numbers.
Definition: sampledSet.H:101
labelList cells_
Cell numbers.
Definition: sampledSet.H:104
label getBoundaryCell(const label) const
Returns cell next to boundary face.
Definition: sampledSet.C:69
void checkDimensions() const
Check for consistent sizing.
Definition: sampledSet.C:47
point pushIn(const point &sample, const label facei) const
Moves sample in direction of -n to it is 'inside' of facei.
Definition: sampledSet.C:218
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:183
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
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 OpenFOAM.
dimensionedScalar sign(const dimensionedScalar &ds)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:680
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
IOerror FatalIOError
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333