Salome HOME
Merge from BR_imps_2013 14/01/2014
[modules/smesh.git] / src / StdMeshers / StdMeshers_Quadrangle_2D.hxx
1 // Copyright (C) 2007-2013  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License.
10 //
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
14 // Lesser General Public License for more details.
15 //
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
19 //
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22 //  File   : StdMeshers_Quadrangle_2D.hxx
23 //           Moved here from SMESH_Quadrangle_2D.hxx
24 //  Author : Paul RASCLE, EDF
25 //  Module : SMESH
26
27 #ifndef _SMESH_QUADRANGLE_2D_HXX_
28 #define _SMESH_QUADRANGLE_2D_HXX_
29
30 #include "SMESH_Algo.hxx"
31 #include "SMESH_ProxyMesh.hxx"
32 #include "SMESH_StdMeshers.hxx"
33 #include "StdMeshers_FaceSide.hxx"
34 #include "StdMeshers_QuadrangleParams.hxx"
35
36 #include <TopoDS_Face.hxx>
37 #include <Bnd_B2d.hxx>
38
39 class SMDS_MeshNode;
40 class SMESH_Mesh;
41 class SMESH_MesherHelper;
42 class SMESH_ProxyMesh;
43 struct uvPtStruct;
44
45
46 enum TSideID { QUAD_BOTTOM_SIDE=0, QUAD_RIGHT_SIDE, QUAD_TOP_SIDE, QUAD_LEFT_SIDE, NB_QUAD_SIDES };
47
48 typedef uvPtStruct UVPtStruct;
49 struct FaceQuadStruct
50 {
51   struct Side // a side of FaceQuadStruct
52   {
53     struct Contact // contact of two sides
54     {
55       int   point; // index of a grid point of this side where two sides meat
56       Side* other_side;
57       int   other_point;
58     };
59     StdMeshers_FaceSidePtr grid;
60     int                    from, to;     // indices of grid points used by the quad
61     std::set<int>          forced_nodes; // indices of forced grid points
62     std::vector<Contact>   contacts;     // contacts with sides of other quads
63     int                    nbNodeOut;    // nb of missing nodes on an opposite shorter side
64
65     Side(StdMeshers_FaceSidePtr theGrid = StdMeshers_FaceSidePtr());
66     Side& operator=(const Side& otherSide);
67     operator StdMeshers_FaceSidePtr() { return grid; }
68     operator const StdMeshers_FaceSidePtr() const { return grid; }
69     void AddContact( int ip, Side* side, int iop );
70     int  ToSideIndex( int quadNodeIndex ) const;
71     int  ToQuadIndex( int sideNodeIndex ) const;
72     bool IsForced( int nodeIndex ) const;
73     bool IsReversed() const { return nbNodeOut ? false : to < from; }
74     int  NbPoints() const { return Abs( to - from ); }
75     double Param( int nodeIndex ) const;
76     double Length( int from=-1, int to=-1) const;
77     gp_XY Value2d( double x ) const;
78     // some sortcuts
79     const vector<UVPtStruct>& GetUVPtStruct(bool isXConst=0, double constValue=0) const
80     { return nbNodeOut ?
81         grid->SimulateUVPtStruct( NbPoints()-nbNodeOut-1, isXConst, constValue ) :
82         grid->GetUVPtStruct( isXConst, constValue );
83     }
84   };
85
86   std::vector< Side >      side;
87   std::vector< UVPtStruct> uv_grid;
88   int                      iSize, jSize;
89   TopoDS_Face              face;
90   Bnd_B2d                  uv_box;
91
92   FaceQuadStruct ( const TopoDS_Face& F = TopoDS_Face() );
93   UVPtStruct& UVPt( int i, int j ) { return uv_grid[ i + j * iSize ]; }
94   void  shift    ( size_t nb, bool keepUnitOri );
95   int & nbNodeOut( int iSide ) { return side[ iSide ].nbNodeOut; }
96   bool  findCell ( const gp_XY& uv, int & i, int & j );
97   bool  isNear   ( const gp_XY& uv, int & i, int & j, int nbLoops=1 );
98   bool  isEqual  ( const gp_XY& uv, int   i, int   j );
99   void  normPa2IJ( double x, double y, int & i, int & j );
100   void  updateUV ( const gp_XY& uv, int   i, int   j, bool isVertical );
101
102   typedef boost::shared_ptr<FaceQuadStruct> Ptr;
103 };
104
105 class STDMESHERS_EXPORT StdMeshers_Quadrangle_2D: public SMESH_2D_Algo
106 {
107 public:
108   StdMeshers_Quadrangle_2D(int hypId, int studyId, SMESH_Gen* gen);
109   virtual ~StdMeshers_Quadrangle_2D();
110
111   virtual bool CheckHypothesis(SMESH_Mesh&         aMesh,
112                                const TopoDS_Shape& aShape,
113                                Hypothesis_Status&  aStatus);
114
115   virtual bool Compute(SMESH_Mesh&         aMesh,
116                        const TopoDS_Shape& aShape);
117
118   virtual bool Evaluate(SMESH_Mesh &         aMesh,
119                         const TopoDS_Shape & aShape,
120                         MapShapeNbElems&     aResMap);
121
122   FaceQuadStruct::Ptr CheckAnd2Dcompute(SMESH_Mesh&         aMesh,
123                                         const TopoDS_Shape& aShape,
124                                         const bool          CreateQuadratic);
125
126   FaceQuadStruct::Ptr CheckNbEdges(SMESH_Mesh&         aMesh,
127                                    const TopoDS_Shape& aShape,
128                                    const bool          considerMesh=false);
129
130 protected:
131
132   bool checkNbEdgesForEvaluate(SMESH_Mesh& aMesh,
133                                const TopoDS_Shape & aShape,
134                                MapShapeNbElems& aResMap,
135                                std::vector<int>& aNbNodes,
136                                bool& IsQuadratic);
137
138   bool setNormalizedGrid(FaceQuadStruct::Ptr quad);
139   
140   void splitQuadFace(SMESHDS_Mesh *       theMeshDS,
141                      const int            theFaceID,
142                      const SMDS_MeshNode* theNode1,
143                      const SMDS_MeshNode* theNode2,
144                      const SMDS_MeshNode* theNode3,
145                      const SMDS_MeshNode* theNode4);
146
147   bool computeQuadDominant(SMESH_Mesh&         aMesh,
148                            const TopoDS_Face&  aFace);
149
150   bool computeQuadDominant(SMESH_Mesh&         aMesh,
151                            const TopoDS_Face&  aFace,
152                            FaceQuadStruct::Ptr quad);
153
154   bool computeQuadPref(SMESH_Mesh&         aMesh,
155                        const TopoDS_Face&  aFace,
156                        FaceQuadStruct::Ptr quad);
157
158   bool computeTriangles(SMESH_Mesh&         aMesh,
159                         const TopoDS_Face&  aFace,
160                         FaceQuadStruct::Ptr quad);
161
162   bool evaluateQuadPref(SMESH_Mesh&         aMesh,
163                         const TopoDS_Shape& aShape,
164                         std::vector<int>&   aNbNodes,
165                         MapShapeNbElems&    aResMap,
166                         bool                isQuadratic);
167
168   bool computeReduced (SMESH_Mesh&         aMesh,
169                        const TopoDS_Face&  aFace,
170                        FaceQuadStruct::Ptr quad);
171
172   void updateDegenUV(FaceQuadStruct::Ptr quad);
173
174   void smooth (FaceQuadStruct::Ptr quad);
175
176   int getCorners(const TopoDS_Face&          theFace,
177                  SMESH_Mesh &                theMesh,
178                  std::list<TopoDS_Edge>&     theWire,
179                  std::vector<TopoDS_Vertex>& theVertices,
180                  int &                       theNbDegenEdges,
181                  const bool                  considerMesh);
182
183   bool getEnforcedUV();
184
185   bool addEnforcedNodes();
186
187   int splitQuad(FaceQuadStruct::Ptr quad, int i, int j);
188
189   typedef std::map< StdMeshers_FaceSidePtr, std::vector< FaceQuadStruct::Ptr > > TQuadsBySide;
190   void updateSideUV( FaceQuadStruct::Side&  side,
191                      int                    iForced,
192                      const TQuadsBySide&    quads,
193                      int *                  iNext=NULL);
194
195
196   // Fields
197
198   bool myQuadranglePreference;
199   bool myTrianglePreference;
200   int  myTriaVertexID;
201   bool myNeedSmooth;
202   const StdMeshers_QuadrangleParams* myParams;
203   StdMeshers_QuadType                myQuadType;
204
205   SMESH_MesherHelper*                myHelper;
206   SMESH_ProxyMesh::Ptr               myProxyMesh;
207   std::list< FaceQuadStruct::Ptr >   myQuadList;
208
209   struct ForcedPoint
210   {
211     gp_XY         uv;
212     gp_XYZ        xyz;
213     TopoDS_Vertex vertex;
214
215     double U() const { return uv.X(); }
216     double V() const { return uv.Y(); }
217     operator const gp_XY& () { return uv; }
218   };
219   std::vector< ForcedPoint >         myForcedPnts;
220 };
221
222 #endif