Salome HOME
bos #20256: [CEA 18523] Porting SMESH to int 64 bits
[modules/smesh.git] / src / StdMeshers / StdMeshers_PolygonPerFace_2D.cxx
1 // Copyright (C) 2007-2021  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, or (at your option) any later version.
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
23 // File      : StdMeshers_PolygonPerFace_2D.cxx
24 // Module    : SMESH
25 // Created   : Fri Oct 20 11:37:07 2006
26 // Author    : Edward AGAPOV (eap)
27 //
28 #include "StdMeshers_PolygonPerFace_2D.hxx"
29
30 #include "SMESH_Comment.hxx"
31 #include "SMESH_Mesh.hxx"
32 #include "SMESH_MesherHelper.hxx"
33 #include "SMESH_ProxyMesh.hxx"
34 #include "SMESH_subMesh.hxx"
35 #include "StdMeshers_FaceSide.hxx"
36 #include "StdMeshers_ViscousLayers2D.hxx"
37
38 #include <TopExp_Explorer.hxx>
39 #include <TopoDS_Face.hxx>
40 #include <TopoDS.hxx>
41
42 #include <vector>
43
44 using namespace std;
45
46 //=======================================================================
47 //function : StdMeshers_PolygonPerFace_2D
48 //purpose  :
49 //=======================================================================
50
51 StdMeshers_PolygonPerFace_2D::StdMeshers_PolygonPerFace_2D(int        hypId,
52                                                            SMESH_Gen* gen)
53   :SMESH_2D_Algo(hypId, gen)
54 {
55   _name = "PolygonPerFace_2D";
56 }
57
58 //=======================================================================
59 //function : CheckHypothesis
60 //purpose  : 
61 //=======================================================================
62
63 bool StdMeshers_PolygonPerFace_2D::CheckHypothesis(SMESH_Mesh&                          /*theMesh*/,
64                                                    const TopoDS_Shape&                  /*theShape*/,
65                                                    SMESH_Hypothesis::Hypothesis_Status& theStatus)
66 {
67   theStatus = HYP_OK;
68   return true;
69 }
70
71 //=======================================================================
72 //function : Compute
73 //purpose  :
74 //=======================================================================
75
76 bool StdMeshers_PolygonPerFace_2D::Compute(SMESH_Mesh&         theMesh,
77                                            const TopoDS_Shape& theShape)
78 {
79   const TopoDS_Face& face = TopoDS::Face( theShape );
80
81   SMESH_MesherHelper helper( theMesh );
82   helper.SetElementsOnShape( true );
83   _quadraticMesh = helper.IsQuadraticSubMesh( face );
84
85   SMESH_ProxyMesh::Ptr proxyMesh = StdMeshers_ViscousLayers2D::Compute( theMesh, face );
86   if ( !proxyMesh )
87     return false;
88
89   TError err;
90   TSideVector wires = StdMeshers_FaceSide::GetFaceWires(face, theMesh,
91                                                         /*skipMediumNodes=*/_quadraticMesh,
92                                                         err, &helper, proxyMesh,
93                                                         /*checkVertexNodes=*/false);
94   if ( wires.size() != 1 )
95     return error( COMPERR_BAD_SHAPE, SMESH_Comment("One wire required, not ") << wires.size() );
96
97   vector<const SMDS_MeshNode*> nodes = wires[0]->GetOrderedNodes();
98   int nbNodes = int( nodes.size() ) - 1; // 1st node is repeated at end
99
100   switch ( nbNodes ) {
101   case 3:
102     helper.AddFace( nodes[0], nodes[1], nodes[2] );
103     break;
104   case 4:
105     helper.AddFace( nodes[0], nodes[1], nodes[2], nodes[3] );
106     break;
107   default:
108     if ( nbNodes < 3 )
109       return error( COMPERR_BAD_INPUT_MESH, "Less that 3 nodes on the wire" );
110     nodes.resize( nodes.size() - 1 );
111     helper.AddPolygonalFace ( nodes );
112   }
113
114   return true;
115 }
116
117 //=======================================================================
118 //function : Evaluate
119 //purpose  : 
120 //=======================================================================
121
122 bool StdMeshers_PolygonPerFace_2D::Evaluate(SMESH_Mesh&         theMesh,
123                                             const TopoDS_Shape& theShape,
124                                             MapShapeNbElems&    theResMap)
125 {
126   // count nb segments
127   int nbLinSegs = 0, nbQuadSegs = 0;
128   TopExp_Explorer edge( theShape, TopAbs_EDGE );
129   for ( ; edge.More(); edge.Next() )
130   {
131     SMESH_subMesh* sm = theMesh.GetSubMesh( edge.Current() );
132     MapShapeNbElems::iterator sm2vec = theResMap.find( sm );
133     if ( sm2vec == theResMap.end() )
134       continue;
135     nbLinSegs  += sm2vec->second.at( SMDSEntity_Edge );
136     nbQuadSegs += sm2vec->second.at( SMDSEntity_Quad_Edge );
137   }
138
139   std::vector<smIdType> aVec( SMDSEntity_Last, 0 );
140   switch ( nbLinSegs + nbQuadSegs ) {
141   case 3:
142     aVec[ nbQuadSegs ? SMDSEntity_Quad_Triangle : SMDSEntity_Triangle ] = 1;
143     break;
144   case 4:
145     aVec[ nbQuadSegs ? SMDSEntity_Quad_Quadrangle : SMDSEntity_Quadrangle ] = 1;
146     break;
147   default:
148     if ( nbLinSegs + nbQuadSegs < 3 )
149       return error( COMPERR_BAD_INPUT_MESH, "Less that 3 nodes on the wire" );
150     aVec[ nbQuadSegs ? SMDSEntity_Quad_Polygon : SMDSEntity_Polygon ] = 1;
151   }
152
153   SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
154   theResMap.insert(std::make_pair(sm,aVec));
155
156   return true;
157 }