Go to the documentation of this file.
33 #include "indexedVertex.H"
57 template<
class Triangulation,
class Type>
60 const Triangulation&
mesh,
72 typename Triangulation::Finite_vertices_iterator vit =
73 mesh.finite_vertices_begin();
74 vit !=
mesh.finite_vertices_end();
80 newField[added++] = field[count];
86 newField.resize(added);
105 typename T::Finite_vertices_iterator vit =
mesh.finite_vertices_begin();
106 vit !=
mesh.finite_vertices_end();
115 std::list<typename T::Vertex_handle> adjVerts;
116 mesh.finite_adjacent_vertices(vit, std::back_inserter(adjVerts));
122 typename std::list<typename T::Vertex_handle>::const_iterator
123 adjVertI = adjVerts.begin();
124 adjVertI != adjVerts.end();
128 typename T::Vertex_handle vh = *adjVertI;
134 globalIndexing.toGlobal(vh->procIndex(), vh->index())
139 pointPoints[vit->index()].
transfer(indices);
167 typename T::Finite_vertices_iterator vit =
mesh.finite_vertices_begin();
168 vit !=
mesh.finite_vertices_end();
177 alignments[vit->index()] = vit->alignment();
195 typename T::Finite_vertices_iterator vit =
mesh.finite_vertices_begin();
196 vit !=
mesh.finite_vertices_end();
216 const label maxRefinementIterations,
217 const scalar defaultCellSize
220 for (
label iter = 0; iter < maxRefinementIterations; ++iter)
226 CellSizeDelaunay::Finite_cells_iterator cit =
227 mesh.finite_cells_begin();
228 cit !=
mesh.finite_cells_end();
232 const point newPoint =
237 cit->vertex(0)->point(),
238 cit->vertex(1)->point(),
239 cit->vertex(2)->point(),
240 cit->vertex(3)->point()
244 if (geometryToConformTo.
inside(newPoint))
246 ptsToInsert.
append(newPoint);
267 int main(
int argc,
char *argv[])
274 label maxRefinementIterations = 2;
275 label maxSmoothingIterations = 200;
276 scalar minResidual = 0;
277 scalar defaultCellSize = 0.001;
278 scalar nearFeatDistSqrCoeff = 1
e-8;
306 "cvSearchableSurfaces",
313 foamyHexMeshDict.subDict(
"geometry"),
314 foamyHexMeshDict.lookupOrDefault(
"singleRegionName",
true)
322 foamyHexMeshDict.subDict(
"surfaceConformation")
335 foamyHexMeshDict.subDict(
"backgroundMeshDecomposition")
354 geometryToConformTo.
geometry()[surfI];
375 nearFeatDistSqrCoeff,
386 geometryToConformTo.
features()[infoFeature];
394 pointAlignment() += norms[nI];
397 pointAlignment().normalize();
398 pointAlignment().orthogonalize();
405 nearFeatDistSqrCoeff,
413 geometryToConformTo.
features()[infoFeature];
421 pointAlignment() += norms[nI];
424 pointAlignment().normalize();
425 pointAlignment().orthogonalize();
433 surface.findNearest(ptField, distField, infoList);
436 surface.getNormal(infoList, normals);
438 pointAlignment.
set(
new triad(normals[0]));
446 CellSizeDelaunay::Vertex_handle vh =
mesh.
insert
451 Vb::vtInternalNearBoundary
457 CellSizeDelaunay::Vertex_handle vh =
mesh.
insert
462 Vb::vtInternalNearBoundary
474 maxRefinementIterations,
500 CellSizeDelaunay::Finite_vertices_iterator vit =
501 mesh.finite_vertices_begin();
502 vit !=
mesh.finite_vertices_end();
506 if (vit->nearBoundary())
508 fixedAlignments[vit->index()] = vit->alignment();
514 for (
label iter = 0; iter < maxSmoothingIterations; iter++)
516 Info<<
"Iteration " << iter;
518 meshDistributor().distribute(
points);
519 meshDistributor().distribute(alignments);
527 const labelList& pPoints = pointPoints[pI];
534 const triad& oldTriad = alignments[pI];
535 triad& newTriad = triadAv[pI];
538 const triad& fixedAlignment = fixedAlignments[pI];
540 forAll(pPoints, adjPointI)
542 const label adjPointIndex = pPoints[adjPointI];
548 triad tmpTriad = alignments[adjPointIndex];
552 if (tmpTriad.
set(dir))
554 tmpTriad[dir] *= (1.0/dist);
558 newTriad += tmpTriad;
567 forAll(fixedAlignment, dirI)
569 if (fixedAlignment.
set(dirI))
577 forAll(fixedAlignment, dirI)
579 if (fixedAlignment.
set(dirI))
581 newTriad.
align(fixedAlignment[dirI]);
585 else if (nFixed == 2)
587 forAll(fixedAlignment, dirI)
589 if (fixedAlignment.
set(dirI))
591 newTriad[dirI] = fixedAlignment[dirI];
601 else if (nFixed == 3)
603 forAll(fixedAlignment, dirI)
605 if (fixedAlignment.
set(dirI))
607 newTriad[dirI] = fixedAlignment[dirI];
621 scalar dotProd = (oldTriad[dir] & newTriad[dir]);
622 scalar
diff =
mag(dotProd) - 1.0;
631 alignments[pI] = triadAv[pI].sortxyz();
636 Info<<
", Residual = " << residual <<
endl;
638 if (residual <= minResidual)
646 OFstream str(runTime.path()/
"alignments.obj");
650 const triad& tri = alignments[pI];
687 filterFarPoints(
mesh, sizes)
700 filterFarPoints(
mesh, alignments)
705 alignmentsIO.write();
707 Info<<
nl <<
"ExecutionTime = " << runTime.elapsedCpuTime() <<
" s"
708 <<
" ClockTime = " << runTime.elapsedClockTime() <<
" s"
vectorField pointField
pointField is a vectorField.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Simple random number generator.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
label index() const
Return index.
const labelListList & featurePointNormals() const
Return the indices of the normals that are adjacent to the.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
static bool & parRun()
Is this a parallel run?
bool hit() const
Is there a hit.
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
void normalize()
Normalize each set axis vector to have a unit magnitude.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Field< triad > triadField
Specialisation of Field<T> for triad.
dimensioned< scalar > mag(const dimensioned< Type > &)
pointFromPoint topoint(const Point &P)
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & printInfo(Ostream &) const
Print information.
PrimitivePatch< face, List, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points)
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
Class containing processor-to-processor mapping information.
void set(T *)
Set pointer to that given.
void orthogonalize()
Orthogonalize this triad so that it is ortho-normal.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void align(const vector &v)
Align this triad with the given vector v.
int main(int argc, char *argv[])
const labelListList & edgeNormals() const
Return the indices of the normals that are adjacent to the.
const word & name() const
Return the name.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
const double e
Elementary charge.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Representation of a 3D Cartesian coordinate system as a Vector of vectors.
void setSize(const label)
Reset size of List.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
PointIndexHit< point > pointIndexHit
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Container for searchableSurfaces.
Store a background polyMesh to use for the decomposition of space and queries for parallel conformalV...
cachedRandom rndGen(label(0), -1)
bool set(const direction d) const
Is the vector in the direction d set.
A list of faces which address into the list of points.