Salome HOME
0020996: EDF 1549 SMESH: Projection 1D/2D algorithms fail
[modules/smesh.git] / src / StdMeshers / StdMeshers_ProjectionUtils.cxx
1 //  Copyright (C) 2007-2010  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
23 //  SMESH SMESH : idl implementation based on 'SMESH' unit's calsses
24 // File      : StdMeshers_ProjectionUtils.cxx
25 // Created   : Fri Oct 27 10:24:28 2006
26 // Author    : Edward AGAPOV (eap)
27 //
28 #include "StdMeshers_ProjectionUtils.hxx"
29
30 #include "StdMeshers_ProjectionSource1D.hxx"
31 #include "StdMeshers_ProjectionSource2D.hxx"
32 #include "StdMeshers_ProjectionSource3D.hxx"
33
34 #include "SMESH_Algo.hxx"
35 #include "SMESH_Block.hxx"
36 #include "SMESH_Gen.hxx"
37 #include "SMESH_Hypothesis.hxx"
38 #include "SMESH_IndexedDataMapOfShapeIndexedMapOfShape.hxx"
39 #include "SMESH_Mesh.hxx"
40 #include "SMESH_MesherHelper.hxx"
41 #include "SMESH_subMesh.hxx"
42 #include "SMESH_subMeshEventListener.hxx"
43 #include "SMDS_EdgePosition.hxx"
44
45 #include "utilities.h"
46
47 #include <BRepAdaptor_Surface.hxx>
48 #include <BRepTools.hxx>
49 #include <BRepTools_WireExplorer.hxx>
50 #include <BRep_Builder.hxx>
51 #include <BRep_Tool.hxx>
52 #include <Bnd_Box.hxx>
53 #include <TopAbs.hxx>
54 #include <TopExp.hxx>
55 #include <TopExp_Explorer.hxx>
56 #include <TopTools_Array1OfShape.hxx>
57 #include <TopTools_DataMapIteratorOfDataMapOfShapeListOfShape.hxx>
58 #include <TopTools_DataMapIteratorOfDataMapOfShapeShape.hxx>
59 #include <TopTools_ListIteratorOfListOfShape.hxx>
60 #include <TopTools_ListOfShape.hxx>
61 #include <TopTools_MapOfShape.hxx>
62 #include <TopoDS.hxx>
63 #include <TopoDS_Compound.hxx>
64 #include <TopoDS_Shape.hxx>
65 #include <gp_Pnt.hxx>
66 #include <gp_Vec.hxx>
67
68 #include <numeric>
69
70 using namespace std;
71
72
73 #define RETURN_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); return false; }
74 #define CONT_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); continue; }
75 #define SHOW_SHAPE(v,msg) \
76 // { \
77 //  if ( (v).IsNull() ) cout << msg << " NULL SHAPE" << endl; \
78 // else if ((v).ShapeType() == TopAbs_VERTEX) {\
79 //   gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( (v) ));\
80 //   cout<<msg<<" "<<shapeIndex((v))<<" ( "<<p.X()<<", "<<p.Y()<<", "<<p.Z()<<" )"<<endl;} \
81 // else {\
82 // cout << msg << " "; TopAbs::Print((v).ShapeType(),cout) <<" "<<shapeIndex((v))<<endl;}\
83 // }
84 #define SHOW_LIST(msg,l) \
85 // { \
86 //     cout << msg << " ";\
87 //     list< TopoDS_Edge >::const_iterator e = l.begin();\
88 //     for ( int i = 0; e != l.end(); ++e, ++i ) {\
89 //       cout << i << "V (" << TopExp::FirstVertex( *e, true ).TShape().operator->() << ") "\
90 //            << i << "E (" << e->TShape().operator->() << "); "; }\
91 //     cout << endl;\
92 //   }
93
94 #define HERE StdMeshers_ProjectionUtils
95
96 namespace {
97
98   static SMESHDS_Mesh* theMeshDS[2] = { 0, 0 }; // used to debug only
99   long shapeIndex(const TopoDS_Shape& S)
100   {
101     if ( theMeshDS[0] && theMeshDS[1] )
102       return max(theMeshDS[0]->ShapeToIndex(S), theMeshDS[1]->ShapeToIndex(S) );
103     return long(S.TShape().operator->());
104   }
105
106   //================================================================================
107   /*!
108    * \brief Write shape for debug purposes
109    */
110   //================================================================================
111
112   bool _StoreBadShape(const TopoDS_Shape& shape)
113   {
114 #ifdef _DEBUG_
115     const char* type[] ={"COMPOUND","COMPSOLID","SOLID","SHELL","FACE","WIRE","EDGE","VERTEX"};
116     BRepTools::Write( shape, SMESH_Comment("/tmp/") << type[shape.ShapeType()] << "_"
117                       << shape.TShape().operator->() << ".brep");
118 #endif
119     return false;
120   }
121   
122   //================================================================================
123   /*!
124    * \brief Reverse order of edges in a list and their orientation
125     * \param edges - list of edges to reverse
126     * \param nbEdges - number of edges to reverse
127    */
128   //================================================================================
129
130   void Reverse( list< TopoDS_Edge > & edges, const int nbEdges, const int firstEdge=0)
131   {
132     SHOW_LIST("BEFORE REVERSE", edges);
133
134     list< TopoDS_Edge >::iterator eIt = edges.begin();
135     std::advance( eIt, firstEdge );
136     list< TopoDS_Edge >::iterator eBackIt = eIt;
137     for ( int i = 0; i < nbEdges; ++i, ++eBackIt )
138       eBackIt->Reverse(); // reverse edge
139     // reverse list
140     --eBackIt;
141     while ( eIt != eBackIt )
142     {
143       std::swap( *eIt, *eBackIt );
144       SHOW_LIST("# AFTER SWAP", edges)
145         if ( (++eIt) != eBackIt )
146           --eBackIt;
147     }
148     SHOW_LIST("ATFER REVERSE", edges)
149   }
150
151   //================================================================================
152   /*!
153    * \brief Check if propagation is possible
154     * \param theMesh1 - source mesh
155     * \param theMesh2 - target mesh
156     * \retval bool - true if possible
157    */
158   //================================================================================
159
160   bool IsPropagationPossible( SMESH_Mesh* theMesh1, SMESH_Mesh* theMesh2 )
161   {
162     if ( theMesh1 != theMesh2 ) {
163       TopoDS_Shape mainShape1 = theMesh1->GetMeshDS()->ShapeToMesh();
164       TopoDS_Shape mainShape2 = theMesh2->GetMeshDS()->ShapeToMesh();
165       return mainShape1.IsSame( mainShape2 );
166     }
167     return true;
168   }
169
170   //================================================================================
171   /*!
172    * \brief Fix up association of edges in faces by possible propagation
173     * \param nbEdges - nb of edges in an outer wire
174     * \param edges1 - edges of one face
175     * \param edges2 - matching edges of another face
176     * \param theMesh1 - mesh 1
177     * \param theMesh2 - mesh 2
178     * \retval bool - true if association was fixed
179    */
180   //================================================================================
181
182   bool FixAssocByPropagation( const int             nbEdges,
183                               list< TopoDS_Edge > & edges1,
184                               list< TopoDS_Edge > & edges2,
185                               SMESH_Mesh*           theMesh1,
186                               SMESH_Mesh*           theMesh2)
187   {
188     if ( nbEdges == 2 && IsPropagationPossible( theMesh1, theMesh2 ) )
189     {
190       list< TopoDS_Edge >::iterator eIt2 = ++edges2.begin(); // 2nd edge of the 2nd face
191       TopoDS_Edge edge2 = HERE::GetPropagationEdge( theMesh1, *eIt2, edges1.front() ).second;
192       if ( !edge2.IsNull() ) { // propagation found for the second edge
193         Reverse( edges2, nbEdges );
194         return true;
195       }
196     }
197     return false;
198   }
199
200   //================================================================================
201   /*!
202    * \brief Look for a group containing a target shape and similar to a source group
203     * \param tgtShape - target edge or face
204     * \param tgtMesh1 - target mesh
205     * \param srcGroup - source group
206     * \retval TopoDS_Shape - found target group
207    */
208   //================================================================================
209
210   TopoDS_Shape FindGroupContaining(const TopoDS_Shape& tgtShape,
211                                    const SMESH_Mesh*   tgtMesh1,
212                                    const TopoDS_Shape& srcGroup)
213   {
214     list<SMESH_subMesh*> subMeshes = tgtMesh1->GetGroupSubMeshesContaining(tgtShape);
215     list<SMESH_subMesh*>::iterator sm = subMeshes.begin();
216     int type, last = TopAbs_SHAPE;
217     StdMeshers_ProjectionUtils util;
218     for ( ; sm != subMeshes.end(); ++sm ) {
219       const TopoDS_Shape & group = (*sm)->GetSubShape();
220       // check if group is similar to srcGroup
221       for ( type = srcGroup.ShapeType(); type < last; ++type)
222         if ( util.Count( srcGroup, (TopAbs_ShapeEnum)type, 0) !=
223              util.Count( group,    (TopAbs_ShapeEnum)type, 0))
224           break;
225       if ( type == last )
226         return group;
227     }
228     return TopoDS_Shape();
229   }
230
231   //================================================================================
232   /*!
233    * \brief Find association of groups at top and bottom of prism
234    */
235   //================================================================================
236
237   bool AssocGroupsByPropagation(const TopoDS_Shape&   theGroup1,
238                                 const TopoDS_Shape&   theGroup2,
239                                 SMESH_Mesh&           theMesh,
240                                 HERE::TShapeShapeMap& theMap)
241   {
242     // If groups are on top and bottom of prism then we can associate
243     // them using "vertical" (or "side") edges and faces of prism since
244     // they connect corresponding vertices and edges of groups.
245
246     TopTools_IndexedMapOfShape subshapes1, subshapes2;
247     TopExp::MapShapes( theGroup1, subshapes1 );
248     TopExp::MapShapes( theGroup2, subshapes2 );
249     TopTools_ListIteratorOfListOfShape ancestIt;
250
251     // Iterate on vertices of group1 to find corresponding vertices in group2
252     // and associate adjacent edges and faces
253
254     TopTools_MapOfShape verticShapes;
255     TopExp_Explorer vExp1( theGroup1, TopAbs_VERTEX );
256     for ( ; vExp1.More(); vExp1.Next() )
257     {
258       const TopoDS_Vertex& v1 = TopoDS::Vertex( vExp1.Current() );
259       if ( theMap.IsBound( v1 )) continue; // already processed
260
261       // Find "vertical" edge ending in v1 and whose other vertex belongs to group2
262       TopoDS_Shape verticEdge, v2;
263       ancestIt.Initialize( theMesh.GetAncestors( v1 ));
264       for ( ; verticEdge.IsNull() && ancestIt.More(); ancestIt.Next() )
265       {
266         if ( ancestIt.Value().ShapeType() != TopAbs_EDGE ) continue;
267         v2 = HERE::GetNextVertex( TopoDS::Edge( ancestIt.Value() ), v1 );
268         if ( subshapes2.Contains( v2 ))
269           verticEdge = ancestIt.Value();
270       }
271       if ( verticEdge.IsNull() )
272         return false;
273
274       HERE::InsertAssociation( v1, v2, theMap);
275
276       // Associate edges by vertical faces sharing the found vertical edge
277       ancestIt.Initialize( theMesh.GetAncestors( verticEdge ) );
278       for ( ; ancestIt.More(); ancestIt.Next() )
279       {
280         if ( ancestIt.Value().ShapeType() != TopAbs_FACE ) continue;
281         if ( !verticShapes.Add( ancestIt.Value() )) continue;
282         const TopoDS_Face& face = TopoDS::Face( ancestIt.Value() );
283
284         // get edges of the face
285         TopoDS_Edge edgeGr1, edgeGr2, verticEdge2;
286         list< TopoDS_Edge > edges;    list< int > nbEdgesInWire;
287         SMESH_Block::GetOrderedEdges( face, v1, edges, nbEdgesInWire);
288         if ( nbEdgesInWire.front() != 4 )
289           return _StoreBadShape( face );
290         list< TopoDS_Edge >::iterator edge = edges.begin();
291         if ( verticEdge.IsSame( *edge )) {
292           edgeGr2     = *(++edge);
293           verticEdge2 = *(++edge);
294           edgeGr1     = *(++edge);
295         } else {
296           edgeGr1     = *(edge++);
297           verticEdge2 = *(edge++);
298           edgeGr2     = *(edge++);
299         }
300
301         HERE::InsertAssociation( edgeGr1, edgeGr2.Reversed(), theMap);
302       }
303     }
304
305     // Associate faces
306     TopoDS_Iterator gr1It( theGroup1 );
307     if ( gr1It.Value().ShapeType() == TopAbs_FACE )
308     {
309       // find a boundary edge of group1 to start from
310       TopoDS_Shape bndEdge;
311       TopExp_Explorer edgeExp1( theGroup1, TopAbs_EDGE );
312       for ( ; bndEdge.IsNull() && edgeExp1.More(); edgeExp1.Next())
313         if ( HERE::IsBoundaryEdge( TopoDS::Edge( edgeExp1.Current()), theGroup1, theMesh ))
314           bndEdge = edgeExp1.Current();
315       if ( bndEdge.IsNull() )
316         return false;
317
318       list< TopoDS_Shape > edges(1, bndEdge);
319       list< TopoDS_Shape >::iterator edge1 = edges.begin();
320       for ( ; edge1 != edges.end(); ++edge1 )
321       {
322         // there must be one or zero not associated faces between ancestors of edge
323         // belonging to theGroup1
324         TopoDS_Shape face1;
325         ancestIt.Initialize( theMesh.GetAncestors( *edge1 ) );
326         for ( ; ancestIt.More() && face1.IsNull(); ancestIt.Next() ) {
327           if ( ancestIt.Value().ShapeType() == TopAbs_FACE &&
328                !theMap.IsBound( ancestIt.Value() ) &&
329                subshapes1.Contains( ancestIt.Value() ))
330             face1 = ancestIt.Value();
331
332           // add edges of face1 to start searching for adjacent faces from
333           for ( TopExp_Explorer e(face1, TopAbs_EDGE); e.More(); e.Next())
334             if ( !edge1->IsSame( e.Current() ))
335               edges.push_back( e.Current() );
336         }
337         if ( !face1.IsNull() ) {
338           // find the corresponding face of theGroup2
339           TopoDS_Shape edge2 = theMap( *edge1 );
340           TopoDS_Shape face2;
341           ancestIt.Initialize( theMesh.GetAncestors( edge2 ) );
342           for ( ; ancestIt.More() && face2.IsNull(); ancestIt.Next() ) {
343             if ( ancestIt.Value().ShapeType() == TopAbs_FACE &&
344                  !theMap.IsBound( ancestIt.Value() ) &&
345                  subshapes2.Contains( ancestIt.Value() ))
346               face2 = ancestIt.Value();
347           }
348           if ( face2.IsNull() )
349             return false;
350
351           HERE::InsertAssociation( face1, face2, theMap);
352         }
353       }
354     }
355     return true;
356   }
357
358   //================================================================================
359   /*!
360    * \brief Return true if uv position of the vIndex-th vertex of edge on face is close
361    * enough to given uv 
362    */
363   //================================================================================
364
365   bool sameVertexUV( const TopoDS_Edge& edge,
366                      const TopoDS_Face& face,
367                      const int&         vIndex,
368                      const gp_Pnt2d&    uv,
369                      const double&      tol2d )
370   {
371     TopoDS_Vertex VV[2];
372     TopExp::Vertices( edge, VV[0], VV[1], true);
373     gp_Pnt2d v1UV = BRep_Tool::Parameters( VV[vIndex], face);
374     double dist2d = v1UV.Distance( uv );
375     return dist2d < tol2d;
376   }
377
378 } // namespace
379
380 //=======================================================================
381 /*!
382  * \brief Looks for association of all subshapes of two shapes
383  * \param theShape1 - shape 1
384  * \param theMesh1 - mesh built on shape 1
385  * \param theShape2 - shape 2
386  * \param theMesh2 - mesh built on shape 2
387  * \param theAssociation - association map to be filled that may
388  *                         contain association of one or two pairs of vertices
389  * \retval bool - true if association found
390  */
391 //=======================================================================
392
393 bool StdMeshers_ProjectionUtils::FindSubShapeAssociation(const TopoDS_Shape& theShape1,
394                                                          SMESH_Mesh*         theMesh1,
395                                                          const TopoDS_Shape& theShape2,
396                                                          SMESH_Mesh*         theMesh2,
397                                                          TShapeShapeMap &    theMap)
398 {
399   // Structure of this long function is following
400   // 1) Group->group projection: theShape1 is a group member,
401   //    theShape2 is a group. We find a group theShape1 is in and recall self.
402   // 2) Accosiate same shapes with different location (partners).
403   // 3) If vertex association is given, perform accosiation according to shape type:
404   //       switch ( ShapeType ) {
405   //         case TopAbs_EDGE:
406   //         case ...:
407   //       }
408   // 4) else try to accosiate in different ways:
409   //       a) accosiate shapes by propagation and other simple cases
410   //            switch ( ShapeType ) {
411   //            case TopAbs_EDGE:
412   //            case ...:
413   //            }
414   //       b) find association of a couple of vertices and recall self.
415   //
416
417   theMeshDS[0] = theMesh1->GetMeshDS(); // debug
418   theMeshDS[1] = theMesh2->GetMeshDS();
419
420   // =================================================================================
421   // 1) Is it the case of associating a group member -> another group? (PAL16202, 16203)
422   // =================================================================================
423   if ( theShape1.ShapeType() != theShape2.ShapeType() ) {
424     TopoDS_Shape group1, group2;
425     if ( theShape1.ShapeType() == TopAbs_COMPOUND ) {
426       group1 = theShape1;
427       group2 = FindGroupContaining( theShape2, theMesh2, group1 );
428     }
429     else if ( theShape2.ShapeType() == TopAbs_COMPOUND ) {
430       group2 = theShape2;
431       group1 = FindGroupContaining( theShape1, theMesh1, group2 );
432     }
433     if ( group1.IsNull() || group2.IsNull() )
434       RETURN_BAD_RESULT("Different shape types");
435     // Associate compounds
436     return FindSubShapeAssociation(group1, theMesh1, group2, theMesh2, theMap );
437   }
438
439   bool bidirect = ( !theShape1.IsSame( theShape2 ));
440
441   // ============
442   // 2) Is partner?
443   // ============
444   bool partner = theShape1.IsPartner( theShape2 );
445   TopTools_DataMapIteratorOfDataMapOfShapeShape vvIt( theMap );
446   for ( ; partner && vvIt.More(); vvIt.Next() )
447     partner = vvIt.Key().IsPartner( vvIt.Value() );
448
449   if ( partner ) // Same shape with different location
450   {
451     // recursively associate all subshapes of theShape1 and theShape2
452     typedef list< pair< TopoDS_Shape, TopoDS_Shape > > TShapePairsList;
453     TShapePairsList shapesQueue( 1, make_pair( theShape1, theShape2 ));
454     TShapePairsList::iterator s1_s2 = shapesQueue.begin();
455     for ( ; s1_s2 != shapesQueue.end(); ++s1_s2 )
456     {
457       InsertAssociation( s1_s2->first, s1_s2->second, theMap, bidirect);
458       TopoDS_Iterator s1It( s1_s2->first), s2It( s1_s2->second );
459       for ( ; s1It.More(); s1It.Next(), s2It.Next() )
460         shapesQueue.push_back( make_pair( s1It.Value(), s2It.Value() ));
461     }
462     return true;
463   }
464
465   if ( !theMap.IsEmpty() )
466   {
467     //======================================================================
468     // 3) HAS initial vertex association
469     //======================================================================
470     switch ( theShape1.ShapeType() ) {
471       // ----------------------------------------------------------------------
472     case TopAbs_EDGE: { // TopAbs_EDGE
473       // ----------------------------------------------------------------------
474       if ( theMap.Extent() != 2 )
475         RETURN_BAD_RESULT("Wrong map extent " << theMap.Extent() );
476       TopoDS_Edge edge1 = TopoDS::Edge( theShape1 );
477       TopoDS_Edge edge2 = TopoDS::Edge( theShape2 );
478       TopoDS_Vertex VV1[2], VV2[2];
479       TopExp::Vertices( edge1, VV1[0], VV1[1] );
480       TopExp::Vertices( edge2, VV2[0], VV2[1] );
481       int i1 = 0, i2 = 0;
482       if ( theMap.IsBound( VV1[ i1 ] )) i1 = 1;
483       if ( theMap.IsBound( VV2[ i2 ] )) i2 = 1;
484       InsertAssociation( VV1[ i1 ], VV2[ i2 ], theMap, bidirect);
485       InsertAssociation( theShape1, theShape2, theMap, bidirect );
486       return true;
487     }
488       // ----------------------------------------------------------------------
489     case TopAbs_FACE: { // TopAbs_FACE
490       // ----------------------------------------------------------------------
491       TopoDS_Face face1 = TopoDS::Face( theShape1 );
492       TopoDS_Face face2 = TopoDS::Face( theShape2 );
493
494       TopoDS_Vertex VV1[2], VV2[2];
495       // find a not closed edge of face1 both vertices of which are associated
496       int nbEdges = 0;
497       TopExp_Explorer exp ( face1, TopAbs_EDGE );
498       for ( ; VV2[ 1 ].IsNull() && exp.More(); exp.Next(), ++nbEdges ) {
499         TopExp::Vertices( TopoDS::Edge( exp.Current() ), VV1[0], VV1[1] );
500         if ( theMap.IsBound( VV1[0] ) ) {
501           VV2[ 0 ] = TopoDS::Vertex( theMap( VV1[0] ));
502           if ( theMap.IsBound( VV1[1] ) && !VV1[0].IsSame( VV1[1] ))
503             VV2[ 1 ] = TopoDS::Vertex( theMap( VV1[1] ));
504         }
505       }
506       if ( VV2[ 1 ].IsNull() ) { // 2 bound vertices not found
507         if ( nbEdges > 1 ) {
508           RETURN_BAD_RESULT("2 bound vertices not found" );
509         } else {
510           VV2[ 1 ] = VV2[ 0 ];
511         }
512       }
513       list< TopoDS_Edge > edges1, edges2;
514       int nbE = FindFaceAssociation( face1, VV1, face2, VV2, edges1, edges2 );
515       if ( !nbE ) RETURN_BAD_RESULT("FindFaceAssociation() failed");
516       FixAssocByPropagation( nbE, edges1, edges2, theMesh1, theMesh2 );
517
518       list< TopoDS_Edge >::iterator eIt1 = edges1.begin();
519       list< TopoDS_Edge >::iterator eIt2 = edges2.begin();
520       for ( ; eIt1 != edges1.end(); ++eIt1, ++eIt2 )
521       {
522         InsertAssociation( *eIt1, *eIt2, theMap, bidirect);
523         VV1[0] = TopExp::FirstVertex( *eIt1, true );
524         VV2[0] = TopExp::FirstVertex( *eIt2, true );
525         InsertAssociation( VV1[0], VV2[0], theMap, bidirect);
526       }
527       InsertAssociation( theShape1, theShape2, theMap, bidirect );
528       return true;
529     }
530       // ----------------------------------------------------------------------
531     case TopAbs_SHELL: // TopAbs_SHELL, TopAbs_SOLID
532     case TopAbs_SOLID: {
533       // ----------------------------------------------------------------------
534       TopoDS_Vertex VV1[2], VV2[2];
535       // try to find a not closed edge of shape1 both vertices of which are associated
536       TopoDS_Edge edge1;
537       TopExp_Explorer exp ( theShape1, TopAbs_EDGE );
538       for ( ; VV2[ 1 ].IsNull() && exp.More(); exp.Next() ) {
539         edge1 = TopoDS::Edge( exp.Current() );
540         TopExp::Vertices( edge1 , VV1[0], VV1[1] );
541         if ( theMap.IsBound( VV1[0] )) {
542           VV2[ 0 ] = TopoDS::Vertex( theMap( VV1[0] ));
543           if ( theMap.IsBound( VV1[1] ) && !VV1[0].IsSame( VV1[1] ))
544             VV2[ 1 ] = TopoDS::Vertex( theMap( VV1[1] ));
545         }
546       }
547       if ( VV2[ 1 ].IsNull() ) // 2 bound vertices not found
548         RETURN_BAD_RESULT("2 bound vertices not found" );
549       // get an edge2 of theShape2 corresponding to edge1
550       TopoDS_Edge edge2 = GetEdgeByVertices( theMesh2, VV2[ 0 ], VV2[ 1 ]);
551       if ( edge2.IsNull() )
552         RETURN_BAD_RESULT("GetEdgeByVertices() failed");
553
554       // build map of edge to faces if shapes are not subshapes of main ones
555       bool isSubOfMain = false;
556       if ( SMESHDS_SubMesh * sm = theMesh1->GetMeshDS()->MeshElements( theShape1 ))
557         isSubOfMain = !sm->IsComplexSubmesh();
558       else
559         isSubOfMain = theMesh1->GetMeshDS()->ShapeToIndex( theShape1 );
560       TAncestorMap e2f1, e2f2;
561       const TAncestorMap& edgeToFace1 = isSubOfMain ? theMesh1->GetAncestorMap() : e2f1;
562       const TAncestorMap& edgeToFace2 = isSubOfMain ? theMesh2->GetAncestorMap() : e2f2;
563       if (!isSubOfMain) {
564         TopExp::MapShapesAndAncestors( theShape1, TopAbs_EDGE, TopAbs_FACE, e2f1 );
565         TopExp::MapShapesAndAncestors( theShape2, TopAbs_EDGE, TopAbs_FACE, e2f2 );
566         if ( !edgeToFace1.Contains( edge1 ))
567           RETURN_BAD_RESULT("edge1 does not belong to theShape1");
568         if ( !edgeToFace2.Contains( edge2 ))
569           RETURN_BAD_RESULT("edge2 does not belong to theShape2");
570       }
571       //
572       // Look for 2 corresponing faces:
573       //
574       TopoDS_Shape F1, F2;
575
576       // get a face sharing edge1 (F1)
577       TopoDS_Shape FF2[2];
578       TopTools_ListIteratorOfListOfShape ancestIt1( edgeToFace1.FindFromKey( edge1 ));
579       for ( ; F1.IsNull() && ancestIt1.More(); ancestIt1.Next() )
580         if ( ancestIt1.Value().ShapeType() == TopAbs_FACE )
581           F1 = ancestIt1.Value().Oriented( TopAbs_FORWARD );
582       if ( F1.IsNull() )
583         RETURN_BAD_RESULT(" Face1 not found");
584
585       // get 2 faces sharing edge2 (one of them is F2)
586       TopTools_ListIteratorOfListOfShape ancestIt2( edgeToFace2.FindFromKey( edge2 ));
587       for ( int i = 0; FF2[1].IsNull() && ancestIt2.More(); ancestIt2.Next() )
588         if ( ancestIt2.Value().ShapeType() == TopAbs_FACE )
589           FF2[ i++ ] = ancestIt2.Value().Oriented( TopAbs_FORWARD );
590
591       // get oriented edge1 and edge2 from F1 and FF2[0]
592       for ( exp.Init( F1, TopAbs_EDGE ); exp.More(); exp.Next() )
593         if ( edge1.IsSame( exp.Current() )) {
594           edge1 = TopoDS::Edge( exp.Current() );
595           break;
596         }
597       for ( exp.Init( FF2[ 0 ], TopAbs_EDGE ); exp.More(); exp.Next() )
598         if ( edge2.IsSame( exp.Current() )) {
599           edge2 = TopoDS::Edge( exp.Current() );
600           break;
601         }
602
603       // compare first vertices of edge1 and edge2
604       TopExp::Vertices( edge1, VV1[0], VV1[1], true );
605       TopExp::Vertices( edge2, VV2[0], VV2[1], true );
606       F2 = FF2[ 0 ]; // (F2 !)
607       if ( !VV1[ 0 ].IsSame( theMap( VV2[ 0 ]))) {
608         edge2.Reverse();
609         if ( FF2[ 1 ].IsNull() )
610           F2.Reverse();
611         else
612           F2 = FF2[ 1 ];
613       }
614
615       TopTools_MapOfShape boundEdges;
616
617       // association of face subshapes and neighbour faces
618       list< pair < TopoDS_Face, TopoDS_Edge > > FE1, FE2;
619       list< pair < TopoDS_Face, TopoDS_Edge > >::iterator fe1, fe2;
620       FE1.push_back( make_pair( TopoDS::Face( F1 ), edge1 ));
621       FE2.push_back( make_pair( TopoDS::Face( F2 ), edge2 ));
622       for ( fe1 = FE1.begin(), fe2 = FE2.begin(); fe1 != FE1.end(); ++fe1, ++fe2 )
623       {
624         const TopoDS_Face& face1 = fe1->first;
625         if ( theMap.IsBound( face1 ) ) continue;
626         const TopoDS_Face& face2 = fe2->first;
627         edge1 = fe1->second;
628         edge2 = fe2->second;
629         TopExp::Vertices( edge1, VV1[0], VV1[1], true );
630         TopExp::Vertices( edge2, VV2[0], VV2[1], true );
631         list< TopoDS_Edge > edges1, edges2;
632         int nbE = FindFaceAssociation( face1, VV1, face2, VV2, edges1, edges2 );
633         if ( !nbE ) RETURN_BAD_RESULT("FindFaceAssociation() failed");
634         InsertAssociation( face1, face2, theMap, bidirect); // assoc faces
635         MESSAGE("Assoc FACE " << theMesh1->GetMeshDS()->ShapeToIndex( face1 )<<
636                 " to "        << theMesh2->GetMeshDS()->ShapeToIndex( face2 ));
637         if ( nbE == 2 && (edge1.IsSame( edges1.front())) != (edge2.IsSame( edges2.front())))
638         {
639           Reverse( edges2, nbE );
640         }
641         list< TopoDS_Edge >::iterator eIt1 = edges1.begin();
642         list< TopoDS_Edge >::iterator eIt2 = edges2.begin();
643         for ( ; eIt1 != edges1.end(); ++eIt1, ++eIt2 )
644         {
645           if ( !boundEdges.Add( *eIt1 )) continue; // already associated
646           InsertAssociation( *eIt1, *eIt2, theMap, bidirect);  // assoc edges
647           MESSAGE("Assoc edge " << theMesh1->GetMeshDS()->ShapeToIndex( *eIt1 )<<
648                   " to "        << theMesh2->GetMeshDS()->ShapeToIndex( *eIt2 ));
649           VV1[0] = TopExp::FirstVertex( *eIt1, true );
650           VV2[0] = TopExp::FirstVertex( *eIt2, true );
651           InsertAssociation( VV1[0], VV2[0], theMap, bidirect); // assoc vertices
652           MESSAGE("Assoc vertex " << theMesh1->GetMeshDS()->ShapeToIndex( VV1[0] )<<
653                   " to "          << theMesh2->GetMeshDS()->ShapeToIndex( VV2[0] ));
654
655           // add adjacent faces to process
656           TopoDS_Face nextFace1 = GetNextFace( edgeToFace1, *eIt1, face1 );
657           TopoDS_Face nextFace2 = GetNextFace( edgeToFace2, *eIt2, face2 );
658           if ( !nextFace1.IsNull() && !nextFace2.IsNull() ) {
659             FE1.push_back( make_pair( nextFace1, *eIt1 ));
660             FE2.push_back( make_pair( nextFace2, *eIt2 ));
661           }
662         }
663       }
664       InsertAssociation( theShape1, theShape2, theMap, bidirect );
665       return true;
666     }
667       // ----------------------------------------------------------------------
668     case TopAbs_COMPOUND: { // GROUP
669       // ----------------------------------------------------------------------
670       // Maybe groups contain only one member
671       TopoDS_Iterator it1( theShape1 ), it2( theShape2 );
672       TopAbs_ShapeEnum memberType = it1.Value().ShapeType();
673       int nbMembers = Count( theShape1, memberType, true );
674       if ( nbMembers == 0 ) return true;
675       if ( nbMembers == 1 ) {
676         return FindSubShapeAssociation( it1.Value(), theMesh1, it2.Value(), theMesh2, theMap );
677       }
678       // Try to make shells of faces
679       //
680       BRep_Builder builder;
681       TopoDS_Shell shell1, shell2;
682       builder.MakeShell(shell1); builder.MakeShell(shell2);
683       if ( memberType == TopAbs_FACE ) {
684         // just add faces of groups to shells
685         for (; it1.More(); it1.Next(), it2.Next() )
686           builder.Add( shell1, it1.Value() ), builder.Add( shell2, it2.Value() );
687       }
688       else if ( memberType == TopAbs_EDGE ) {
689         // Try to add faces sharing more than one edge of a group or
690         // sharing all its vertices with the group
691         TopTools_IndexedMapOfShape groupVertices[2];
692         TopExp::MapShapes( theShape1, TopAbs_VERTEX, groupVertices[0]);
693         TopExp::MapShapes( theShape2, TopAbs_VERTEX, groupVertices[1]);
694         //
695         TopTools_MapOfShape groupEdges[2], addedFaces[2];
696         bool hasInitAssoc = (!theMap.IsEmpty()), initAssocOK = !hasInitAssoc;
697         for (; it1.More(); it1.Next(), it2.Next() ) {
698           groupEdges[0].Add( it1.Value() );
699           groupEdges[1].Add( it2.Value() );
700           if ( !initAssocOK ) {
701             // for shell association there must be an edge with both vertices bound
702             TopoDS_Vertex v1, v2;
703             TopExp::Vertices( TopoDS::Edge( it1.Value()), v1, v2 );
704             initAssocOK = ( theMap.IsBound( v1 ) && theMap.IsBound( v2 ));
705           }
706         }
707         for (int is2ndGroup = 0; initAssocOK && is2ndGroup < 2; ++is2ndGroup) {
708           const TopoDS_Shape& group = is2ndGroup ? theShape2: theShape1;
709           SMESH_Mesh*         mesh  = is2ndGroup ? theMesh2 : theMesh1;
710           TopoDS_Shell&       shell = is2ndGroup ? shell2   : shell1;
711           for ( TopoDS_Iterator it( group ); it.More(); it.Next() ) {
712             const TopoDS_Edge& edge = TopoDS::Edge( it.Value() );
713             TopoDS_Face face;
714             for ( int iF = 0; iF < 2; ++iF ) { // loop on 2 faces sharing edge
715               face = GetNextFace(mesh->GetAncestorMap(), edge, face);
716               if ( !face.IsNull() ) {
717                 int nbGroupEdges = 0;
718                 for ( TopExp_Explorer f( face, TopAbs_EDGE ); f.More(); f.Next())
719                   if ( groupEdges[ is2ndGroup ].Contains( f.Current() ))
720                     if ( ++nbGroupEdges > 1 )
721                       break;
722                 bool add = (nbGroupEdges > 1 || Count( face, TopAbs_EDGE, true ) == 1 );
723                 if ( !add ) {
724                   add = true;
725                   for ( TopExp_Explorer v( face, TopAbs_VERTEX ); add && v.More(); v.Next())
726                     add = groupVertices[ is2ndGroup ].Contains( v.Current() );
727                 }
728                 if ( add && addedFaces[ is2ndGroup ].Add( face ))
729                   builder.Add( shell, face );
730               }
731             }
732           }
733         }
734       } else {
735         RETURN_BAD_RESULT("Unexpected group type");
736       }
737       // Associate shells
738       //
739       int nbFaces1 = Count( shell1, TopAbs_FACE, 0 );
740       int nbFaces2 = Count( shell2, TopAbs_FACE, 0 );
741       if ( nbFaces1 != nbFaces2 )
742         RETURN_BAD_RESULT("Different nb of faces found for shells");
743       if ( nbFaces1 > 0 ) {
744         bool ok = false;
745         if ( nbFaces1 == 1 ) {
746           TopoDS_Shape F1 = TopoDS_Iterator( shell1 ).Value();
747           TopoDS_Shape F2 = TopoDS_Iterator( shell2 ).Value();
748           ok = FindSubShapeAssociation( F1, theMesh1, F2, theMesh2, theMap );
749         }
750         else {
751           ok = FindSubShapeAssociation(shell1, theMesh1, shell2, theMesh2, theMap );
752         }
753         // Check if all members are mapped 
754         if ( ok ) {
755           TopTools_MapOfShape boundMembers[2];
756           TopoDS_Iterator mIt;
757           for ( mIt.Initialize( theShape1 ); mIt.More(); mIt.Next())
758             if ( theMap.IsBound( mIt.Value() )) {
759               boundMembers[0].Add( mIt.Value() );
760               boundMembers[1].Add( theMap( mIt.Value() ));
761             }
762           if ( boundMembers[0].Extent() != nbMembers ) {
763             // make compounds of not bound members
764             TopoDS_Compound comp[2];
765             for ( int is2ndGroup = 0; is2ndGroup < 2; ++is2ndGroup ) {
766               builder.MakeCompound( comp[is2ndGroup] );
767               for ( mIt.Initialize( is2ndGroup ? theShape2:theShape1 ); mIt.More(); mIt.Next())
768                 if ( ! boundMembers[ is2ndGroup ].Contains( mIt.Value() ))
769                   builder.Add( comp[ is2ndGroup ], mIt.Value() );
770             }
771             // check if theMap contains initial association for the comp's
772             bool hasInitialAssoc = false;
773             if ( memberType == TopAbs_EDGE ) {
774               for ( TopExp_Explorer v( comp[0], TopAbs_VERTEX ); v.More(); v.Next())
775                 if ( theMap.IsBound( v.Current() )) {
776                   hasInitialAssoc = true;
777                   break;
778                 }
779             }
780             if ( hasInitialAssoc == bool( !theMap.IsEmpty() ))
781               ok = FindSubShapeAssociation( comp[0], theMesh1, comp[1], theMesh2, theMap );
782             else {
783               TShapeShapeMap tmpMap;
784               ok = FindSubShapeAssociation( comp[0], theMesh1, comp[1], theMesh2, tmpMap );
785               if ( ok ) {
786                 TopTools_DataMapIteratorOfDataMapOfShapeShape mapIt( tmpMap );
787                 for ( ; mapIt.More(); mapIt.Next() )
788                   theMap.Bind( mapIt.Key(), mapIt.Value());
789               }
790             }
791           }
792         }
793         return ok;
794       }
795       // Each edge of an edge group is shared by own faces
796       // ------------------------------------------------------------------
797       //
798       // map vertices to edges sharing them, avoid doubling edges in lists
799       TopTools_DataMapOfShapeListOfShape v2e[2];
800       for (int isFirst = 0; isFirst < 2; ++isFirst ) {
801         const TopoDS_Shape& group = isFirst ? theShape1 : theShape2;
802         TopTools_DataMapOfShapeListOfShape& veMap = v2e[ isFirst ? 0 : 1 ];
803         TopTools_MapOfShape addedEdges;
804         for ( TopExp_Explorer e( group, TopAbs_EDGE ); e.More(); e.Next() ) {
805           const TopoDS_Shape& edge = e.Current();
806           if ( addedEdges.Add( edge )) {
807             for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next()) {
808               const TopoDS_Shape& vertex = v.Current();
809               if ( !veMap.IsBound( vertex )) {
810                 TopTools_ListOfShape l;
811                 veMap.Bind( vertex, l );
812               }
813               veMap( vertex ).Append( edge );
814             }
815           }
816         }   
817       }
818       while ( !v2e[0].IsEmpty() )
819       {
820         // find a bound vertex
821         TopoDS_Vertex V[2];
822         TopTools_DataMapIteratorOfDataMapOfShapeListOfShape v2eIt( v2e[0] );
823         for ( ; v2eIt.More(); v2eIt.Next())
824           if ( theMap.IsBound( v2eIt.Key() )) {
825             V[0] = TopoDS::Vertex( v2eIt.Key() );
826             V[1] = TopoDS::Vertex( theMap( V[0] ));
827             break;
828           }
829         if ( V[0].IsNull() )
830           RETURN_BAD_RESULT("No more bound vertices");
831
832         while ( !V[0].IsNull() && v2e[0].IsBound( V[0] )) {
833           TopTools_ListOfShape& edges0 = v2e[0]( V[0] );
834           TopTools_ListOfShape& edges1 = v2e[1]( V[1] );
835           int nbE0 = edges0.Extent(), nbE1 = edges1.Extent();
836           if ( nbE0 != nbE1 )
837             RETURN_BAD_RESULT("Different nb of edges: "<< nbE0 << " != " << nbE1);
838
839           if ( nbE0 == 1 )
840           {
841             TopoDS_Edge e0 = TopoDS::Edge( edges0.First() );
842             TopoDS_Edge e1 = TopoDS::Edge( edges1.First() );
843             v2e[0].UnBind( V[0] );
844             v2e[1].UnBind( V[1] );
845             InsertAssociation( e0, e1, theMap, bidirect );
846             MESSAGE("Assoc edge " << theMesh1->GetMeshDS()->ShapeToIndex( e0 )<<
847                     " to "        << theMesh2->GetMeshDS()->ShapeToIndex( e1 ));
848             V[0] = GetNextVertex( e0, V[0] );
849             V[1] = GetNextVertex( e1, V[1] );
850             if ( !V[0].IsNull() ) {
851               InsertAssociation( V[0], V[1], theMap, bidirect );
852               MESSAGE("Assoc vertex " << theMesh1->GetMeshDS()->ShapeToIndex( V[0] )<<
853                       " to "          << theMesh2->GetMeshDS()->ShapeToIndex( V[1] ));
854             }
855           }
856           else if ( nbE0 == 2 )
857           {
858             // one of edges must have both ends bound
859             TopoDS_Vertex v0e0 = GetNextVertex( TopoDS::Edge( edges0.First() ), V[0] );
860             TopoDS_Vertex v1e0 = GetNextVertex( TopoDS::Edge( edges0.Last() ),  V[0] );
861             TopoDS_Vertex v0e1 = GetNextVertex( TopoDS::Edge( edges1.First() ), V[1] );
862             TopoDS_Vertex v1e1 = GetNextVertex( TopoDS::Edge( edges1.Last() ),  V[1] );
863             TopoDS_Shape e0b, e1b, e0n, e1n, v1b; // bound and not-bound
864             TopoDS_Vertex v0n, v1n;
865             if ( theMap.IsBound( v0e0 )) {
866               v0n = v1e0; e0b = edges0.First(); e0n = edges0.Last(); v1b = theMap( v0e0 );
867             } else if ( theMap.IsBound( v1e0 )) {
868               v0n = v0e0; e0n = edges0.First(); e0b = edges0.Last(); v1b = theMap( v1e0 );
869             } else {
870               RETURN_BAD_RESULT("None of vertices bound");
871             }
872             if ( v1b.IsSame( v1e1 )) {
873               v1n = v0e1; e1n = edges1.First(); e1b = edges1.Last();
874             } else {
875               v1n = v1e1; e1b = edges1.First(); e1n = edges1.Last();
876             }
877             InsertAssociation( e0b, e1b, theMap, bidirect );
878             InsertAssociation( e0n, e1n, theMap, bidirect );
879             InsertAssociation( v0n, v1n, theMap, bidirect );
880             MESSAGE("Assoc edge " << theMesh1->GetMeshDS()->ShapeToIndex( e0b )<<
881                     " to "        << theMesh2->GetMeshDS()->ShapeToIndex( e1b ));
882             MESSAGE("Assoc edge " << theMesh1->GetMeshDS()->ShapeToIndex( e0n )<<
883                     " to "        << theMesh2->GetMeshDS()->ShapeToIndex( e1n ));
884             MESSAGE("Assoc vertex " << theMesh1->GetMeshDS()->ShapeToIndex( v0n )<<
885                     " to "          << theMesh2->GetMeshDS()->ShapeToIndex( v1n ));
886             v2e[0].UnBind( V[0] );
887             v2e[1].UnBind( V[1] );
888             V[0] = v0n;
889             V[1] = v1n;
890           }
891           else {
892             RETURN_BAD_RESULT("Not implemented");
893           }
894         }
895       } //while ( !v2e[0].IsEmpty() )
896       return true;
897     }
898
899     default:
900       RETURN_BAD_RESULT("Unexpected shape type");
901
902     } // end switch by shape type
903   } // end case of available initial vertex association
904
905   //======================================================================
906   // 4) NO INITIAL VERTEX ASSOCIATION
907   //======================================================================
908
909   switch ( theShape1.ShapeType() ) {
910
911   case TopAbs_EDGE: {
912     // ----------------------------------------------------------------------
913     TopoDS_Edge edge1 = TopoDS::Edge( theShape1 );
914     TopoDS_Edge edge2 = TopoDS::Edge( theShape2 );
915     if ( IsPropagationPossible( theMesh1, theMesh2 ))
916     {
917       TopoDS_Edge prpEdge = GetPropagationEdge( theMesh1, edge2, edge1 ).second;
918       if ( !prpEdge.IsNull() )
919       {
920         TopoDS_Vertex VV1[2], VV2[2];
921         TopExp::Vertices( edge1,   VV1[0], VV1[1], true );
922         TopExp::Vertices( prpEdge, VV2[0], VV2[1], true );
923         InsertAssociation( VV1[ 0 ], VV2[ 0 ], theMap, bidirect);
924         InsertAssociation( VV1[ 1 ], VV2[ 1 ], theMap, bidirect);
925         if ( VV1[0].IsSame( VV1[1] ) || // one of edges is closed
926              VV2[0].IsSame( VV2[1] ) )
927         {
928           InsertAssociation( edge1, prpEdge, theMap, bidirect); // insert with a proper orientation
929         }
930         InsertAssociation( theShape1, theShape2, theMap, bidirect );
931         return true; // done
932       }
933     }
934     if ( IsClosedEdge( edge1 ) && IsClosedEdge( edge2 ))
935     {
936       // TODO: find out a proper orientation (is it possible?)
937       InsertAssociation( edge1, edge2, theMap, bidirect); // insert with a proper orientation
938       InsertAssociation( TopExp::FirstVertex(edge1), TopExp::FirstVertex(edge2),
939                          theMap, bidirect);
940       InsertAssociation( theShape1, theShape2, theMap, bidirect );
941       return true; // done
942     }
943     break; // try by vertex closeness
944   }
945
946   case TopAbs_FACE: {
947     // ----------------------------------------------------------------------
948     if ( IsPropagationPossible( theMesh1, theMesh2 )) // try by propagation in one mesh
949     {
950       TopoDS_Face face1 = TopoDS::Face(theShape1);
951       TopoDS_Face face2 = TopoDS::Face(theShape2);
952       TopoDS_Edge edge1, edge2;
953       // get outer edge of theShape1
954       edge1 = TopoDS::Edge( OuterShape( face1, TopAbs_EDGE ));
955       // find out if any edge of face2 is a propagation edge of outer edge1
956       map<int,TopoDS_Edge> propag_edges; // use map to find the closest propagation edge
957       for ( TopExp_Explorer exp( face2, TopAbs_EDGE ); exp.More(); exp.Next() ) {
958         edge2 = TopoDS::Edge( exp.Current() );
959         pair<int,TopoDS_Edge> step_edge = GetPropagationEdge( theMesh1, edge2, edge1 );
960         if ( !step_edge.second.IsNull() ) { // propagation found
961           propag_edges.insert( step_edge );
962           if ( step_edge.first == 1 ) break; // most close found
963         }
964       }
965       if ( !propag_edges.empty() ) // propagation found
966       {
967         edge2 = propag_edges.begin()->second;
968         TopoDS_Vertex VV1[2], VV2[2];
969         TopExp::Vertices( edge1, VV1[0], VV1[1], true );
970         TopExp::Vertices( edge2, VV2[0], VV2[1], true );
971         list< TopoDS_Edge > edges1, edges2;
972         int nbE = FindFaceAssociation( face1, VV1, face2, VV2, edges1, edges2 );
973         if ( !nbE ) RETURN_BAD_RESULT("FindFaceAssociation() failed");
974         if ( nbE == 2 ) // only 2 edges
975         {
976           // take care of proper association of propagated edges
977           bool same1 = edge1.IsSame( edges1.front() );
978           bool same2 = edge2.IsSame( edges2.front() );
979           if ( same1 != same2 )
980             Reverse(edges2, nbE);
981         }
982         // store association
983         list< TopoDS_Edge >::iterator eIt1 = edges1.begin();
984         list< TopoDS_Edge >::iterator eIt2 = edges2.begin();
985         for ( ; eIt1 != edges1.end(); ++eIt1, ++eIt2 )
986         {
987           InsertAssociation( *eIt1, *eIt2, theMap, bidirect);
988           VV1[0] = TopExp::FirstVertex( *eIt1, true );
989           VV2[0] = TopExp::FirstVertex( *eIt2, true );
990           InsertAssociation( VV1[0], VV2[0], theMap, bidirect);
991         }
992         InsertAssociation( theShape1, theShape2, theMap, bidirect );
993         return true;
994       }
995     }
996     break; // try by vertex closeness
997   }
998   case TopAbs_COMPOUND: {
999     // ----------------------------------------------------------------------
1000     if ( IsPropagationPossible( theMesh1, theMesh2 )) {
1001
1002       // try to accosiate all using propagation
1003       if ( AssocGroupsByPropagation( theShape1, theShape2, *theMesh1, theMap ))
1004         return true;
1005
1006       // find a boundary edge for theShape1
1007       TopoDS_Edge E;
1008       for(TopExp_Explorer exp(theShape1, TopAbs_EDGE); exp.More(); exp.Next() ) {
1009         E = TopoDS::Edge( exp.Current() );
1010         if ( IsBoundaryEdge( E, theShape1, *theMesh1 ))
1011           break;
1012         else
1013           E.Nullify();
1014       }
1015       if ( E.IsNull() )
1016         break; // try by vertex closeness
1017
1018       // find association for vertices of edge E
1019       TopoDS_Vertex VV1[2], VV2[2];
1020       for(TopExp_Explorer eexp(E, TopAbs_VERTEX); eexp.More(); eexp.Next()) {
1021         TopoDS_Vertex V1 = TopoDS::Vertex( eexp.Current() );
1022         // look for an edge ending in E whose one vertex is in theShape1
1023         // and the other, in theShape2
1024         const TopTools_ListOfShape& Ancestors = theMesh1->GetAncestors(V1);
1025         TopTools_ListIteratorOfListOfShape ita(Ancestors);
1026         for(; ita.More(); ita.Next()) {
1027           if( ita.Value().ShapeType() != TopAbs_EDGE ) continue;
1028           TopoDS_Edge edge = TopoDS::Edge(ita.Value());
1029           bool FromShape1 = false;
1030           for(TopExp_Explorer expe(theShape1, TopAbs_EDGE); expe.More(); expe.Next() ) {
1031             if(edge.IsSame(expe.Current())) {
1032               FromShape1 = true;
1033               break;
1034             }
1035           }
1036           if(!FromShape1) {
1037             // is it an edge between theShape1 and theShape2?
1038             TopExp_Explorer expv(edge, TopAbs_VERTEX);
1039             TopoDS_Vertex V2 = TopoDS::Vertex( expv.Current() );
1040             if(V2.IsSame(V1)) {
1041               expv.Next();
1042               V2 = TopoDS::Vertex( expv.Current() );
1043             }
1044             bool FromShape2 = false;
1045             for ( expv.Init( theShape2, TopAbs_VERTEX ); expv.More(); expv.Next()) {
1046               if ( V2.IsSame( expv.Current() )) {
1047                 FromShape2 = true;
1048                 break;
1049               }
1050             }
1051             if ( FromShape2 ) {
1052               if ( VV1[0].IsNull() )
1053                 VV1[0] = V1, VV2[0] = V2;
1054               else
1055                 VV1[1] = V1, VV2[1] = V2;
1056               break; // from loop on ancestors of V1
1057             }
1058           }
1059         }
1060       }
1061       if ( !VV1[1].IsNull() ) {
1062         InsertAssociation( VV1[0], VV2[0], theMap, bidirect);
1063         InsertAssociation( VV1[1], VV2[1], theMap, bidirect);
1064         return FindSubShapeAssociation( theShape1, theMesh1, theShape2, theMesh2, theMap);
1065       }
1066     }
1067     break; // try by vertex closeness
1068   }
1069   default:;
1070   }
1071
1072   // Find association by closeness of vertices
1073   // ------------------------------------------
1074
1075   TopTools_IndexedMapOfShape vMap1, vMap2;
1076   TopExp::MapShapes( theShape1, TopAbs_VERTEX, vMap1 );
1077   TopExp::MapShapes( theShape2, TopAbs_VERTEX, vMap2 );
1078
1079   if ( vMap1.Extent() != vMap2.Extent() )
1080     RETURN_BAD_RESULT("Different nb of vertices");
1081
1082   if ( vMap1.Extent() == 1 ) {
1083     InsertAssociation( vMap1(1), vMap2(1), theMap, bidirect);
1084     if ( theShape1.ShapeType() == TopAbs_EDGE ) {
1085       InsertAssociation( theShape1, theShape2, theMap, bidirect );
1086       return true;
1087     }
1088     return FindSubShapeAssociation( theShape1, theMesh1, theShape2, theMesh2, theMap);
1089   }
1090
1091   // Find transformation to make the shapes be of similar size at same location
1092
1093   Bnd_Box box[2];
1094   for ( int i = 1; i <= vMap1.Extent(); ++i ) {
1095     box[ 0 ].Add( BRep_Tool::Pnt ( TopoDS::Vertex( vMap1( i ))));
1096     box[ 1 ].Add( BRep_Tool::Pnt ( TopoDS::Vertex( vMap2( i ))));
1097   }
1098
1099   gp_Pnt gc[2]; // box center
1100   double x0,y0,z0, x1,y1,z1;
1101   box[0].Get( x0,y0,z0, x1,y1,z1 );
1102   gc[0] = 0.5 * ( gp_XYZ( x0,y0,z0 ) + gp_XYZ( x1,y1,z1 ));
1103   box[1].Get( x0,y0,z0, x1,y1,z1 );
1104   gc[1] = 0.5 * ( gp_XYZ( x0,y0,z0 ) + gp_XYZ( x1,y1,z1 ));
1105
1106   // 1 -> 2
1107   gp_Vec vec01( gc[0], gc[1] );
1108   double scale = sqrt( box[1].SquareExtent() / box[0].SquareExtent() );
1109
1110   // Find 2 closest vertices
1111
1112   TopoDS_Vertex VV1[2], VV2[2];
1113   // get 2 linked vertices of shape 1 not belonging to an inner wire of a face
1114   TopoDS_Shape edge = theShape1;
1115   TopExp_Explorer expF( theShape1, TopAbs_FACE ), expE;
1116   if ( expF.More() ) {
1117     for ( ; expF.More(); expF.Next() ) {
1118       edge.Nullify();
1119       TopoDS_Shape wire = OuterShape( TopoDS::Face( expF.Current() ), TopAbs_WIRE );
1120       for ( expE.Init( wire, TopAbs_EDGE ); edge.IsNull() && expE.More(); expE.Next() )
1121         if ( !IsClosedEdge( TopoDS::Edge( expE.Current() )))
1122           edge = expE.Current();
1123       if ( !edge.IsNull() )
1124         break;
1125     }
1126   } else if (edge.ShapeType() != TopAbs_EDGE) { // no faces
1127     edge.Nullify();
1128     for ( expE.Init( theShape1, TopAbs_EDGE ); edge.IsNull() && expE.More(); expE.Next() )
1129       if ( !IsClosedEdge( TopoDS::Edge( expE.Current() )))
1130         edge = expE.Current();
1131   }
1132   if ( edge.IsNull() || edge.ShapeType() != TopAbs_EDGE )
1133     RETURN_BAD_RESULT("Edge not found");
1134
1135   TopExp::Vertices( TopoDS::Edge( edge ), VV1[0], VV1[1]);
1136   if ( VV1[0].IsSame( VV1[1] ))
1137     RETURN_BAD_RESULT("Only closed edges");
1138
1139   // find vertices closest to 2 linked vertices of shape 1
1140   for ( int i1 = 0; i1 < 2; ++i1 )
1141   {
1142     double dist2 = DBL_MAX;
1143     gp_Pnt p1 = BRep_Tool::Pnt( VV1[ i1 ]);
1144     p1.Translate( vec01 );
1145     p1.Scale( gc[1], scale );
1146     for ( int i2 = 1; i2 <= vMap2.Extent(); ++i2 )
1147     {
1148       TopoDS_Vertex V2 = TopoDS::Vertex( vMap2( i2 ));
1149       gp_Pnt p2 = BRep_Tool::Pnt ( V2 );
1150       double d2 = p1.SquareDistance( p2 );
1151       if ( d2 < dist2 && !V2.IsSame( VV2[ 0 ])) {
1152         VV2[ i1 ] = V2; dist2 = d2;
1153       }
1154     }
1155   }
1156
1157   InsertAssociation( VV1[ 0 ], VV2[ 0 ], theMap, bidirect);
1158   InsertAssociation( VV1[ 1 ], VV2[ 1 ], theMap, bidirect);
1159   MESSAGE("Initial assoc VERT " << theMesh1->GetMeshDS()->ShapeToIndex( VV1[ 0 ] )<<
1160           " to "                << theMesh2->GetMeshDS()->ShapeToIndex( VV2[ 0 ] )<<
1161           "\nand         VERT " << theMesh1->GetMeshDS()->ShapeToIndex( VV1[ 1 ] )<<
1162           " to "                << theMesh2->GetMeshDS()->ShapeToIndex( VV2[ 1 ] ));
1163   if ( theShape1.ShapeType() == TopAbs_EDGE ) {
1164     InsertAssociation( theShape1, theShape2, theMap, bidirect );
1165     return true;
1166   }
1167
1168   return FindSubShapeAssociation( theShape1, theMesh1, theShape2, theMesh2, theMap );
1169 }
1170
1171 //================================================================================
1172 /*!
1173  * \brief Find association of edges of faces
1174  * \param face1 - face 1
1175  * \param VV1 - vertices of face 1
1176  * \param face2 - face 2
1177  * \param VV2 - vertices of face 2 associated with ones of face 1
1178  * \param edges1 - out list of edges of face 1
1179  * \param edges2 - out list of edges of face 2
1180  * \retval int - nb of edges in an outer wire in a success case, else zero
1181  */
1182 //================================================================================
1183
1184 int StdMeshers_ProjectionUtils::FindFaceAssociation(const TopoDS_Face&    face1,
1185                                                     TopoDS_Vertex         VV1[2],
1186                                                     const TopoDS_Face&    face2,
1187                                                     TopoDS_Vertex         VV2[2],
1188                                                     list< TopoDS_Edge > & edges1,
1189                                                     list< TopoDS_Edge > & edges2)
1190 {
1191   list< int > nbEInW1, nbEInW2;
1192   for ( int outer_wire_algo = 0; outer_wire_algo < 2; ++outer_wire_algo )
1193   {
1194     edges1.clear();
1195     edges2.clear();
1196
1197     if ( SMESH_Block::GetOrderedEdges( face1, VV1[0], edges1, nbEInW1, outer_wire_algo) !=
1198          SMESH_Block::GetOrderedEdges( face2, VV2[0], edges2, nbEInW2, outer_wire_algo) )
1199       CONT_BAD_RESULT("Different number of wires in faces ");
1200
1201     if ( nbEInW1.front() != nbEInW2.front() )
1202       CONT_BAD_RESULT("Different number of edges in faces: " <<
1203                       nbEInW1.front() << " != " << nbEInW2.front());
1204
1205     // Define if we need to reverse one of wires to make edges in lists match each other
1206
1207     bool reverse = false;
1208
1209     list< TopoDS_Edge >::iterator edgeIt;
1210     if ( !VV1[1].IsSame( TopExp::LastVertex( edges1.front(), true ))) {
1211       reverse = true;
1212       edgeIt = --edges1.end();
1213       // check if the second vertex belongs to the first or last edge in the wire
1214       if ( !VV1[1].IsSame( TopExp::FirstVertex( *edgeIt, true ))) {
1215         bool KO = true; // belongs to none
1216         if ( nbEInW1.size() > 1 ) { // several wires
1217           edgeIt = edges1.begin();
1218           std::advance( edgeIt, nbEInW1.front()-1 );
1219           KO = !VV1[1].IsSame( TopExp::FirstVertex( *edgeIt, true ));
1220         }
1221         if ( KO )
1222           CONT_BAD_RESULT("GetOrderedEdges() failed");
1223       }
1224     }
1225     edgeIt = --edges2.end();
1226     if ( !VV2[1].IsSame( TopExp::LastVertex( edges2.front(), true ))) {
1227       reverse = !reverse;
1228       // check if the second vertex belongs to the first or last edge in the wire
1229       if ( !VV2[1].IsSame( TopExp::FirstVertex( *edgeIt, true ))) {
1230         bool KO = true; // belongs to none
1231         if ( nbEInW2.size() > 1 ) { // several wires
1232           edgeIt = edges2.begin();
1233           std::advance( edgeIt, nbEInW2.front()-1 );
1234           KO = !VV2[1].IsSame( TopExp::FirstVertex( *edgeIt, true ));
1235         }
1236         if ( KO )
1237           CONT_BAD_RESULT("GetOrderedEdges() failed");
1238       }
1239     }
1240     if ( reverse )
1241     {
1242       Reverse( edges2 , nbEInW2.front());
1243       if (( VV1[1].IsSame( TopExp::LastVertex( edges1.front(), true ))) !=
1244           ( VV2[1].IsSame( TopExp::LastVertex( edges2.front(), true ))))
1245         CONT_BAD_RESULT("GetOrderedEdges() failed");
1246     }
1247
1248     // Try to orient internal wires (0020996)
1249     if ( nbEInW1.size() > 1 )
1250     {
1251       // Try by UV similarity. Check that Vec(VV1[0],VV1[1]) in 2D on face1 is the same
1252       // as Vec(VV2[0],VV2[1]) on face2
1253       double vTol = BRep_Tool::Tolerance( VV1[0] );
1254       BRepAdaptor_Surface surface1( face1, false );
1255       double vTolUV =
1256         surface1.UResolution( vTol ) + surface1.VResolution( vTol ); // let's be tolerant
1257       gp_Pnt2d v1UV = BRep_Tool::Parameters( VV1[0], face1 );
1258       gp_Pnt2d v2UV = BRep_Tool::Parameters( VV2[0], face2 );
1259       gp_Pnt2d v1UV1 = BRep_Tool::Parameters( VV1[1], face1 );
1260       gp_Pnt2d v2UV1 = BRep_Tool::Parameters( VV2[1], face2 );
1261       gp_Vec2d e1Vec( v1UV, v1UV1 );
1262       gp_Vec2d e2Vec( v2UV, v2UV1 );
1263       if ( Abs( e1Vec.X()-e2Vec.X()) < vTolUV && Abs( e1Vec.Y()-e2Vec.Y()) < vTolUV )
1264       {
1265         gp_XY dUV = v2UV.XY() - v1UV.XY();
1266         // skip edges of the outer wire
1267         list< int >::iterator nbEInW = nbEInW1.begin();
1268         list< TopoDS_Edge >::iterator edge1Beg = edges1.begin(), edge2Beg = edges2.begin();
1269         for ( int i = 0; i < *nbEInW; ++i )
1270           ++edge1Beg, ++edge2Beg;
1271         for ( ++nbEInW; nbEInW != nbEInW1.end(); ++nbEInW ) // loop on inner wires
1272         {
1273           // reach an end of edges of a current innner wire
1274           list< TopoDS_Edge >::iterator edge1End = edge1Beg, edge2End = edge2Beg;
1275           for ( int i = 0; i < *nbEInW; ++i )
1276             ++edge1End, ++edge2End;
1277           // rotate edges2 untill coincident with edges1 in 2D
1278           v1UV = BRep_Tool::Parameters( TopExp::FirstVertex(*edge1Beg,true), face1 );
1279           v2UV = BRep_Tool::Parameters( TopExp::LastVertex (*edge1Beg,true), face1 );
1280           v1UV.ChangeCoord() += dUV;
1281           v2UV.ChangeCoord() += dUV;
1282           int i = *nbEInW;
1283           while ( --i > 0 && !sameVertexUV( *edge2Beg, face2, 0, v1UV, vTolUV ))
1284             edges2.splice( edge2End, edges2, edge2Beg++ ); // move edge2Beg to place before edge2End
1285           // reverse edges2 if needed
1286           if ( sameVertexUV( *edge2Beg, face2, 0, v1UV, vTolUV ) &&
1287                !sameVertexUV( *edge2Beg, face2, 1, v2UV, vTolUV ))
1288           {
1289             Reverse( edges2 , *nbEInW, distance( edges2.begin(),edge2Beg ));
1290             edge2End = edges2.begin();
1291             std::advance( edge2End, std::accumulate( nbEInW1.begin(), nbEInW, *nbEInW));
1292           }
1293           // prepare to the next wire loop
1294           edge1Beg = edge1End, edge2Beg = edge2End;
1295         }
1296       }
1297     }
1298     break; // OK
1299
1300   } // loop algos getting an outer wire
1301
1302   return nbEInW2.front();
1303 }
1304
1305 //=======================================================================
1306 //function : InitVertexAssociation
1307 //purpose  : 
1308 //=======================================================================
1309
1310 void StdMeshers_ProjectionUtils::InitVertexAssociation( const SMESH_Hypothesis* theHyp,
1311                                                         TShapeShapeMap &        theAssociationMap,
1312                                                         const TopoDS_Shape&     theTargetShape)
1313 {
1314   string hypName = theHyp->GetName();
1315   if ( hypName == "ProjectionSource1D" ) {
1316     const StdMeshers_ProjectionSource1D * hyp =
1317       static_cast<const StdMeshers_ProjectionSource1D*>( theHyp );
1318     if ( hyp->HasVertexAssociation() )
1319       InsertAssociation( hyp->GetSourceVertex(),hyp->GetTargetVertex(),theAssociationMap);
1320   }
1321   else if ( hypName == "ProjectionSource2D" ) {
1322     const StdMeshers_ProjectionSource2D * hyp =
1323       static_cast<const StdMeshers_ProjectionSource2D*>( theHyp );
1324     if ( hyp->HasVertexAssociation() ) {
1325       InsertAssociation( hyp->GetSourceVertex(1),hyp->GetTargetVertex(1),theAssociationMap);
1326       InsertAssociation( hyp->GetSourceVertex(2),hyp->GetTargetVertex(2),theAssociationMap);
1327     }
1328   }
1329   else if ( hypName == "ProjectionSource3D" ) {
1330     const StdMeshers_ProjectionSource3D * hyp =
1331       static_cast<const StdMeshers_ProjectionSource3D*>( theHyp );
1332     if ( hyp->HasVertexAssociation() ) {
1333       InsertAssociation( hyp->GetSourceVertex(1),hyp->GetTargetVertex(1),theAssociationMap);
1334       InsertAssociation( hyp->GetSourceVertex(2),hyp->GetTargetVertex(2),theAssociationMap);
1335     }
1336   }
1337 }
1338
1339 //=======================================================================
1340 /*!
1341  * \brief Inserts association theShape1 <-> theShape2 to TShapeShapeMap
1342  * \param theShape1 - shape 1
1343  * \param theShape2 - shape 2
1344  * \param theAssociationMap - association map 
1345  * \retval bool - true if there was no association for these shapes before
1346  */
1347 //=======================================================================
1348
1349 bool StdMeshers_ProjectionUtils::InsertAssociation( const TopoDS_Shape& theShape1,
1350                                                     const TopoDS_Shape& theShape2,
1351                                                     TShapeShapeMap &    theAssociationMap,
1352                                                     const bool          theBidirectional)
1353 {
1354   if ( !theShape1.IsNull() && !theShape2.IsNull() ) {
1355     SHOW_SHAPE(theShape1,"Assoc ");
1356     SHOW_SHAPE(theShape2," to ");
1357     bool isNew = ( theAssociationMap.Bind( theShape1, theShape2 ));
1358     if ( theBidirectional )
1359       theAssociationMap.Bind( theShape2, theShape1 );
1360     return isNew;
1361   }
1362   else {
1363     throw SALOME_Exception("StdMeshers_ProjectionUtils: attempt to associate NULL shape");
1364   }
1365   return false;
1366 }
1367
1368 //=======================================================================
1369 /*!
1370  * \brief Finds an edge by its vertices in a main shape of the mesh
1371  * \param aMesh - the mesh
1372  * \param V1 - vertex 1
1373  * \param V2 - vertex 2
1374  * \retval TopoDS_Edge - found edge
1375  */
1376 //=======================================================================
1377
1378 TopoDS_Edge StdMeshers_ProjectionUtils::GetEdgeByVertices( SMESH_Mesh*          theMesh,
1379                                                            const TopoDS_Vertex& theV1,
1380                                                            const TopoDS_Vertex& theV2)
1381 {
1382   if ( theMesh && !theV1.IsNull() && !theV2.IsNull() )
1383   {
1384     TopTools_ListIteratorOfListOfShape ancestorIt( theMesh->GetAncestors( theV1 ));
1385     for ( ; ancestorIt.More(); ancestorIt.Next() )
1386       if ( ancestorIt.Value().ShapeType() == TopAbs_EDGE )
1387         for ( TopExp_Explorer expV ( ancestorIt.Value(), TopAbs_VERTEX );
1388               expV.More();
1389               expV.Next() )
1390           if ( theV2.IsSame( expV.Current() ))
1391             return TopoDS::Edge( ancestorIt.Value() );
1392   }
1393   return TopoDS_Edge();
1394 }
1395
1396 //================================================================================
1397 /*!
1398  * \brief Return another face sharing an edge
1399  * \param edgeToFaces - data map of descendants to ancestors
1400  * \param edge - edge
1401  * \param face - face
1402  * \retval TopoDS_Face - found face
1403  */
1404 //================================================================================
1405
1406 TopoDS_Face StdMeshers_ProjectionUtils::GetNextFace( const TAncestorMap& edgeToFaces,
1407                                                      const TopoDS_Edge&  edge,
1408                                                      const TopoDS_Face&  face)
1409 {
1410 //   if ( !edge.IsNull() && !face.IsNull() && edgeToFaces.Contains( edge ))
1411   if ( !edge.IsNull() && edgeToFaces.Contains( edge )) // PAL16202
1412   {
1413     TopTools_ListIteratorOfListOfShape ancestorIt( edgeToFaces.FindFromKey( edge ));
1414     for ( ; ancestorIt.More(); ancestorIt.Next() )
1415       if ( ancestorIt.Value().ShapeType() == TopAbs_FACE &&
1416            !face.IsSame( ancestorIt.Value() ))
1417         return TopoDS::Face( ancestorIt.Value() );
1418   }
1419   return TopoDS_Face();
1420 }
1421
1422 //================================================================================
1423 /*!
1424  * \brief Return other vertex of an edge
1425  */
1426 //================================================================================
1427
1428 TopoDS_Vertex StdMeshers_ProjectionUtils::GetNextVertex(const TopoDS_Edge&   edge,
1429                                                         const TopoDS_Vertex& vertex)
1430 {
1431   TopoDS_Vertex vF,vL;
1432   TopExp::Vertices(edge,vF,vL);
1433   if ( vF.IsSame( vL ))
1434     return TopoDS_Vertex();
1435   return vertex.IsSame( vF ) ? vL : vF; 
1436 }
1437
1438 //================================================================================
1439 /*!
1440  * \brief Return a propagation edge
1441  * \param aMesh - mesh
1442  * \param theEdge - edge to find by propagation
1443  * \param fromEdge - start edge for propagation
1444  * \retval pair<int,TopoDS_Edge> - propagation step and found edge
1445  */
1446 //================================================================================
1447
1448 pair<int,TopoDS_Edge>
1449 StdMeshers_ProjectionUtils::GetPropagationEdge( SMESH_Mesh*        aMesh,
1450                                                 const TopoDS_Edge& theEdge,
1451                                                 const TopoDS_Edge& fromEdge)
1452 {
1453   SMESH_IndexedMapOfShape aChain;
1454   int step = 0;
1455
1456   // List of edges, added to chain on the previous cycle pass
1457   TopTools_ListOfShape listPrevEdges;
1458   listPrevEdges.Append(fromEdge);
1459
1460   // Collect all edges pass by pass
1461   while (listPrevEdges.Extent() > 0) {
1462     step++;
1463     // List of edges, added to chain on this cycle pass
1464     TopTools_ListOfShape listCurEdges;
1465
1466     // Find the next portion of edges
1467     TopTools_ListIteratorOfListOfShape itE (listPrevEdges);
1468     for (; itE.More(); itE.Next()) {
1469       TopoDS_Shape anE = itE.Value();
1470
1471       // Iterate on faces, having edge <anE>
1472       TopTools_ListIteratorOfListOfShape itA (aMesh->GetAncestors(anE));
1473       for (; itA.More(); itA.Next()) {
1474         TopoDS_Shape aW = itA.Value();
1475
1476         // There are objects of different type among the ancestors of edge
1477         if (aW.ShapeType() == TopAbs_WIRE) {
1478           TopoDS_Shape anOppE;
1479
1480           BRepTools_WireExplorer aWE (TopoDS::Wire(aW));
1481           Standard_Integer nb = 1, found = 0;
1482           TopTools_Array1OfShape anEdges (1,4);
1483           for (; aWE.More(); aWE.Next(), nb++) {
1484             if (nb > 4) {
1485               found = 0;
1486               break;
1487             }
1488             anEdges(nb) = aWE.Current();
1489             if (anEdges(nb).IsSame(anE)) found = nb;
1490           }
1491
1492           if (nb == 5 && found > 0) {
1493             // Quadrangle face found, get an opposite edge
1494             Standard_Integer opp = found + 2;
1495             if (opp > 4) opp -= 4;
1496             anOppE = anEdges(opp);
1497
1498             // add anOppE to aChain if ...
1499             if (!aChain.Contains(anOppE)) { // ... anOppE is not in aChain
1500               // Add found edge to the chain oriented so that to
1501               // have it co-directed with a forward MainEdge
1502               TopAbs_Orientation ori = anE.Orientation();
1503               if ( anEdges(opp).Orientation() == anEdges(found).Orientation() )
1504                 ori = TopAbs::Reverse( ori );
1505               anOppE.Orientation( ori );
1506               if ( anOppE.IsSame( theEdge ))
1507                 return make_pair( step, TopoDS::Edge( anOppE ));
1508               aChain.Add(anOppE);
1509               listCurEdges.Append(anOppE);
1510             }
1511           } // if (nb == 5 && found > 0)
1512         } // if (aF.ShapeType() == TopAbs_WIRE)
1513       } // for (; itF.More(); itF.Next())
1514     } // for (; itE.More(); itE.Next())
1515
1516     listPrevEdges = listCurEdges;
1517   } // while (listPrevEdges.Extent() > 0)
1518
1519   return make_pair( INT_MAX, TopoDS_Edge());
1520 }
1521
1522 //================================================================================
1523   /*!
1524    * \brief Find corresponding nodes on two faces
1525     * \param face1 - the first face
1526     * \param mesh1 - mesh containing elements on the first face
1527     * \param face2 - the second face
1528     * \param mesh2 - mesh containing elements on the second face
1529     * \param assocMap - map associating subshapes of the faces
1530     * \param node1To2Map - map containing found matching nodes
1531     * \retval bool - is a success
1532    */
1533 //================================================================================
1534
1535 bool StdMeshers_ProjectionUtils::
1536 FindMatchingNodesOnFaces( const TopoDS_Face&     face1,
1537                           SMESH_Mesh*            mesh1,
1538                           const TopoDS_Face&     face2,
1539                           SMESH_Mesh*            mesh2,
1540                           const TShapeShapeMap & assocMap,
1541                           TNodeNodeMap &         node1To2Map)
1542 {
1543   SMESHDS_Mesh* meshDS1 = mesh1->GetMeshDS();
1544   SMESHDS_Mesh* meshDS2 = mesh2->GetMeshDS();
1545   
1546   SMESH_MesherHelper helper1( *mesh1 );
1547   SMESH_MesherHelper helper2( *mesh2 );
1548
1549   // Get corresponding submeshes and roughly check match of meshes
1550
1551   SMESHDS_SubMesh * SM2 = meshDS2->MeshElements( face2 );
1552   SMESHDS_SubMesh * SM1 = meshDS1->MeshElements( face1 );
1553   if ( !SM2 || !SM1 )
1554     RETURN_BAD_RESULT("Empty submeshes");
1555   if ( SM2->NbNodes()    != SM1->NbNodes() ||
1556        SM2->NbElements() != SM1->NbElements() )
1557     RETURN_BAD_RESULT("Different meshes on corresponding faces "
1558                       << meshDS1->ShapeToIndex( face1 ) << " and "
1559                       << meshDS2->ShapeToIndex( face2 ));
1560   if ( SM2->NbElements() == 0 )
1561     RETURN_BAD_RESULT("Empty submeshes");
1562
1563   helper1.SetSubShape( face1 );
1564   helper2.SetSubShape( face2 );
1565   if ( helper1.HasSeam() != helper2.HasSeam() )
1566     RETURN_BAD_RESULT("Different faces' geometry");
1567
1568   // Data to call SMESH_MeshEditor::FindMatchingNodes():
1569
1570   // 1. Nodes of corresponding links:
1571
1572   // get 2 matching edges, try to find not seam ones
1573   TopoDS_Edge edge1, edge2, seam1, seam2, anyEdge1, anyEdge2;
1574   TopExp_Explorer eE( OuterShape( face2, TopAbs_WIRE ), TopAbs_EDGE );
1575   do {
1576     // edge 2
1577     TopoDS_Edge e2 = TopoDS::Edge( eE.Current() );
1578     eE.Next();
1579     // edge 1
1580     if ( !assocMap.IsBound( e2 ))
1581       RETURN_BAD_RESULT("Association not found for edge " << meshDS2->ShapeToIndex( e2 ));
1582     TopoDS_Edge e1 = TopoDS::Edge( assocMap( e2 ));
1583     if ( !helper1.IsSubShape( e1, face1 ))
1584       RETURN_BAD_RESULT("Wrong association, edge " << meshDS1->ShapeToIndex( e1 ) <<
1585                         " isn't a subshape of face " << meshDS1->ShapeToIndex( face1 ));
1586     // check that there are nodes on edges
1587     SMESHDS_SubMesh * eSM1 = meshDS1->MeshElements( e1 );
1588     SMESHDS_SubMesh * eSM2 = meshDS2->MeshElements( e2 );
1589     bool nodesOnEdges = ( eSM1 && eSM2 && eSM1->NbNodes() && eSM2->NbNodes() );
1590     // check that the nodes on edges belong to faces
1591     // (as NETGEN ignores nodes on the degenerated geom edge)
1592     bool nodesOfFaces = false;
1593     if ( nodesOnEdges ) {
1594       const SMDS_MeshNode* n1 = eSM1->GetNodes()->next();
1595       const SMDS_MeshNode* n2 = eSM2->GetNodes()->next();
1596       nodesOfFaces = ( n1->GetInverseElementIterator(SMDSAbs_Face)->more() &&
1597                        n2->GetInverseElementIterator(SMDSAbs_Face)->more() );
1598     }
1599     if ( nodesOfFaces )
1600     {
1601       if ( helper2.IsRealSeam( e2 )) {
1602         seam1 = e1; seam2 = e2;
1603       }
1604       else {
1605         edge1 = e1; edge2 = e2;
1606       }
1607     }
1608     else {
1609       anyEdge1 = e1; anyEdge2 = e2;
1610     }
1611   } while ( edge2.IsNull() && eE.More() );
1612   //
1613   if ( edge2.IsNull() ) {
1614     edge1 = seam1; edge2 = seam2;
1615   }
1616   bool hasNodesOnEdge = (! edge2.IsNull() );
1617   if ( !hasNodesOnEdge ) {
1618     // 0020338 - nb segments == 1
1619     edge1 = anyEdge1; edge2 = anyEdge2;
1620   }
1621
1622   // get 2 matching vertices
1623   TopoDS_Vertex V2 = TopExp::FirstVertex( TopoDS::Edge( edge2 ));
1624   if ( !assocMap.IsBound( V2 ))
1625     RETURN_BAD_RESULT("Association not found for vertex " << meshDS2->ShapeToIndex( V2 ));
1626   TopoDS_Vertex V1 = TopoDS::Vertex( assocMap( V2 ));
1627
1628   // nodes on vertices
1629   const SMDS_MeshNode* vNode1 = SMESH_Algo::VertexNode( V1, meshDS1 );
1630   const SMDS_MeshNode* vNode2 = SMESH_Algo::VertexNode( V2, meshDS2 );
1631   if ( !vNode1 ) RETURN_BAD_RESULT("No node on vertex #" << meshDS1->ShapeToIndex( V1 ));
1632   if ( !vNode2 ) RETURN_BAD_RESULT("No node on vertex #" << meshDS2->ShapeToIndex( V2 ));
1633
1634   // nodes on edges linked with nodes on vertices
1635   const SMDS_MeshNode* nullNode = 0;
1636   vector< const SMDS_MeshNode*> eNode1( 2, nullNode );
1637   vector< const SMDS_MeshNode*> eNode2( 2, nullNode );
1638   if ( hasNodesOnEdge )
1639   {
1640     int nbNodeToGet = 1;
1641     if ( IsClosedEdge( edge1 ) || IsClosedEdge( edge2 ) )
1642       nbNodeToGet = 2;
1643     for ( int is2 = 0; is2 < 2; ++is2 )
1644     {
1645       TopoDS_Edge &     edge  = is2 ? edge2 : edge1;
1646       SMESHDS_Mesh *    smDS  = is2 ? meshDS2 : meshDS1;
1647       SMESHDS_SubMesh* edgeSM = smDS->MeshElements( edge );
1648       // nodes linked with ones on vertices
1649       const SMDS_MeshNode*           vNode = is2 ? vNode2 : vNode1;
1650       vector< const SMDS_MeshNode*>& eNode = is2 ? eNode2 : eNode1;
1651       int nbGotNode = 0;
1652       SMDS_ElemIteratorPtr vElem = vNode->GetInverseElementIterator(SMDSAbs_Edge);
1653       while ( vElem->more() && nbGotNode != nbNodeToGet ) {
1654         const SMDS_MeshElement* elem = vElem->next();
1655         if ( edgeSM->Contains( elem ))
1656           eNode[ nbGotNode++ ] = 
1657             ( elem->GetNode(0) == vNode ) ? elem->GetNode(1) : elem->GetNode(0);
1658       }
1659       if ( nbGotNode > 1 ) // sort found nodes by param on edge
1660       {
1661         SMESH_MesherHelper* helper = is2 ? &helper2 : &helper1;
1662         double u0 = helper->GetNodeU( edge, eNode[ 0 ]);
1663         double u1 = helper->GetNodeU( edge, eNode[ 1 ]);
1664         if ( u0 > u1 ) std::swap( eNode[ 0 ], eNode[ 1 ]);
1665       }
1666       if ( nbGotNode == 0 )
1667         RETURN_BAD_RESULT("Found no nodes on edge " << smDS->ShapeToIndex( edge ) <<
1668                           " linked to " << vNode );
1669     }
1670   }
1671   else // 0020338 - nb segments == 1
1672   {
1673     // get 2 other matching vertices
1674     V2 = TopExp::LastVertex( TopoDS::Edge( edge2 ));
1675     if ( !assocMap.IsBound( V2 ))
1676       RETURN_BAD_RESULT("Association not found for vertex " << meshDS2->ShapeToIndex( V2 ));
1677     V1 = TopoDS::Vertex( assocMap( V2 ));
1678
1679     // nodes on vertices
1680     eNode1[0] = SMESH_Algo::VertexNode( V1, meshDS1 );
1681     eNode2[0] = SMESH_Algo::VertexNode( V2, meshDS2 );
1682     if ( !eNode1[0] ) RETURN_BAD_RESULT("No node on vertex #" << meshDS1->ShapeToIndex( V1 ));
1683     if ( !eNode2[0] ) RETURN_BAD_RESULT("No node on vertex #" << meshDS2->ShapeToIndex( V2 ));
1684   }
1685
1686   // 2. face sets
1687
1688   set<const SMDS_MeshElement*> Elems1, Elems2;
1689   for ( int is2 = 0; is2 < 2; ++is2 )
1690   {
1691     set<const SMDS_MeshElement*> & elems = is2 ? Elems2 : Elems1;
1692     SMESHDS_SubMesh*                  sm = is2 ? SM2 : SM1;
1693     SMESH_MesherHelper*           helper = is2 ? &helper2 : &helper1;
1694     const TopoDS_Face &             face = is2 ? face2 : face1;
1695     SMDS_ElemIteratorPtr eIt = sm->GetElements();
1696
1697     if ( !helper->IsRealSeam( is2 ? edge2 : edge1 ))
1698     {
1699       while ( eIt->more() ) elems.insert( eIt->next() );
1700     }
1701     else
1702     {
1703       // the only suitable edge is seam, i.e. it is a sphere.
1704       // FindMatchingNodes() will not know which way to go from any edge.
1705       // So we ignore all faces having nodes on edges or vertices except
1706       // one of faces sharing current start nodes
1707
1708       // find a face to keep
1709       const SMDS_MeshElement* faceToKeep = 0;
1710       const SMDS_MeshNode* vNode = is2 ? vNode2 : vNode1;
1711       const SMDS_MeshNode* eNode = is2 ? eNode2[0] : eNode1[0];
1712       TIDSortedElemSet inSet, notInSet;
1713
1714       const SMDS_MeshElement* f1 =
1715         SMESH_MeshEditor::FindFaceInSet( vNode, eNode, inSet, notInSet );
1716       if ( !f1 ) RETURN_BAD_RESULT("The first face on seam not found");
1717       notInSet.insert( f1 );
1718
1719       const SMDS_MeshElement* f2 =
1720         SMESH_MeshEditor::FindFaceInSet( vNode, eNode, inSet, notInSet );
1721       if ( !f2 ) RETURN_BAD_RESULT("The second face on seam not found");
1722
1723       // select a face with less UV of vNode
1724       const SMDS_MeshNode* notSeamNode[2] = {0, 0};
1725       for ( int iF = 0; iF < 2; ++iF ) {
1726         const SMDS_MeshElement* f = ( iF ? f2 : f1 );
1727         for ( int i = 0; !notSeamNode[ iF ] && i < f->NbNodes(); ++i ) {
1728           const SMDS_MeshNode* node = f->GetNode( i );
1729           if ( !helper->IsSeamShape( node->GetPosition()->GetShapeId() ))
1730             notSeamNode[ iF ] = node;
1731         }
1732       }
1733       gp_Pnt2d uv1 = helper->GetNodeUV( face, vNode, notSeamNode[0] );
1734       gp_Pnt2d uv2 = helper->GetNodeUV( face, vNode, notSeamNode[1] );
1735       if ( uv1.X() + uv1.Y() > uv2.X() + uv2.Y() )
1736         faceToKeep = f2;
1737       else
1738         faceToKeep = f1;
1739
1740       // fill elem set
1741       elems.insert( faceToKeep );
1742       while ( eIt->more() ) {
1743         const SMDS_MeshElement* f = eIt->next();
1744         int nbNodes = f->NbNodes();
1745         if ( f->IsQuadratic() )
1746           nbNodes /= 2;
1747         bool onBnd = false;
1748         for ( int i = 0; !onBnd && i < nbNodes; ++i ) {
1749           const SMDS_MeshNode* node = f->GetNode( i );
1750           onBnd = ( node->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE);
1751         }
1752         if ( !onBnd )
1753           elems.insert( f );
1754       }
1755       // add also faces adjacent to faceToKeep
1756       int nbNodes = faceToKeep->NbNodes();
1757       if ( faceToKeep->IsQuadratic() ) nbNodes /= 2;
1758       notInSet.insert( f1 );
1759       notInSet.insert( f2 );
1760       for ( int i = 0; i < nbNodes; ++i ) {
1761         const SMDS_MeshNode* n1 = faceToKeep->GetNode( i );
1762         const SMDS_MeshNode* n2 = faceToKeep->GetNode( i+1 % nbNodes );
1763         f1 = SMESH_MeshEditor::FindFaceInSet( n1, n2, inSet, notInSet );
1764         if ( f1 )
1765           elems.insert( f1 );
1766       }
1767     } // case on a sphere
1768   } // loop on 2 faces
1769
1770   //  int quadFactor = (*Elems1.begin())->IsQuadratic() ? 2 : 1;
1771
1772   node1To2Map.clear();
1773   int res = SMESH_MeshEditor::FindMatchingNodes( Elems1, Elems2,
1774                                                  vNode1, vNode2,
1775                                                  eNode1[0], eNode2[0],
1776                                                  node1To2Map);
1777   if ( res != SMESH_MeshEditor::SEW_OK )
1778     RETURN_BAD_RESULT("FindMatchingNodes() result " << res );
1779
1780   // On a sphere, add matching nodes on the edge
1781
1782   if ( helper1.IsRealSeam( edge1 ))
1783   {
1784     // sort nodes on edges by param on edge
1785     map< double, const SMDS_MeshNode* > u2nodesMaps[2];
1786     for ( int is2 = 0; is2 < 2; ++is2 )
1787     {
1788       TopoDS_Edge &     edge  = is2 ? edge2 : edge1;
1789       SMESHDS_Mesh *    smDS  = is2 ? meshDS2 : meshDS1;
1790       SMESHDS_SubMesh* edgeSM = smDS->MeshElements( edge );
1791       map< double, const SMDS_MeshNode* > & pos2nodes = u2nodesMaps[ is2 ];
1792
1793       SMDS_NodeIteratorPtr nIt = edgeSM->GetNodes();
1794       while ( nIt->more() ) {
1795         const SMDS_MeshNode* node = nIt->next();
1796         const SMDS_EdgePosition* pos =
1797           static_cast<const SMDS_EdgePosition*>(node->GetPosition().get());
1798         pos2nodes.insert( make_pair( pos->GetUParameter(), node ));
1799       }
1800       if ( pos2nodes.size() != edgeSM->NbNodes() )
1801         RETURN_BAD_RESULT("Equal params of nodes on edge "
1802                           << smDS->ShapeToIndex( edge ) << " of face " << is2 );
1803     }
1804     if ( u2nodesMaps[0].size() != u2nodesMaps[1].size() )
1805       RETURN_BAD_RESULT("Different nb of new nodes on edges or wrong params");
1806
1807     // compare edge orientation
1808     double u1 = helper1.GetNodeU( edge1, vNode1 );
1809     double u2 = helper2.GetNodeU( edge2, vNode2 );
1810     bool isFirst1 = ( u1 < u2nodesMaps[0].begin()->first );
1811     bool isFirst2 = ( u2 < u2nodesMaps[1].begin()->first );
1812     bool reverse ( isFirst1 != isFirst2 );
1813
1814     // associate matching nodes
1815     map< double, const SMDS_MeshNode* >::iterator u_Node1, u_Node2, end1;
1816     map< double, const SMDS_MeshNode* >::reverse_iterator uR_Node2;
1817     u_Node1 = u2nodesMaps[0].begin();
1818     u_Node2 = u2nodesMaps[1].begin();
1819     uR_Node2 = u2nodesMaps[1].rbegin();
1820     end1 = u2nodesMaps[0].end();
1821     for ( ; u_Node1 != end1; ++u_Node1 ) {
1822       const SMDS_MeshNode* n1 = u_Node1->second;
1823       const SMDS_MeshNode* n2 = ( reverse ? (uR_Node2++)->second : (u_Node2++)->second );
1824       node1To2Map.insert( make_pair( n1, n2 ));
1825     }
1826
1827     // associate matching nodes on the last vertices
1828     V2 = TopExp::LastVertex( TopoDS::Edge( edge2 ));
1829     if ( !assocMap.IsBound( V2 ))
1830       RETURN_BAD_RESULT("Association not found for vertex " << meshDS2->ShapeToIndex( V2 ));
1831     V1 = TopoDS::Vertex( assocMap( V2 ));
1832     vNode1 = SMESH_Algo::VertexNode( V1, meshDS1 );
1833     vNode2 = SMESH_Algo::VertexNode( V2, meshDS2 );
1834     if ( !vNode1 ) RETURN_BAD_RESULT("No node on vertex #" << meshDS1->ShapeToIndex( V1 ));
1835     if ( !vNode2 ) RETURN_BAD_RESULT("No node on vertex #" << meshDS2->ShapeToIndex( V2 ));
1836     node1To2Map.insert( make_pair( vNode1, vNode2 ));
1837   }
1838
1839 // don't know why this condition is usually true :(
1840 //   if ( node1To2Map.size() * quadFactor < SM1->NbNodes() )
1841 //     MESSAGE("FindMatchingNodes() found too few node pairs starting from nodes ("
1842 //             << vNode1->GetID() << " - " << eNode1[0]->GetID() << ") ("
1843 //             << vNode2->GetID() << " - " << eNode2[0]->GetID() << "):"
1844 //             << node1To2Map.size() * quadFactor << " < " << SM1->NbNodes());
1845   
1846   return true;
1847 }
1848
1849 //================================================================================
1850 /*!
1851  * \brief Check if the first and last vertices of an edge are the same
1852  * \param anEdge - the edge to check
1853  * \retval bool - true if same
1854  */
1855 //================================================================================
1856
1857 bool StdMeshers_ProjectionUtils::IsClosedEdge( const TopoDS_Edge& anEdge )
1858 {
1859   return TopExp::FirstVertex( anEdge ).IsSame( TopExp::LastVertex( anEdge ));
1860 }
1861
1862 //================================================================================
1863   /*!
1864    * \brief Return any subshape of a face belonging to the outer wire
1865     * \param face - the face
1866     * \param type - type of subshape to return
1867     * \retval TopoDS_Shape - the found subshape
1868    */
1869 //================================================================================
1870
1871 TopoDS_Shape StdMeshers_ProjectionUtils::OuterShape( const TopoDS_Face& face,
1872                                                      TopAbs_ShapeEnum   type)
1873 {
1874   TopExp_Explorer exp( BRepTools::OuterWire( face ), type );
1875   if ( exp.More() )
1876     return exp.Current();
1877   return TopoDS_Shape();
1878 }
1879
1880 //================================================================================
1881   /*!
1882    * \brief Check that submesh is computed and try to compute it if is not
1883     * \param sm - submesh to compute
1884     * \param iterationNb - int used to stop infinite recursive call
1885     * \retval bool - true if computed
1886    */
1887 //================================================================================
1888
1889 bool StdMeshers_ProjectionUtils::MakeComputed(SMESH_subMesh * sm, const int iterationNb)
1890 {
1891   if ( iterationNb > 10 )
1892     RETURN_BAD_RESULT("Infinite recursive projection");
1893   if ( !sm )
1894     RETURN_BAD_RESULT("NULL submesh");
1895   if ( sm->IsMeshComputed() )
1896     return true;
1897
1898   SMESH_Mesh* mesh = sm->GetFather();
1899   SMESH_Gen* gen   = mesh->GetGen();
1900   SMESH_Algo* algo = gen->GetAlgo( *mesh, sm->GetSubShape() );
1901   if ( !algo )
1902   {
1903     if ( sm->GetSubShape().ShapeType() != TopAbs_COMPOUND )
1904       RETURN_BAD_RESULT("No algo assigned to submesh " << sm->GetId());
1905     // group
1906     bool computed = true;
1907     for ( TopoDS_Iterator grMember( sm->GetSubShape() ); grMember.More(); grMember.Next())
1908       if ( SMESH_subMesh* grSub = mesh->GetSubMesh( grMember.Value() ))
1909         if ( !MakeComputed( grSub, iterationNb + 1 ))
1910           computed = false;
1911     return computed;
1912   }
1913
1914   string algoType = algo->GetName();
1915   if ( algoType.substr(0, 11) != "Projection_")
1916     return gen->Compute( *mesh, sm->GetSubShape() );
1917
1918   // try to compute source mesh
1919
1920   const list <const SMESHDS_Hypothesis *> & hyps =
1921     algo->GetUsedHypothesis( *mesh, sm->GetSubShape() );
1922
1923   TopoDS_Shape srcShape;
1924   SMESH_Mesh* srcMesh = 0;
1925   list <const SMESHDS_Hypothesis*>::const_iterator hIt = hyps.begin();
1926   for ( ; srcShape.IsNull() && hIt != hyps.end(); ++hIt ) {
1927     string hypName = (*hIt)->GetName();
1928     if ( hypName == "ProjectionSource1D" ) {
1929       const StdMeshers_ProjectionSource1D * hyp =
1930         static_cast<const StdMeshers_ProjectionSource1D*>( *hIt );
1931       srcShape = hyp->GetSourceEdge();
1932       srcMesh = hyp->GetSourceMesh();
1933     }
1934     else if ( hypName == "ProjectionSource2D" ) {
1935       const StdMeshers_ProjectionSource2D * hyp =
1936         static_cast<const StdMeshers_ProjectionSource2D*>( *hIt );
1937       srcShape = hyp->GetSourceFace();
1938       srcMesh = hyp->GetSourceMesh();
1939     }
1940     else if ( hypName == "ProjectionSource3D" ) {
1941       const StdMeshers_ProjectionSource3D * hyp =
1942         static_cast<const StdMeshers_ProjectionSource3D*>( *hIt );
1943       srcShape = hyp->GetSource3DShape();
1944       srcMesh = hyp->GetSourceMesh();
1945     }
1946   }
1947   if ( srcShape.IsNull() ) // no projection source defined
1948     return gen->Compute( *mesh, sm->GetSubShape() );
1949
1950   if ( srcShape.IsSame( sm->GetSubShape() ))
1951     RETURN_BAD_RESULT("Projection from self");
1952     
1953   if ( !srcMesh )
1954     srcMesh = mesh;
1955
1956   if ( MakeComputed( srcMesh->GetSubMesh( srcShape ), iterationNb + 1 ))
1957     return gen->Compute( *mesh, sm->GetSubShape() );
1958
1959   return false;
1960 }
1961
1962 //================================================================================
1963 /*!
1964  * \brief Count nb of subshapes
1965  * \param shape - the shape
1966  * \param type - the type of subshapes to count
1967  * \retval int - the calculated number
1968  */
1969 //================================================================================
1970
1971 int StdMeshers_ProjectionUtils::Count(const TopoDS_Shape&    shape,
1972                                       const TopAbs_ShapeEnum type,
1973                                       const bool             ignoreSame)
1974 {
1975   if ( ignoreSame ) {
1976     TopTools_IndexedMapOfShape map;
1977     TopExp::MapShapes( shape, type, map );
1978     return map.Extent();
1979   }
1980   else {
1981     int nb = 0;
1982     for ( TopExp_Explorer exp( shape, type ); exp.More(); exp.Next() )
1983       ++nb;
1984     return nb;
1985   }
1986 }
1987
1988 //================================================================================
1989 /*!
1990  * \brief Return true if edge is a boundary of edgeContainer
1991  */
1992 //================================================================================
1993
1994 bool StdMeshers_ProjectionUtils::IsBoundaryEdge(const TopoDS_Edge&  edge,
1995                                                 const TopoDS_Shape& edgeContainer,
1996                                                 SMESH_Mesh&         mesh)
1997 {
1998   TopTools_IndexedMapOfShape facesOfEdgeContainer, facesNearEdge;
1999   TopExp::MapShapes( edgeContainer, TopAbs_FACE, facesOfEdgeContainer );
2000
2001   const TopTools_ListOfShape& EAncestors = mesh.GetAncestors(edge);
2002   TopTools_ListIteratorOfListOfShape itea(EAncestors);
2003   for(; itea.More(); itea.Next()) {
2004     if( itea.Value().ShapeType() == TopAbs_FACE &&
2005         facesOfEdgeContainer.Contains( itea.Value() ))
2006     {
2007       facesNearEdge.Add( itea.Value() );
2008       if ( facesNearEdge.Extent() > 1 )
2009         return false;
2010     }
2011   }
2012   return ( facesNearEdge.Extent() == 1 );
2013 }
2014
2015
2016 namespace {
2017
2018   SMESH_subMeshEventListener* GetSrcSubMeshListener();
2019
2020   //================================================================================
2021   /*!
2022    * \brief Listener that resets an event listener on source submesh when 
2023    * "ProjectionSource*D" hypothesis is modified
2024    */
2025   //================================================================================
2026
2027   struct HypModifWaiter: SMESH_subMeshEventListener
2028   {
2029     HypModifWaiter():SMESH_subMeshEventListener(0){} // won't be deleted by submesh
2030
2031     void ProcessEvent(const int event, const int eventType, SMESH_subMesh* subMesh,
2032                       EventListenerData*, const SMESH_Hypothesis*)
2033     {
2034       if ( event     == SMESH_subMesh::MODIF_HYP &&
2035            eventType == SMESH_subMesh::ALGO_EVENT)
2036       {
2037         // delete current source listener
2038         subMesh->DeleteEventListener( GetSrcSubMeshListener() );
2039         // let algo set a new one
2040         SMESH_Gen* gen = subMesh->GetFather()->GetGen();
2041         if ( SMESH_Algo* algo = gen->GetAlgo( *subMesh->GetFather(),
2042                                               subMesh->GetSubShape() ))
2043           algo->SetEventListener( subMesh );
2044       }
2045     }
2046   };
2047   //================================================================================
2048   /*!
2049    * \brief return static HypModifWaiter
2050    */
2051   //================================================================================
2052
2053   SMESH_subMeshEventListener* GetHypModifWaiter() {
2054     static HypModifWaiter aHypModifWaiter;
2055     return &aHypModifWaiter;
2056   }
2057   //================================================================================
2058   /*!
2059    * \brief return static listener for source shape submeshes
2060    */
2061   //================================================================================
2062
2063   SMESH_subMeshEventListener* GetSrcSubMeshListener() {
2064     static SMESH_subMeshEventListener srcListener(0); // won't be deleted by submesh
2065     return &srcListener;
2066   }
2067 }
2068
2069 //================================================================================
2070 /*!
2071  * \brief Set event listeners to submesh with projection algo
2072  * \param subMesh - submesh with projection algo
2073  * \param srcShape - source shape
2074  * \param srcMesh - source mesh
2075  */
2076 //================================================================================
2077
2078 void StdMeshers_ProjectionUtils::SetEventListener(SMESH_subMesh* subMesh,
2079                                                   TopoDS_Shape   srcShape,
2080                                                   SMESH_Mesh*    srcMesh)
2081 {
2082   // Set listener that resets an event listener on source submesh when
2083   // "ProjectionSource*D" hypothesis is modified since source shape can be changed
2084   subMesh->SetEventListener( GetHypModifWaiter(),0,subMesh);
2085
2086   // Set an event listener to submesh of the source shape
2087   if ( !srcShape.IsNull() )
2088   {
2089     if ( !srcMesh )
2090       srcMesh = subMesh->GetFather();
2091
2092     SMESH_subMesh* srcShapeSM = srcMesh->GetSubMesh( srcShape );
2093
2094     if ( srcShapeSM != subMesh ) {
2095       if ( srcShapeSM->GetSubMeshDS() &&
2096            srcShapeSM->GetSubMeshDS()->IsComplexSubmesh() )
2097       {  // source shape is a group
2098         TopExp_Explorer it(srcShapeSM->GetSubShape(), // explore the group into subshapes...
2099                            subMesh->GetSubShape().ShapeType()); // ...of target shape type
2100         for (; it.More(); it.Next())
2101         {
2102           SMESH_subMesh* srcSM = srcMesh->GetSubMesh( it.Current() );
2103           if ( srcSM != subMesh )
2104           {
2105             SMESH_subMeshEventListenerData* data =
2106               srcSM->GetEventListenerData(GetSrcSubMeshListener());
2107             if ( data )
2108               data->mySubMeshes.push_back( subMesh );
2109             else
2110               data = SMESH_subMeshEventListenerData::MakeData( subMesh );
2111             subMesh->SetEventListener ( GetSrcSubMeshListener(), data, srcSM );
2112           }
2113         }
2114       }
2115       else
2116       {
2117         subMesh->SetEventListener( GetSrcSubMeshListener(),
2118                                    SMESH_subMeshEventListenerData::MakeData( subMesh ),
2119                                    srcShapeSM );
2120       }
2121     }
2122   }
2123 }