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