Loading...
Searching...
No Matches
distanceSurface.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) 2016-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 "distanceSurface.H"
30#include "dictionary.H"
31#include "volFields.H"
34#include "fvMesh.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
41}
42
43
44// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45
46const Foam::Enum
47<
48 Foam::distanceSurface::topologyFilterType
49>
50Foam::distanceSurface::topoFilterNames_
51({
52 { topologyFilterType::NONE, "none" },
53 { topologyFilterType::LARGEST_REGION, "largestRegion" },
54 { topologyFilterType::NEAREST_POINTS, "nearestPoints" },
55 { topologyFilterType::PROXIMITY_REGIONS, "proximityRegions" },
56 { topologyFilterType::PROXIMITY_FACES, "proximityFaces" },
57 { topologyFilterType::PROXIMITY_FACES, "proximity" },
58});
59
61// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
62
63namespace Foam
64{
65
66// Check that all point hits are valid
67static inline void checkAllHits(const UList<pointIndexHit>& nearest)
68{
69 label notHit = 0;
70 for (const pointIndexHit& pHit : nearest)
71 {
72 if (!pHit.hit())
73 {
74 ++notHit;
75 }
76 }
77
78 if (notHit)
79 {
81 << "Had " << notHit << " faces/cells from "
82 << nearest.size() << " without a point hit." << nl
83 << "May be caused by a severely degenerate input surface" << nl
84 << abort(FatalError);
85 }
86}
87
88
89// Normal distance from surface hit point to a point in the mesh
90static inline scalar normalDistance_zero
91(
92 const point& pt,
93 const pointIndexHit& pHit,
94 const vector& norm
95)
96{
97 const vector diff(pt - pHit.point());
99 return (diff & norm);
100}
101
102
103// Signed distance from surface hit point to a point in the mesh,
104// the sign is dictated by the normal
105static inline scalar normalDistance_nonzero
106(
107 const point& pt,
108 const pointIndexHit& pHit,
109 const vector& norm
110)
111{
112 const vector diff(pt - pHit.point());
113 const scalar normDist = (diff & norm);
114
115 return Foam::sign(normDist) * Foam::mag(diff);
116}
117
118
119// Normal distance from surface hit point to a point in the mesh
120static inline void calcNormalDistance_zero
121(
123 const pointField& points,
124 const List<pointIndexHit>& nearest,
125 const vectorField& normals
126)
127{
128 forAll(nearest, i)
129 {
130 distance[i] =
131 normalDistance_zero(points[i], nearest[i], normals[i]);
132 }
133}
134
135
136// Signed distance from surface hit point -> point in the mesh,
137// the sign is dictated by the normal
138static inline void calcNormalDistance_nonzero
139(
141 const pointField& points,
142 const List<pointIndexHit>& nearest,
143 const vectorField& normals
144)
145{
146 forAll(nearest, i)
147 {
148 distance[i] =
149 normalDistance_nonzero(points[i], nearest[i], normals[i]);
150 }
151}
152
153
154// Close to the surface: normal distance from surface hit point
155// Far from surface: distance from surface hit point
156//
157// Note
158// This switch may be helpful when working directly with
159// distance/gradient fields. Has low overhead otherwise.
160// May be replaced in the future (2020-11)
161static inline void calcNormalDistance_filtered
162(
164 const bitSet& ignoreLocation,
165 const pointField& points,
166 const List<pointIndexHit>& nearest,
167 const vectorField& normals
168)
169{
170 forAll(nearest, i)
171 {
172 if (ignoreLocation.test(i))
173 {
174 distance[i] =
175 normalDistance_nonzero(points[i], nearest[i], normals[i]);
176 }
177 else
178 {
179 distance[i] =
180 normalDistance_zero(points[i], nearest[i], normals[i]);
181 }
182 }
183}
184
185
186// Flat surfaces (eg, a plane) have an extreme change in
187// the normal at the edge, which creates a zero-crossing
188// extending to infinity.
190// Ad hoc treatment: require that the surface hit point is within a
191// somewhat generous bounding box for the cell (+10%)
192//
193// Provisioning for filtering based on the cell points,
194// but its usefulness remains to be seen (2020-12-09)
195template<bool WantPointFilter = false>
197(
198 bitSet& ignoreCells,
199 const List<pointIndexHit>& nearest,
200 const polyMesh& mesh,
201 const scalar boundBoxInflate = 0.1 // 10% to catch corners
202)
203{
204 // A deny filter. Initially false (accept everything)
205 ignoreCells.resize(mesh.nCells());
206
207 bitSet ignorePoints;
208 if (WantPointFilter)
209 {
210 // Create deny filter
211 ignorePoints.resize(mesh.nPoints(), true);
212 }
213
214 forAll(nearest, celli)
215 {
216 const point& pt = nearest[celli].point();
217
218 boundBox cellBb(mesh.cellBb(celli));
219 cellBb.inflate(boundBoxInflate);
220
221 if (!cellBb.contains(pt))
222 {
223 ignoreCells.set(celli);
224 }
225 else if (WantPointFilter)
226 {
227 // Good cell candidate, do not ignore its points
228 ignorePoints.unset(mesh.cellPoints(celli));
229 }
230 }
231
232 return ignorePoints;
233}
235
236} // End namespace Foam
237
238
239// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
240
242(
243 const word& defaultSurfaceName,
244 const polyMesh& mesh,
245 const dictionary& dict
246)
247:
248 mesh_(mesh),
249 geometryPtr_
250 (
252 (
253 dict.get<word>("surfaceType"),
255 (
256 dict.getOrDefault("surfaceName", defaultSurfaceName),
257 mesh.time().constant(), // directory
258 "triSurface", // instance
259 mesh.time(), // registry
260 IOobject::MUST_READ,
261 IOobject::NO_WRITE
262 ),
263 dict
264 )
265 ),
266 distance_(dict.getOrDefault<scalar>("distance", 0)),
267 withZeroDistance_(equal(distance_, 0)),
268 withSignDistance_
269 (
270 withZeroDistance_
271 || (distance_ < 0)
272 || dict.getOrDefault<bool>("signed", true)
273 ),
274
275 isoParams_
276 (
277 dict,
278 isoSurfaceParams::ALGO_TOPO,
279 isoSurfaceParams::filterType::DIAGCELL
280 ),
281 topoFilter_
282 (
283 topoFilterNames_.getOrDefault
284 (
285 "topology",
286 dict,
287 topologyFilterType::NONE
288 )
289 ),
290 nearestPoints_(),
291 maxDistanceSqr_(Foam::sqr(GREAT)),
292 absProximity_(dict.getOrDefault<scalar>("absProximity", 1e-5)),
293 cellDistancePtr_(nullptr),
294 pointDistance_(),
295 surface_(),
296 meshCells_(),
297 isoSurfacePtr_(nullptr)
298{
299 if (topologyFilterType::NEAREST_POINTS == topoFilter_)
300 {
301 dict.readEntry("nearestPoints", nearestPoints_);
302 }
303
304 if (dict.readIfPresent("maxDistance", maxDistanceSqr_))
305 {
306 maxDistanceSqr_ = Foam::sqr(maxDistanceSqr_);
307 }
308}
309
310
312(
313 const polyMesh& mesh,
314 const word& surfaceType,
315 const word& surfaceName,
316 const isoSurfaceParams& params,
317 const bool interpolate
318)
319:
321 (
322 mesh,
324 surfaceType,
325 surfaceName,
326 scalar(0),
327 true, // redundant - must be signed
328 params
329 )
330{}
331
332
334(
335 const polyMesh& mesh,
336 const bool interpolate,
337 const word& surfaceType,
338 const word& surfaceName,
339 const scalar distance,
340 const bool useSignedDistance,
341 const isoSurfaceParams& params
342)
343:
345 (
346 mesh,
349 (
350 surfaceType,
352 (
353 surfaceName, // name
354 mesh.time().constant(), // directory
355 "triSurface", // instance
356 mesh.time(), // registry
357 IOobject::MUST_READ,
358 IOobject::NO_WRITE
359 ),
360 dictionary()
361 ),
363 useSignedDistance,
364 params
365 )
366{}
367
368
370(
371 const polyMesh& mesh,
372 const bool interpolate,
374 const scalar distance,
375 const bool useSignedDistance,
376 const isoSurfaceParams& params
377)
378:
379 mesh_(mesh),
380 geometryPtr_(surface),
381 distance_(distance),
382 withZeroDistance_(equal(distance_, 0)),
383 withSignDistance_
384 (
385 withZeroDistance_
386 || (distance_ < 0)
387 || useSignedDistance
388 ),
389
390 isoParams_(params),
391 topoFilter_(topologyFilterType::NONE),
392 nearestPoints_(),
393 maxDistanceSqr_(Foam::sqr(GREAT)),
394 absProximity_(1e-5),
395 cellDistancePtr_(nullptr),
396 pointDistance_(),
397 surface_(),
398 meshCells_(),
399 isoSurfacePtr_(nullptr)
400{}
401
402
403// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
404
406{
407 if (debug)
408 {
409 Pout<< "distanceSurface::createGeometry updating geometry." << endl;
410 }
411
412 // Clear any previously stored topologies
413 isoSurfacePtr_.reset(nullptr);
414 surface_.clear();
415 meshCells_.clear();
416
417 // Doing searches on this surface
418 const searchableSurface& geom = geometryPtr_();
419
420 const fvMesh& fvmesh = static_cast<const fvMesh&>(mesh_);
421
422 // Distance to cell centres
423 // ~~~~~~~~~~~~~~~~~~~~~~~~
424
425 cellDistancePtr_.reset
426 (
428 (
429 IOobject
430 (
431 "distanceSurface.cellDistance",
432 fvmesh.time().timeName(),
433 fvmesh.time(),
437 ),
438 fvmesh,
440 )
441 );
442 auto& cellDistance = *cellDistancePtr_;
443
444
445 // For distance = 0 we apply additional geometric filtering
446 // to limit the extent of open edges.
447 //
448 // Does not work with ALGO_POINT
449
450 bitSet ignoreCells, ignoreCellPoints;
451
452 const bool filterCells =
453 (
454 withZeroDistance_
455 && isoParams_.algorithm() != isoSurfaceParams::ALGO_POINT
456 );
457
458
459 // Internal field
460 {
461 const pointField& cc = fvmesh.C();
462 scalarField& fld = cellDistance.primitiveFieldRef();
463
464 List<pointIndexHit> nearest;
465 geom.findNearest
466 (
467 cc,
468 // Use initialized field (GREAT) to limit search too
469 fld,
470 nearest
471 );
472 checkAllHits(nearest);
473
474 // Geometric pre-filtering when distance == 0
475
476 // NOTE (2021-05-31)
477 // Can skip the prefilter if we use proximity-regions filter anyhow
478 // but it makes the iso algorithm more expensive and doesn't help
479 // unless we start relying on area-based weighting for rejecting regions.
480
481 if (filterCells)
482 {
483 ignoreCellPoints = simpleGeometricFilter<false>
484 (
485 ignoreCells,
486 nearest,
487 fvmesh,
488
489 // Inflate bound box.
490 // - To catch corners: approx. 10%
491 // - Extra generous for PROXIMITY_REGIONS
492 // (extra weighting for 'bad' faces)
493 (
494 topologyFilterType::PROXIMITY_REGIONS == topoFilter_
495 ? 1
496 : 0.1
497 )
498 );
499 }
500
501 if (withSignDistance_)
502 {
503 vectorField norms;
504 geom.getNormal(nearest, norms);
505
506 if (filterCells)
507 {
508 // With inside/outside switching (see note above)
510 (
511 fld,
512 ignoreCells,
513 cc,
514 nearest,
515 norms
516 );
517 }
518 else if (withZeroDistance_)
519 {
520 calcNormalDistance_zero(fld, cc, nearest, norms);
521 }
522 else
523 {
524 calcNormalDistance_nonzero(fld, cc, nearest, norms);
525 }
526 }
527 else
528 {
529 calcAbsoluteDistance(fld, cc, nearest);
530 }
531 }
532
533
534 // Patch fields
535 {
536 forAll(fvmesh.C().boundaryField(), patchi)
537 {
538 const pointField& cc = fvmesh.C().boundaryField()[patchi];
539 scalarField& fld = cellDistance.boundaryFieldRef()[patchi];
540
541 List<pointIndexHit> nearest;
542 geom.findNearest
543 (
544 cc,
545 scalarField(cc.size(), GREAT),
546 nearest
547 );
548 checkAllHits(nearest);
549
550 if (withSignDistance_)
551 {
552 vectorField norms;
553 geom.getNormal(nearest, norms);
554
555 if (withZeroDistance_)
556 {
557 // Slight inconsistency in boundary vs interior when
558 // cells are filtered, but the patch fields are only
559 // used by isoSurfacePoint, and filtering is disabled
560 // for that anyhow.
561
562 calcNormalDistance_zero(fld, cc, nearest, norms);
563 }
564 else
565 {
566 calcNormalDistance_nonzero(fld, cc, nearest, norms);
567 }
568 }
569 else
570 {
571 calcAbsoluteDistance(fld, cc, nearest);
572 }
573 }
574 }
575
576
577 // On processor patches the mesh.C() will already be the cell centre
578 // on the opposite side so no need to swap cellDistance.
579
580
581 // Distance to points
582 pointDistance_.resize(fvmesh.nPoints());
583 pointDistance_ = GREAT;
584 {
585 const pointField& pts = fvmesh.points();
586 scalarField& fld = pointDistance_;
587
588 List<pointIndexHit> nearest;
589 geom.findNearest
590 (
591 pts,
592 // Use initialized field (GREAT) to limit search too
593 pointDistance_,
594 nearest
595 );
596 checkAllHits(nearest);
597
598 if (withSignDistance_)
599 {
600 vectorField norms;
601 geom.getNormal(nearest, norms);
602
603 if (filterCells)
604 {
605 // With inside/outside switching (see note above)
607 (
608 fld,
609 ignoreCellPoints,
610 pts,
611 nearest,
612 norms
613 );
614 }
615 else if (withZeroDistance_)
616 {
617 calcNormalDistance_zero(fld, pts, nearest, norms);
618 }
619 else
620 {
621 calcNormalDistance_nonzero(fld, pts, nearest, norms);
622 }
623 }
624 else
625 {
626 calcAbsoluteDistance(fld, pts, nearest);
627 }
628 }
629
630
631 // Don't need ignoreCells if there is nothing to ignore.
632 if (ignoreCells.none())
633 {
634 ignoreCells.clearStorage();
635 }
636 else if (filterCells && topologyFilterType::NONE != topoFilter_)
637 {
638 // For refine blocked cells (eg, checking actual cells cut)
639 isoSurfaceBase isoCutter
640 (
641 mesh_,
642 cellDistance,
643 pointDistance_,
644 distance_
645 );
646
647 if (topologyFilterType::LARGEST_REGION == topoFilter_)
648 {
649 refineBlockedCells(ignoreCells, isoCutter);
650 filterKeepLargestRegion(ignoreCells);
651 }
652 else if (topologyFilterType::NEAREST_POINTS == topoFilter_)
653 {
654 refineBlockedCells(ignoreCells, isoCutter);
655 filterKeepNearestRegions(ignoreCells);
656 }
657
658 // Note: apply similar filtering for PROXIMITY_REGIONS later instead
659 }
660
661 // Don't need point filter beyond this point
662 ignoreCellPoints.clearStorage();
663
664
665 if (debug)
666 {
667 Pout<< "Writing cell distance:" << cellDistance.objectPath() << endl;
668 cellDistance.write();
669
670 pointScalarField pDist
671 (
673 (
674 "distanceSurface.pointDistance",
675 fvmesh.time().timeName(),
676 fvmesh.time(),
680 ),
681 pointMesh::New(fvmesh),
683 );
684 pDist.primitiveFieldRef() = pointDistance_;
685
686 Pout<< "Writing point distance:" << pDist.objectPath() << endl;
687 pDist.write();
688 }
689
690 isoSurfacePtr_.reset
691 (
693 (
694 isoParams_,
695 cellDistance,
696 pointDistance_,
697 distance_,
698 ignoreCells
699 )
700 );
701
702
703 // Restrict ignored cells to those actually cut
704 if (filterCells && topologyFilterType::PROXIMITY_REGIONS == topoFilter_)
705 {
706 isoSurfaceBase isoCutter
707 (
708 mesh_,
709 cellDistance,
710 pointDistance_,
711 distance_
712 );
713
714 refineBlockedCells(ignoreCells, isoCutter);
715 }
716
717
718 // ALGO_POINT still needs cell, point fields (for interpolate)
719 // The others can do straight transfer
720
721 // But also flatten into a straight transfer for proximity filtering
722 if
723 (
724 isoParams_.algorithm() != isoSurfaceParams::ALGO_POINT
725 || topologyFilterType::PROXIMITY_FACES == topoFilter_
726 || topologyFilterType::PROXIMITY_REGIONS == topoFilter_
727 )
728 {
729 surface_.transfer(static_cast<meshedSurface&>(*isoSurfacePtr_));
730 meshCells_.transfer(isoSurfacePtr_->meshCells());
731
732 isoSurfacePtr_.reset(nullptr);
733 cellDistancePtr_.reset(nullptr);
734 pointDistance_.clear();
735 }
736
737 if (topologyFilterType::PROXIMITY_FACES == topoFilter_)
738 {
739 filterFaceProximity();
740 }
741 else if (topologyFilterType::PROXIMITY_REGIONS == topoFilter_)
742 {
743 filterRegionProximity(ignoreCells);
744 }
745
746 if (debug)
748 print(Pout, debug);
749 Pout<< endl;
750 }
751}
752
753
754void Foam::distanceSurface::print(Ostream& os, int level) const
755{
756 os << " surface:" << surfaceName()
757 << " distance:" << distance()
758 << " topology:" << topoFilterNames_[topoFilter_];
759
760 isoParams_.print(os);
761
762 if (level)
763 {
764 os << " faces:" << surface().surfFaces().size()
765 << " points:" << surface().points().size();
766 }
767}
768
769
770// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
fileName objectPath() const
The complete path + object name.
Definition IOobjectI.H:313
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
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
void resize(const label numElem, const unsigned int val=0u)
Reset addressable list size, does not shrink the allocated size.
void clearStorage()
Clear the list and delete storage.
const point_type & point() const noexcept
Return point, no checks.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
bool none() const
True if no bits in this bitset are set.
Definition bitSetI.H:414
bool test(const label pos) const
Test for True value at specified position, never auto-vivify entries.
Definition bitSet.H:334
bitSet & unset(const bitSet &other)
Unset (subtract) the bits specified in the other bitset, which is a set difference corresponds to the...
Definition bitSetI.H:540
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
bool contains(const point &pt) const
Contains point? (inside or on edge).
Definition boundBoxI.H:455
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Definition boundBoxI.H:381
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A surface defined by a distance from an input searchable surface. Uses an iso-surface algorithm (cell...
void filterFaceProximity()
Adjust extracted iso-surface to remove far faces.
void createGeometry()
Create/recreate the distance surface.
void filterRegionProximity(bitSet &ignoreCells) const
Remove region(s) that have far faces.
scalar distance() const noexcept
The distance to the underlying searchableSurface.
distanceSurface(const word &defaultSurfaceName, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
void print(Ostream &os, int level=0) const
Print information.
void filterKeepNearestRegions(bitSet &ignoreCells) const
Keep region(s) closest to the nearest points.
const meshedSurface & surface() const
The underlying surface.
bool refineBlockedCells(bitSet &ignoreCells, const isoSurfaceBase &isoContext) const
Re-filter the blocked cells based on the anticipated cuts.
void filterKeepLargestRegion(bitSet &ignoreCells) const
Keep region with the most cuts (after region split).
const word & surfaceName() const
The name of the underlying searchableSurface.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const volVectorField & C() const
Return cell centres as volVectorField.
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
Low-level components common to various iso-surface algorithms.
static autoPtr< isoSurfaceBase > New(const isoSurfaceParams &params, const volScalarField &cellValues, const scalarField &pointValues, const scalar iso, const bitSet &ignoreCells=bitSet())
Create for specified algorithm type.
Preferences for controlling iso-surface algorithms.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
label nPoints() const noexcept
Number of mesh points.
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const =0
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const =0
From a set of points and indices get the normal.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
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
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
Different types of constants.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
static void calcNormalDistance_zero(scalarField &distance, const pointField &points, const List< pointIndexHit > &nearest, const vectorField &normals)
scalar distance(const vector &p1, const vector &p2)
Definition curveTools.C:12
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
static scalar normalDistance_zero(const point &pt, const pointIndexHit &pHit, const vector &norm)
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
bool equal(const T &a, const T &b)
Compare two values for equality.
Definition label.H:180
static void calcNormalDistance_filtered(scalarField &distance, const bitSet &ignoreLocation, const pointField &points, const List< pointIndexHit > &nearest, const vectorField &normals)
static scalar normalDistance_nonzero(const point &pt, const pointIndexHit &pHit, const vector &norm)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
MeshedSurface< face > meshedSurface
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition triad.C:373
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
static void calcNormalDistance_nonzero(scalarField &distance, const pointField &points, const List< pointIndexHit > &nearest, const vectorField &normals)
static void checkAllHits(const UList< pointIndexHit > &nearest)
static bitSet simpleGeometricFilter(bitSet &ignoreCells, const List< pointIndexHit > &nearest, const polyMesh &mesh, const scalar boundBoxInflate=0.1)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
const pointField & pts
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299