smesh/src/StdMeshers/StdMeshers_Quadrangle_2D.hxx

253 lines
9.2 KiB
C++
Raw Normal View History

2013-04-01 19:05:47 +06:00
// Copyright (C) 2007-2013 CEA/DEN, EDF R&D, OPEN CASCADE
2009-02-17 10:27:49 +05:00
//
2012-08-09 16:03:55 +06:00
// Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
// CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
2004-06-18 14:34:31 +06:00
//
2012-08-09 16:03:55 +06:00
// This library is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 2.1 of the License.
2004-06-18 14:34:31 +06:00
//
2012-08-09 16:03:55 +06:00
// This library is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
// Lesser General Public License for more details.
2004-06-18 14:34:31 +06:00
//
2012-08-09 16:03:55 +06:00
// You should have received a copy of the GNU Lesser General Public
// License along with this library; if not, write to the Free Software
// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
2004-06-18 14:34:31 +06:00
//
2012-08-09 16:03:55 +06:00
// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
2009-02-17 10:27:49 +05:00
//
2004-06-18 14:34:31 +06:00
// File : StdMeshers_Quadrangle_2D.hxx
// Moved here from SMESH_Quadrangle_2D.hxx
// Author : Paul RASCLE, EDF
// Module : SMESH
2012-08-09 16:03:55 +06:00
2004-06-18 14:34:31 +06:00
#ifndef _SMESH_QUADRANGLE_2D_HXX_
#define _SMESH_QUADRANGLE_2D_HXX_
2012-12-13 17:41:29 +06:00
#include "SMESH_Algo.hxx"
#include "SMESH_ProxyMesh.hxx"
#include "SMESH_StdMeshers.hxx"
2014-01-15 15:41:17 +06:00
#include "StdMeshers_FaceSide.hxx"
2012-08-09 16:03:55 +06:00
#include "StdMeshers_QuadrangleParams.hxx"
#include <TopoDS_Face.hxx>
2014-01-15 15:41:17 +06:00
#include <Bnd_B2d.hxx>
2012-08-09 16:03:55 +06:00
2012-12-13 17:41:29 +06:00
class SMDS_MeshNode;
class SMESH_Mesh;
class SMESH_MesherHelper;
2012-12-13 17:41:29 +06:00
class SMESH_ProxyMesh;
struct uvPtStruct;
2004-12-01 15:48:31 +05:00
2013-02-12 20:37:44 +06:00
enum TSideID { QUAD_BOTTOM_SIDE=0, QUAD_RIGHT_SIDE, QUAD_TOP_SIDE, QUAD_LEFT_SIDE, NB_QUAD_SIDES };
2004-06-18 14:34:31 +06:00
typedef uvPtStruct UVPtStruct;
2014-01-15 15:41:17 +06:00
struct FaceQuadStruct
2004-06-18 14:34:31 +06:00
{
2014-01-15 15:41:17 +06:00
struct Side // a side of FaceQuadStruct
{
struct Contact // contact of two sides
{
int point; // index of a grid point of this side where two sides meat
Side* other_side;
int other_point;
};
StdMeshers_FaceSidePtr grid;
int from, to; // indices of grid points used by the quad
int di; // +1 or -1 depending on IsReversed()
2014-01-15 15:41:17 +06:00
std::set<int> forced_nodes; // indices of forced grid points
std::vector<Contact> contacts; // contacts with sides of other quads
int nbNodeOut; // nb of missing nodes on an opposite shorter side
Side(StdMeshers_FaceSidePtr theGrid = StdMeshers_FaceSidePtr());
Side& operator=(const Side& otherSide);
operator StdMeshers_FaceSidePtr() { return grid; }
operator const StdMeshers_FaceSidePtr() const { return grid; }
void AddContact( int ip, Side* side, int iop );
int ToSideIndex( int quadNodeIndex ) const;
int ToQuadIndex( int sideNodeIndex ) const;
bool IsForced( int nodeIndex ) const;
bool IsReversed() const { return nbNodeOut ? false : to < from; }
bool Reverse(bool keepGrid);
2014-01-15 15:41:17 +06:00
int NbPoints() const { return Abs( to - from ); }
double Param( int nodeIndex ) const;
double Length( int from=-1, int to=-1) const;
gp_XY Value2d( double x ) const;
const UVPtStruct& First() const { return GetUVPtStruct()[ from ]; }
const UVPtStruct& Last() const {
return GetUVPtStruct()[ to-nbNodeOut-(IsReversed() ? -1 : +1)];
}
2014-01-15 15:41:17 +06:00
// some sortcuts
const vector<UVPtStruct>& GetUVPtStruct(bool isXConst=0, double constValue=0) const
{ return nbNodeOut ?
grid->SimulateUVPtStruct( NbPoints()-nbNodeOut-1, isXConst, constValue ) :
grid->GetUVPtStruct( isXConst, constValue );
}
};
struct SideIterator // iterator on UVPtStruct of a Side
{
const UVPtStruct *uvPtr, *uvEnd;
int dPtr, counter;
SideIterator(): uvPtr(0), uvEnd(0), dPtr(0), counter(0) {}
void Init( const Side& side ) {
dPtr = counter = 0;
uvPtr = uvEnd = 0;
if ( side.NbPoints() > 0 ) {
uvPtr = & side.First();
uvEnd = & side.Last();
dPtr = ( uvEnd > uvPtr ) ? +1 : -1;
uvEnd += dPtr;
}
}
bool More() const { return uvPtr != uvEnd; }
void Next() { uvPtr += dPtr; ++counter; }
UVPtStruct& UVPt() const { return (UVPtStruct&) *uvPtr; }
UVPtStruct& operator[](int i) { return (UVPtStruct&) uvPtr[ i*dPtr]; }
int Count() const { return counter; }
};
2014-01-15 15:41:17 +06:00
std::vector< Side > side;
std::vector< UVPtStruct> uv_grid;
int iSize, jSize;
TopoDS_Face face;
Bnd_B2d uv_box;
std::string name; // to ease debugging
2014-01-15 15:41:17 +06:00
FaceQuadStruct ( const TopoDS_Face& F = TopoDS_Face(), const std::string& nm="main" );
2014-01-15 15:41:17 +06:00
UVPtStruct& UVPt( int i, int j ) { return uv_grid[ i + j * iSize ]; }
void shift ( size_t nb, bool keepUnitOri, bool keepGrid=false );
2014-01-15 15:41:17 +06:00
int & nbNodeOut( int iSide ) { return side[ iSide ].nbNodeOut; }
bool findCell ( const gp_XY& uv, int & i, int & j );
bool isNear ( const gp_XY& uv, int & i, int & j, int nbLoops=1 );
bool isEqual ( const gp_XY& uv, int i, int j );
void normPa2IJ( double x, double y, int & i, int & j );
void updateUV ( const gp_XY& uv, int i, int j, bool isVertical );
typedef boost::shared_ptr<FaceQuadStruct> Ptr;
};
2004-06-18 14:34:31 +06:00
class STDMESHERS_EXPORT StdMeshers_Quadrangle_2D: public SMESH_2D_Algo
2004-06-18 14:34:31 +06:00
{
public:
StdMeshers_Quadrangle_2D(int hypId, int studyId, SMESH_Gen* gen);
virtual ~StdMeshers_Quadrangle_2D();
2013-07-31 14:12:54 +06:00
virtual bool CheckHypothesis(SMESH_Mesh& aMesh,
2004-06-18 14:34:31 +06:00
const TopoDS_Shape& aShape,
2013-07-31 14:12:54 +06:00
Hypothesis_Status& aStatus);
2004-06-18 14:34:31 +06:00
2013-07-31 14:12:54 +06:00
virtual bool Compute(SMESH_Mesh& aMesh,
2012-08-09 16:03:55 +06:00
const TopoDS_Shape& aShape);
2013-07-31 14:12:54 +06:00
virtual bool Evaluate(SMESH_Mesh & aMesh,
const TopoDS_Shape & aShape,
MapShapeNbElems& aResMap);
2004-06-18 14:34:31 +06:00
2013-07-31 14:12:54 +06:00
FaceQuadStruct::Ptr CheckAnd2Dcompute(SMESH_Mesh& aMesh,
2013-02-12 20:37:44 +06:00
const TopoDS_Shape& aShape,
2013-07-31 14:12:54 +06:00
const bool CreateQuadratic);
2004-06-18 14:34:31 +06:00
2013-07-31 14:12:54 +06:00
FaceQuadStruct::Ptr CheckNbEdges(SMESH_Mesh& aMesh,
const TopoDS_Shape& aShape,
const bool considerMesh=false);
2012-08-09 16:03:55 +06:00
protected:
bool checkNbEdgesForEvaluate(SMESH_Mesh& aMesh,
2012-08-09 16:03:55 +06:00
const TopoDS_Shape & aShape,
MapShapeNbElems& aResMap,
std::vector<int>& aNbNodes,
bool& IsQuadratic);
2014-01-15 15:41:17 +06:00
bool setNormalizedGrid(FaceQuadStruct::Ptr quad);
2009-02-17 10:27:49 +05:00
2014-01-15 15:41:17 +06:00
void splitQuadFace(SMESHDS_Mesh * theMeshDS,
const int theFaceID,
const SMDS_MeshNode* theNode1,
const SMDS_MeshNode* theNode2,
const SMDS_MeshNode* theNode3,
const SMDS_MeshNode* theNode4);
bool computeQuadDominant(SMESH_Mesh& aMesh,
const TopoDS_Face& aFace);
2004-06-18 14:34:31 +06:00
bool computeQuadDominant(SMESH_Mesh& aMesh,
const TopoDS_Face& aFace,
FaceQuadStruct::Ptr quad);
bool computeQuadPref(SMESH_Mesh& aMesh,
const TopoDS_Face& aFace,
2013-02-12 20:37:44 +06:00
FaceQuadStruct::Ptr quad);
2014-01-15 15:41:17 +06:00
bool computeTriangles(SMESH_Mesh& aMesh,
const TopoDS_Face& aFace,
FaceQuadStruct::Ptr quad);
bool evaluateQuadPref(SMESH_Mesh& aMesh,
2012-08-09 16:03:55 +06:00
const TopoDS_Shape& aShape,
2013-02-12 20:37:44 +06:00
std::vector<int>& aNbNodes,
MapShapeNbElems& aResMap,
bool isQuadratic);
2012-08-09 16:03:55 +06:00
bool computeReduced (SMESH_Mesh& aMesh,
const TopoDS_Face& aFace,
2013-02-12 20:37:44 +06:00
FaceQuadStruct::Ptr quad);
2012-08-09 16:03:55 +06:00
void updateDegenUV(FaceQuadStruct::Ptr quad);
void smooth (FaceQuadStruct::Ptr quad);
2004-06-18 14:34:31 +06:00
int getCorners(const TopoDS_Face& theFace,
SMESH_Mesh & theMesh,
std::list<TopoDS_Edge>& theWire,
std::vector<TopoDS_Vertex>& theVertices,
int & theNbDegenEdges,
const bool considerMesh);
2014-01-15 15:41:17 +06:00
bool getEnforcedUV();
bool addEnforcedNodes();
int splitQuad(FaceQuadStruct::Ptr quad, int i, int j);
void shiftQuad(FaceQuadStruct::Ptr& quad, const int num );
2014-01-15 15:41:17 +06:00
typedef std::map< StdMeshers_FaceSidePtr, std::vector< FaceQuadStruct::Ptr > > TQuadsBySide;
void updateSideUV( FaceQuadStruct::Side& side,
int iForced,
const TQuadsBySide& quads,
int * iNext=NULL);
// Fields
bool myQuadranglePreference;
2009-02-17 10:27:49 +05:00
bool myTrianglePreference;
2012-12-13 17:41:29 +06:00
int myTriaVertexID;
bool myNeedSmooth;
2014-01-15 15:41:17 +06:00
const StdMeshers_QuadrangleParams* myParams;
StdMeshers_QuadType myQuadType;
SMESH_MesherHelper* myHelper;
SMESH_ProxyMesh::Ptr myProxyMesh;
std::list< FaceQuadStruct::Ptr > myQuadList;
struct ForcedPoint
{
gp_XY uv;
gp_XYZ xyz;
TopoDS_Vertex vertex;
double U() const { return uv.X(); }
double V() const { return uv.Y(); }
operator const gp_XY& () { return uv; }
};
std::vector< ForcedPoint > myForcedPnts;
2004-06-18 14:34:31 +06:00
};
#endif