Basic infrastructure for face motion controller - cell alignment and size

storage and lookup.
This commit is contained in:
graham
2009-07-10 20:23:33 +01:00
parent 02d4810d54
commit 8415cc64bc
5 changed files with 188 additions and 21 deletions

View File

@ -1989,18 +1989,20 @@ void Foam::conformalVoronoiMesh::conformToSurface()
void Foam::conformalVoronoiMesh::move()
{
timeCheck();
scalar relaxation = relaxationModel_->relaxation();
Info<< nl << " Relaxation = " << relaxation << endl;
timeCheck();
pointField dualVertices(number_of_cells());
pointField displacementAccumulator(startOfSurfacePointPairs_, vector::zero);
List<bool> pointToBeRetained(startOfSurfacePointPairs_, true);
DynamicList<point> newPointsToInsert;
label dualVerti = 0;
// Find the dual point of each tetrahedron and assign it an index.
@ -2029,10 +2031,146 @@ void Foam::conformalVoronoiMesh::move()
}
}
// setVertexAlignmentDirections();
dualVertices.setSize(dualVerti);
timeCheck();
Info<< nl << " Calculating target cell alignment and size" << endl;
for
(
Triangulation::Finite_vertices_iterator vit = finite_vertices_begin();
vit != finite_vertices_end();
vit++
)
{
if (vit->internalOrBoundaryPoint())
{
point pt(topoint(vit->point()));
vit->alignment() = requiredAlignment(pt);
vit->targetCellSize() = targetCellSize(pt);
}
}
timeCheck();
Info<< nl << " Looping over all dual faces" << endl;
vectorField cartesianDirections(3);
cartesianDirections[0] = vector(0,0,1);
cartesianDirections[1] = vector(0,1,0);
cartesianDirections[2] = vector(1,0,0);
for
(
Triangulation::Finite_edges_iterator eit = finite_edges_begin();
eit != finite_edges_end();
++eit
)
{
if
(
eit->first->vertex(eit->second)->internalOrBoundaryPoint()
&& eit->first->vertex(eit->third)->internalOrBoundaryPoint()
)
{
Cell_circulator ccStart = incident_cells(*eit);
Cell_circulator cc = ccStart;
DynamicList<label> verticesOnFace;
do
{
if (!is_infinite(cc))
{
if (cc->cellIndex() < 0)
{
FatalErrorIn("conformalVoronoiMesh::move")
<< "Dual face uses circumcenter defined by a "
<< " Delaunay tetrahedron with no internal "
<< "or boundary points."
<< exit(FatalError);
}
verticesOnFace.append(cc->cellIndex());
}
} while (++cc != ccStart);
verticesOnFace.shrink();
}
Cell_handle c = eit->first;
Vertex_handle vA = c->vertex(eit->second);
Vertex_handle vB = c->vertex(eit->third);
point dVA = topoint(vA->point());
point dVB = topoint(vB->point());
Field<vector> alignmentDirsA = vA->alignment() & cartesianDirections;
Field<vector> alignmentDirsB = vB->alignment() & cartesianDirections;
Field<vector> alignmentDirs(3);
forAll(alignmentDirsA, aA)
{
const vector& a(alignmentDirsA[aA]);
scalar maxDotProduct = 0.0;
forAll(alignmentDirsB, aB)
{
const vector& b(alignmentDirsB[aB]);
scalar dotProduct = a & b;
if (mag(dotProduct) > maxDotProduct)
{
maxDotProduct = mag(dotProduct);
alignmentDirs[aA] = a + sign(dotProduct)*b;
alignmentDirs[aA] /= mag(alignmentDirs[aA]);
}
}
}
vector rAB = dVA - dVB;
scalar rABMag = mag(rAB);
forAll(alignmentDirs, aD)
{
vector& alignmentDir = alignmentDirs[aD];
if ((rAB & alignmentDir) < 0)
{
// swap the direction of the alignment so that has the
// same sense as rAB
alignmentDir *= -1;
}
scalar alignmentDotProd = ((rAB/rABMag) & alignmentDir);
scalar targetCellSize =
0.5*(vA->targetCellSize() + vB->targetCellSize());
scalar targetFaceArea = sqr(targetCellSize);
if
(
alignmentDotProd
> cvMeshControls().cosAlignmentAcceptanceAngle()
)
{
alignmentDir *= 0.5*targetCellSize;
vector delta = alignmentDir - 0.5*rAB;
}
}
}
}

