Salome HOME
86edbacfda44b4c686139bcbf8486f20d54ecd49
[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       if ( edge1.Orientation() >= TopAbs_INTERNAL ) edge1.Orientation( TopAbs_FORWARD );
479       if ( edge2.Orientation() >= TopAbs_INTERNAL ) edge2.Orientation( TopAbs_FORWARD );
480       TopoDS_Vertex VV1[2], VV2[2];
481       TopExp::Vertices( edge1, VV1[0], VV1[1] );
482       TopExp::Vertices( edge2, VV2[0], VV2[1] );
483       int i1 = 0, i2 = 0;
484       if ( theMap.IsBound( VV1[ i1 ] )) i1 = 1;
485       if ( theMap.IsBound( VV2[ i2 ] )) i2 = 1;
486       InsertAssociation( VV1[ i1 ], VV2[ i2 ], theMap, bidirect);
487       InsertAssociation( theShape1, theShape2, theMap, bidirect );
488       return true;
489     }
490       // ----------------------------------------------------------------------
491     case TopAbs_FACE: { // TopAbs_FACE
492       // ----------------------------------------------------------------------
493       TopoDS_Face face1 = TopoDS::Face( theShape1 );
494       TopoDS_Face face2 = TopoDS::Face( theShape2 );
495       if ( face1.Orientation() >= TopAbs_INTERNAL ) face1.Orientation( TopAbs_FORWARD );
496       if ( face2.Orientation() >= TopAbs_INTERNAL ) face2.Orientation( TopAbs_FORWARD );
497
498       TopoDS_Vertex VV1[2], VV2[2];
499       // find a not closed edge of face1 both vertices of which are associated
500       int nbEdges = 0;
501       TopExp_Explorer exp ( face1, TopAbs_EDGE );
502       for ( ; VV2[ 1 ].IsNull() && exp.More(); exp.Next(), ++nbEdges ) {
503         TopExp::Vertices( TopoDS::Edge( exp.Current() ), VV1[0], VV1[1] );
504         if ( theMap.IsBound( VV1[0] ) ) {
505           VV2[ 0 ] = TopoDS::Vertex( theMap( VV1[0] ));
506           if ( theMap.IsBound( VV1[1] ) && !VV1[0].IsSame( VV1[1] ))
507             VV2[ 1 ] = TopoDS::Vertex( theMap( VV1[1] ));
508         }
509       }
510       if ( VV2[ 1 ].IsNull() ) { // 2 bound vertices not found
511         if ( nbEdges > 1 ) {
512           RETURN_BAD_RESULT("2 bound vertices not found" );
513         } else {
514           VV2[ 1 ] = VV2[ 0 ];
515         }
516       }
517       list< TopoDS_Edge > edges1, edges2;
518       int nbE = FindFaceAssociation( face1, VV1, face2, VV2, edges1, edges2 );
519       if ( !nbE ) RETURN_BAD_RESULT("FindFaceAssociation() failed");
520       FixAssocByPropagation( nbE, edges1, edges2, theMesh1, theMesh2 );
521
522       list< TopoDS_Edge >::iterator eIt1 = edges1.begin();
523       list< TopoDS_Edge >::iterator eIt2 = edges2.begin();
524       for ( ; eIt1 != edges1.end(); ++eIt1, ++eIt2 )
525       {
526         InsertAssociation( *eIt1, *eIt2, theMap, bidirect);
527         VV1[0] = TopExp::FirstVertex( *eIt1, true );
528         VV2[0] = TopExp::FirstVertex( *eIt2, true );
529         InsertAssociation( VV1[0], VV2[0], theMap, bidirect);
530       }
531       InsertAssociation( theShape1, theShape2, theMap, bidirect );
532       return true;
533     }
534       // ----------------------------------------------------------------------
535     case TopAbs_SHELL: // TopAbs_SHELL, TopAbs_SOLID
536     case TopAbs_SOLID: {
537       // ----------------------------------------------------------------------
538       TopoDS_Vertex VV1[2], VV2[2];
539       // try to find a not closed edge of shape1 both vertices of which are associated
540       TopoDS_Edge edge1;
541       TopExp_Explorer exp ( theShape1, TopAbs_EDGE );
542       for ( ; VV2[ 1 ].IsNull() && exp.More(); exp.Next() ) {
543         edge1 = TopoDS::Edge( exp.Current() );
544         if ( edge1.Orientation() >= TopAbs_INTERNAL ) edge1.Orientation( TopAbs_FORWARD );
545         TopExp::Vertices( edge1 , VV1[0], VV1[1] );
546         if ( theMap.IsBound( VV1[0] )) {
547           VV2[ 0 ] = TopoDS::Vertex( theMap( VV1[0] ));
548           if ( theMap.IsBound( VV1[1] ) && !VV1[0].IsSame( VV1[1] ))
549             VV2[ 1 ] = TopoDS::Vertex( theMap( VV1[1] ));
550         }
551       }
552       if ( VV2[ 1 ].IsNull() ) // 2 bound vertices not found
553         RETURN_BAD_RESULT("2 bound vertices not found" );
554       // get an edge2 of theShape2 corresponding to edge1
555       TopoDS_Edge edge2 = GetEdgeByVertices( theMesh2, VV2[ 0 ], VV2[ 1 ]);
556       if ( edge2.IsNull() )
557         RETURN_BAD_RESULT("GetEdgeByVertices() failed");
558
559       // build map of edge to faces if shapes are not subshapes of main ones
560       bool isSubOfMain = false;
561       if ( SMESHDS_SubMesh * sm = theMesh1->GetMeshDS()->MeshElements( theShape1 ))
562         isSubOfMain = !sm->IsComplexSubmesh();
563       else
564         isSubOfMain = theMesh1->GetMeshDS()->ShapeToIndex( theShape1 );
565       TAncestorMap e2f1, e2f2;
566       const TAncestorMap& edgeToFace1 = isSubOfMain ? theMesh1->GetAncestorMap() : e2f1;
567       const TAncestorMap& edgeToFace2 = isSubOfMain ? theMesh2->GetAncestorMap() : e2f2;
568       if (!isSubOfMain) {
569         TopExp::MapShapesAndAncestors( theShape1, TopAbs_EDGE, TopAbs_FACE, e2f1 );
570         TopExp::MapShapesAndAncestors( theShape2, TopAbs_EDGE, TopAbs_FACE, e2f2 );
571         if ( !edgeToFace1.Contains( edge1 ))
572           RETURN_BAD_RESULT("edge1 does not belong to theShape1");
573         if ( !edgeToFace2.Contains( edge2 ))
574           RETURN_BAD_RESULT("edge2 does not belong to theShape2");
575       }
576       //
577       // Look for 2 corresponing faces:
578       //
579       TopoDS_Shape F1, F2;
580
581       // get a face sharing edge1 (F1)
582       TopoDS_Shape FF2[2];
583       TopTools_ListIteratorOfListOfShape ancestIt1( edgeToFace1.FindFromKey( edge1 ));
584       for ( ; F1.IsNull() && ancestIt1.More(); ancestIt1.Next() )
585         if ( ancestIt1.Value().ShapeType() == TopAbs_FACE )
586           F1 = ancestIt1.Value().Oriented( TopAbs_FORWARD );
587       if ( F1.IsNull() )
588         RETURN_BAD_RESULT(" Face1 not found");
589
590       // get 2 faces sharing edge2 (one of them is F2)
591       TopTools_ListIteratorOfListOfShape ancestIt2( edgeToFace2.FindFromKey( edge2 ));
592       for ( int i = 0; FF2[1].IsNull() && ancestIt2.More(); ancestIt2.Next() )
593         if ( ancestIt2.Value().ShapeType() == TopAbs_FACE )
594           FF2[ i++ ] = ancestIt2.Value().Oriented( TopAbs_FORWARD );
595
596       // get oriented edge1 and edge2 from F1 and FF2[0]
597       for ( exp.Init( F1, TopAbs_EDGE ); exp.More(); exp.Next() )
598         if ( edge1.IsSame( exp.Current() )) {
599           edge1 = TopoDS::Edge( exp.Current() );
600           break;
601         }
602       for ( exp.Init( FF2[ 0 ], TopAbs_EDGE ); exp.More(); exp.Next() )
603         if ( edge2.IsSame( exp.Current() )) {
604           edge2 = TopoDS::Edge( exp.Current() );
605           break;
606         }
607
608       // compare first vertices of edge1 and edge2
609       TopExp::Vertices( edge1, VV1[0], VV1[1], true );
610       TopExp::Vertices( edge2, VV2[0], VV2[1], true );
611       F2 = FF2[ 0 ]; // (F2 !)
612       if ( !VV1[ 0 ].IsSame( theMap( VV2[ 0 ]))) {
613         edge2.Reverse();
614         if ( FF2[ 1 ].IsNull() )
615           F2.Reverse();
616         else
617           F2 = FF2[ 1 ];
618       }
619
620       TopTools_MapOfShape boundEdges;
621
622       // association of face subshapes and neighbour faces
623       list< pair < TopoDS_Face, TopoDS_Edge > > FE1, FE2;
624       list< pair < TopoDS_Face, TopoDS_Edge > >::iterator fe1, fe2;
625       FE1.push_back( make_pair( TopoDS::Face( F1 ), edge1 ));
626       FE2.push_back( make_pair( TopoDS::Face( F2 ), edge2 ));
627       for ( fe1 = FE1.begin(), fe2 = FE2.begin(); fe1 != FE1.end(); ++fe1, ++fe2 )
628       {
629         const TopoDS_Face& face1 = fe1->first;
630         if ( theMap.IsBound( face1 ) ) continue;
631         const TopoDS_Face& face2 = fe2->first;
632         edge1 = fe1->second;
633         edge2 = fe2->second;
634         TopExp::Vertices( edge1, VV1[0], VV1[1], true );
635         TopExp::Vertices( edge2, VV2[0], VV2[1], true );
636         list< TopoDS_Edge > edges1, edges2;
637         int nbE = FindFaceAssociation( face1, VV1, face2, VV2, edges1, edges2 );
638         if ( !nbE ) RETURN_BAD_RESULT("FindFaceAssociation() failed");
639         InsertAssociation( face1, face2, theMap, bidirect); // assoc faces
640         MESSAGE("Assoc FACE " << theMesh1->GetMeshDS()->ShapeToIndex( face1 )<<
641                 " to "        << theMesh2->GetMeshDS()->ShapeToIndex( face2 ));
642         if ( nbE == 2 && (edge1.IsSame( edges1.front())) != (edge2.IsSame( edges2.front())))
643         {
644           Reverse( edges2, nbE );
645         }
646         list< TopoDS_Edge >::iterator eIt1 = edges1.begin();
647         list< TopoDS_Edge >::iterator eIt2 = edges2.begin();
648         for ( ; eIt1 != edges1.end(); ++eIt1, ++eIt2 )
649         {
650           if ( !boundEdges.Add( *eIt1 )) continue; // already associated
651           InsertAssociation( *eIt1, *eIt2, theMap, bidirect);  // assoc edges
652           MESSAGE("Assoc edge " << theMesh1->GetMeshDS()->ShapeToIndex( *eIt1 )<<
653                   " to "        << theMesh2->GetMeshDS()->ShapeToIndex( *eIt2 ));
654           VV1[0] = TopExp::FirstVertex( *eIt1, true );
655           VV2[0] = TopExp::FirstVertex( *eIt2, true );
656           InsertAssociation( VV1[0], VV2[0], theMap, bidirect); // assoc vertices
657           MESSAGE("Assoc vertex " << theMesh1->GetMeshDS()->ShapeToIndex( VV1[0] )<<
658                   " to "          << theMesh2->GetMeshDS()->ShapeToIndex( VV2[0] ));
659
660           // add adjacent faces to process
661           TopoDS_Face nextFace1 = GetNextFace( edgeToFace1, *eIt1, face1 );
662           TopoDS_Face nextFace2 = GetNextFace( edgeToFace2, *eIt2, face2 );
663           if ( !nextFace1.IsNull() && !nextFace2.IsNull() ) {
664             FE1.push_back( make_pair( nextFace1, *eIt1 ));
665             FE2.push_back( make_pair( nextFace2, *eIt2 ));
666           }
667         }
668       }
669       InsertAssociation( theShape1, theShape2, theMap, bidirect );
670       return true;
671     }
672       // ----------------------------------------------------------------------
673     case TopAbs_COMPOUND: { // GROUP
674       // ----------------------------------------------------------------------
675       // Maybe groups contain only one member
676       TopoDS_Iterator it1( theShape1 ), it2( theShape2 );
677       TopAbs_ShapeEnum memberType = it1.Value().ShapeType();
678       int nbMembers = Count( theShape1, memberType, true );
679       if ( nbMembers == 0 ) return true;
680       if ( nbMembers == 1 ) {
681         return FindSubShapeAssociation( it1.Value(), theMesh1, it2.Value(), theMesh2, theMap );
682       }
683       // Try to make shells of faces
684       //
685       BRep_Builder builder;
686       TopoDS_Shell shell1, shell2;
687       builder.MakeShell(shell1); builder.MakeShell(shell2);
688       if ( memberType == TopAbs_FACE ) {
689         // just add faces of groups to shells
690         for (; it1.More(); it1.Next(), it2.Next() )
691           builder.Add( shell1, it1.Value() ), builder.Add( shell2, it2.Value() );
692       }
693       else if ( memberType == TopAbs_EDGE ) {
694         // Try to add faces sharing more than one edge of a group or
695         // sharing all its vertices with the group
696         TopTools_IndexedMapOfShape groupVertices[2];
697         TopExp::MapShapes( theShape1, TopAbs_VERTEX, groupVertices[0]);
698         TopExp::MapShapes( theShape2, TopAbs_VERTEX, groupVertices[1]);
699         //
700         TopTools_MapOfShape groupEdges[2], addedFaces[2];
701         bool hasInitAssoc = (!theMap.IsEmpty()), initAssocOK = !hasInitAssoc;
702         for (; it1.More(); it1.Next(), it2.Next() ) {
703           groupEdges[0].Add( it1.Value() );
704           groupEdges[1].Add( it2.Value() );
705           if ( !initAssocOK ) {
706             // for shell association there must be an edge with both vertices bound
707             TopoDS_Vertex v1, v2;
708             TopExp::Vertices( TopoDS::Edge( it1.Value().Oriented(TopAbs_FORWARD)), v1, v2 );
709             initAssocOK = ( theMap.IsBound( v1 ) && theMap.IsBound( v2 ));
710           }
711         }
712         for (int is2ndGroup = 0; initAssocOK && is2ndGroup < 2; ++is2ndGroup) {
713           const TopoDS_Shape& group = is2ndGroup ? theShape2: theShape1;
714           SMESH_Mesh*         mesh  = is2ndGroup ? theMesh2 : theMesh1;
715           TopoDS_Shell&       shell = is2ndGroup ? shell2   : shell1;
716           for ( TopoDS_Iterator it( group ); it.More(); it.Next() ) {
717             const TopoDS_Edge& edge = TopoDS::Edge( it.Value() );
718             TopoDS_Face face;
719             for ( int iF = 0; iF < 2; ++iF ) { // loop on 2 faces sharing edge
720               face = GetNextFace(mesh->GetAncestorMap(), edge, face);
721               if ( !face.IsNull() ) {
722                 int nbGroupEdges = 0;
723                 for ( TopExp_Explorer f( face, TopAbs_EDGE ); f.More(); f.Next())
724                   if ( groupEdges[ is2ndGroup ].Contains( f.Current() ))
725                     if ( ++nbGroupEdges > 1 )
726                       break;
727                 bool add = (nbGroupEdges > 1 || Count( face, TopAbs_EDGE, true ) == 1 );
728                 if ( !add ) {
729                   add = true;
730                   for ( TopExp_Explorer v( face, TopAbs_VERTEX ); add && v.More(); v.Next())
731                     add = groupVertices[ is2ndGroup ].Contains( v.Current() );
732                 }
733                 if ( add && addedFaces[ is2ndGroup ].Add( face ))
734                   builder.Add( shell, face );
735               }
736             }
737           }
738         }
739       } else {
740         RETURN_BAD_RESULT("Unexpected group type");
741       }
742       // Associate shells
743       //
744       int nbFaces1 = Count( shell1, TopAbs_FACE, 0 );
745       int nbFaces2 = Count( shell2, TopAbs_FACE, 0 );
746       if ( nbFaces1 != nbFaces2 )
747         RETURN_BAD_RESULT("Different nb of faces found for shells");
748       if ( nbFaces1 > 0 ) {
749         bool ok = false;
750         if ( nbFaces1 == 1 ) {
751           TopoDS_Shape F1 = TopoDS_Iterator( shell1 ).Value();
752           TopoDS_Shape F2 = TopoDS_Iterator( shell2 ).Value();
753           ok = FindSubShapeAssociation( F1, theMesh1, F2, theMesh2, theMap );
754         }
755         else {
756           ok = FindSubShapeAssociation(shell1, theMesh1, shell2, theMesh2, theMap );
757         }
758         // Check if all members are mapped 
759         if ( ok ) {
760           TopTools_MapOfShape boundMembers[2];
761           TopoDS_Iterator mIt;
762           for ( mIt.Initialize( theShape1 ); mIt.More(); mIt.Next())
763             if ( theMap.IsBound( mIt.Value() )) {
764               boundMembers[0].Add( mIt.Value() );
765               boundMembers[1].Add( theMap( mIt.Value() ));
766             }
767           if ( boundMembers[0].Extent() != nbMembers ) {
768             // make compounds of not bound members
769             TopoDS_Compound comp[2];
770             for ( int is2ndGroup = 0; is2ndGroup < 2; ++is2ndGroup ) {
771               builder.MakeCompound( comp[is2ndGroup] );
772               for ( mIt.Initialize( is2ndGroup ? theShape2:theShape1 ); mIt.More(); mIt.Next())
773                 if ( ! boundMembers[ is2ndGroup ].Contains( mIt.Value() ))
774                   builder.Add( comp[ is2ndGroup ], mIt.Value() );
775             }
776             // check if theMap contains initial association for the comp's
777             bool hasInitialAssoc = false;
778             if ( memberType == TopAbs_EDGE ) {
779               for ( TopExp_Explorer v( comp[0], TopAbs_VERTEX ); v.More(); v.Next())
780                 if ( theMap.IsBound( v.Current() )) {
781                   hasInitialAssoc = true;
782                   break;
783                 }
784             }
785             if ( hasInitialAssoc == bool( !theMap.IsEmpty() ))
786               ok = FindSubShapeAssociation( comp[0], theMesh1, comp[1], theMesh2, theMap );
787             else {
788               TShapeShapeMap tmpMap;
789               ok = FindSubShapeAssociation( comp[0], theMesh1, comp[1], theMesh2, tmpMap );
790               if ( ok ) {
791                 TopTools_DataMapIteratorOfDataMapOfShapeShape mapIt( tmpMap );
792                 for ( ; mapIt.More(); mapIt.Next() )
793                   theMap.Bind( mapIt.Key(), mapIt.Value());
794               }
795             }
796           }
797         }
798         return ok;
799       }
800       // Each edge of an edge group is shared by own faces
801       // ------------------------------------------------------------------
802       //
803       // map vertices to edges sharing them, avoid doubling edges in lists
804       TopTools_DataMapOfShapeListOfShape v2e[2];
805       for (int isFirst = 0; isFirst < 2; ++isFirst ) {
806         const TopoDS_Shape& group = isFirst ? theShape1 : theShape2;
807         TopTools_DataMapOfShapeListOfShape& veMap = v2e[ isFirst ? 0 : 1 ];
808         TopTools_MapOfShape addedEdges;
809         for ( TopExp_Explorer e( group, TopAbs_EDGE ); e.More(); e.Next() ) {
810           const TopoDS_Shape& edge = e.Current();
811           if ( addedEdges.Add( edge )) {
812             for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next()) {
813               const TopoDS_Shape& vertex = v.Current();
814               if ( !veMap.IsBound( vertex )) {
815                 TopTools_ListOfShape l;
816                 veMap.Bind( vertex, l );
817               }
818               veMap( vertex ).Append( edge );
819             }
820           }
821         }   
822       }
823       while ( !v2e[0].IsEmpty() )
824       {
825         // find a bound vertex
826         TopoDS_Vertex V[2];
827         TopTools_DataMapIteratorOfDataMapOfShapeListOfShape v2eIt( v2e[0] );
828         for ( ; v2eIt.More(); v2eIt.Next())
829           if ( theMap.IsBound( v2eIt.Key() )) {
830             V[0] = TopoDS::Vertex( v2eIt.Key() );
831             V[1] = TopoDS::Vertex( theMap( V[0] ));
832             break;
833           }
834         if ( V[0].IsNull() )
835           RETURN_BAD_RESULT("No more bound vertices");
836
837         while ( !V[0].IsNull() && v2e[0].IsBound( V[0] )) {
838           TopTools_ListOfShape& edges0 = v2e[0]( V[0] );
839           TopTools_ListOfShape& edges1 = v2e[1]( V[1] );
840           int nbE0 = edges0.Extent(), nbE1 = edges1.Extent();
841           if ( nbE0 != nbE1 )
842             RETURN_BAD_RESULT("Different nb of edges: "<< nbE0 << " != " << nbE1);
843
844           if ( nbE0 == 1 )
845           {
846             TopoDS_Edge e0 = TopoDS::Edge( edges0.First() );
847             TopoDS_Edge e1 = TopoDS::Edge( edges1.First() );
848             v2e[0].UnBind( V[0] );
849             v2e[1].UnBind( V[1] );
850             InsertAssociation( e0, e1, theMap, bidirect );
851             MESSAGE("Assoc edge " << theMesh1->GetMeshDS()->ShapeToIndex( e0 )<<
852                     " to "        << theMesh2->GetMeshDS()->ShapeToIndex( e1 ));
853             V[0] = GetNextVertex( e0, V[0] );
854             V[1] = GetNextVertex( e1, V[1] );
855             if ( !V[0].IsNull() ) {
856               InsertAssociation( V[0], V[1], theMap, bidirect );
857               MESSAGE("Assoc vertex " << theMesh1->GetMeshDS()->ShapeToIndex( V[0] )<<
858                       " to "          << theMesh2->GetMeshDS()->ShapeToIndex( V[1] ));
859             }
860           }
861           else if ( nbE0 == 2 )
862           {
863             // one of edges must have both ends bound
864             TopoDS_Vertex v0e0 = GetNextVertex( TopoDS::Edge( edges0.First() ), V[0] );
865             TopoDS_Vertex v1e0 = GetNextVertex( TopoDS::Edge( edges0.Last() ),  V[0] );
866             TopoDS_Vertex v0e1 = GetNextVertex( TopoDS::Edge( edges1.First() ), V[1] );
867             TopoDS_Vertex v1e1 = GetNextVertex( TopoDS::Edge( edges1.Last() ),  V[1] );
868             TopoDS_Shape e0b, e1b, e0n, e1n, v1b; // bound and not-bound
869             TopoDS_Vertex v0n, v1n;
870             if ( theMap.IsBound( v0e0 )) {
871               v0n = v1e0; e0b = edges0.First(); e0n = edges0.Last(); v1b = theMap( v0e0 );
872             } else if ( theMap.IsBound( v1e0 )) {
873               v0n = v0e0; e0n = edges0.First(); e0b = edges0.Last(); v1b = theMap( v1e0 );
874             } else {
875               RETURN_BAD_RESULT("None of vertices bound");
876             }
877             if ( v1b.IsSame( v1e1 )) {
878               v1n = v0e1; e1n = edges1.First(); e1b = edges1.Last();
879             } else {
880               v1n = v1e1; e1b = edges1.First(); e1n = edges1.Last();
881             }
882             InsertAssociation( e0b, e1b, theMap, bidirect );
883             InsertAssociation( e0n, e1n, theMap, bidirect );
884             InsertAssociation( v0n, v1n, theMap, bidirect );
885             MESSAGE("Assoc edge " << theMesh1->GetMeshDS()->ShapeToIndex( e0b )<<
886                     " to "        << theMesh2->GetMeshDS()->ShapeToIndex( e1b ));
887             MESSAGE("Assoc edge " << theMesh1->GetMeshDS()->ShapeToIndex( e0n )<<
888                     " to "        << theMesh2->GetMeshDS()->ShapeToIndex( e1n ));
889             MESSAGE("Assoc vertex " << theMesh1->GetMeshDS()->ShapeToIndex( v0n )<<
890                     " to "          << theMesh2->GetMeshDS()->ShapeToIndex( v1n ));
891             v2e[0].UnBind( V[0] );
892             v2e[1].UnBind( V[1] );
893             V[0] = v0n;
894             V[1] = v1n;
895           }
896           else {
897             RETURN_BAD_RESULT("Not implemented");
898           }
899         }
900       } //while ( !v2e[0].IsEmpty() )
901       return true;
902     }
903
904     default:
905       RETURN_BAD_RESULT("Unexpected shape type");
906
907     } // end switch by shape type
908   } // end case of available initial vertex association
909
910   //======================================================================
911   // 4) NO INITIAL VERTEX ASSOCIATION
912   //======================================================================
913
914   switch ( theShape1.ShapeType() ) {
915
916   case TopAbs_EDGE: {
917     // ----------------------------------------------------------------------
918     TopoDS_Edge edge1 = TopoDS::Edge( theShape1 );
919     TopoDS_Edge edge2 = TopoDS::Edge( theShape2 );
920     if ( IsPropagationPossible( theMesh1, theMesh2 ))
921     {
922       TopoDS_Edge prpEdge = GetPropagationEdge( theMesh1, edge2, edge1 ).second;
923       if ( !prpEdge.IsNull() )
924       {
925         TopoDS_Vertex VV1[2], VV2[2];
926         TopExp::Vertices( edge1,   VV1[0], VV1[1], true );
927         TopExp::Vertices( prpEdge, VV2[0], VV2[1], true );
928         InsertAssociation( VV1[ 0 ], VV2[ 0 ], theMap, bidirect);
929         InsertAssociation( VV1[ 1 ], VV2[ 1 ], theMap, bidirect);
930         if ( VV1[0].IsSame( VV1[1] ) || // one of edges is closed
931              VV2[0].IsSame( VV2[1] ) )
932         {
933           InsertAssociation( edge1, prpEdge, theMap, bidirect); // insert with a proper orientation
934         }
935         InsertAssociation( theShape1, theShape2, theMap, bidirect );
936         return true; // done
937       }
938     }
939     if ( IsClosedEdge( edge1 ) && IsClosedEdge( edge2 ))
940     {
941       // TODO: find out a proper orientation (is it possible?)
942       InsertAssociation( edge1, edge2, theMap, bidirect); // insert with a proper orientation
943       InsertAssociation( TopExp::FirstVertex(edge1), TopExp::FirstVertex(edge2),
944                          theMap, bidirect);
945       InsertAssociation( theShape1, theShape2, theMap, bidirect );
946       return true; // done
947     }
948     break; // try by vertex closeness
949   }
950
951   case TopAbs_FACE: {
952     // ----------------------------------------------------------------------
953     if ( IsPropagationPossible( theMesh1, theMesh2 )) // try by propagation in one mesh
954     {
955       TopoDS_Face face1 = TopoDS::Face(theShape1);
956       TopoDS_Face face2 = TopoDS::Face(theShape2);
957       if ( face1.Orientation() >= TopAbs_INTERNAL ) face1.Orientation( TopAbs_FORWARD );
958       if ( face2.Orientation() >= TopAbs_INTERNAL ) face2.Orientation( TopAbs_FORWARD );
959       TopoDS_Edge edge1, edge2;
960       // get outer edge of theShape1
961       edge1 = TopoDS::Edge( OuterShape( face1, TopAbs_EDGE ));
962       // find out if any edge of face2 is a propagation edge of outer edge1
963       map<int,TopoDS_Edge> propag_edges; // use map to find the closest propagation edge
964       for ( TopExp_Explorer exp( face2, TopAbs_EDGE ); exp.More(); exp.Next() ) {
965         edge2 = TopoDS::Edge( exp.Current() );
966         pair<int,TopoDS_Edge> step_edge = GetPropagationEdge( theMesh1, edge2, edge1 );
967         if ( !step_edge.second.IsNull() ) { // propagation found
968           propag_edges.insert( step_edge );
969           if ( step_edge.first == 1 ) break; // most close found
970         }
971       }
972       if ( !propag_edges.empty() ) // propagation found
973       {
974         edge2 = propag_edges.begin()->second;
975         TopoDS_Vertex VV1[2], VV2[2];
976         TopExp::Vertices( edge1, VV1[0], VV1[1], true );
977         TopExp::Vertices( edge2, VV2[0], VV2[1], true );
978         list< TopoDS_Edge > edges1, edges2;
979         int nbE = FindFaceAssociation( face1, VV1, face2, VV2, edges1, edges2 );
980         if ( !nbE ) RETURN_BAD_RESULT("FindFaceAssociation() failed");
981         if ( nbE == 2 ) // only 2 edges
982         {
983           // take care of proper association of propagated edges
984           bool same1 = edge1.IsSame( edges1.front() );
985           bool same2 = edge2.IsSame( edges2.front() );
986           if ( same1 != same2 )
987             Reverse(edges2, nbE);
988         }
989         // store association
990         list< TopoDS_Edge >::iterator eIt1 = edges1.begin();
991         list< TopoDS_Edge >::iterator eIt2 = edges2.begin();
992         for ( ; eIt1 != edges1.end(); ++eIt1, ++eIt2 )
993         {
994           InsertAssociation( *eIt1, *eIt2, theMap, bidirect);
995           VV1[0] = TopExp::FirstVertex( *eIt1, true );
996           VV2[0] = TopExp::FirstVertex( *eIt2, true );
997           InsertAssociation( VV1[0], VV2[0], theMap, bidirect);
998         }
999         InsertAssociation( theShape1, theShape2, theMap, bidirect );
1000         return true;
1001       }
1002     }
1003     break; // try by vertex closeness
1004   }
1005   case TopAbs_COMPOUND: {
1006     // ----------------------------------------------------------------------
1007     if ( IsPropagationPossible( theMesh1, theMesh2 )) {
1008
1009       // try to accosiate all using propagation
1010       if ( AssocGroupsByPropagation( theShape1, theShape2, *theMesh1, theMap ))
1011         return true;
1012
1013       // find a boundary edge for theShape1
1014       TopoDS_Edge E;
1015       for(TopExp_Explorer exp(theShape1, TopAbs_EDGE); exp.More(); exp.Next() ) {
1016         E = TopoDS::Edge( exp.Current() );
1017         if ( IsBoundaryEdge( E, theShape1, *theMesh1 ))
1018           break;
1019         else
1020           E.Nullify();
1021       }
1022       if ( E.IsNull() )
1023         break; // try by vertex closeness
1024
1025       // find association for vertices of edge E
1026       TopoDS_Vertex VV1[2], VV2[2];
1027       for(TopExp_Explorer eexp(E, TopAbs_VERTEX); eexp.More(); eexp.Next()) {
1028         TopoDS_Vertex V1 = TopoDS::Vertex( eexp.Current() );
1029         // look for an edge ending in E whose one vertex is in theShape1
1030         // and the other, in theShape2
1031         const TopTools_ListOfShape& Ancestors = theMesh1->GetAncestors(V1);
1032         TopTools_ListIteratorOfListOfShape ita(Ancestors);
1033         for(; ita.More(); ita.Next()) {
1034           if( ita.Value().ShapeType() != TopAbs_EDGE ) continue;
1035           TopoDS_Edge edge = TopoDS::Edge(ita.Value());
1036           bool FromShape1 = false;
1037           for(TopExp_Explorer expe(theShape1, TopAbs_EDGE); expe.More(); expe.Next() ) {
1038             if(edge.IsSame(expe.Current())) {
1039               FromShape1 = true;
1040               break;
1041             }
1042           }
1043           if(!FromShape1) {
1044             // is it an edge between theShape1 and theShape2?
1045             TopExp_Explorer expv(edge, TopAbs_VERTEX);
1046             TopoDS_Vertex V2 = TopoDS::Vertex( expv.Current() );
1047             if(V2.IsSame(V1)) {
1048               expv.Next();
1049               V2 = TopoDS::Vertex( expv.Current() );
1050             }
1051             bool FromShape2 = false;
1052             for ( expv.Init( theShape2, TopAbs_VERTEX ); expv.More(); expv.Next()) {
1053               if ( V2.IsSame( expv.Current() )) {
1054                 FromShape2 = true;
1055                 break;
1056               }
1057             }
1058             if ( FromShape2 ) {
1059               if ( VV1[0].IsNull() )
1060                 VV1[0] = V1, VV2[0] = V2;
1061               else
1062                 VV1[1] = V1, VV2[1] = V2;
1063               break; // from loop on ancestors of V1
1064             }
1065           }
1066         }
1067       }
1068       if ( !VV1[1].IsNull() ) {
1069         InsertAssociation( VV1[0], VV2[0], theMap, bidirect);
1070         InsertAssociation( VV1[1], VV2[1], theMap, bidirect);
1071         return FindSubShapeAssociation( theShape1, theMesh1, theShape2, theMesh2, theMap);
1072       }
1073     }
1074     break; // try by vertex closeness
1075   }
1076   default:;
1077   }
1078
1079   // Find association by closeness of vertices
1080   // ------------------------------------------
1081
1082   TopTools_IndexedMapOfShape vMap1, vMap2;
1083   TopExp::MapShapes( theShape1, TopAbs_VERTEX, vMap1 );
1084   TopExp::MapShapes( theShape2, TopAbs_VERTEX, vMap2 );
1085
1086   if ( vMap1.Extent() != vMap2.Extent() )
1087     RETURN_BAD_RESULT("Different nb of vertices");
1088
1089   if ( vMap1.Extent() == 1 ) {
1090     InsertAssociation( vMap1(1), vMap2(1), theMap, bidirect);
1091     if ( theShape1.ShapeType() == TopAbs_EDGE ) {
1092       InsertAssociation( theShape1, theShape2, theMap, bidirect );
1093       return true;
1094     }
1095     return FindSubShapeAssociation( theShape1, theMesh1, theShape2, theMesh2, theMap);
1096   }
1097
1098   // Find transformation to make the shapes be of similar size at same location
1099
1100   Bnd_Box box[2];
1101   for ( int i = 1; i <= vMap1.Extent(); ++i ) {
1102     box[ 0 ].Add( BRep_Tool::Pnt ( TopoDS::Vertex( vMap1( i ))));
1103     box[ 1 ].Add( BRep_Tool::Pnt ( TopoDS::Vertex( vMap2( i ))));
1104   }
1105
1106   gp_Pnt gc[2]; // box center
1107   double x0,y0,z0, x1,y1,z1;
1108   box[0].Get( x0,y0,z0, x1,y1,z1 );
1109   gc[0] = 0.5 * ( gp_XYZ( x0,y0,z0 ) + gp_XYZ( x1,y1,z1 ));
1110   box[1].Get( x0,y0,z0, x1,y1,z1 );
1111   gc[1] = 0.5 * ( gp_XYZ( x0,y0,z0 ) + gp_XYZ( x1,y1,z1 ));
1112
1113   // 1 -> 2
1114   gp_Vec vec01( gc[0], gc[1] );
1115   double scale = sqrt( box[1].SquareExtent() / box[0].SquareExtent() );
1116
1117   // Find 2 closest vertices
1118
1119   TopoDS_Vertex VV1[2], VV2[2];
1120   // get 2 linked vertices of shape 1 not belonging to an inner wire of a face
1121   TopoDS_Shape edge = theShape1;
1122   TopExp_Explorer expF( theShape1, TopAbs_FACE ), expE;
1123   if ( expF.More() ) {
1124     for ( ; expF.More(); expF.Next() ) {
1125       edge.Nullify();
1126       TopoDS_Shape wire = OuterShape( TopoDS::Face( expF.Current() ), TopAbs_WIRE );
1127       for ( expE.Init( wire, TopAbs_EDGE ); edge.IsNull() && expE.More(); expE.Next() )
1128         if ( !IsClosedEdge( TopoDS::Edge( expE.Current() )))
1129           edge = expE.Current();
1130       if ( !edge.IsNull() )
1131         break;
1132     }
1133   } else if (edge.ShapeType() != TopAbs_EDGE) { // no faces
1134     edge.Nullify();
1135     for ( expE.Init( theShape1, TopAbs_EDGE ); edge.IsNull() && expE.More(); expE.Next() )
1136       if ( !IsClosedEdge( TopoDS::Edge( expE.Current() )))
1137         edge = expE.Current();
1138   }
1139   if ( edge.IsNull() || edge.ShapeType() != TopAbs_EDGE )
1140     RETURN_BAD_RESULT("Edge not found");
1141
1142   TopExp::Vertices( TopoDS::Edge( edge.Oriented(TopAbs_FORWARD)), VV1[0], VV1[1]);
1143   if ( VV1[0].IsSame( VV1[1] ))
1144     RETURN_BAD_RESULT("Only closed edges");
1145
1146   // find vertices closest to 2 linked vertices of shape 1
1147   for ( int i1 = 0; i1 < 2; ++i1 )
1148   {
1149     double dist2 = DBL_MAX;
1150     gp_Pnt p1 = BRep_Tool::Pnt( VV1[ i1 ]);
1151     p1.Translate( vec01 );
1152     p1.Scale( gc[1], scale );
1153     for ( int i2 = 1; i2 <= vMap2.Extent(); ++i2 )
1154     {
1155       TopoDS_Vertex V2 = TopoDS::Vertex( vMap2( i2 ));
1156       gp_Pnt p2 = BRep_Tool::Pnt ( V2 );
1157       double d2 = p1.SquareDistance( p2 );
1158       if ( d2 < dist2 && !V2.IsSame( VV2[ 0 ])) {
1159         VV2[ i1 ] = V2; dist2 = d2;
1160       }
1161     }
1162   }
1163
1164   InsertAssociation( VV1[ 0 ], VV2[ 0 ], theMap, bidirect);
1165   InsertAssociation( VV1[ 1 ], VV2[ 1 ], theMap, bidirect);
1166   MESSAGE("Initial assoc VERT " << theMesh1->GetMeshDS()->ShapeToIndex( VV1[ 0 ] )<<
1167           " to "                << theMesh2->GetMeshDS()->ShapeToIndex( VV2[ 0 ] )<<
1168           "\nand         VERT " << theMesh1->GetMeshDS()->ShapeToIndex( VV1[ 1 ] )<<
1169           " to "                << theMesh2->GetMeshDS()->ShapeToIndex( VV2[ 1 ] ));
1170   if ( theShape1.ShapeType() == TopAbs_EDGE ) {
1171     InsertAssociation( theShape1, theShape2, theMap, bidirect );
1172     return true;
1173   }
1174
1175   return FindSubShapeAssociation( theShape1, theMesh1, theShape2, theMesh2, theMap );
1176 }
1177
1178 //================================================================================
1179 /*!
1180  * \brief Find association of edges of faces
1181  * \param face1 - face 1
1182  * \param VV1 - vertices of face 1
1183  * \param face2 - face 2
1184  * \param VV2 - vertices of face 2 associated with ones of face 1
1185  * \param edges1 - out list of edges of face 1
1186  * \param edges2 - out list of edges of face 2
1187  * \retval int - nb of edges in an outer wire in a success case, else zero
1188  */
1189 //================================================================================
1190
1191 int StdMeshers_ProjectionUtils::FindFaceAssociation(const TopoDS_Face&    face1,
1192                                                     TopoDS_Vertex         VV1[2],
1193                                                     const TopoDS_Face&    face2,
1194                                                     TopoDS_Vertex         VV2[2],
1195                                                     list< TopoDS_Edge > & edges1,
1196                                                     list< TopoDS_Edge > & edges2)
1197 {
1198   bool OK = false;
1199   list< int > nbEInW1, nbEInW2;
1200   int i_ok_wire_algo = -1;
1201   for ( int outer_wire_algo = 0; outer_wire_algo < 2 && !OK; ++outer_wire_algo )
1202   {
1203     edges1.clear();
1204     edges2.clear();
1205
1206     if ( SMESH_Block::GetOrderedEdges( face1, VV1[0], edges1, nbEInW1, outer_wire_algo) !=
1207          SMESH_Block::GetOrderedEdges( face2, VV2[0], edges2, nbEInW2, outer_wire_algo) )
1208       CONT_BAD_RESULT("Different number of wires in faces ");
1209
1210     if ( nbEInW1 != nbEInW2 )
1211       CONT_BAD_RESULT("Different number of edges in faces: " <<
1212                       nbEInW1.front() << " != " << nbEInW2.front());
1213
1214     i_ok_wire_algo = outer_wire_algo;
1215
1216     // Define if we need to reverse one of wires to make edges in lists match each other
1217
1218     bool reverse = false;
1219
1220     list< TopoDS_Edge >::iterator edgeIt;
1221     if ( !VV1[1].IsSame( TopExp::LastVertex( edges1.front(), true ))) {
1222       reverse = true;
1223       edgeIt = --edges1.end();
1224       // check if the second vertex belongs to the first or last edge in the wire
1225       if ( !VV1[1].IsSame( TopExp::FirstVertex( *edgeIt, true ))) {
1226         bool KO = true; // belongs to none
1227         if ( nbEInW1.size() > 1 ) { // several wires
1228           edgeIt = edges1.begin();
1229           std::advance( edgeIt, nbEInW1.front()-1 );
1230           KO = !VV1[1].IsSame( TopExp::FirstVertex( *edgeIt, true ));
1231         }
1232         if ( KO )
1233           CONT_BAD_RESULT("GetOrderedEdges() failed");
1234       }
1235     }
1236     edgeIt = --edges2.end();
1237     if ( !VV2[1].IsSame( TopExp::LastVertex( edges2.front(), true ))) {
1238       reverse = !reverse;
1239       // check if the second vertex belongs to the first or last edge in the wire
1240       if ( !VV2[1].IsSame( TopExp::FirstVertex( *edgeIt, true ))) {
1241         bool KO = true; // belongs to none
1242         if ( nbEInW2.size() > 1 ) { // several wires
1243           edgeIt = edges2.begin();
1244           std::advance( edgeIt, nbEInW2.front()-1 );
1245           KO = !VV2[1].IsSame( TopExp::FirstVertex( *edgeIt, true ));
1246         }
1247         if ( KO )
1248           CONT_BAD_RESULT("GetOrderedEdges() failed");
1249       }
1250     }
1251     if ( reverse )
1252     {
1253       Reverse( edges2 , nbEInW2.front());
1254       if (( VV1[1].IsSame( TopExp::LastVertex( edges1.front(), true ))) !=
1255           ( VV2[1].IsSame( TopExp::LastVertex( edges2.front(), true ))))
1256         CONT_BAD_RESULT("GetOrderedEdges() failed");
1257     }
1258     OK = true;
1259
1260   } // loop algos getting an outer wire
1261   
1262   // Try to orient all (if !OK) or only internal wires (issue 0020996) by UV similarity
1263   if (( !OK || nbEInW1.size() > 1 ) && i_ok_wire_algo > -1 )
1264   {
1265     // Check that Vec(VV1[0],VV1[1]) in 2D on face1 is the same
1266     // as Vec(VV2[0],VV2[1]) on face2
1267     double vTol = BRep_Tool::Tolerance( VV1[0] );
1268     BRepAdaptor_Surface surface1( face1, false );
1269     double vTolUV =
1270       surface1.UResolution( vTol ) + surface1.VResolution( vTol ); // let's be tolerant
1271     gp_Pnt2d v0f1UV = BRep_Tool::Parameters( VV1[0], face1 );
1272     gp_Pnt2d v0f2UV = BRep_Tool::Parameters( VV2[0], face2 );
1273     gp_Pnt2d v1f1UV = BRep_Tool::Parameters( VV1[1], face1 );
1274     gp_Pnt2d v1f2UV = BRep_Tool::Parameters( VV2[1], face2 );
1275     gp_Vec2d v01f1Vec( v0f1UV, v1f1UV );
1276     gp_Vec2d v01f2Vec( v0f2UV, v1f2UV );
1277     if ( Abs( v01f1Vec.X()-v01f2Vec.X()) < vTolUV && Abs( v01f1Vec.Y()-v01f2Vec.Y()) < vTolUV )
1278     {
1279       if ( i_ok_wire_algo != 1 )
1280       {
1281         edges1.clear();
1282         edges2.clear();
1283         SMESH_Block::GetOrderedEdges( face1, VV1[0], edges1, nbEInW1, i_ok_wire_algo);
1284         SMESH_Block::GetOrderedEdges( face2, VV2[0], edges2, nbEInW2, i_ok_wire_algo);
1285       }
1286       gp_XY dUV = v0f2UV.XY() - v0f1UV.XY(); // UV shift between 2 faces
1287       // skip edges of the outer wire (if the outer wire is OK)
1288       list< int >::iterator nbEInW = nbEInW1.begin();
1289       list< TopoDS_Edge >::iterator edge1Beg = edges1.begin(), edge2Beg = edges2.begin();
1290       if ( OK )
1291       {
1292         for ( int i = 0; i < *nbEInW; ++i )
1293           ++edge1Beg, ++edge2Beg;
1294         ++nbEInW;
1295       }
1296       for ( ; nbEInW != nbEInW1.end(); ++nbEInW ) // loop on wires
1297       {
1298         // reach an end of edges of a current wire
1299         list< TopoDS_Edge >::iterator edge1End = edge1Beg, edge2End = edge2Beg;
1300         for ( int i = 0; i < *nbEInW; ++i )
1301           ++edge1End, ++edge2End;
1302         // rotate edges2 untill coincident with edges1 in 2D
1303         v0f1UV = BRep_Tool::Parameters( TopExp::FirstVertex(*edge1Beg,true), face1 );
1304         v1f1UV = BRep_Tool::Parameters( TopExp::LastVertex (*edge1Beg,true), face1 );
1305         v0f1UV.ChangeCoord() += dUV;
1306         v1f1UV.ChangeCoord() += dUV;
1307         int i = *nbEInW;
1308         while ( --i > 0 && !sameVertexUV( *edge2Beg, face2, 0, v0f1UV, vTolUV ))
1309           edges2.splice( edge2End, edges2, edge2Beg++ ); // move edge2Beg to place before edge2End
1310         if ( sameVertexUV( *edge2Beg, face2, 0, v0f1UV, vTolUV ))
1311         {
1312           if ( nbEInW == nbEInW1.begin() )
1313             OK = true; // OK is for the first wire
1314           // reverse edges2 if needed
1315           if ( !sameVertexUV( *edge2Beg, face2, 1, v1f1UV, vTolUV ))
1316           {
1317             Reverse( edges2 , *nbEInW, distance( edges2.begin(),edge2Beg ));
1318             // set correct edge2End
1319             edge2End = edges2.begin();
1320             std::advance( edge2End, std::accumulate( nbEInW1.begin(), nbEInW, *nbEInW));
1321           }
1322         }
1323         // prepare to the next wire loop
1324         edge1Beg = edge1End, edge2Beg = edge2End;
1325       }
1326     }
1327   }
1328
1329   return OK ? nbEInW1.front() : 0;
1330 }
1331
1332 //=======================================================================
1333 //function : InitVertexAssociation
1334 //purpose  : 
1335 //=======================================================================
1336
1337 void StdMeshers_ProjectionUtils::InitVertexAssociation( const SMESH_Hypothesis* theHyp,
1338                                                         TShapeShapeMap &        theAssociationMap,
1339                                                         const TopoDS_Shape&     theTargetShape)
1340 {
1341   string hypName = theHyp->GetName();
1342   if ( hypName == "ProjectionSource1D" ) {
1343     const StdMeshers_ProjectionSource1D * hyp =
1344       static_cast<const StdMeshers_ProjectionSource1D*>( theHyp );
1345     if ( hyp->HasVertexAssociation() )
1346       InsertAssociation( hyp->GetSourceVertex(),hyp->GetTargetVertex(),theAssociationMap);
1347   }
1348   else if ( hypName == "ProjectionSource2D" ) {
1349     const StdMeshers_ProjectionSource2D * hyp =
1350       static_cast<const StdMeshers_ProjectionSource2D*>( theHyp );
1351     if ( hyp->HasVertexAssociation() ) {
1352       InsertAssociation( hyp->GetSourceVertex(1),hyp->GetTargetVertex(1),theAssociationMap);
1353       InsertAssociation( hyp->GetSourceVertex(2),hyp->GetTargetVertex(2),theAssociationMap);
1354     }
1355   }
1356   else if ( hypName == "ProjectionSource3D" ) {
1357     const StdMeshers_ProjectionSource3D * hyp =
1358       static_cast<const StdMeshers_ProjectionSource3D*>( theHyp );
1359     if ( hyp->HasVertexAssociation() ) {
1360       InsertAssociation( hyp->GetSourceVertex(1),hyp->GetTargetVertex(1),theAssociationMap);
1361       InsertAssociation( hyp->GetSourceVertex(2),hyp->GetTargetVertex(2),theAssociationMap);
1362     }
1363   }
1364 }
1365
1366 //=======================================================================
1367 /*!
1368  * \brief Inserts association theShape1 <-> theShape2 to TShapeShapeMap
1369  * \param theShape1 - shape 1
1370  * \param theShape2 - shape 2
1371  * \param theAssociationMap - association map 
1372  * \retval bool - true if there was no association for these shapes before
1373  */
1374 //=======================================================================
1375
1376 bool StdMeshers_ProjectionUtils::InsertAssociation( const TopoDS_Shape& theShape1,
1377                                                     const TopoDS_Shape& theShape2,
1378                                                     TShapeShapeMap &    theAssociationMap,
1379                                                     const bool          theBidirectional)
1380 {
1381   if ( !theShape1.IsNull() && !theShape2.IsNull() ) {
1382     SHOW_SHAPE(theShape1,"Assoc ");
1383     SHOW_SHAPE(theShape2," to ");
1384     bool isNew = ( theAssociationMap.Bind( theShape1, theShape2 ));
1385     if ( theBidirectional )
1386       theAssociationMap.Bind( theShape2, theShape1 );
1387     return isNew;
1388   }
1389   else {
1390     throw SALOME_Exception("StdMeshers_ProjectionUtils: attempt to associate NULL shape");
1391   }
1392   return false;
1393 }
1394
1395 //=======================================================================
1396 /*!
1397  * \brief Finds an edge by its vertices in a main shape of the mesh
1398  * \param aMesh - the mesh
1399  * \param V1 - vertex 1
1400  * \param V2 - vertex 2
1401  * \retval TopoDS_Edge - found edge
1402  */
1403 //=======================================================================
1404
1405 TopoDS_Edge StdMeshers_ProjectionUtils::GetEdgeByVertices( SMESH_Mesh*          theMesh,
1406                                                            const TopoDS_Vertex& theV1,
1407                                                            const TopoDS_Vertex& theV2)
1408 {
1409   if ( theMesh && !theV1.IsNull() && !theV2.IsNull() )
1410   {
1411     TopTools_ListIteratorOfListOfShape ancestorIt( theMesh->GetAncestors( theV1 ));
1412     for ( ; ancestorIt.More(); ancestorIt.Next() )
1413       if ( ancestorIt.Value().ShapeType() == TopAbs_EDGE )
1414         for ( TopExp_Explorer expV ( ancestorIt.Value(), TopAbs_VERTEX );
1415               expV.More();
1416               expV.Next() )
1417           if ( theV2.IsSame( expV.Current() ))
1418             return TopoDS::Edge( ancestorIt.Value() );
1419   }
1420   return TopoDS_Edge();
1421 }
1422
1423 //================================================================================
1424 /*!
1425  * \brief Return another face sharing an edge
1426  * \param edgeToFaces - data map of descendants to ancestors
1427  * \param edge - edge
1428  * \param face - face
1429  * \retval TopoDS_Face - found face
1430  */
1431 //================================================================================
1432
1433 TopoDS_Face StdMeshers_ProjectionUtils::GetNextFace( const TAncestorMap& edgeToFaces,
1434                                                      const TopoDS_Edge&  edge,
1435                                                      const TopoDS_Face&  face)
1436 {
1437 //   if ( !edge.IsNull() && !face.IsNull() && edgeToFaces.Contains( edge ))
1438   if ( !edge.IsNull() && edgeToFaces.Contains( edge )) // PAL16202
1439   {
1440     TopTools_ListIteratorOfListOfShape ancestorIt( edgeToFaces.FindFromKey( edge ));
1441     for ( ; ancestorIt.More(); ancestorIt.Next() )
1442       if ( ancestorIt.Value().ShapeType() == TopAbs_FACE &&
1443            !face.IsSame( ancestorIt.Value() ))
1444         return TopoDS::Face( ancestorIt.Value() );
1445   }
1446   return TopoDS_Face();
1447 }
1448
1449 //================================================================================
1450 /*!
1451  * \brief Return other vertex of an edge
1452  */
1453 //================================================================================
1454
1455 TopoDS_Vertex StdMeshers_ProjectionUtils::GetNextVertex(const TopoDS_Edge&   edge,
1456                                                         const TopoDS_Vertex& vertex)
1457 {
1458   TopoDS_Vertex vF,vL;
1459   TopExp::Vertices(edge,vF,vL);
1460   if ( vF.IsSame( vL ))
1461     return TopoDS_Vertex();
1462   return vertex.IsSame( vF ) ? vL : vF; 
1463 }
1464
1465 //================================================================================
1466 /*!
1467  * \brief Return a propagation edge
1468  * \param aMesh - mesh
1469  * \param theEdge - edge to find by propagation
1470  * \param fromEdge - start edge for propagation
1471  * \retval pair<int,TopoDS_Edge> - propagation step and found edge
1472  */
1473 //================================================================================
1474
1475 pair<int,TopoDS_Edge>
1476 StdMeshers_ProjectionUtils::GetPropagationEdge( SMESH_Mesh*        aMesh,
1477                                                 const TopoDS_Edge& theEdge,
1478                                                 const TopoDS_Edge& fromEdge)
1479 {
1480   SMESH_IndexedMapOfShape aChain;
1481   int step = 0;
1482
1483   // List of edges, added to chain on the previous cycle pass
1484   TopTools_ListOfShape listPrevEdges;
1485   listPrevEdges.Append(fromEdge);
1486
1487   // Collect all edges pass by pass
1488   while (listPrevEdges.Extent() > 0) {
1489     step++;
1490     // List of edges, added to chain on this cycle pass
1491     TopTools_ListOfShape listCurEdges;
1492
1493     // Find the next portion of edges
1494     TopTools_ListIteratorOfListOfShape itE (listPrevEdges);
1495     for (; itE.More(); itE.Next()) {
1496       TopoDS_Shape anE = itE.Value();
1497
1498       // Iterate on faces, having edge <anE>
1499       TopTools_ListIteratorOfListOfShape itA (aMesh->GetAncestors(anE));
1500       for (; itA.More(); itA.Next()) {
1501         TopoDS_Shape aW = itA.Value();
1502
1503         // There are objects of different type among the ancestors of edge
1504         if (aW.ShapeType() == TopAbs_WIRE) {
1505           TopoDS_Shape anOppE;
1506
1507           BRepTools_WireExplorer aWE (TopoDS::Wire(aW));
1508           Standard_Integer nb = 1, found = 0;
1509           TopTools_Array1OfShape anEdges (1,4);
1510           for (; aWE.More(); aWE.Next(), nb++) {
1511             if (nb > 4) {
1512               found = 0;
1513               break;
1514             }
1515             anEdges(nb) = aWE.Current();
1516             if (anEdges(nb).IsSame(anE)) found = nb;
1517           }
1518
1519           if (nb == 5 && found > 0) {
1520             // Quadrangle face found, get an opposite edge
1521             Standard_Integer opp = found + 2;
1522             if (opp > 4) opp -= 4;
1523             anOppE = anEdges(opp);
1524
1525             // add anOppE to aChain if ...
1526             if (!aChain.Contains(anOppE)) { // ... anOppE is not in aChain
1527               // Add found edge to the chain oriented so that to
1528               // have it co-directed with a forward MainEdge
1529               TopAbs_Orientation ori = anE.Orientation();
1530               if ( anEdges(opp).Orientation() == anEdges(found).Orientation() )
1531                 ori = TopAbs::Reverse( ori );
1532               anOppE.Orientation( ori );
1533               if ( anOppE.IsSame( theEdge ))
1534                 return make_pair( step, TopoDS::Edge( anOppE ));
1535               aChain.Add(anOppE);
1536               listCurEdges.Append(anOppE);
1537             }
1538           } // if (nb == 5 && found > 0)
1539         } // if (aF.ShapeType() == TopAbs_WIRE)
1540       } // for (; itF.More(); itF.Next())
1541     } // for (; itE.More(); itE.Next())
1542
1543     listPrevEdges = listCurEdges;
1544   } // while (listPrevEdges.Extent() > 0)
1545
1546   return make_pair( INT_MAX, TopoDS_Edge());
1547 }
1548
1549 //================================================================================
1550   /*!
1551    * \brief Find corresponding nodes on two faces
1552     * \param face1 - the first face
1553     * \param mesh1 - mesh containing elements on the first face
1554     * \param face2 - the second face
1555     * \param mesh2 - mesh containing elements on the second face
1556     * \param assocMap - map associating subshapes of the faces
1557     * \param node1To2Map - map containing found matching nodes
1558     * \retval bool - is a success
1559    */
1560 //================================================================================
1561
1562 bool StdMeshers_ProjectionUtils::
1563 FindMatchingNodesOnFaces( const TopoDS_Face&     face1,
1564                           SMESH_Mesh*            mesh1,
1565                           const TopoDS_Face&     face2,
1566                           SMESH_Mesh*            mesh2,
1567                           const TShapeShapeMap & assocMap,
1568                           TNodeNodeMap &         node1To2Map)
1569 {
1570   SMESHDS_Mesh* meshDS1 = mesh1->GetMeshDS();
1571   SMESHDS_Mesh* meshDS2 = mesh2->GetMeshDS();
1572   
1573   SMESH_MesherHelper helper1( *mesh1 );
1574   SMESH_MesherHelper helper2( *mesh2 );
1575
1576   // Get corresponding submeshes and roughly check match of meshes
1577
1578   SMESHDS_SubMesh * SM2 = meshDS2->MeshElements( face2 );
1579   SMESHDS_SubMesh * SM1 = meshDS1->MeshElements( face1 );
1580   if ( !SM2 || !SM1 )
1581     RETURN_BAD_RESULT("Empty submeshes");
1582   if ( SM2->NbNodes()    != SM1->NbNodes() ||
1583        SM2->NbElements() != SM1->NbElements() )
1584     RETURN_BAD_RESULT("Different meshes on corresponding faces "
1585                       << meshDS1->ShapeToIndex( face1 ) << " and "
1586                       << meshDS2->ShapeToIndex( face2 ));
1587   if ( SM2->NbElements() == 0 )
1588     RETURN_BAD_RESULT("Empty submeshes");
1589
1590   helper1.SetSubShape( face1 );
1591   helper2.SetSubShape( face2 );
1592   if ( helper1.HasSeam() != helper2.HasSeam() )
1593     RETURN_BAD_RESULT("Different faces' geometry");
1594
1595   // Data to call SMESH_MeshEditor::FindMatchingNodes():
1596
1597   // 1. Nodes of corresponding links:
1598
1599   // get 2 matching edges, try to find not seam ones
1600   TopoDS_Edge edge1, edge2, seam1, seam2, anyEdge1, anyEdge2;
1601   TopExp_Explorer eE( OuterShape( face2, TopAbs_WIRE ), TopAbs_EDGE );
1602   do {
1603     // edge 2
1604     TopoDS_Edge e2 = TopoDS::Edge( eE.Current() );
1605     eE.Next();
1606     // edge 1
1607     if ( !assocMap.IsBound( e2 ))
1608       RETURN_BAD_RESULT("Association not found for edge " << meshDS2->ShapeToIndex( e2 ));
1609     TopoDS_Edge e1 = TopoDS::Edge( assocMap( e2 ));
1610     if ( !helper1.IsSubShape( e1, face1 ))
1611       RETURN_BAD_RESULT("Wrong association, edge " << meshDS1->ShapeToIndex( e1 ) <<
1612                         " isn't a subshape of face " << meshDS1->ShapeToIndex( face1 ));
1613     // check that there are nodes on edges
1614     SMESHDS_SubMesh * eSM1 = meshDS1->MeshElements( e1 );
1615     SMESHDS_SubMesh * eSM2 = meshDS2->MeshElements( e2 );
1616     bool nodesOnEdges = ( eSM1 && eSM2 && eSM1->NbNodes() && eSM2->NbNodes() );
1617     // check that the nodes on edges belong to faces
1618     // (as NETGEN ignores nodes on the degenerated geom edge)
1619     bool nodesOfFaces = false;
1620     if ( nodesOnEdges ) {
1621       const SMDS_MeshNode* n1 = eSM1->GetNodes()->next();
1622       const SMDS_MeshNode* n2 = eSM2->GetNodes()->next();
1623       nodesOfFaces = ( n1->GetInverseElementIterator(SMDSAbs_Face)->more() &&
1624                        n2->GetInverseElementIterator(SMDSAbs_Face)->more() );
1625     }
1626     if ( nodesOfFaces )
1627     {
1628       if ( helper2.IsRealSeam( e2 )) {
1629         seam1 = e1; seam2 = e2;
1630       }
1631       else {
1632         edge1 = e1; edge2 = e2;
1633       }
1634     }
1635     else {
1636       anyEdge1 = e1; anyEdge2 = e2;
1637     }
1638   } while ( edge2.IsNull() && eE.More() );
1639   //
1640   if ( edge2.IsNull() ) {
1641     edge1 = seam1; edge2 = seam2;
1642   }
1643   bool hasNodesOnEdge = (! edge2.IsNull() );
1644   if ( !hasNodesOnEdge ) {
1645     // 0020338 - nb segments == 1
1646     edge1 = anyEdge1; edge2 = anyEdge2;
1647   }
1648
1649   // get 2 matching vertices
1650   TopoDS_Vertex V2 = TopExp::FirstVertex( TopoDS::Edge( edge2 ));
1651   if ( !assocMap.IsBound( V2 ))
1652     RETURN_BAD_RESULT("Association not found for vertex " << meshDS2->ShapeToIndex( V2 ));
1653   TopoDS_Vertex V1 = TopoDS::Vertex( assocMap( V2 ));
1654
1655   // nodes on vertices
1656   const SMDS_MeshNode* vNode1 = SMESH_Algo::VertexNode( V1, meshDS1 );
1657   const SMDS_MeshNode* vNode2 = SMESH_Algo::VertexNode( V2, meshDS2 );
1658   if ( !vNode1 ) RETURN_BAD_RESULT("No node on vertex #" << meshDS1->ShapeToIndex( V1 ));
1659   if ( !vNode2 ) RETURN_BAD_RESULT("No node on vertex #" << meshDS2->ShapeToIndex( V2 ));
1660
1661   // nodes on edges linked with nodes on vertices
1662   const SMDS_MeshNode* nullNode = 0;
1663   vector< const SMDS_MeshNode*> eNode1( 2, nullNode );
1664   vector< const SMDS_MeshNode*> eNode2( 2, nullNode );
1665   if ( hasNodesOnEdge )
1666   {
1667     int nbNodeToGet = 1;
1668     if ( IsClosedEdge( edge1 ) || IsClosedEdge( edge2 ) )
1669       nbNodeToGet = 2;
1670     for ( int is2 = 0; is2 < 2; ++is2 )
1671     {
1672       TopoDS_Edge &     edge  = is2 ? edge2 : edge1;
1673       SMESHDS_Mesh *    smDS  = is2 ? meshDS2 : meshDS1;
1674       SMESHDS_SubMesh* edgeSM = smDS->MeshElements( edge );
1675       // nodes linked with ones on vertices
1676       const SMDS_MeshNode*           vNode = is2 ? vNode2 : vNode1;
1677       vector< const SMDS_MeshNode*>& eNode = is2 ? eNode2 : eNode1;
1678       int nbGotNode = 0;
1679       SMDS_ElemIteratorPtr vElem = vNode->GetInverseElementIterator(SMDSAbs_Edge);
1680       while ( vElem->more() && nbGotNode != nbNodeToGet ) {
1681         const SMDS_MeshElement* elem = vElem->next();
1682         if ( edgeSM->Contains( elem ))
1683           eNode[ nbGotNode++ ] = 
1684             ( elem->GetNode(0) == vNode ) ? elem->GetNode(1) : elem->GetNode(0);
1685       }
1686       if ( nbGotNode > 1 ) // sort found nodes by param on edge
1687       {
1688         SMESH_MesherHelper* helper = is2 ? &helper2 : &helper1;
1689         double u0 = helper->GetNodeU( edge, eNode[ 0 ]);
1690         double u1 = helper->GetNodeU( edge, eNode[ 1 ]);
1691         if ( u0 > u1 ) std::swap( eNode[ 0 ], eNode[ 1 ]);
1692       }
1693       if ( nbGotNode == 0 )
1694         RETURN_BAD_RESULT("Found no nodes on edge " << smDS->ShapeToIndex( edge ) <<
1695                           " linked to " << vNode );
1696     }
1697   }
1698   else // 0020338 - nb segments == 1
1699   {
1700     // get 2 other matching vertices
1701     V2 = TopExp::LastVertex( TopoDS::Edge( edge2 ));
1702     if ( !assocMap.IsBound( V2 ))
1703       RETURN_BAD_RESULT("Association not found for vertex " << meshDS2->ShapeToIndex( V2 ));
1704     V1 = TopoDS::Vertex( assocMap( V2 ));
1705
1706     // nodes on vertices
1707     eNode1[0] = SMESH_Algo::VertexNode( V1, meshDS1 );
1708     eNode2[0] = SMESH_Algo::VertexNode( V2, meshDS2 );
1709     if ( !eNode1[0] ) RETURN_BAD_RESULT("No node on vertex #" << meshDS1->ShapeToIndex( V1 ));
1710     if ( !eNode2[0] ) RETURN_BAD_RESULT("No node on vertex #" << meshDS2->ShapeToIndex( V2 ));
1711   }
1712
1713   // 2. face sets
1714
1715   set<const SMDS_MeshElement*> Elems1, Elems2;
1716   for ( int is2 = 0; is2 < 2; ++is2 )
1717   {
1718     set<const SMDS_MeshElement*> & elems = is2 ? Elems2 : Elems1;
1719     SMESHDS_SubMesh*                  sm = is2 ? SM2 : SM1;
1720     SMESH_MesherHelper*           helper = is2 ? &helper2 : &helper1;
1721     const TopoDS_Face &             face = is2 ? face2 : face1;
1722     SMDS_ElemIteratorPtr eIt = sm->GetElements();
1723
1724     if ( !helper->IsRealSeam( is2 ? edge2 : edge1 ))
1725     {
1726       while ( eIt->more() ) elems.insert( eIt->next() );
1727     }
1728     else
1729     {
1730       // the only suitable edge is seam, i.e. it is a sphere.
1731       // FindMatchingNodes() will not know which way to go from any edge.
1732       // So we ignore all faces having nodes on edges or vertices except
1733       // one of faces sharing current start nodes
1734
1735       // find a face to keep
1736       const SMDS_MeshElement* faceToKeep = 0;
1737       const SMDS_MeshNode* vNode = is2 ? vNode2 : vNode1;
1738       const SMDS_MeshNode* eNode = is2 ? eNode2[0] : eNode1[0];
1739       TIDSortedElemSet inSet, notInSet;
1740
1741       const SMDS_MeshElement* f1 =
1742         SMESH_MeshEditor::FindFaceInSet( vNode, eNode, inSet, notInSet );
1743       if ( !f1 ) RETURN_BAD_RESULT("The first face on seam not found");
1744       notInSet.insert( f1 );
1745
1746       const SMDS_MeshElement* f2 =
1747         SMESH_MeshEditor::FindFaceInSet( vNode, eNode, inSet, notInSet );
1748       if ( !f2 ) RETURN_BAD_RESULT("The second face on seam not found");
1749
1750       // select a face with less UV of vNode
1751       const SMDS_MeshNode* notSeamNode[2] = {0, 0};
1752       for ( int iF = 0; iF < 2; ++iF ) {
1753         const SMDS_MeshElement* f = ( iF ? f2 : f1 );
1754         for ( int i = 0; !notSeamNode[ iF ] && i < f->NbNodes(); ++i ) {
1755           const SMDS_MeshNode* node = f->GetNode( i );
1756           if ( !helper->IsSeamShape( node->getshapeId() ))
1757             notSeamNode[ iF ] = node;
1758         }
1759       }
1760       gp_Pnt2d uv1 = helper->GetNodeUV( face, vNode, notSeamNode[0] );
1761       gp_Pnt2d uv2 = helper->GetNodeUV( face, vNode, notSeamNode[1] );
1762       if ( uv1.X() + uv1.Y() > uv2.X() + uv2.Y() )
1763         faceToKeep = f2;
1764       else
1765         faceToKeep = f1;
1766
1767       // fill elem set
1768       elems.insert( faceToKeep );
1769       while ( eIt->more() ) {
1770         const SMDS_MeshElement* f = eIt->next();
1771         int nbNodes = f->NbNodes();
1772         if ( f->IsQuadratic() )
1773           nbNodes /= 2;
1774         bool onBnd = false;
1775         for ( int i = 0; !onBnd && i < nbNodes; ++i ) {
1776           const SMDS_MeshNode* node = f->GetNode( i );
1777           onBnd = ( node->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE);
1778         }
1779         if ( !onBnd )
1780           elems.insert( f );
1781       }
1782       // add also faces adjacent to faceToKeep
1783       int nbNodes = faceToKeep->NbNodes();
1784       if ( faceToKeep->IsQuadratic() ) nbNodes /= 2;
1785       notInSet.insert( f1 );
1786       notInSet.insert( f2 );
1787       for ( int i = 0; i < nbNodes; ++i ) {
1788         const SMDS_MeshNode* n1 = faceToKeep->GetNode( i );
1789         const SMDS_MeshNode* n2 = faceToKeep->GetNode( i+1 % nbNodes );
1790         f1 = SMESH_MeshEditor::FindFaceInSet( n1, n2, inSet, notInSet );
1791         if ( f1 )
1792           elems.insert( f1 );
1793       }
1794     } // case on a sphere
1795   } // loop on 2 faces
1796
1797   //  int quadFactor = (*Elems1.begin())->IsQuadratic() ? 2 : 1;
1798
1799   node1To2Map.clear();
1800   int res = SMESH_MeshEditor::FindMatchingNodes( Elems1, Elems2,
1801                                                  vNode1, vNode2,
1802                                                  eNode1[0], eNode2[0],
1803                                                  node1To2Map);
1804   if ( res != SMESH_MeshEditor::SEW_OK )
1805     RETURN_BAD_RESULT("FindMatchingNodes() result " << res );
1806
1807   // On a sphere, add matching nodes on the edge
1808
1809   if ( helper1.IsRealSeam( edge1 ))
1810   {
1811     // sort nodes on edges by param on edge
1812     map< double, const SMDS_MeshNode* > u2nodesMaps[2];
1813     for ( int is2 = 0; is2 < 2; ++is2 )
1814     {
1815       TopoDS_Edge &     edge  = is2 ? edge2 : edge1;
1816       SMESHDS_Mesh *    smDS  = is2 ? meshDS2 : meshDS1;
1817       SMESHDS_SubMesh* edgeSM = smDS->MeshElements( edge );
1818       map< double, const SMDS_MeshNode* > & pos2nodes = u2nodesMaps[ is2 ];
1819
1820       SMDS_NodeIteratorPtr nIt = edgeSM->GetNodes();
1821       while ( nIt->more() ) {
1822         const SMDS_MeshNode* node = nIt->next();
1823         const SMDS_EdgePosition* pos =
1824           static_cast<const SMDS_EdgePosition*>(node->GetPosition());
1825         pos2nodes.insert( make_pair( pos->GetUParameter(), node ));
1826       }
1827       if ( pos2nodes.size() != edgeSM->NbNodes() )
1828         RETURN_BAD_RESULT("Equal params of nodes on edge "
1829                           << smDS->ShapeToIndex( edge ) << " of face " << is2 );
1830     }
1831     if ( u2nodesMaps[0].size() != u2nodesMaps[1].size() )
1832       RETURN_BAD_RESULT("Different nb of new nodes on edges or wrong params");
1833
1834     // compare edge orientation
1835     double u1 = helper1.GetNodeU( edge1, vNode1 );
1836     double u2 = helper2.GetNodeU( edge2, vNode2 );
1837     bool isFirst1 = ( u1 < u2nodesMaps[0].begin()->first );
1838     bool isFirst2 = ( u2 < u2nodesMaps[1].begin()->first );
1839     bool reverse ( isFirst1 != isFirst2 );
1840
1841     // associate matching nodes
1842     map< double, const SMDS_MeshNode* >::iterator u_Node1, u_Node2, end1;
1843     map< double, const SMDS_MeshNode* >::reverse_iterator uR_Node2;
1844     u_Node1 = u2nodesMaps[0].begin();
1845     u_Node2 = u2nodesMaps[1].begin();
1846     uR_Node2 = u2nodesMaps[1].rbegin();
1847     end1 = u2nodesMaps[0].end();
1848     for ( ; u_Node1 != end1; ++u_Node1 ) {
1849       const SMDS_MeshNode* n1 = u_Node1->second;
1850       const SMDS_MeshNode* n2 = ( reverse ? (uR_Node2++)->second : (u_Node2++)->second );
1851       node1To2Map.insert( make_pair( n1, n2 ));
1852     }
1853
1854     // associate matching nodes on the last vertices
1855     V2 = TopExp::LastVertex( TopoDS::Edge( edge2 ));
1856     if ( !assocMap.IsBound( V2 ))
1857       RETURN_BAD_RESULT("Association not found for vertex " << meshDS2->ShapeToIndex( V2 ));
1858     V1 = TopoDS::Vertex( assocMap( V2 ));
1859     vNode1 = SMESH_Algo::VertexNode( V1, meshDS1 );
1860     vNode2 = SMESH_Algo::VertexNode( V2, meshDS2 );
1861     if ( !vNode1 ) RETURN_BAD_RESULT("No node on vertex #" << meshDS1->ShapeToIndex( V1 ));
1862     if ( !vNode2 ) RETURN_BAD_RESULT("No node on vertex #" << meshDS2->ShapeToIndex( V2 ));
1863     node1To2Map.insert( make_pair( vNode1, vNode2 ));
1864   }
1865
1866 // don't know why this condition is usually true :(
1867 //   if ( node1To2Map.size() * quadFactor < SM1->NbNodes() )
1868 //     MESSAGE("FindMatchingNodes() found too few node pairs starting from nodes ("
1869 //             << vNode1->GetID() << " - " << eNode1[0]->GetID() << ") ("
1870 //             << vNode2->GetID() << " - " << eNode2[0]->GetID() << "):"
1871 //             << node1To2Map.size() * quadFactor << " < " << SM1->NbNodes());
1872   
1873   return true;
1874 }
1875
1876 //================================================================================
1877 /*!
1878  * \brief Check if the first and last vertices of an edge are the same
1879  * \param anEdge - the edge to check
1880  * \retval bool - true if same
1881  */
1882 //================================================================================
1883
1884 bool StdMeshers_ProjectionUtils::IsClosedEdge( const TopoDS_Edge& anEdge )
1885 {
1886   return TopExp::FirstVertex( anEdge ).IsSame( TopExp::LastVertex( anEdge ));
1887 }
1888
1889 //================================================================================
1890   /*!
1891    * \brief Return any subshape of a face belonging to the outer wire
1892     * \param face - the face
1893     * \param type - type of subshape to return
1894     * \retval TopoDS_Shape - the found subshape
1895    */
1896 //================================================================================
1897
1898 TopoDS_Shape StdMeshers_ProjectionUtils::OuterShape( const TopoDS_Face& face,
1899                                                      TopAbs_ShapeEnum   type)
1900 {
1901   TopExp_Explorer exp( BRepTools::OuterWire( face ), type );
1902   if ( exp.More() )
1903     return exp.Current();
1904   return TopoDS_Shape();
1905 }
1906
1907 //================================================================================
1908   /*!
1909    * \brief Check that submesh is computed and try to compute it if is not
1910     * \param sm - submesh to compute
1911     * \param iterationNb - int used to stop infinite recursive call
1912     * \retval bool - true if computed
1913    */
1914 //================================================================================
1915
1916 bool StdMeshers_ProjectionUtils::MakeComputed(SMESH_subMesh * sm, const int iterationNb)
1917 {
1918   if ( iterationNb > 10 )
1919     RETURN_BAD_RESULT("Infinite recursive projection");
1920   if ( !sm )
1921     RETURN_BAD_RESULT("NULL submesh");
1922   if ( sm->IsMeshComputed() )
1923     return true;
1924
1925   SMESH_Mesh* mesh = sm->GetFather();
1926   SMESH_Gen* gen   = mesh->GetGen();
1927   SMESH_Algo* algo = gen->GetAlgo( *mesh, sm->GetSubShape() );
1928   if ( !algo )
1929   {
1930     if ( sm->GetSubShape().ShapeType() != TopAbs_COMPOUND )
1931       RETURN_BAD_RESULT("No algo assigned to submesh " << sm->GetId());
1932     // group
1933     bool computed = true;
1934     for ( TopoDS_Iterator grMember( sm->GetSubShape() ); grMember.More(); grMember.Next())
1935       if ( SMESH_subMesh* grSub = mesh->GetSubMesh( grMember.Value() ))
1936         if ( !MakeComputed( grSub, iterationNb + 1 ))
1937           computed = false;
1938     return computed;
1939   }
1940
1941   string algoType = algo->GetName();
1942   if ( algoType.substr(0, 11) != "Projection_")
1943     return gen->Compute( *mesh, sm->GetSubShape() );
1944
1945   // try to compute source mesh
1946
1947   const list <const SMESHDS_Hypothesis *> & hyps =
1948     algo->GetUsedHypothesis( *mesh, sm->GetSubShape() );
1949
1950   TopoDS_Shape srcShape;
1951   SMESH_Mesh* srcMesh = 0;
1952   list <const SMESHDS_Hypothesis*>::const_iterator hIt = hyps.begin();
1953   for ( ; srcShape.IsNull() && hIt != hyps.end(); ++hIt ) {
1954     string hypName = (*hIt)->GetName();
1955     if ( hypName == "ProjectionSource1D" ) {
1956       const StdMeshers_ProjectionSource1D * hyp =
1957         static_cast<const StdMeshers_ProjectionSource1D*>( *hIt );
1958       srcShape = hyp->GetSourceEdge();
1959       srcMesh = hyp->GetSourceMesh();
1960     }
1961     else if ( hypName == "ProjectionSource2D" ) {
1962       const StdMeshers_ProjectionSource2D * hyp =
1963         static_cast<const StdMeshers_ProjectionSource2D*>( *hIt );
1964       srcShape = hyp->GetSourceFace();
1965       srcMesh = hyp->GetSourceMesh();
1966     }
1967     else if ( hypName == "ProjectionSource3D" ) {
1968       const StdMeshers_ProjectionSource3D * hyp =
1969         static_cast<const StdMeshers_ProjectionSource3D*>( *hIt );
1970       srcShape = hyp->GetSource3DShape();
1971       srcMesh = hyp->GetSourceMesh();
1972     }
1973   }
1974   if ( srcShape.IsNull() ) // no projection source defined
1975     return gen->Compute( *mesh, sm->GetSubShape() );
1976
1977   if ( srcShape.IsSame( sm->GetSubShape() ))
1978     RETURN_BAD_RESULT("Projection from self");
1979     
1980   if ( !srcMesh )
1981     srcMesh = mesh;
1982
1983   if ( MakeComputed( srcMesh->GetSubMesh( srcShape ), iterationNb + 1 ))
1984     return gen->Compute( *mesh, sm->GetSubShape() );
1985
1986   return false;
1987 }
1988
1989 //================================================================================
1990 /*!
1991  * \brief Count nb of subshapes
1992  * \param shape - the shape
1993  * \param type - the type of subshapes to count
1994  * \retval int - the calculated number
1995  */
1996 //================================================================================
1997
1998 int StdMeshers_ProjectionUtils::Count(const TopoDS_Shape&    shape,
1999                                       const TopAbs_ShapeEnum type,
2000                                       const bool             ignoreSame)
2001 {
2002   if ( ignoreSame ) {
2003     TopTools_IndexedMapOfShape map;
2004     TopExp::MapShapes( shape, type, map );
2005     return map.Extent();
2006   }
2007   else {
2008     int nb = 0;
2009     for ( TopExp_Explorer exp( shape, type ); exp.More(); exp.Next() )
2010       ++nb;
2011     return nb;
2012   }
2013 }
2014
2015 //================================================================================
2016 /*!
2017  * \brief Return true if edge is a boundary of edgeContainer
2018  */
2019 //================================================================================
2020
2021 bool StdMeshers_ProjectionUtils::IsBoundaryEdge(const TopoDS_Edge&  edge,
2022                                                 const TopoDS_Shape& edgeContainer,
2023                                                 SMESH_Mesh&         mesh)
2024 {
2025   TopTools_IndexedMapOfShape facesOfEdgeContainer, facesNearEdge;
2026   TopExp::MapShapes( edgeContainer, TopAbs_FACE, facesOfEdgeContainer );
2027
2028   const TopTools_ListOfShape& EAncestors = mesh.GetAncestors(edge);
2029   TopTools_ListIteratorOfListOfShape itea(EAncestors);
2030   for(; itea.More(); itea.Next()) {
2031     if( itea.Value().ShapeType() == TopAbs_FACE &&
2032         facesOfEdgeContainer.Contains( itea.Value() ))
2033     {
2034       facesNearEdge.Add( itea.Value() );
2035       if ( facesNearEdge.Extent() > 1 )
2036         return false;
2037     }
2038   }
2039   return ( facesNearEdge.Extent() == 1 );
2040 }
2041
2042
2043 namespace {
2044
2045   SMESH_subMeshEventListener* GetSrcSubMeshListener();
2046
2047   //================================================================================
2048   /*!
2049    * \brief Listener that resets an event listener on source submesh when 
2050    * "ProjectionSource*D" hypothesis is modified
2051    */
2052   //================================================================================
2053
2054   struct HypModifWaiter: SMESH_subMeshEventListener
2055   {
2056     HypModifWaiter():SMESH_subMeshEventListener(0){} // won't be deleted by submesh
2057
2058     void ProcessEvent(const int event, const int eventType, SMESH_subMesh* subMesh,
2059                       EventListenerData*, const SMESH_Hypothesis*)
2060     {
2061       if ( event     == SMESH_subMesh::MODIF_HYP &&
2062            eventType == SMESH_subMesh::ALGO_EVENT)
2063       {
2064         // delete current source listener
2065         subMesh->DeleteEventListener( GetSrcSubMeshListener() );
2066         // let algo set a new one
2067         SMESH_Gen* gen = subMesh->GetFather()->GetGen();
2068         if ( SMESH_Algo* algo = gen->GetAlgo( *subMesh->GetFather(),
2069                                               subMesh->GetSubShape() ))
2070           algo->SetEventListener( subMesh );
2071       }
2072     }
2073   };
2074   //================================================================================
2075   /*!
2076    * \brief return static HypModifWaiter
2077    */
2078   //================================================================================
2079
2080   SMESH_subMeshEventListener* GetHypModifWaiter() {
2081     static HypModifWaiter aHypModifWaiter;
2082     return &aHypModifWaiter;
2083   }
2084   //================================================================================
2085   /*!
2086    * \brief return static listener for source shape submeshes
2087    */
2088   //================================================================================
2089
2090   SMESH_subMeshEventListener* GetSrcSubMeshListener() {
2091     static SMESH_subMeshEventListener srcListener(0); // won't be deleted by submesh
2092     return &srcListener;
2093   }
2094 }
2095
2096 //================================================================================
2097 /*!
2098  * \brief Set event listeners to submesh with projection algo
2099  * \param subMesh - submesh with projection algo
2100  * \param srcShape - source shape
2101  * \param srcMesh - source mesh
2102  */
2103 //================================================================================
2104
2105 void StdMeshers_ProjectionUtils::SetEventListener(SMESH_subMesh* subMesh,
2106                                                   TopoDS_Shape   srcShape,
2107                                                   SMESH_Mesh*    srcMesh)
2108 {
2109   // Set listener that resets an event listener on source submesh when
2110   // "ProjectionSource*D" hypothesis is modified since source shape can be changed
2111   subMesh->SetEventListener( GetHypModifWaiter(),0,subMesh);
2112
2113   // Set an event listener to submesh of the source shape
2114   if ( !srcShape.IsNull() )
2115   {
2116     if ( !srcMesh )
2117       srcMesh = subMesh->GetFather();
2118
2119     SMESH_subMesh* srcShapeSM = srcMesh->GetSubMesh( srcShape );
2120
2121     if ( srcShapeSM != subMesh ) {
2122       if ( srcShapeSM->GetSubMeshDS() &&
2123            srcShapeSM->GetSubMeshDS()->IsComplexSubmesh() )
2124       {  // source shape is a group
2125         TopExp_Explorer it(srcShapeSM->GetSubShape(), // explore the group into subshapes...
2126                            subMesh->GetSubShape().ShapeType()); // ...of target shape type
2127         for (; it.More(); it.Next())
2128         {
2129           SMESH_subMesh* srcSM = srcMesh->GetSubMesh( it.Current() );
2130           if ( srcSM != subMesh )
2131           {
2132             SMESH_subMeshEventListenerData* data =
2133               srcSM->GetEventListenerData(GetSrcSubMeshListener());
2134             if ( data )
2135               data->mySubMeshes.push_back( subMesh );
2136             else
2137               data = SMESH_subMeshEventListenerData::MakeData( subMesh );
2138             subMesh->SetEventListener ( GetSrcSubMeshListener(), data, srcSM );
2139           }
2140         }
2141       }
2142       else
2143       {
2144         subMesh->SetEventListener( GetSrcSubMeshListener(),
2145                                    SMESH_subMeshEventListenerData::MakeData( subMesh ),
2146                                    srcShapeSM );
2147       }
2148     }
2149   }
2150 }