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