added solutionD and geometricD

This commit is contained in:
mattijs
2009-02-24 19:20:55 +00:00
parent c8944ce200
commit c49b302aa3
14 changed files with 411 additions and 236 deletions

View File

@ -23,16 +23,13 @@ Foam::label Foam::checkGeometry(const polyMesh& mesh, const bool allGeometry)
scalar minDistSqr = magSqr(1e-6 * globalBb.span()); scalar minDistSqr = magSqr(1e-6 * globalBb.span());
// Non-empty directions // Non-empty directions
const Vector<label> validDirs = (mesh.directions() + Vector<label>::one)/2; const Vector<label> validDirs = (mesh.geometricD() + Vector<label>::one)/2;
Info<< " Mesh (non-empty, non-wedge) directions " << validDirs << endl;
Info<< " Mesh (non-empty) directions " << validDirs << endl; const Vector<label> solDirs = (mesh.solutionD() + Vector<label>::one)/2;
Info<< " Mesh (non-empty) directions " << solDirs << endl;
scalar nGeomDims = mesh.nGeometricD(); if (mesh.nGeometricD() < 3)
Info<< " Mesh (non-empty, non-wedge) dimensions "
<< nGeomDims << endl;
if (nGeomDims < 3)
{ {
pointSet nonAlignedPoints(mesh, "nonAlignedEdges", mesh.nPoints()/100); pointSet nonAlignedPoints(mesh, "nonAlignedEdges", mesh.nPoints()/100);

View File

@ -54,40 +54,79 @@ void Foam::polyMesh::calcDirections() const
{ {
for (direction cmpt=0; cmpt<vector::nComponents; cmpt++) for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
{ {
directions_[cmpt] = 1; solutionD_[cmpt] = 1;
} }
label nEmptyPatches = 0; // Knock out empty and wedge directions. Note:they will be present on all
// domains.
vector dirVec = vector::zero; label nEmptyPatches = 0;
label nWedgePatches = 0;
vector emptyDirVec = vector::zero;
vector wedgeDirVec = vector::zero;
forAll(boundaryMesh(), patchi) forAll(boundaryMesh(), patchi)
{
if (isA<emptyPolyPatch>(boundaryMesh()[patchi]))
{ {
if (boundaryMesh()[patchi].size()) if (boundaryMesh()[patchi].size())
{
if (isA<emptyPolyPatch>(boundaryMesh()[patchi]))
{ {
nEmptyPatches++; nEmptyPatches++;
dirVec += sum(cmptMag(boundaryMesh()[patchi].faceAreas())); emptyDirVec += sum(cmptMag(boundaryMesh()[patchi].faceAreas()));
}
else if (isA<wedgePolyPatch>(boundaryMesh()[patchi]))
{
const wedgePolyPatch& wpp = refCast<const wedgePolyPatch>
(
boundaryMesh()[patchi]
);
nWedgePatches++;
wedgeDirVec += cmptMag(wpp.centreNormal());
} }
} }
} }
if (nEmptyPatches) if (nEmptyPatches)
{ {
reduce(dirVec, sumOp<vector>()); reduce(emptyDirVec, sumOp<vector>());
dirVec /= mag(dirVec); emptyDirVec /= mag(emptyDirVec);
for (direction cmpt=0; cmpt<vector::nComponents; cmpt++) for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
{ {
if (dirVec[cmpt] > 1e-6) if (emptyDirVec[cmpt] > 1e-6)
{ {
directions_[cmpt] = -1; solutionD_[cmpt] = -1;
} }
else else
{ {
directions_[cmpt] = 1; solutionD_[cmpt] = 1;
}
}
}
// Knock out wedge directions
geometricD_ = solutionD_;
if (nWedgePatches)
{
reduce(wedgeDirVec, sumOp<vector>());
wedgeDirVec /= mag(wedgeDirVec);
for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
{
if (wedgeDirVec[cmpt] > 1e-6)
{
geometricD_[cmpt] = -1;
}
else
{
geometricD_[cmpt] = 1;
} }
} }
} }
@ -163,7 +202,8 @@ Foam::polyMesh::polyMesh(const IOobject& io)
*this *this
), ),
bounds_(points_), bounds_(points_),
directions_(Vector<label>::zero), geometricD_(Vector<label>::zero),
solutionD_(Vector<label>::zero),
pointZones_ pointZones_
( (
IOobject IOobject
@ -350,7 +390,8 @@ Foam::polyMesh::polyMesh
0 0
), ),
bounds_(points_, syncPar), bounds_(points_, syncPar),
directions_(Vector<label>::zero), geometricD_(Vector<label>::zero),
solutionD_(Vector<label>::zero),
pointZones_ pointZones_
( (
IOobject IOobject
@ -505,7 +546,8 @@ Foam::polyMesh::polyMesh
0 0
), ),
bounds_(points_, syncPar), bounds_(points_, syncPar),
directions_(Vector<label>::zero), geometricD_(Vector<label>::zero),
solutionD_(Vector<label>::zero),
pointZones_ pointZones_
( (
IOobject IOobject
@ -766,44 +808,37 @@ const Foam::fileName& Foam::polyMesh::facesInstance() const
} }
const Foam::Vector<Foam::label>& Foam::polyMesh::directions() const const Foam::Vector<Foam::label>& Foam::polyMesh::geometricD() const
{ {
if (directions_.x() == 0) if (geometricD_.x() == 0)
{ {
calcDirections(); calcDirections();
} }
return directions_; return geometricD_;
} }
Foam::label Foam::polyMesh::nGeometricD() const Foam::label Foam::polyMesh::nGeometricD() const
{ {
label nWedges = 0; return cmptSum(geometricD() + Vector<label>::one)/2;
forAll(boundary_, patchi)
{
if (isA<wedgePolyPatch>(boundary_[patchi]))
{
nWedges++;
}
} }
if (nWedges != 0 && nWedges != 2 && nWedges != 4)
const Foam::Vector<Foam::label>& Foam::polyMesh::solutionD() const
{ {
FatalErrorIn("label polyMesh::nGeometricD() const") if (solutionD_.x() == 0)
<< "Number of wedge patches " << nWedges << " is incorrect, " {
"should be 0, 2 or 4" calcDirections();
<< exit(FatalError);
} }
return nSolutionD() - nWedges/2; return solutionD_;
} }
Foam::label Foam::polyMesh::nSolutionD() const Foam::label Foam::polyMesh::nSolutionD() const
{ {
return cmptSum(directions() + Vector<label>::one)/2; return cmptSum(solutionD() + Vector<label>::one)/2;
} }
@ -823,6 +858,10 @@ void Foam::polyMesh::addPatches
<< abort(FatalError); << abort(FatalError);
} }
// Reset valid directions
geometricD_ = Vector<label>::zero;
solutionD_ = Vector<label>::zero;
boundary_.setSize(p.size()); boundary_.setSize(p.size());
// Copy the patch pointers // Copy the patch pointers
@ -1037,6 +1076,10 @@ Foam::tmp<Foam::scalarField> Foam::polyMesh::movePoints
faceZones_.movePoints(points_); faceZones_.movePoints(points_);
cellZones_.movePoints(points_); cellZones_.movePoints(points_);
// Reset valid directions (could change with rotation)
geometricD_ = Vector<label>::zero;
solutionD_ = Vector<label>::zero;
// Hack until proper callbacks. Below are all the polyMeh MeshObjects with a // Hack until proper callbacks. Below are all the polyMeh MeshObjects with a
// movePoints function. // movePoints function.

