mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
Storing and reinserting surface conformation instead of recalculating at each
iteration. Rebuilding every 10 steps as a hard-coded experiment.
This commit is contained in:
@ -63,7 +63,7 @@ int main(int argc, char *argv[])
|
|||||||
|
|
||||||
mesh.move();
|
mesh.move();
|
||||||
|
|
||||||
mesh.conformToSurface();
|
// mesh.conformToSurface();
|
||||||
}
|
}
|
||||||
|
|
||||||
mesh.writeMesh();
|
mesh.writeMesh();
|
||||||
|
|||||||
@ -74,6 +74,7 @@ Foam::conformalVoronoiMesh::conformalVoronoiMesh
|
|||||||
storedSizes_(),
|
storedSizes_(),
|
||||||
storedAlignments_(),
|
storedAlignments_(),
|
||||||
sizeAndAlignmentTree_(),
|
sizeAndAlignmentTree_(),
|
||||||
|
surfaceConformationVertices_(),
|
||||||
initialPointsMethod_
|
initialPointsMethod_
|
||||||
(
|
(
|
||||||
initialPointsMethod::New
|
initialPointsMethod::New
|
||||||
@ -105,6 +106,8 @@ Foam::conformalVoronoiMesh::conformalVoronoiMesh
|
|||||||
|
|
||||||
conformToSurface();
|
conformToSurface();
|
||||||
|
|
||||||
|
storeSurfaceConformation();
|
||||||
|
|
||||||
if(cvMeshControls().objOutput())
|
if(cvMeshControls().objOutput())
|
||||||
{
|
{
|
||||||
writePoints("allInitialPoints.obj", false);
|
writePoints("allInitialPoints.obj", false);
|
||||||
@ -1314,6 +1317,61 @@ void Foam::conformalVoronoiMesh::addSurfaceAndEdgeHits
|
|||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
void Foam::conformalVoronoiMesh::storeSurfaceConformation()
|
||||||
|
{
|
||||||
|
Info<< nl << " Storing surface conformation." << endl;
|
||||||
|
|
||||||
|
surfaceConformationVertices_.setSize
|
||||||
|
(
|
||||||
|
number_of_vertices() - startOfSurfacePointPairs_
|
||||||
|
);
|
||||||
|
|
||||||
|
label surfPtI = 0;
|
||||||
|
|
||||||
|
for
|
||||||
|
(
|
||||||
|
Triangulation::Finite_vertices_iterator vit = finite_vertices_begin();
|
||||||
|
vit != finite_vertices_end();
|
||||||
|
vit++
|
||||||
|
)
|
||||||
|
{
|
||||||
|
if (vit->index() >= startOfSurfacePointPairs_)
|
||||||
|
{
|
||||||
|
if (!vit->pairPoint())
|
||||||
|
{
|
||||||
|
FatalErrorIn("storeSurfaceConformation()")
|
||||||
|
<< "Trying to store a vertex that is not a surface point"
|
||||||
|
<< exit(FatalError);
|
||||||
|
}
|
||||||
|
|
||||||
|
surfaceConformationVertices_[surfPtI] = Vb(vit->point());
|
||||||
|
|
||||||
|
surfaceConformationVertices_[surfPtI].index() =
|
||||||
|
vit->index() - startOfSurfacePointPairs_;
|
||||||
|
|
||||||
|
surfaceConformationVertices_[surfPtI].type() =
|
||||||
|
vit->type() - startOfSurfacePointPairs_;
|
||||||
|
|
||||||
|
surfPtI++;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
Info<< " Stored " << surfaceConformationVertices_.size()
|
||||||
|
<< " vertices" << endl;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
void Foam::conformalVoronoiMesh::reinsertSurfaceConformation()
|
||||||
|
{
|
||||||
|
Info<< nl << " Reinserting stored surface conformation" << endl;
|
||||||
|
|
||||||
|
forAll(surfaceConformationVertices_, v)
|
||||||
|
{
|
||||||
|
insertVb(surfaceConformationVertices_[v], startOfSurfacePointPairs_);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
void Foam::conformalVoronoiMesh::calcDualMesh
|
void Foam::conformalVoronoiMesh::calcDualMesh
|
||||||
(
|
(
|
||||||
pointField& points,
|
pointField& points,
|
||||||
@ -2559,8 +2617,28 @@ void Foam::conformalVoronoiMesh::move()
|
|||||||
reinsertFeaturePoints();
|
reinsertFeaturePoints();
|
||||||
startOfInternalPoints_ = number_of_vertices();
|
startOfInternalPoints_ = number_of_vertices();
|
||||||
|
|
||||||
|
timeCheck();
|
||||||
|
|
||||||
Info<< nl << " Reinserting entire tessellation" << endl;
|
Info<< nl << " Reinserting entire tessellation" << endl;
|
||||||
insertPoints(pointsToInsert);
|
insertPoints(pointsToInsert);
|
||||||
|
|
||||||
|
timeCheck();
|
||||||
|
|
||||||
|
if (runTime_.timeIndex() % 10 == 0)
|
||||||
|
{
|
||||||
|
Info<< nl << " Rebuilding surface conformation "
|
||||||
|
<< "HARD CODED TO EVERY 10 STEPS" << endl;
|
||||||
|
|
||||||
|
conformToSurface();
|
||||||
|
storeSurfaceConformation();
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
startOfSurfacePointPairs_ = number_of_vertices();
|
||||||
|
reinsertSurfaceConformation();
|
||||||
|
}
|
||||||
|
|
||||||
|
timeCheck();
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
|||||||
@ -132,6 +132,10 @@ class conformalVoronoiMesh
|
|||||||
//- Search tree for size and alignment lookup points
|
//- Search tree for size and alignment lookup points
|
||||||
mutable autoPtr<indexedOctree<treeDataPoint> > sizeAndAlignmentTree_;
|
mutable autoPtr<indexedOctree<treeDataPoint> > sizeAndAlignmentTree_;
|
||||||
|
|
||||||
|
//- Store the surface and feature edge conformation locations to be
|
||||||
|
// reinserted
|
||||||
|
List<Vb> surfaceConformationVertices_;
|
||||||
|
|
||||||
//- Method for inserting initial points. Runtime selectable.
|
//- Method for inserting initial points. Runtime selectable.
|
||||||
autoPtr<initialPointsMethod> initialPointsMethod_;
|
autoPtr<initialPointsMethod> initialPointsMethod_;
|
||||||
|
|
||||||
@ -200,8 +204,9 @@ class conformalVoronoiMesh
|
|||||||
const vector& n
|
const vector& n
|
||||||
);
|
);
|
||||||
|
|
||||||
//- Insert a Vb (a typedef of CGAL::indexedVertex<K>)
|
//- Insert a Vb (a typedef of CGAL::indexedVertex<K>) with the
|
||||||
inline void insertVb(const Vb& v);
|
// possibility of an offset for the index and the type
|
||||||
|
inline void insertVb(const Vb& v, label offset = 0);
|
||||||
|
|
||||||
//- Insert pairs of points on the surface with the given normals, at the
|
//- Insert pairs of points on the surface with the given normals, at the
|
||||||
// specified spacing
|
// specified spacing
|
||||||
@ -362,6 +367,14 @@ class conformalVoronoiMesh
|
|||||||
autoPtr<indexedOctree<treeDataPoint> >& edgeLocationTree
|
autoPtr<indexedOctree<treeDataPoint> >& edgeLocationTree
|
||||||
) const;
|
) const;
|
||||||
|
|
||||||
|
//- Store the surface conformation with the indices offset to be
|
||||||
|
// relative to zero
|
||||||
|
void storeSurfaceConformation();
|
||||||
|
|
||||||
|
//- Reinsert the surface conformation re-offsetting indices to be
|
||||||
|
// relative to new number of internal vertices
|
||||||
|
void reinsertSurfaceConformation();
|
||||||
|
|
||||||
//- Dual calculation
|
//- Dual calculation
|
||||||
void calcDualMesh
|
void calcDualMesh
|
||||||
(
|
(
|
||||||
|
|||||||
@ -202,7 +202,7 @@ inline void Foam::conformalVoronoiMesh::insertPointPair
|
|||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
inline void Foam::conformalVoronoiMesh::insertVb(const Vb& v)
|
inline void Foam::conformalVoronoiMesh::insertVb(const Vb& v, label offset)
|
||||||
{
|
{
|
||||||
const Point& Pt(v.point());
|
const Point& Pt(v.point());
|
||||||
|
|
||||||
@ -218,8 +218,16 @@ inline void Foam::conformalVoronoiMesh::insertVb(const Vb& v)
|
|||||||
}
|
}
|
||||||
else
|
else
|
||||||
{
|
{
|
||||||
vh->index() = v.index();
|
vh->index() = v.index() + offset;
|
||||||
vh->type() = v.type();
|
|
||||||
|
if (v.pairPoint())
|
||||||
|
{
|
||||||
|
vh->type() = v.type() + offset;
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
vh->type() = v.type();
|
||||||
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|||||||
Reference in New Issue
Block a user