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