47const Foam::scalar Foam::featurePointConformer::tolParallel = 1
e-3;
52Foam::vector Foam::featurePointConformer::sharedFaceNormal
54 const extendedFeatureEdgeMesh& feMesh,
59 const labelList& edgeInormals = feMesh.edgeNormals()[edgeI];
60 const labelList& nextEdgeInormals = feMesh.edgeNormals()[nextEdgeI];
62 const vector& A1 = feMesh.normals()[edgeInormals[0]];
63 const vector& A2 = feMesh.normals()[edgeInormals[1]];
65 const vector& B1 = feMesh.normals()[nextEdgeInormals[0]];
66 const vector& B2 = feMesh.normals()[nextEdgeInormals[1]];
73 const scalar A1B1 =
mag((A1 & B1) - 1.0);
74 const scalar A1B2 =
mag((A1 & B2) - 1.0);
75 const scalar A2B1 =
mag((A2 & B1) - 1.0);
76 const scalar A2B2 =
mag((A2 & B2) - 1.0);
83 if (A1B1 < A1B2 && A1B1 < A2B1 && A1B1 < A2B2)
87 else if (A1B2 < A1B1 && A1B2 < A2B1 && A1B2 < A2B2)
91 else if (A2B1 < A1B1 && A2B1 < A1B2 && A2B1 < A2B2)
102Foam::label Foam::featurePointConformer::getSign
120void Foam::featurePointConformer::addMasterAndSlavePoints
137 const vertexType masterType = masterPointsTypes[pI];
144 foamyHexMesh_.vertexCount() +
pts.
size(),
150 const label masterIndex =
pts.
last().index();
154 const planeDynList& masterPointPlanes = masterPointReflections[pI];
156 forAll(masterPointPlanes, planeI)
160 const plane& reflPlane = masterPointPlanes[planeI]();
162 const Foam::point slavePt = reflPlane.mirror(masterPt);
164 const vertexType slaveType =
176 foamyHexMesh_.vertexCount() +
pts.
size(),
182 ftPtPairs_.addPointPair
194void Foam::featurePointConformer::createMasterAndSlavePoints
211 && !foamyHexMesh_.decomposition().positionOnThisProcessor(featPt)
213 || geometryToConformTo_.outside(featPt)
219 const scalar ppDist = foamyHexMesh_.pointPairDistance(featPt);
226 const labelList& featPtEdges = feMesh.featurePointEdges()[ptI];
231 pointEdgeTypes.calcPointFeatureEdgesTypes();
241 const edgeStatus eStatusCurr = feMesh.getEdgeStatus(circ.curr());
254 label
sign = getSign(eStatusCurr);
256 const vector n = sharedFaceNormal(feMesh, circ.curr(), circ.next());
258 const vector pointMotionDirection =
sign*0.5*ppDist*
n;
264 if (masterPoints.empty())
269 planeDynList firstPlane;
272 masterPoints.append(pt);
274 masterPointsTypes.append
281 masterPointReflections.insert
283 masterPoints.size() - 1,
284 std::move(firstPlane)
326 masterPoints.last() += pointMotionDirection;
328 masterPointReflections[masterPoints.size() - 1].append
336 addMasterAndSlavePoints
340 masterPointReflections,
347void Foam::featurePointConformer::createMixedFeaturePoints
352 if (foamyHexMeshControls_.mixedFeaturePointPPDistanceCoeff() < 0)
355 <<
"Skipping specialised handling for mixed feature points" <<
endl;
362 geometryToConformTo_.features()
373 label ptI = feMesh.mixedStart();
374 ptI < feMesh.nonFeatureStart();
383 && !foamyHexMesh_.decomposition().positionOnThisProcessor(featPt))
388 const labelList& pEds = pointsEdges[ptI];
393 pFEdgeTypes.calcPointFeatureEdgesTypes();
395 bool specialisedSuccess =
false;
397 if (foamyHexMeshControls_.specialiseFeaturePoints())
400 createSpecialisedFeaturePoint
402 feMesh, pEds, pFEdgeTypes, allEdStat, ptI,
pts
406 if (!specialisedSuccess && foamyHexMeshControls_.edgeAiming())
448 const scalar edgeGroupDistance =
449 foamyHexMesh_.mixedFeaturePointDistance(pt);
453 const label edgeI = pEds[
e];
456 pt + edgeGroupDistance*feMesh.edgeDirection(edgeI, ptI);
460 foamyHexMesh_.createEdgePointGroup(feMesh, edgeHit,
pts);
472 geometryToConformTo_.features()
481 label ptI = feMesh.convexStart();
482 ptI < feMesh.mixedStart();
487 createMasterAndSlavePoints(feMesh, ptI,
pts);
490 if (foamyHexMeshControls_.guardFeaturePoints())
495 label ptI = feMesh.mixedStart();
496 ptI < feMesh.nonFeatureStart();
504 feMesh.points()[ptI],
516Foam::featurePointConformer::featurePointConformer
521 foamyHexMesh_(foamyHexMesh),
522 foamyHexMeshControls_(foamyHexMesh.foamyHexMeshControls()),
523 geometryToConformTo_(foamyHexMesh.geometryToConformTo()),
524 featurePointVertices_(),
525 ftPtPairs_(foamyHexMesh)
528 <<
"Conforming to feature points" <<
endl;
531 <<
indent <<
"Circulating edges is: "
532 << foamyHexMeshControls_.circulateEdges().c_str() << nl
533 << indent <<
"Guarding feature points is: "
534 << foamyHexMeshControls_.guardFeaturePoints().c_str() << nl
535 << indent <<
"Snapping to feature points is: "
536 << foamyHexMeshControls_.snapFeaturePoints().c_str() << nl
537 << indent <<
"Specialising feature points is: "
538 << foamyHexMeshControls_.specialiseFeaturePoints().c_str()
544 createFeaturePoints(
pts);
546 createMixedFeaturePoints(
pts);
568 if (foamyHexMeshControls_.objOutput())
570 DelaunayMeshTools::writeOBJ(
"featureVertices.obj",
pts);
573 featurePointVertices_.transfer(
pts);
591 decomposition.distributePoints(featurePointVertices_);
594 forAll(featurePointVertices_, vI)
608 forAll(featurePointVertices_, vI)
610 const label currentIndex = featurePointVertices_[vI].index();
612 const auto newIndexIter = oldToNewIndices.cfind(currentIndex);
614 if (newIndexIter.good())
616 featurePointVertices_[vI].index() = *newIndexIter;
620 ftPtPairs_.reIndex(oldToNewIndices);
CGAL::indexedVertex< K > Vb
bool externalBoundaryPoint() const
bool internalBoundaryPoint() const
Like Foam::Circulator, but with const-access iterators.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void append(const T &val)
Append an element at the end of the list.
A HashTable to objects of type <T> with a label key.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
void size(const label n)
Older name for setAddressableSize.
T & last()
Access last element of the list, position [size()-1].
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Store a background polyMesh to use for the decomposition of space and queries for parallel conformalV...
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Hold the types of feature edges attached to the point.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
List< labelList > labelListList
List of labelList.
dimensionedScalar sign(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Ostream & indent(Ostream &os)
Indent stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
vector point
Point is a vector.
vectorField pointField
pointField is a vectorField.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.