View File

@ -36,7 +36,7 @@ Description
#include <CGAL/Triangulation_3.h>
#include "List.H"
#include "vector.H"
#include "tensor.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -65,7 +65,9 @@ class indexedVertex
// Lowest numbered is inside one (master).
int type_;
Foam::List<Foam::vector> alignmentDirections_;
Foam::tensor alignment_;
Foam::scalar targetCellSize_;
public:
@ -94,7 +96,8 @@ public:
Vb(),
index_(INTERNAL_POINT),
type_(INTERNAL_POINT),
alignmentDirections_(0)
alignment_(),
targetCellSize_()
{}
indexedVertex(const Point& p)
@ -102,7 +105,8 @@ public:
Vb(p),
index_(INTERNAL_POINT),
type_(INTERNAL_POINT),
alignmentDirections_(0)
alignment_(),
targetCellSize_()
{}
indexedVertex(const Point& p, Cell_handle f)
@ -110,7 +114,8 @@ public:
Vb(f, p),
index_(INTERNAL_POINT),
type_(INTERNAL_POINT),
alignmentDirections_(0)
alignment_(),
targetCellSize_()
{}
indexedVertex(Cell_handle f)
@ -118,7 +123,8 @@ public:
Vb(f),
index_(INTERNAL_POINT),
type_(INTERNAL_POINT),
alignmentDirections_(0)
alignment_(),
targetCellSize_()
{}
@ -142,14 +148,24 @@ public:
return type_;
}
inline Foam::List<Foam::vector>& alignmentDirections()
inline Foam::tensor& alignment()
{
return alignmentDirections_;
return alignment_;
}
inline const Foam::List<Foam::vector>& alignmentDirections() const
inline const Foam::tensor& alignment() const
{
return alignmentDirections_;
return alignment_;
}
inline Foam::scalar& targetCellSize()
{
return targetCellSize_;
}
inline Foam::scalar targetCellSize() const
{
return targetCellSize_;
}
inline bool uninitialised() const

View File

@ -94,6 +94,14 @@ Foam::cvControls::cvControls
motionDict.lookup("alignmentSearchSpokes")
);
cosAlignmentAcceptanceAngle_ = cos
(
readScalar
(
motionDict.lookup("alignmentAcceptanceAngle")
)
);
const dictionary& insertionDict
(
motionDict.subDict("pointInsertionCriteria")
@ -109,11 +117,6 @@ Foam::cvControls::cvControls
insertionDict.lookup("faceAreaRatioCoeff")
);
alignmentAcceptanceAngle_ = readScalar
(
insertionDict.lookup("alignmentAcceptanceAngle")
);
const dictionary& removalDict
(
motionDict.subDict("pointRemovalCriteria")

View File

@ -108,14 +108,16 @@ class cvControls
// direction
label alignmentSearchSpokes_;
//- cosine of angle of alignment with required direction where a face
// will be accepted for rotation
scalar cosAlignmentAcceptanceAngle_;
// Point insertion criteria
scalar cellCentreInsertionDistCoeff_;
scalar faceAreaRatioCoeff_;
scalar alignmentAcceptanceAngle_;
// Point removal criteria
@ -184,6 +186,8 @@ public:
//- Return the number of alignmentSearchSpokes to use
inline label alignmentSearchSpokes() const;
//- Return the cosAlignmentAcceptanceAngle
inline scalar cosAlignmentAcceptanceAngle() const;
};

View File

@ -98,4 +98,10 @@ inline Foam::label Foam::cvControls::alignmentSearchSpokes() const
}
inline Foam::scalar Foam::cvControls::cosAlignmentAcceptanceAngle() const
{
return cosAlignmentAcceptanceAngle_;
}
// ************************************************************************* //