View File

@ -120,9 +120,13 @@ private:
// Created from points on construction, updated when the mesh moves // Created from points on construction, updated when the mesh moves
boundBox bounds_; boundBox bounds_;
//- vector of non-constrained directions in mesh
// defined according to the presence of empty and wedge patches
mutable Vector<label> geometricD_;
//- vector of valid directions in mesh //- vector of valid directions in mesh
// defined according to the presence of empty patches // defined according to the presence of empty patches
mutable Vector<label> directions_; mutable Vector<label> solutionD_;
// Zoning information // Zoning information
@ -309,17 +313,22 @@ public:
return bounds_; return bounds_;
} }
//- Return the vector of valid directions in mesh. //- Return the vector of geometric directions in mesh.
// Defined according to the presence of empty patches. // Defined according to the presence of empty and wedge patches.
// 1 indicates valid direction and -1 an invalid direction. // 1 indicates unconstrained direction and -1 a constrained
const Vector<label>& directions() const; // direction.
const Vector<label>& geometricD() const;
//- Return the number of valid geometric dimensions in the mesh //- Return the number of valid geometric dimensions in the mesh
label nGeometricD() const; label nGeometricD() const;
//- Return the number of valid solution dimensions in the mesh. //- Return the vector of solved-for directions in mesh.
// For wedge cases this includes the circumferential direction // Differs from geometricD in that it includes for wedge cases
// in case of swirl. // the circumferential direction in case of swirl.
// 1 indicates valid direction and -1 an invalid direction.
const Vector<label>& solutionD() const;
//- Return the number of valid solved-for dimensions in the mesh
label nSolutionD() const; label nSolutionD() const;
//- Return point zone mesh //- Return point zone mesh

