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