Salome HOME
Merge from V6_main (04/10/2012)
[modules/smesh.git] / src / StdMeshers / StdMeshers_Projection_1D2D.cxx
1 // Copyright (C) 2007-2012  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_Projection_1D2D.cxx
23 // Module    : SMESH
24 // Author    : Edward AGAPOV (eap)
25 //
26 #include "StdMeshers_Projection_1D2D.hxx"
27
28 #include "SMESH_Gen.hxx"
29 #include "SMESH_MesherHelper.hxx"
30 #include "SMESH_subMesh.hxx"
31 #include "SMESH_subMeshEventListener.hxx"
32 #include "StdMeshers_FaceSide.hxx"
33 #include "StdMeshers_ProjectionSource2D.hxx"
34 #include "StdMeshers_ProjectionUtils.hxx"
35
36 #include <TopoDS.hxx>
37 #include <TopExp_Explorer.hxx>
38
39 using namespace std;
40
41 namespace
42 {
43   // --------------------------------------------------------------------------------
44   /*!
45    * \brief an event listener updating submehses of EDGEs according to
46    *        events on the target FACE submesh
47    */
48   struct EventProparatorToEdges : public SMESH_subMeshEventListener
49   {
50     EventProparatorToEdges(): SMESH_subMeshEventListener(/*isDeletable=*/false,
51                                                          "Projection_1D2D::EventProparatorToEdges")
52     {}
53     static EventProparatorToEdges* Instance() { static EventProparatorToEdges E; return &E; }
54
55     static void Set(SMESH_subMesh* faceSubMesh)
56     {
57       SMESH_subMeshEventListenerData* edgeSubMeshes =
58         new SMESH_subMeshEventListenerData(/*isDeletable=*/true);
59       SMESH_Mesh* mesh = faceSubMesh->GetFather();
60       TopExp_Explorer eExp( faceSubMesh->GetSubShape(), TopAbs_EDGE );
61       for ( ; eExp.More(); eExp.Next() )
62         edgeSubMeshes->mySubMeshes.push_back( mesh->GetSubMesh( eExp.Current() ));
63
64       // set a listener
65       faceSubMesh->SetEventListener( Instance(), edgeSubMeshes, faceSubMesh );
66     }
67   };
68   // --------------------------------------------------------------------------------
69   /*!
70    * \brief Structure used to temporary remove EventProparatorToEdges from faceSubMesh
71    *  in order to prevent propagation of CLEAN event from FACE to EDGEs during 
72    *  StdMeshers_Projection_1D2D::Compute(). The CLEAN event is emmited by Pattern mapper
73    * and causes removal of faces generated on adjacent FACEs.
74    */
75   struct UnsetterOfEventProparatorToEdges
76   {
77     SMESH_subMesh* _faceSubMesh;
78     UnsetterOfEventProparatorToEdges( SMESH_subMesh* faceSubMesh ):_faceSubMesh(faceSubMesh)
79     {
80       faceSubMesh->DeleteEventListener( EventProparatorToEdges::Instance() );
81     }
82     ~UnsetterOfEventProparatorToEdges()
83     {
84       EventProparatorToEdges::Set(_faceSubMesh);
85     }
86   };
87 }
88
89 //=======================================================================
90 //function : StdMeshers_Projection_1D2D
91 //purpose  : 
92 //=======================================================================
93
94 StdMeshers_Projection_1D2D::StdMeshers_Projection_1D2D(int hypId, int studyId, SMESH_Gen* gen)
95   :StdMeshers_Projection_2D(hypId, studyId, gen)
96 {
97   _name = "Projection_1D2D";
98   _requireDiscreteBoundary = false;
99   _supportSubmeshes = true;
100 }
101
102 //=======================================================================
103 //function : Compute
104 //purpose  : 
105 //=======================================================================
106
107 bool StdMeshers_Projection_1D2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theShape)
108 {
109   UnsetterOfEventProparatorToEdges eventBarrier( theMesh.GetSubMesh( theShape ));
110
111   if ( !StdMeshers_Projection_2D::Compute(theMesh, theShape))
112     return false;
113
114   SMESHDS_Mesh * meshDS = theMesh.GetMeshDS();
115
116   SMESHDS_SubMesh * faceSubMesh = meshDS->MeshElements( theShape );
117   if ( !faceSubMesh || faceSubMesh->NbElements() == 0 ) return false;
118   _quadraticMesh = faceSubMesh->GetElements()->next()->IsQuadratic();
119
120   SMESH_MesherHelper helper( theMesh );
121   helper.SetSubShape( theShape );
122
123   TopoDS_Face F = TopoDS::Face( theShape );
124   TError err;
125   TSideVector wires = StdMeshers_FaceSide::GetFaceWires( F, theMesh,
126                                                          /*ignoreMediumNodes=*/false, err);
127   if ( err && !err->IsOK() )
128     return error( err );
129
130   for ( size_t iWire = 0; iWire < wires.size(); ++iWire )
131   {
132     vector<const SMDS_MeshNode*> nodes = wires[ iWire ]->GetOrderedNodes();
133     if ( nodes.empty() )
134       return error("Wrong nodes on a wire");
135
136     // check that all nodes are shared by faces generated on F
137     for ( size_t i = 0; i < nodes.size(); ++i )
138     {
139       SMDS_ElemIteratorPtr fIt = nodes[i]->GetInverseElementIterator(SMDSAbs_Face);
140       bool faceFound = false;
141       while ( !faceFound && fIt->more() )
142         faceFound = ( helper.GetSubShapeID() == fIt->next()->getshapeId() );
143       if ( !faceFound )
144         return error("The existing 1D mesh mismatches the generated 2D mesh");
145     }
146
147     const bool checkExisting = ( wires[ iWire ]->NbSegments() || helper.HasSeam() );
148
149     if ( _quadraticMesh )
150     {
151       for ( size_t i = 2; i < nodes.size(); i += 2 )
152       {
153         if ( checkExisting && meshDS->FindEdge( nodes[i-2], nodes[i], nodes[i-1]))
154           continue;
155         SMDS_MeshElement* e = meshDS->AddEdge( nodes[i-2], nodes[i], nodes[i-1] );
156         meshDS->SetMeshElementOnShape( e, nodes[i-1]->getshapeId() );
157       }
158     }
159     else
160     {
161       int edgeID = meshDS->ShapeToIndex( wires[ iWire ]->Edge(0) );
162       for ( size_t i = 1; i < nodes.size(); ++i )
163       {
164         if ( checkExisting && meshDS->FindEdge( nodes[i-1], nodes[i]))
165           continue;
166         SMDS_MeshElement* e = meshDS->AddEdge( nodes[i-1], nodes[i] );
167         if ( nodes[i-1]->getshapeId() != edgeID &&
168              nodes[i  ]->getshapeId() != edgeID )
169         {
170           edgeID = helper.GetMediumPos( nodes[i-1], nodes[i] ).first;
171           if ( edgeID < 1 ) edgeID = helper.GetSubShapeID();
172         }
173         meshDS->SetMeshElementOnShape( e, edgeID );
174       }
175     }
176   }   
177
178   return true;
179 }
180
181 //=======================================================================
182 //function : Evaluate
183 //purpose  : 
184 //=======================================================================
185
186 bool StdMeshers_Projection_1D2D::Evaluate(SMESH_Mesh&         theMesh,
187                                           const TopoDS_Shape& theShape,
188                                           MapShapeNbElems&    aResMap)
189 {
190   if ( !StdMeshers_Projection_2D::Evaluate(theMesh,theShape,aResMap))
191     return false;
192
193   TopoDS_Shape srcFace = _sourceHypo->GetSourceFace();
194   SMESH_Mesh * srcMesh = _sourceHypo->GetSourceMesh();
195   if ( !srcMesh ) srcMesh = & theMesh;
196   SMESH_subMesh* srcFaceSM = srcMesh->GetSubMesh( srcFace );
197
198   typedef StdMeshers_ProjectionUtils SPU;
199   SPU::TShapeShapeMap shape2ShapeMap;
200   SPU::InitVertexAssociation( _sourceHypo, shape2ShapeMap );
201   if ( !SPU::FindSubShapeAssociation( theShape, &theMesh, srcFace, srcMesh, shape2ShapeMap))
202     return error(COMPERR_BAD_SHAPE,"Topology of source and target faces seems different" );
203
204   MapShapeNbElems srcResMap;
205   if ( !srcFaceSM->IsMeshComputed() )
206     _gen->Evaluate( *srcMesh, srcFace, srcResMap);
207
208   SMESH_subMeshIteratorPtr smIt = srcFaceSM->getDependsOnIterator(/*includeSelf=*/false,
209                                                                   /*complexShapeFirst=*/true);
210   while ( smIt->more() )
211   {
212     SMESH_subMesh* srcSM = smIt->next();
213     TopAbs_ShapeEnum shapeType = srcSM->GetSubShape().ShapeType();
214     if ( shapeType == TopAbs_EDGE )
215     {
216       std::vector<int> aVec;
217       SMESHDS_SubMesh* srcSubMeshDS = srcSM->GetSubMeshDS();
218       if ( srcSubMeshDS && srcSubMeshDS->NbElements() )
219       {
220         aVec.resize(SMDSEntity_Last, 0);
221         SMDS_ElemIteratorPtr eIt = srcSubMeshDS->GetElements();
222         _quadraticMesh = ( eIt->more() && eIt->next()->IsQuadratic() );
223
224         aVec[SMDSEntity_Node] = srcSubMeshDS->NbNodes();
225         aVec[_quadraticMesh ? SMDSEntity_Quad_Edge : SMDSEntity_Edge] = srcSubMeshDS->NbElements();
226       }
227       else
228       {
229         if ( srcResMap.empty() )
230           if ( !_gen->Evaluate( *srcMesh, srcSM->GetSubShape(), srcResMap ))
231             return error(COMPERR_BAD_INPUT_MESH,"Source mesh not evaluatable");
232         aVec = srcResMap[ srcSM ];
233         if ( aVec.empty() )
234           return error(COMPERR_BAD_INPUT_MESH,"Source mesh is wrongly evaluated");
235       }
236       TopoDS_Shape tgtEdge = shape2ShapeMap( srcSM->GetSubShape(), /*isSrc=*/true  );
237       SMESH_subMesh* tgtSM = theMesh.GetSubMesh( tgtEdge );
238       aResMap.insert(std::make_pair(tgtSM,aVec));
239     }
240     if ( shapeType == TopAbs_VERTEX ) break;
241   }
242
243   return true;
244 }
245
246 //=======================================================================
247 //function : SetEventListener
248 //purpose  : Sets a default event listener to submesh of the source face.
249 //           faceSubMesh - submesh where algo is set
250 // After being set, event listener is notified on each event of a submesh.
251 // This method is called when a submesh gets HYP_OK algo_state.
252 // Arranges that CLEAN event is translated from source submesh to
253 // the faceSubMesh submesh.
254 //=======================================================================
255
256 void StdMeshers_Projection_1D2D::SetEventListener(SMESH_subMesh* faceSubMesh)
257 {
258   // set a listener of events on a source submesh
259   StdMeshers_Projection_2D::SetEventListener(faceSubMesh);
260
261   // set a listener to the target FACE submesh in order to update submehses
262   // of EDGEs according to events on the target FACE submesh
263   EventProparatorToEdges::Set( faceSubMesh );
264 }
265