Salome HOME
b686fe1e0d2a42636b4e6fd8d3af27f5455fe3f6
[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     bool Reverse();
75     int  NbPoints() const { return Abs( to - from ); }
76     double Param( int nodeIndex ) const;
77     double Length( int from=-1, int to=-1) const;
78     gp_XY Value2d( double x ) const;
79     const UVPtStruct& First() const { return GetUVPtStruct()[ from ]; }
80     const UVPtStruct& Last()  const {
81       return GetUVPtStruct()[ to-nbNodeOut-(IsReversed() ? -1 : +1)];
82     }
83     // some sortcuts
84     const vector<UVPtStruct>& GetUVPtStruct(bool isXConst=0, double constValue=0) const
85     { return nbNodeOut ?
86         grid->SimulateUVPtStruct( NbPoints()-nbNodeOut-1, isXConst, constValue ) :
87         grid->GetUVPtStruct( isXConst, constValue );
88     }
89   };
90   struct SideIterator // iterator on UVPtStruct of a Side
91   {
92     const UVPtStruct *uvPtr, *uvEnd;
93     int               dPtr, counter;
94     SideIterator(): uvPtr(0), uvEnd(0), dPtr(0), counter(0) {}
95     void Init( const Side& side ) {
96       dPtr  = counter = 0;
97       uvPtr = uvEnd = 0;
98       if ( side.NbPoints() > 0 ) {
99         uvPtr = & side.First();
100         uvEnd = & side.Last();
101         dPtr  = ( uvEnd > uvPtr ) ? +1 : -1;
102         uvEnd += dPtr;
103       }
104     }
105     bool More() const { return uvPtr != uvEnd; }
106     void Next() { uvPtr += dPtr; ++counter; }
107     UVPtStruct& UVPt() const { return (UVPtStruct&) *uvPtr; }
108     int  Count() const { return counter; }
109   };
110
111   std::vector< Side >      side;
112   std::vector< UVPtStruct> uv_grid;
113   int                      iSize, jSize;
114   TopoDS_Face              face;
115   Bnd_B2d                  uv_box;
116   std::string              name; // to ease debugging
117
118   FaceQuadStruct ( const TopoDS_Face& F = TopoDS_Face(), const std::string& nm="main" );
119   UVPtStruct& UVPt( int i, int j ) { return uv_grid[ i + j * iSize ]; }
120   void  shift    ( size_t nb, bool keepUnitOri );
121   int & nbNodeOut( int iSide ) { return side[ iSide ].nbNodeOut; }
122   bool  findCell ( const gp_XY& uv, int & i, int & j );
123   bool  isNear   ( const gp_XY& uv, int & i, int & j, int nbLoops=1 );
124   bool  isEqual  ( const gp_XY& uv, int   i, int   j );
125   void  normPa2IJ( double x, double y, int & i, int & j );
126   void  updateUV ( const gp_XY& uv, int   i, int   j, bool isVertical );
127
128   typedef boost::shared_ptr<FaceQuadStruct> Ptr;
129 };
130
131 class STDMESHERS_EXPORT StdMeshers_Quadrangle_2D: public SMESH_2D_Algo
132 {
133 public:
134   StdMeshers_Quadrangle_2D(int hypId, int studyId, SMESH_Gen* gen);
135   virtual ~StdMeshers_Quadrangle_2D();
136
137   virtual bool CheckHypothesis(SMESH_Mesh&         aMesh,
138                                const TopoDS_Shape& aShape,
139                                Hypothesis_Status&  aStatus);
140
141   virtual bool Compute(SMESH_Mesh&         aMesh,
142                        const TopoDS_Shape& aShape);
143
144   virtual bool Evaluate(SMESH_Mesh &         aMesh,
145                         const TopoDS_Shape & aShape,
146                         MapShapeNbElems&     aResMap);
147
148   FaceQuadStruct::Ptr CheckAnd2Dcompute(SMESH_Mesh&         aMesh,
149                                         const TopoDS_Shape& aShape,
150                                         const bool          CreateQuadratic);
151
152   FaceQuadStruct::Ptr CheckNbEdges(SMESH_Mesh&         aMesh,
153                                    const TopoDS_Shape& aShape,
154                                    const bool          considerMesh=false);
155
156 protected:
157
158   bool checkNbEdgesForEvaluate(SMESH_Mesh& aMesh,
159                                const TopoDS_Shape & aShape,
160                                MapShapeNbElems& aResMap,
161                                std::vector<int>& aNbNodes,
162                                bool& IsQuadratic);
163
164   bool setNormalizedGrid(FaceQuadStruct::Ptr quad);
165   
166   void splitQuadFace(SMESHDS_Mesh *       theMeshDS,
167                      const int            theFaceID,
168                      const SMDS_MeshNode* theNode1,
169                      const SMDS_MeshNode* theNode2,
170                      const SMDS_MeshNode* theNode3,
171                      const SMDS_MeshNode* theNode4);
172
173   bool computeQuadDominant(SMESH_Mesh&         aMesh,
174                            const TopoDS_Face&  aFace);
175
176   bool computeQuadDominant(SMESH_Mesh&         aMesh,
177                            const TopoDS_Face&  aFace,
178                            FaceQuadStruct::Ptr quad);
179
180   bool computeQuadPref(SMESH_Mesh&         aMesh,
181                        const TopoDS_Face&  aFace,
182                        FaceQuadStruct::Ptr quad);
183
184   bool computeTriangles(SMESH_Mesh&         aMesh,
185                         const TopoDS_Face&  aFace,
186                         FaceQuadStruct::Ptr quad);
187
188   bool evaluateQuadPref(SMESH_Mesh&         aMesh,
189                         const TopoDS_Shape& aShape,
190                         std::vector<int>&   aNbNodes,
191                         MapShapeNbElems&    aResMap,
192                         bool                isQuadratic);
193
194   bool computeReduced (SMESH_Mesh&         aMesh,
195                        const TopoDS_Face&  aFace,
196                        FaceQuadStruct::Ptr quad);
197
198   void updateDegenUV(FaceQuadStruct::Ptr quad);
199
200   void smooth (FaceQuadStruct::Ptr quad);
201
202   int getCorners(const TopoDS_Face&          theFace,
203                  SMESH_Mesh &                theMesh,
204                  std::list<TopoDS_Edge>&     theWire,
205                  std::vector<TopoDS_Vertex>& theVertices,
206                  int &                       theNbDegenEdges,
207                  const bool                  considerMesh);
208
209   bool getEnforcedUV();
210
211   bool addEnforcedNodes();
212
213   int splitQuad(FaceQuadStruct::Ptr quad, int i, int j);
214
215   typedef std::map< StdMeshers_FaceSidePtr, std::vector< FaceQuadStruct::Ptr > > TQuadsBySide;
216   void updateSideUV( FaceQuadStruct::Side&  side,
217                      int                    iForced,
218                      const TQuadsBySide&    quads,
219                      int *                  iNext=NULL);
220
221
222   // Fields
223
224   bool myQuadranglePreference;
225   bool myTrianglePreference;
226   int  myTriaVertexID;
227   bool myNeedSmooth;
228   const StdMeshers_QuadrangleParams* myParams;
229   StdMeshers_QuadType                myQuadType;
230
231   SMESH_MesherHelper*                myHelper;
232   SMESH_ProxyMesh::Ptr               myProxyMesh;
233   std::list< FaceQuadStruct::Ptr >   myQuadList;
234
235   struct ForcedPoint
236   {
237     gp_XY         uv;
238     gp_XYZ        xyz;
239     TopoDS_Vertex vertex;
240
241     double U() const { return uv.X(); }
242     double V() const { return uv.Y(); }
243     operator const gp_XY& () { return uv; }
244   };
245   std::vector< ForcedPoint >         myForcedPnts;
246 };
247
248 #endif