View File

@ -68,6 +68,10 @@ void Foam::polyMesh::clearGeom()
boundary_[patchI].clearGeom(); boundary_[patchI].clearGeom();
} }
// Reset valid directions (could change with rotation)
geometricD_ = Vector<label>::zero;
solutionD_ = Vector<label>::zero;
pointMesh::Delete(*this); pointMesh::Delete(*this);
} }
@ -87,6 +91,10 @@ void Foam::polyMesh::clearAddressing()
// recalculation // recalculation
deleteDemandDrivenData(globalMeshDataPtr_); deleteDemandDrivenData(globalMeshDataPtr_);
// Reset valid directions
geometricD_ = Vector<label>::zero;
solutionD_ = Vector<label>::zero;
pointMesh::Delete(*this); pointMesh::Delete(*this);
} }

View File

@ -217,7 +217,8 @@ Foam::polyMesh::polyMesh
boundaryFaces.size() + 1 // add room for a default patch boundaryFaces.size() + 1 // add room for a default patch
), ),
bounds_(points_, syncPar), bounds_(points_, syncPar),
directions_(Vector<label>::zero), geometricD_(Vector<label>::zero),
solutionD_(Vector<label>::zero),
pointZones_ pointZones_
( (
IOobject IOobject

View File

@ -68,6 +68,11 @@ void Foam::polyMesh::updateMesh(const mapPolyMesh& mpm)
newMotionPoints.map(oldMotionPoints, mpm.pointMap()); newMotionPoints.map(oldMotionPoints, mpm.pointMap());
} }
// Reset valid directions (could change by faces put into empty patches)
geometricD_ = Vector<label>::zero;
solutionD_ = Vector<label>::zero;
// Hack until proper callbacks. Below are all the polyMesh-MeshObjects. // Hack until proper callbacks. Below are all the polyMesh-MeshObjects.
// pointMesh // pointMesh

View File

