1 // Copyright (C) 2007-2016 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
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.
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.
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
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
23 // File : StdMeshers_PolygonPerFace_2D.cxx
25 // Created : Fri Oct 20 11:37:07 2006
26 // Author : Edward AGAPOV (eap)
28 #include "StdMeshers_PolygonPerFace_2D.hxx"
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"
38 #include <TopExp_Explorer.hxx>
39 #include <TopoDS_Face.hxx>
46 //=======================================================================
47 //function : StdMeshers_PolygonPerFace_2D
49 //=======================================================================
51 StdMeshers_PolygonPerFace_2D::StdMeshers_PolygonPerFace_2D(int hypId,
54 :SMESH_2D_Algo(hypId, studyId, gen)
56 _name = "PolygonPerFace_2D";
59 //=======================================================================
60 //function : CheckHypothesis
62 //=======================================================================
64 bool StdMeshers_PolygonPerFace_2D::CheckHypothesis(SMESH_Mesh& theMesh,
65 const TopoDS_Shape& theShape,
66 SMESH_Hypothesis::Hypothesis_Status& theStatus)
72 //=======================================================================
75 //=======================================================================
77 bool StdMeshers_PolygonPerFace_2D::Compute(SMESH_Mesh& theMesh,
78 const TopoDS_Shape& theShape)
80 const TopoDS_Face& face = TopoDS::Face( theShape );
82 SMESH_MesherHelper helper( theMesh );
83 helper.SetElementsOnShape( true );
84 _quadraticMesh = helper.IsQuadraticSubMesh( face );
86 SMESH_ProxyMesh::Ptr proxyMesh = StdMeshers_ViscousLayers2D::Compute( theMesh, face );
91 TSideVector wires = StdMeshers_FaceSide::GetFaceWires(face, theMesh,
92 /*skipMediumNodes=*/_quadraticMesh,
94 /*checkVertexNodes=*/false);
95 if ( wires.size() != 1 )
96 return error( COMPERR_BAD_SHAPE, SMESH_Comment("One wire required, not ") << wires.size() );
98 vector<const SMDS_MeshNode*> nodes = wires[0]->GetOrderedNodes();
99 int nbNodes = int( nodes.size() ) - 1; // 1st node is repeated at end
103 helper.AddFace( nodes[0], nodes[1], nodes[2] );
106 helper.AddFace( nodes[0], nodes[1], nodes[2], nodes[3] );
110 return error( COMPERR_BAD_INPUT_MESH, "Less that 3 nodes on the wire" );
111 nodes.resize( nodes.size() - 1 );
112 helper.AddPolygonalFace ( nodes );
118 //=======================================================================
119 //function : Evaluate
121 //=======================================================================
123 bool StdMeshers_PolygonPerFace_2D::Evaluate(SMESH_Mesh& theMesh,
124 const TopoDS_Shape& theShape,
125 MapShapeNbElems& theResMap)
128 int nbLinSegs = 0, nbQuadSegs = 0;
129 TopExp_Explorer edge( theShape, TopAbs_EDGE );
130 for ( ; edge.More(); edge.Next() )
132 SMESH_subMesh* sm = theMesh.GetSubMesh( edge.Current() );
133 MapShapeNbElems::iterator sm2vec = theResMap.find( sm );
134 if ( sm2vec == theResMap.end() )
136 nbLinSegs += sm2vec->second.at( SMDSEntity_Edge );
137 nbQuadSegs += sm2vec->second.at( SMDSEntity_Quad_Edge );
140 std::vector<int> aVec( SMDSEntity_Last, 0 );
141 switch ( nbLinSegs + nbQuadSegs ) {
143 aVec[ nbQuadSegs ? SMDSEntity_Quad_Triangle : SMDSEntity_Triangle ] = 1;
146 aVec[ nbQuadSegs ? SMDSEntity_Quad_Quadrangle : SMDSEntity_Quadrangle ] = 1;
149 if ( nbLinSegs + nbQuadSegs < 3 )
150 return error( COMPERR_BAD_INPUT_MESH, "Less that 3 nodes on the wire" );
151 aVec[ nbQuadSegs ? SMDSEntity_Quad_Polygon : SMDSEntity_Polygon ] = 1;
154 SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
155 theResMap.insert(std::make_pair(sm,aVec));