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