@ -140,11 +140,11 @@ void Foam::quadraticFitSnGradData::findFaceDirs
#ifndef SPHERICAL_GEOMETRY #ifndef SPHERICAL_GEOMETRY
if (mesh.nGeometricD() <= 2) // find the normal direcion if (mesh.nGeometricD() <= 2) // find the normal direcion
{ {
if (mesh.directions()[0] == -1) if (mesh.geometricD()[0] == -1)
{ {
kdir = vector(1, 0, 0); kdir = vector(1, 0, 0);
} }
else if (mesh.directions()[1] == -1) else if (mesh.geometricD()[1] == -1)
{ {
kdir = vector(0, 1, 0); kdir = vector(0, 1, 0);
} }
@ -153,7 +153,7 @@ void Foam::quadraticFitSnGradData::findFaceDirs
kdir = vector(0, 0, 1); kdir = vector(0, 0, 1);
} }
} }
else // 3D so find a direction in the place of the face else // 3D so find a direction in the plane of the face
{ {
const face& f = mesh.faces()[faci]; const face& f = mesh.faces()[faci];
kdir = mesh.points()[f[0]] - mesh.points()[f[1]]; kdir = mesh.points()[f[0]] - mesh.points()[f[1]];

View File

@ -713,7 +713,7 @@ Foam::fvMatrix<Type>::H() const
( (
pow pow
( (
psi_.mesh().directions(), psi_.mesh().solutionD(),
pTraits<typename powProduct<Vector<label>, Type::rank>::type>::zero pTraits<typename powProduct<Vector<label>, Type::rank>::type>::zero
) )
); );

View File

@ -82,7 +82,7 @@ Foam::lduMatrix::solverPerformance Foam::fvMatrix<Type>::solve
( (
pow pow
( (
psi_.mesh().directions(), psi_.mesh().solutionD(),
pTraits<typename powProduct<Vector<label>, Type::rank>::type>::zero pTraits<typename powProduct<Vector<label>, Type::rank>::type>::zero
) )
); );

View File

@ -83,11 +83,11 @@ void Foam::FitData<FitDataType, ExtendedStencil, Polynomial>::findFaceDirs
# ifndef SPHERICAL_GEOMETRY # ifndef SPHERICAL_GEOMETRY
if (mesh.nGeometricD() <= 2) // find the normal direction if (mesh.nGeometricD() <= 2) // find the normal direction
{ {
if (mesh.directions()[0] == -1) if (mesh.geometricD()[0] == -1)
{ {
kdir = vector(1, 0, 0); kdir = vector(1, 0, 0);
} }
else if (mesh.directions()[1] == -1) else if (mesh.geometricD()[1] == -1)
{ {
kdir = vector(0, 1, 0); kdir = vector(0, 1, 0);
} }
@ -115,7 +115,7 @@ void Foam::FitData<FitDataType, ExtendedStencil, Polynomial>::findFaceDirs
if (magk < SMALL) if (magk < SMALL)
{ {
FatalErrorIn("findFaceDirs") << " calculated kdir = zero" FatalErrorIn("findFaceDirs(..)") << " calculated kdir = zero"
<< exit(FatalError); << exit(FatalError);
} }
else else

View File

@ -25,7 +25,7 @@ License
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
#include "meshTools.H" #include "meshTools.H"
#include "primitiveMesh.H" #include "polyMesh.H"
#include "hexMatcher.H" #include "hexMatcher.H"
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
@ -635,6 +635,103 @@ Foam::label Foam::meshTools::walkFace
} }
void Foam::meshTools::constrainToMeshCentre(const polyMesh& mesh, point& pt)
{
const Vector<label>& dirs = mesh.geometricD();
const point& min = mesh.bounds().min();
const point& max = mesh.bounds().max();
for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
{
if (dirs[cmpt] == -1)
{
pt[cmpt] = 0.5*(min[cmpt]+max[cmpt]);
}
}
}
void Foam::meshTools::constrainToMeshCentre
(
const polyMesh& mesh,
pointField& pts
)
{
const Vector<label>& dirs = mesh.geometricD();
const point& min = mesh.bounds().min();
const point& max = mesh.bounds().max();
bool isConstrained = false;
for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
{
if (dirs[cmpt] == -1)
{
isConstrained = true;
break;
}
}
if (isConstrained)
{
forAll(pts, i)
{
for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
{
if (dirs[cmpt] == -1)
{
pts[i][cmpt] = 0.5*(min[cmpt]+max[cmpt]);
}
}
}
}
}
//- Set the constrained components of directions/velocity to zero
void Foam::meshTools::constrainDirection(const polyMesh& mesh, vector& d)
{
const Vector<label>& dirs = mesh.geometricD();
for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
{
if (dirs[cmpt] == -1)
{
d[cmpt] = 0.0;
}
}
}
void Foam::meshTools::constrainDirection(const polyMesh& mesh, vectorField& d)
{
const Vector<label>& dirs = mesh.geometricD();
bool isConstrained = false;
for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
{
if (dirs[cmpt] == -1)
{
isConstrained = true;
break;
}
}
if (isConstrained)
{
forAll(d, i)
{
for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
{
if (dirs[cmpt] == -1)
{
d[i][cmpt] = 0.0;
}
}
}
}
}
void Foam::meshTools::getParallelEdges void Foam::meshTools::getParallelEdges
( (
const primitiveMesh& mesh, const primitiveMesh& mesh,

View File

@ -50,7 +50,7 @@ namespace Foam
{ {
class primitiveMesh; class primitiveMesh;
class polyMesh;
/*---------------------------------------------------------------------------*\ /*---------------------------------------------------------------------------*\
Namespace meshTools Declaration Namespace meshTools Declaration
@ -81,6 +81,8 @@ namespace meshTools
static const label pXpYpZMask = 1 << pXpYpZ; static const label pXpYpZMask = 1 << pXpYpZ;
// Normal handling
//- Check if n is in same direction as normals of all faceLabels //- Check if n is in same direction as normals of all faceLabels
bool visNormal bool visNormal
( (
@ -96,6 +98,9 @@ namespace meshTools
//- Normalized edge vector //- Normalized edge vector
vector normEdgeVec(const primitiveMesh&, const label edgeI); vector normEdgeVec(const primitiveMesh&, const label edgeI);
// OBJ writing
//- Write obj representation of point //- Write obj representation of point
void writeOBJ void writeOBJ
( (
@ -103,7 +108,6 @@ namespace meshTools
const point& pt const point& pt
); );
//- Write obj representation of faces subset //- Write obj representation of faces subset
void writeOBJ void writeOBJ
( (
@ -131,6 +135,9 @@ namespace meshTools
const labelList& cellLabels const labelList& cellLabels
); );
// Cell/face/edge walking
//- Is edge used by cell //- Is edge used by cell
bool edgeOnCell bool edgeOnCell
( (
@ -239,9 +246,18 @@ namespace meshTools
); );
// // Constraints on position
//- Set the constrained components of position to mesh centre
void constrainToMeshCentre(const polyMesh&, point&);
void constrainToMeshCentre(const polyMesh&, pointField&);
//- Set the constrained components of directions/velocity to zero
void constrainDirection(const polyMesh&, vector&);
void constrainDirection(const polyMesh&, vectorField&);
// Hex only functionality. // Hex only functionality.
//
//- Given edge on hex find other 'parallel', non-connected edges. //- Given edge on hex find other 'parallel', non-connected edges.
void getParallelEdges void getParallelEdges

View File

@ -42,8 +42,7 @@ addToRunTimeSelectionTable(LESdelta, cubeRootVolDelta, dictionary);
void cubeRootVolDelta::calcDelta() void cubeRootVolDelta::calcDelta()
{ {
const Vector<label>& directions = mesh().directions(); label nD = mesh().nGeometricD();
label nD = (directions.nComponents + cmptSum(directions))/2;
if (nD == 3) if (nD == 3)
{ {
@ -55,14 +54,15 @@ void cubeRootVolDelta::calcDelta()
<< "Case is 2D, LES is not strictly applicable\n" << "Case is 2D, LES is not strictly applicable\n"
<< endl; << endl;
const Vector<label>& directions = mesh().geometricD();
scalar thickness = 0.0; scalar thickness = 0.0;
for (direction dir=0; dir<directions.nComponents; dir++) for (direction dir=0; dir<directions.nComponents; dir++)
{ {
if (directions[dir] == -1) if (directions[dir] == -1)
{ {
boundBox bb(mesh().points(), false); thickness = mesh().bounds().span()[dir];
break;
thickness = bb.span()[dir];
} }
} }

View File

@ -41,8 +41,7 @@ namespace Foam
void Foam::IDDESDelta::calcDelta() void Foam::IDDESDelta::calcDelta()
{ {
const Vector<label>& directions = mesh().directions(); label nD = mesh().nGeometricD();
label nD = (directions.nComponents + cmptSum(directions))/2;
// initialise hwn as wall distance // initialise hwn as wall distance
volScalarField hwn = wallDist(mesh()).y(); volScalarField hwn = wallDist(mesh()).y();