68scalar calcVertexNormalWeight
76 label index =
f.find(pI);
97 auto& pointNormals = tpointNormals.ref();
109 const label faceI =
pFaces[fI];
114 scalar weight = calcVertexNormalWeight
122 pointNormals[pI] += weight * areaNorm;
128 return tpointNormals;
135 const bitSet& isFeaturePoint,
136 const List<surfaceFeatures::edgeStatus>& edgeStat
153 const edge&
e =
s.edges()[edgeI];
156 if (!isFeaturePoint.
test(
e[i]))
158 pointNormals[
e[i]] =
Zero;
168 const labelList& eFaces = edgeFaces[edgeI];
172 for (
const label facei : eFaces)
174 n +=
s.faceNormals()[facei];
180 const edge&
e =
s.edges()[edgeI];
183 if (!isFeaturePoint.
test(
e[i]))
185 pointNormals[
e[i]] +=
n;
194 if (nNormals[pointI] > 0)
202 forAll(pointNormals, pointI)
204 if (
mag(
mag(pointNormals[pointI])-1) > SMALL)
212 return tpointNormals;
216void detectSelfIntersections
219 bitSet& isEdgeIntersecting
229 isEdgeIntersecting =
false;
233 const edge&
e = edges[edgeI];
241 treeDataTriSurface::findSelfIntersectOp(
tree, edgeI)
247 isEdgeIntersecting.
set(edgeI);
253label detectIntersectionPoints
255 const scalar tolerance,
258 const scalar extendFactor,
262 const bool checkSelfIntersect,
263 const bitSet& initialIsEdgeIntersecting,
265 bitSet& isPointOnHitEdge,
269 const pointField initialLocalPoints(initialPoints,
s.meshPoints());
270 const List<labelledTri>& localFaces =
s.localFaces();
274 isPointOnHitEdge.
setSize(
s.nPoints());
275 isPointOnHitEdge =
false;
280 scalar tol =
Foam::max(tolerance, 10*
s.tolerance());
286 pointField start(initialLocalPoints+tol*displacement);
287 pointField end(initialLocalPoints+extendFactor*displacement);
289 List<pointIndexHit> hits;
290 s.findLineAny(start, end, hits);
297 && !localFaces[hits[pointI].index()].
found(pointI)
300 scale[pointI] =
Foam::max(0.0, scale[pointI]-0.2);
302 isPointOnHitEdge.
set(pointI);
310 if (checkSelfIntersect)
312 bitSet isEdgeIntersecting;
313 detectSelfIntersections(
s, isEdgeIntersecting);
321 const edge&
e = edges[edgeI];
323 if (isEdgeIntersecting[edgeI] && !initialIsEdgeIntersecting[edgeI])
325 if (isPointOnHitEdge.
set(
e[0]))
327 label start =
s.meshPoints()[
e[0]];
330 Pout<<
"Additional self intersection at "
337 if (isPointOnHitEdge.
set(
e[1]))
339 label
end =
s.meshPoints()[
e[1]];
342 Pout<<
"Additional self intersection at "
365 auto& res = tres.ref();
373 const edge&
e = edges[edgeI];
374 const scalar w = edgeWeights[edgeI];
376 res[
e[0]] += w*
fld[
e[1]];
377 sumWeight[
e[0]] += w;
379 res[
e[1]] += w*
fld[
e[0]];
380 sumWeight[
e[1]] += w;
388 if (
mag(sumWeight[pointI]) < VSMALL)
391 res[pointI] =
fld[pointI];
395 res[pointI] /= sumWeight[pointI];
406 const bitSet& isAffectedPoint,
419 const edge&
e = edges[edgeI];
422 edgeWeights[edgeI] = 1.0/(
Foam::max(w, SMALL));
431 if (isAffectedPoint.
test(pointI))
436 0.5*
fld[pointI] + 0.5*avgFld[pointI]
448 const bitSet& isFeaturePoint,
449 const List<surfaceFeatures::edgeStatus>& edgeStat,
450 const bitSet& isAffectedPoint
457 bitSet isSmoothPoint(isAffectedPoint);
460 bitSet newIsSmoothPoint(isSmoothPoint);
463 const edge&
e = edges[edgeI];
464 if (isSmoothPoint.test(
e[0]))
466 newIsSmoothPoint.set(
e[1]);
468 if (isSmoothPoint.test(
e[1]))
470 newIsSmoothPoint.set(
e[0]);
473 Info<<
"From points-to-smooth " << isSmoothPoint.count()
474 <<
" to " << newIsSmoothPoint.count()
476 isSmoothPoint.transfer(newIsSmoothPoint);
480 for (label i = 0; i < nSmooth; i++)
486 forAll(isSmoothPoint, pointI)
488 if (isSmoothPoint[pointI] && !isFeaturePoint[pointI])
495 avg += faceCentres[
pFaces[pFaceI]];
498 newPoints[meshPoints[pointI]] = avg;
514 const edge&
e = edges[edgeI];
518 pointSum[
e[0]] += eMid;
520 pointSum[
e[1]] += eMid;
529 isSmoothPoint[pointI]
530 && isFeaturePoint[pointI]
531 && nPointSum[pointI] >= 2
534 newPoints[meshPoints[pointI]] =
535 pointSum[pointI]/nPointSum[pointI];
540 s.movePoints(newPoints);
545 bitSet newIsSmoothPoint(isSmoothPoint);
548 const edge&
e = edges[edgeI];
549 if (isSmoothPoint[
e[0]])
551 newIsSmoothPoint.set(
e[1]);
553 if (isSmoothPoint[
e[1]])
555 newIsSmoothPoint.set(
e[0]);
558 Info<<
"From points-to-smooth " << isSmoothPoint.count()
559 <<
" to " << newIsSmoothPoint.count()
561 isSmoothPoint.transfer(newIsSmoothPoint);
570int main(
int argc,
char *argv[])
574 "Inflates surface according to point normals."
580 "Creates inflated version of surface using point normals. "
581 "Takes surface, distance to inflate and additional safety factor"
585 "checkSelfIntersection",
586 "Also check for self-intersection"
592 "Number of smoothing iterations (default 20)"
603 "Switch on additional debug information"
615 const auto inputName =
args.get<
word>(1);
617 const auto extendFactor =
args.get<scalar>(3);
618 const bool checkSelfIntersect =
args.found(
"checkSelfIntersection");
619 const auto nSmooth =
args.getOrDefault<label>(
"nSmooth", 10);
620 const auto featureAngle =
args.getOrDefault<scalar>(
"featureAngle", 180);
624 Info<<
"Inflating surface " << inputName
625 <<
" according to point normals" <<
nl
627 <<
" safety factor : " << extendFactor <<
nl
628 <<
" self intersection test : " << checkSelfIntersect <<
nl
629 <<
" smoothing iterations : " << nSmooth <<
nl
630 <<
" feature angle : " << featureAngle <<
nl
634 if (extendFactor < 1 || extendFactor > 10)
637 <<
"Illegal safety factor " << extendFactor
638 <<
". It is usually 1..2"
661 Info<<
"Detected features:" <<
nl
662 <<
" feature points : " << features.featurePoints().size()
663 <<
" out of " <<
s.nPoints() <<
nl
664 <<
" feature edges : " << features.featureEdges().size()
665 <<
" out of " <<
s.nEdges() <<
nl
668 bitSet isFeaturePoint(
s.nPoints(), features.featurePoints());
670 const List<surfaceFeatures::edgeStatus> edgeStat(features.toStatus());
677 Info<<
"Constructing field scale\n" <<
endl;
696 Info<<
"Calculating vertex normals\n" <<
endl;
709 Info<<
"Constructing field pointDisplacement\n" <<
endl;
723 (
distance*scale.field()*pointNormals)
731 bitSet isScaledPoint(
s.nPoints());
735 bitSet initialIsEdgeIntersecting;
736 if (checkSelfIntersect)
738 detectSelfIntersections(
s, initialIsEdgeIntersecting);
744 initialIsEdgeIntersecting.
setSize(
s.nEdges(),
true);
757 newPoints[meshPoints[i]] += pointDisplacement[i];
759 s.movePoints(newPoints);
760 Info<<
"Bounding box : " <<
s.bounds() <<
endl;
766 if (
s.pointFaces()[pointI].size())
768 scale[pointI] =
mag(pointDisplacement[pointI])/
distance;
783 bitSet isAffectedPoint;
784 label nIntersections = detectIntersectionPoints
793 initialIsEdgeIntersecting,
798 Info<<
"Detected " << nIntersections <<
" intersections" <<
nl
801 if (nIntersections == 0)
809 forAll(isAffectedPoint, pointI)
811 if (isAffectedPoint.
test(pointI))
813 isScaledPoint.set(pointI);
819 for (label i = 0; i < nSmooth; i++)
822 oldScale.rename(
"oldScale");
826 bitSet(
s.nPoints(),
true),
835 pointDisplacement.normalise();
836 pointDisplacement.field() *=
distance*scale.field();
840 lloydsSmoothing(nSmooth,
s, isFeaturePoint, edgeStat, isAffectedPoint);
848 label meshPointI = meshPoints[i];
849 pointDisplacement[i] =
pts[meshPointI]-initialPoints[meshPointI];
860 Info<<
"Detected overall intersecting points " << isScaledPoint.count()
861 <<
" out of " <<
s.nPoints() <<
nl <<
endl;
868 for (
const label pointi : isScaledPoint)
870 str.write(initialPoints[meshPoints[pointi]]);
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))
void normalise()
Inplace normalise this field. Generally a no-op except for vector fields.
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
void resize_nocopy(const label numElem)
Currently identical to resize. Subject to future change (Oct-2021).
void setSize(const label n, unsigned int val=0u)
Alias for resize().
label nPoints() const
Number of points supporting patch faces.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & points() const noexcept
Return reference to global points.
const labelListList & pointFaces() const
Return point-face addressing.
void size(const label n)
Older name for setAddressableSize.
static void noFunctionObjects(bool addWithOption=false)
Remove '-noFunctionObjects' option and ignore any occurrences.
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
static void noParallel()
Remove the parallel options.
static void addOption(const word &optName, const string ¶m="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
static void addNote(const string ¬e)
Add extra notes for the usage information.
void set(const bitSet &bitset)
Set specified bits from another bitset.
bool test(const label pos) const
Test for True value at specified position, never auto-vivify entries.
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Non-pointer based hierarchical recursive searching.
Holds feature edges/points of surface.
@ NONE
Not a classified feature edge.
A class for managing temporary objects.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
IOoject and searching on triSurface.
Triangulated surface description with patch information.
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.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const dimensionedScalar mp
Proton mass.
Namespace for handling debugging switches.
List< edge > edgeList
List of edge.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
scalar distance(const vector &p1, const vector &p2)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
List< labelList > labelListList
List of labelList.
List< label > labelList
A List of labels.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
messageStream Info
Information stream (stdout output on master, null elsewhere).
DimensionedField< vector, triSurfacePointGeoMesh > triSurfacePointVectorField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
errorManip< error > abort(error &err)
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
static constexpr const zero Zero
Global zero (0).
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.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a).
DimensionedField< scalar, triSurfacePointGeoMesh > triSurfacePointScalarField
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} pointMap[start]=pointMap[end]=Foam::min(start, end);} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};constexpr label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< DynamicList< face > > pFaces[nBCs]
Tree tree(triangles.begin(), triangles.end())
Foam::argList args(argc, argv)
#define forAll(list, i)
Loop across all elements in list.