Salome HOME
OCCT dev version porting (6.7.2)
[modules/smesh.git] / src / StdMeshers / StdMeshers_Projection_1D2D.cxx
1 // Copyright (C) 2007-2014  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 // 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   // 1) Project faces
112
113   if ( !StdMeshers_Projection_2D::Compute(theMesh, theShape))
114     return false;
115
116   // 2) Create segments
117
118   SMESHDS_Mesh * meshDS = theMesh.GetMeshDS();
119
120   SMESHDS_SubMesh * faceSubMesh = meshDS->MeshElements( theShape );
121   if ( !faceSubMesh || faceSubMesh->NbElements() == 0 ) return false;
122   _quadraticMesh = faceSubMesh->GetElements()->next()->IsQuadratic();
123
124   SMESH_MesherHelper helper( theMesh );
125   helper.SetSubShape( theShape );
126
127   if ( _quadraticMesh )
128   {
129     // 2a) Move some medium nodes from FACE to EDGES. They are on FACE because
130     // EDGEs are discreteized later than FACE, in this case.
131
132     SMESH_MesherHelper posFixer( theMesh );
133     posFixer.ToFixNodeParameters( true );
134     SMDS_ElemIteratorPtr fIt = faceSubMesh->GetElements();
135     vector< const SMDS_MeshNode* > nodes;
136     double dummyU, tol = 1e-7;
137     while ( fIt->more() ) // loop on mesh faces created by StdMeshers_Projection_2D
138     {
139       const SMDS_MeshElement* f = fIt->next();
140       //if ( !f->IsQuadratic() ) continue;
141       nodes.assign( SMDS_MeshElement::iterator( f->interlacedNodesElemIterator() ),
142                     SMDS_MeshElement::iterator() );
143       nodes.push_back( nodes[0] );
144       for ( size_t i = 2; i < nodes.size(); i += 2 )
145       {
146         pair<int, TopAbs_ShapeEnum> idType = helper.GetMediumPos( nodes[i], nodes[i-2] );
147         if ( idType.second == TopAbs_EDGE &&
148              idType.first  != nodes[i-1]->getshapeId() )
149         {
150           faceSubMesh->RemoveNode( nodes[i-1], /*isDeleted=*/false );
151           meshDS->SetNodeOnEdge( (SMDS_MeshNode*) nodes[i-1], idType.first );
152           posFixer.SetSubShape( idType.first );
153           posFixer.CheckNodeU( TopoDS::Edge( posFixer.GetSubShape() ),
154                                nodes[i-1], dummyU=0., tol, /*force=*/true );
155         }
156       }
157     }
158   }
159   TopoDS_Face F = TopoDS::Face( theShape );
160   TError err;
161   TSideVector wires = StdMeshers_FaceSide::GetFaceWires( F, theMesh,
162                                                          /*ignoreMediumNodes=*/false, err);
163   if ( err && !err->IsOK() )
164     return error( err );
165
166   for ( size_t iWire = 0; iWire < wires.size(); ++iWire )
167   {
168     vector<const SMDS_MeshNode*> nodes = wires[ iWire ]->GetOrderedNodes();
169     if ( nodes.empty() )
170       return error("Wrong nodes on a wire");
171
172     // check that all nodes are shared by faces generated on F
173     for ( size_t i = 0; i < nodes.size(); ++i )
174     {
175       SMDS_ElemIteratorPtr fIt = nodes[i]->GetInverseElementIterator(SMDSAbs_Face);
176       bool faceFound = false;
177       while ( !faceFound && fIt->more() )
178         faceFound = ( helper.GetSubShapeID() == fIt->next()->getshapeId() );
179       if ( !faceFound )
180         return error("The existing 1D mesh mismatches the generated 2D mesh");
181     }
182
183     const bool checkExisting = ( wires[ iWire ]->NbSegments() || helper.HasSeam() );
184
185     if ( _quadraticMesh )
186     {
187       for ( size_t i = 2; i < nodes.size(); i += 2 )
188       {
189         if ( checkExisting && meshDS->FindEdge( nodes[i-2], nodes[i], nodes[i-1]))
190           continue;
191         SMDS_MeshElement* e = meshDS->AddEdge( nodes[i-2], nodes[i], nodes[i-1] );
192         meshDS->SetMeshElementOnShape( e, nodes[i-1]->getshapeId() );
193       }
194     }
195     else
196     {
197       int edgeID = meshDS->ShapeToIndex( wires[ iWire ]->Edge(0) );
198       for ( size_t i = 1; i < nodes.size(); ++i )
199       {
200         if ( checkExisting && meshDS->FindEdge( nodes[i-1], nodes[i]))
201           continue;
202         SMDS_MeshElement* e = meshDS->AddEdge( nodes[i-1], nodes[i] );
203         if ( nodes[i-1]->getshapeId() != edgeID &&
204              nodes[i  ]->getshapeId() != edgeID )
205         {
206           edgeID = helper.GetMediumPos( nodes[i-1], nodes[i] ).first;
207           if ( edgeID < 1 ) edgeID = helper.GetSubShapeID();
208         }
209         meshDS->SetMeshElementOnShape( e, edgeID );
210       }
211     }
212   }
213
214   return true;
215 }
216
217 //=======================================================================
218 //function : Evaluate
219 //purpose  : 
220 //=======================================================================
221
222 bool StdMeshers_Projection_1D2D::Evaluate(SMESH_Mesh&         theMesh,
223                                           const TopoDS_Shape& theShape,
224                                           MapShapeNbElems&    aResMap)
225 {
226   if ( !StdMeshers_Projection_2D::Evaluate(theMesh,theShape,aResMap))
227     return false;
228
229   TopoDS_Shape srcFace = _sourceHypo->GetSourceFace();
230   SMESH_Mesh * srcMesh = _sourceHypo->GetSourceMesh();
231   if ( !srcMesh ) srcMesh = & theMesh;
232   SMESH_subMesh* srcFaceSM = srcMesh->GetSubMesh( srcFace );
233
234   namespace SPU = StdMeshers_ProjectionUtils;
235   SPU::TShapeShapeMap shape2ShapeMap;
236   SPU::InitVertexAssociation( _sourceHypo, shape2ShapeMap );
237   if ( !SPU::FindSubShapeAssociation( theShape, &theMesh, srcFace, srcMesh, shape2ShapeMap))
238     return error(COMPERR_BAD_SHAPE,"Topology of source and target faces seems different" );
239
240   MapShapeNbElems srcResMap;
241   if ( !srcFaceSM->IsMeshComputed() )
242     _gen->Evaluate( *srcMesh, srcFace, srcResMap);
243
244   SMESH_subMeshIteratorPtr smIt = srcFaceSM->getDependsOnIterator(/*includeSelf=*/false,
245                                                                   /*complexShapeFirst=*/true);
246   while ( smIt->more() )
247   {
248     SMESH_subMesh* srcSM = smIt->next();
249     TopAbs_ShapeEnum shapeType = srcSM->GetSubShape().ShapeType();
250     if ( shapeType == TopAbs_EDGE )
251     {
252       std::vector<int> aVec;
253       SMESHDS_SubMesh* srcSubMeshDS = srcSM->GetSubMeshDS();
254       if ( srcSubMeshDS && srcSubMeshDS->NbElements() )
255       {
256         aVec.resize(SMDSEntity_Last, 0);
257         SMDS_ElemIteratorPtr eIt = srcSubMeshDS->GetElements();
258         _quadraticMesh = ( eIt->more() && eIt->next()->IsQuadratic() );
259
260         aVec[SMDSEntity_Node] = srcSubMeshDS->NbNodes();
261         aVec[_quadraticMesh ? SMDSEntity_Quad_Edge : SMDSEntity_Edge] = srcSubMeshDS->NbElements();
262       }
263       else
264       {
265         if ( srcResMap.empty() )
266           if ( !_gen->Evaluate( *srcMesh, srcSM->GetSubShape(), srcResMap ))
267             return error(COMPERR_BAD_INPUT_MESH,"Source mesh not evaluatable");
268         aVec = srcResMap[ srcSM ];
269         if ( aVec.empty() )
270           return error(COMPERR_BAD_INPUT_MESH,"Source mesh is wrongly evaluated");
271       }
272       TopoDS_Shape tgtEdge = shape2ShapeMap( srcSM->GetSubShape(), /*isSrc=*/true  );
273       SMESH_subMesh* tgtSM = theMesh.GetSubMesh( tgtEdge );
274       aResMap.insert(std::make_pair(tgtSM,aVec));
275     }
276     if ( shapeType == TopAbs_VERTEX ) break;
277   }
278
279   return true;
280 }
281
282 //=======================================================================
283 //function : SetEventListener
284 //purpose  : Sets a default event listener to submesh of the source face.
285 //           faceSubMesh - submesh where algo is set
286 // After being set, event listener is notified on each event of a submesh.
287 // This method is called when a submesh gets HYP_OK algo_state.
288 // Arranges that CLEAN event is translated from source submesh to
289 // the faceSubMesh submesh.
290 //=======================================================================
291
292 void StdMeshers_Projection_1D2D::SetEventListener(SMESH_subMesh* faceSubMesh)
293 {
294   // set a listener of events on a source submesh
295   StdMeshers_Projection_2D::SetEventListener(faceSubMesh);
296
297   // set a listener to the target FACE submesh in order to update submehses
298   // of EDGEs according to events on the target FACE submesh
299   EventProparatorToEdges::Set( faceSubMesh );
300 }
301