Back to index

salome-smesh  6.5.0
StdMeshers_ProjectionUtils.cxx
Go to the documentation of this file.
00001 // Copyright (C) 2007-2012  CEA/DEN, EDF R&D, OPEN CASCADE
00002 //
00003 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
00004 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
00005 //
00006 // This library is free software; you can redistribute it and/or
00007 // modify it under the terms of the GNU Lesser General Public
00008 // License as published by the Free Software Foundation; either
00009 // version 2.1 of the License.
00010 //
00011 // This library is distributed in the hope that it will be useful,
00012 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00013 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00014 // Lesser General Public License for more details.
00015 //
00016 // You should have received a copy of the GNU Lesser General Public
00017 // License along with this library; if not, write to the Free Software
00018 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
00019 //
00020 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
00021 //
00022 
00023 //  SMESH SMESH : idl implementation based on 'SMESH' unit's calsses
00024 // File      : StdMeshers_ProjectionUtils.cxx
00025 // Created   : Fri Oct 27 10:24:28 2006
00026 // Author    : Edward AGAPOV (eap)
00027 //
00028 #include "StdMeshers_ProjectionUtils.hxx"
00029 
00030 #include "StdMeshers_ProjectionSource1D.hxx"
00031 #include "StdMeshers_ProjectionSource2D.hxx"
00032 #include "StdMeshers_ProjectionSource3D.hxx"
00033 
00034 #include "SMESH_Algo.hxx"
00035 #include "SMESH_Block.hxx"
00036 #include "SMESH_Gen.hxx"
00037 #include "SMESH_Hypothesis.hxx"
00038 #include "SMESH_Mesh.hxx"
00039 #include "SMESH_MesherHelper.hxx"
00040 #include "SMESH_subMesh.hxx"
00041 #include "SMESH_subMeshEventListener.hxx"
00042 #include "SMDS_EdgePosition.hxx"
00043 
00044 #include "utilities.h"
00045 
00046 #include <BRepAdaptor_Surface.hxx>
00047 #include <BRepTools.hxx>
00048 #include <BRepTools_WireExplorer.hxx>
00049 #include <BRep_Builder.hxx>
00050 #include <BRep_Tool.hxx>
00051 #include <Bnd_Box.hxx>
00052 #include <TopAbs.hxx>
00053 #include <TopExp.hxx>
00054 #include <TopExp_Explorer.hxx>
00055 #include <TopTools_Array1OfShape.hxx>
00056 #include <TopTools_DataMapIteratorOfDataMapOfShapeListOfShape.hxx>
00057 #include <TopTools_DataMapIteratorOfDataMapOfShapeShape.hxx>
00058 #include <TopTools_IndexedMapOfShape.hxx>
00059 #include <TopTools_ListIteratorOfListOfShape.hxx>
00060 #include <TopTools_ListOfShape.hxx>
00061 #include <TopTools_MapOfShape.hxx>
00062 #include <TopoDS.hxx>
00063 #include <TopoDS_Compound.hxx>
00064 #include <TopoDS_Shape.hxx>
00065 #include <gp_Pnt.hxx>
00066 #include <gp_Vec.hxx>
00067 
00068 #include <numeric>
00069 
00070 using namespace std;
00071 
00072 
00073 #define RETURN_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); return false; }
00074 #define CONT_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); continue; }
00075 #define SHOW_SHAPE(v,msg) \
00076 // { \
00077 //  if ( (v).IsNull() ) cout << msg << " NULL SHAPE" << endl; \
00078 // else if ((v).ShapeType() == TopAbs_VERTEX) {\
00079 //   gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( (v) ));\
00080 //   cout<<msg<<" "<<shapeIndex((v))<<" ( "<<p.X()<<", "<<p.Y()<<", "<<p.Z()<<" )"<<endl;} \
00081 // else {\
00082 // cout << msg << " "; TopAbs::Print((v).ShapeType(),cout) <<" "<<shapeIndex((v))<<endl;}\
00083 // }
00084 #define SHOW_LIST(msg,l) \
00085 // { \
00086 //     cout << msg << " ";\
00087 //     list< TopoDS_Edge >::const_iterator e = l.begin();\
00088 //     for ( int i = 0; e != l.end(); ++e, ++i ) {\
00089 //       cout << i << "V (" << TopExp::FirstVertex( *e, true ).TShape().operator->() << ") "\
00090 //            << i << "E (" << e->TShape().operator->() << "); "; }\
00091 //     cout << endl;\
00092 //   }
00093 
00094 #define HERE StdMeshers_ProjectionUtils
00095 
00096 namespace {
00097 
00098   static SMESHDS_Mesh* theMeshDS[2] = { 0, 0 }; // used to debug only
00099   long shapeIndex(const TopoDS_Shape& S)
00100   {
00101     if ( theMeshDS[0] && theMeshDS[1] )
00102       return max(theMeshDS[0]->ShapeToIndex(S), theMeshDS[1]->ShapeToIndex(S) );
00103     return long(S.TShape().operator->());
00104   }
00105 
00106   //================================================================================
00110   //================================================================================
00111 
00112   bool _StoreBadShape(const TopoDS_Shape& shape)
00113   {
00114 #ifdef _DEBUG_
00115     const char* type[] ={"COMPOUND","COMPSOLID","SOLID","SHELL","FACE","WIRE","EDGE","VERTEX"};
00116     BRepTools::Write( shape, SMESH_Comment("/tmp/") << type[shape.ShapeType()] << "_"
00117                       << shape.TShape().operator->() << ".brep");
00118 #endif
00119     return false;
00120   }
00121   
00122   //================================================================================
00128   //================================================================================
00129 
00130   void Reverse( list< TopoDS_Edge > & edges, const int nbEdges, const int firstEdge=0)
00131   {
00132     SHOW_LIST("BEFORE REVERSE", edges);
00133 
00134     list< TopoDS_Edge >::iterator eIt = edges.begin();
00135     std::advance( eIt, firstEdge );
00136     list< TopoDS_Edge >::iterator eBackIt = eIt;
00137     for ( int i = 0; i < nbEdges; ++i, ++eBackIt )
00138       eBackIt->Reverse(); // reverse edge
00139     // reverse list
00140     --eBackIt;
00141     while ( eIt != eBackIt )
00142     {
00143       std::swap( *eIt, *eBackIt );
00144       SHOW_LIST("# AFTER SWAP", edges)
00145         if ( (++eIt) != eBackIt )
00146           --eBackIt;
00147     }
00148     SHOW_LIST("ATFER REVERSE", edges)
00149   }
00150 
00151   //================================================================================
00158   //================================================================================
00159 
00160   bool IsPropagationPossible( SMESH_Mesh* theMesh1, SMESH_Mesh* theMesh2 )
00161   {
00162     if ( theMesh1 != theMesh2 ) {
00163       TopoDS_Shape mainShape1 = theMesh1->GetMeshDS()->ShapeToMesh();
00164       TopoDS_Shape mainShape2 = theMesh2->GetMeshDS()->ShapeToMesh();
00165       return mainShape1.IsSame( mainShape2 );
00166     }
00167     return true;
00168   }
00169 
00170   //================================================================================
00180   //================================================================================
00181 
00182   bool FixAssocByPropagation( const int             nbEdges,
00183                               list< TopoDS_Edge > & edges1,
00184                               list< TopoDS_Edge > & edges2,
00185                               SMESH_Mesh*           theMesh1,
00186                               SMESH_Mesh*           theMesh2)
00187   {
00188     if ( nbEdges == 2 && IsPropagationPossible( theMesh1, theMesh2 ) )
00189     {
00190       list< TopoDS_Edge >::iterator eIt2 = ++edges2.begin(); // 2nd edge of the 2nd face
00191       TopoDS_Edge edge2 = HERE::GetPropagationEdge( theMesh1, *eIt2, edges1.front() ).second;
00192       if ( !edge2.IsNull() ) { // propagation found for the second edge
00193         Reverse( edges2, nbEdges );
00194         return true;
00195       }
00196     }
00197     return false;
00198   }
00199 
00200   //================================================================================
00208   //================================================================================
00209 
00210   TopoDS_Shape FindGroupContaining(const TopoDS_Shape& tgtShape,
00211                                    const SMESH_Mesh*   tgtMesh1,
00212                                    const TopoDS_Shape& srcGroup)
00213   {
00214     list<SMESH_subMesh*> subMeshes = tgtMesh1->GetGroupSubMeshesContaining(tgtShape);
00215     list<SMESH_subMesh*>::iterator sm = subMeshes.begin();
00216     int type, last = TopAbs_SHAPE;
00217     StdMeshers_ProjectionUtils util;
00218     for ( ; sm != subMeshes.end(); ++sm ) {
00219       const TopoDS_Shape & group = (*sm)->GetSubShape();
00220       // check if group is similar to srcGroup
00221       for ( type = srcGroup.ShapeType(); type < last; ++type)
00222         if ( util.Count( srcGroup, (TopAbs_ShapeEnum)type, 0) !=
00223              util.Count( group,    (TopAbs_ShapeEnum)type, 0))
00224           break;
00225       if ( type == last )
00226         return group;
00227     }
00228     return TopoDS_Shape();
00229   }
00230 
00231   //================================================================================
00235   //================================================================================
00236 
00237   bool AssocGroupsByPropagation(const TopoDS_Shape&   theGroup1,
00238                                 const TopoDS_Shape&   theGroup2,
00239                                 SMESH_Mesh&           theMesh,
00240                                 HERE::TShapeShapeMap& theMap)
00241   {
00242     // If groups are on top and bottom of prism then we can associate
00243     // them using "vertical" (or "side") edges and faces of prism since
00244     // they connect corresponding vertices and edges of groups.
00245 
00246     TopTools_IndexedMapOfShape subshapes1, subshapes2;
00247     TopExp::MapShapes( theGroup1, subshapes1 );
00248     TopExp::MapShapes( theGroup2, subshapes2 );
00249     TopTools_ListIteratorOfListOfShape ancestIt;
00250 
00251     // Iterate on vertices of group1 to find corresponding vertices in group2
00252     // and associate adjacent edges and faces
00253 
00254     TopTools_MapOfShape verticShapes;
00255     TopExp_Explorer vExp1( theGroup1, TopAbs_VERTEX );
00256     for ( ; vExp1.More(); vExp1.Next() )
00257     {
00258       const TopoDS_Vertex& v1 = TopoDS::Vertex( vExp1.Current() );
00259       if ( theMap.IsBound( v1 )) continue; // already processed
00260 
00261       // Find "vertical" edge ending in v1 and whose other vertex belongs to group2
00262       TopoDS_Shape verticEdge, v2;
00263       ancestIt.Initialize( theMesh.GetAncestors( v1 ));
00264       for ( ; verticEdge.IsNull() && ancestIt.More(); ancestIt.Next() )
00265       {
00266         if ( ancestIt.Value().ShapeType() != TopAbs_EDGE ) continue;
00267         v2 = HERE::GetNextVertex( TopoDS::Edge( ancestIt.Value() ), v1 );
00268         if ( subshapes2.Contains( v2 ))
00269           verticEdge = ancestIt.Value();
00270       }
00271       if ( verticEdge.IsNull() )
00272         return false;
00273 
00274       HERE::InsertAssociation( v1, v2, theMap);
00275 
00276       // Associate edges by vertical faces sharing the found vertical edge
00277       ancestIt.Initialize( theMesh.GetAncestors( verticEdge ) );
00278       for ( ; ancestIt.More(); ancestIt.Next() )
00279       {
00280         if ( ancestIt.Value().ShapeType() != TopAbs_FACE ) continue;
00281         if ( !verticShapes.Add( ancestIt.Value() )) continue;
00282         const TopoDS_Face& face = TopoDS::Face( ancestIt.Value() );
00283 
00284         // get edges of the face
00285         TopoDS_Edge edgeGr1, edgeGr2, verticEdge2;
00286         list< TopoDS_Edge > edges;    list< int > nbEdgesInWire;
00287         SMESH_Block::GetOrderedEdges( face, v1, edges, nbEdgesInWire);
00288         if ( nbEdgesInWire.front() != 4 )
00289           return _StoreBadShape( face );
00290         list< TopoDS_Edge >::iterator edge = edges.begin();
00291         if ( verticEdge.IsSame( *edge )) {
00292           edgeGr2     = *(++edge);
00293           verticEdge2 = *(++edge);
00294           edgeGr1     = *(++edge);
00295         } else {
00296           edgeGr1     = *(edge++);
00297           verticEdge2 = *(edge++);
00298           edgeGr2     = *(edge++);
00299         }
00300 
00301         HERE::InsertAssociation( edgeGr1, edgeGr2.Reversed(), theMap);
00302       }
00303     }
00304 
00305     // Associate faces
00306     TopoDS_Iterator gr1It( theGroup1 );
00307     if ( gr1It.Value().ShapeType() == TopAbs_FACE )
00308     {
00309       // find a boundary edge of group1 to start from
00310       TopoDS_Shape bndEdge = StdMeshers_ProjectionUtils::GetBoundaryEdge( theGroup1, theMesh );
00311       if ( bndEdge.IsNull() )
00312         return false;
00313 
00314       list< TopoDS_Shape > edges(1, bndEdge);
00315       list< TopoDS_Shape >::iterator edge1 = edges.begin();
00316       for ( ; edge1 != edges.end(); ++edge1 )
00317       {
00318         // there must be one or zero not associated faces between ancestors of edge
00319         // belonging to theGroup1
00320         TopoDS_Shape face1;
00321         ancestIt.Initialize( theMesh.GetAncestors( *edge1 ) );
00322         for ( ; ancestIt.More() && face1.IsNull(); ancestIt.Next() ) {
00323           if ( ancestIt.Value().ShapeType() == TopAbs_FACE &&
00324                !theMap.IsBound( ancestIt.Value() ) &&
00325                subshapes1.Contains( ancestIt.Value() ))
00326             face1 = ancestIt.Value();
00327 
00328           // add edges of face1 to start searching for adjacent faces from
00329           for ( TopExp_Explorer e(face1, TopAbs_EDGE); e.More(); e.Next())
00330             if ( !edge1->IsSame( e.Current() ))
00331               edges.push_back( e.Current() );
00332         }
00333         if ( !face1.IsNull() ) {
00334           // find the corresponding face of theGroup2
00335           TopoDS_Shape edge2 = theMap( *edge1 );
00336           TopoDS_Shape face2;
00337           ancestIt.Initialize( theMesh.GetAncestors( edge2 ) );
00338           for ( ; ancestIt.More() && face2.IsNull(); ancestIt.Next() ) {
00339             if ( ancestIt.Value().ShapeType() == TopAbs_FACE &&
00340                  !theMap.IsBound( ancestIt.Value(), /*is2nd=*/true ) &&
00341                  subshapes2.Contains( ancestIt.Value() ))
00342               face2 = ancestIt.Value();
00343           }
00344           if ( face2.IsNull() )
00345             return false;
00346 
00347           HERE::InsertAssociation( face1, face2, theMap);
00348         }
00349       }
00350     }
00351     return true;
00352   }
00353 
00354   //================================================================================
00359   //================================================================================
00360 
00361   bool sameVertexUV( const TopoDS_Edge& edge,
00362                      const TopoDS_Face& face,
00363                      const int&         vIndex,
00364                      const gp_Pnt2d&    uv,
00365                      const double&      tol2d )
00366   {
00367     TopoDS_Vertex VV[2];
00368     TopExp::Vertices( edge, VV[0], VV[1], true);
00369     gp_Pnt2d v1UV = BRep_Tool::Parameters( VV[vIndex], face);
00370     double dist2d = v1UV.Distance( uv );
00371     return dist2d < tol2d;
00372   }
00373 
00374   //================================================================================
00378   //================================================================================
00379 
00380   TopoDS_Shape getOuterEdge( const TopoDS_Shape theShape1, SMESH_Mesh& mesh )
00381   {
00382     TopoDS_Shape edge;
00383     if ( theShape1.ShapeType() == TopAbs_COMPOUND )
00384     {
00385       TopoDS_Iterator it( theShape1 );
00386       if ( it.Value().ShapeType() == TopAbs_FACE ) // group of FACEs
00387       {
00388         // look for a boundary EDGE of a group
00389         edge = StdMeshers_ProjectionUtils::GetBoundaryEdge( theShape1, mesh );
00390         if ( !edge.IsNull() )
00391           return edge;
00392       }
00393     }
00394     edge = theShape1;
00395     TopExp_Explorer expF( theShape1, TopAbs_FACE ), expE;
00396     if ( expF.More() ) {
00397       for ( ; expF.More(); expF.Next() ) {
00398         edge.Nullify();
00399         TopoDS_Shape wire =
00400           StdMeshers_ProjectionUtils::OuterShape( TopoDS::Face( expF.Current() ), TopAbs_WIRE );
00401         for ( expE.Init( wire, TopAbs_EDGE ); edge.IsNull() && expE.More(); expE.Next() )
00402           if ( !SMESH_MesherHelper::IsClosedEdge( TopoDS::Edge( expE.Current() )))
00403             edge = expE.Current();
00404         if ( !edge.IsNull() )
00405           break;
00406       }
00407     } else if (edge.ShapeType() != TopAbs_EDGE) { // no faces
00408       edge.Nullify();
00409       for ( expE.Init( theShape1, TopAbs_EDGE ); edge.IsNull() && expE.More(); expE.Next() )
00410         if ( !SMESH_MesherHelper::IsClosedEdge( TopoDS::Edge( expE.Current() )))
00411           edge = expE.Current();
00412     }
00413     return edge;
00414   }
00415 
00416 } // namespace
00417 
00418 //=======================================================================
00429 //=======================================================================
00430 
00431 bool StdMeshers_ProjectionUtils::FindSubShapeAssociation(const TopoDS_Shape& theShape1,
00432                                                          SMESH_Mesh*         theMesh1,
00433                                                          const TopoDS_Shape& theShape2,
00434                                                          SMESH_Mesh*         theMesh2,
00435                                                          TShapeShapeMap &    theMap)
00436 {
00437   // Structure of this long function is following
00438   // 1) Group->group projection: theShape1 is a group member,
00439   //    theShape2 is a group. We find a group theShape1 is in and recall self.
00440   // 2) Accosiate same shapes with different location (partners).
00441   // 3) If vertex association is given, perform accosiation according to shape type:
00442   //       switch ( ShapeType ) {
00443   //         case TopAbs_EDGE:
00444   //         case ...:
00445   //       }
00446   // 4) else try to accosiate in different ways:
00447   //       a) accosiate shapes by propagation and other simple cases
00448   //            switch ( ShapeType ) {
00449   //            case TopAbs_EDGE:
00450   //            case ...:
00451   //            }
00452   //       b) find association of a couple of vertices and recall self.
00453   //
00454 
00455   theMeshDS[0] = theMesh1->GetMeshDS(); // debug
00456   theMeshDS[1] = theMesh2->GetMeshDS();
00457 
00458   // =================================================================================
00459   // 1) Is it the case of associating a group member -> another group? (PAL16202, 16203)
00460   // =================================================================================
00461   if ( theShape1.ShapeType() != theShape2.ShapeType() ) {
00462     TopoDS_Shape group1, group2;
00463     if ( theShape1.ShapeType() == TopAbs_COMPOUND ) {
00464       group1 = theShape1;
00465       group2 = FindGroupContaining( theShape2, theMesh2, group1 );
00466     }
00467     else if ( theShape2.ShapeType() == TopAbs_COMPOUND ) {
00468       group2 = theShape2;
00469       group1 = FindGroupContaining( theShape1, theMesh1, group2 );
00470     }
00471     if ( group1.IsNull() || group2.IsNull() )
00472       RETURN_BAD_RESULT("Different shape types");
00473     // Associate compounds
00474     return FindSubShapeAssociation(group1, theMesh1, group2, theMesh2, theMap );
00475   }
00476 
00477   // ============
00478   // 2) Is partner?
00479   // ============
00480   bool partner = theShape1.IsPartner( theShape2 );
00481   TopTools_DataMapIteratorOfDataMapOfShapeShape vvIt( theMap._map1to2 );
00482   for ( ; partner && vvIt.More(); vvIt.Next() )
00483     partner = vvIt.Key().IsPartner( vvIt.Value() );
00484 
00485   if ( partner ) // Same shape with different location
00486   {
00487     // recursively associate all sub-shapes of theShape1 and theShape2
00488     typedef list< pair< TopoDS_Shape, TopoDS_Shape > > TShapePairsList;
00489     TShapePairsList shapesQueue( 1, make_pair( theShape1, theShape2 ));
00490     TShapePairsList::iterator s1_s2 = shapesQueue.begin();
00491     for ( ; s1_s2 != shapesQueue.end(); ++s1_s2 )
00492     {
00493       InsertAssociation( s1_s2->first, s1_s2->second, theMap );
00494       TopoDS_Iterator s1It( s1_s2->first), s2It( s1_s2->second );
00495       for ( ; s1It.More(); s1It.Next(), s2It.Next() )
00496         shapesQueue.push_back( make_pair( s1It.Value(), s2It.Value() ));
00497     }
00498     return true;
00499   }
00500 
00501   if ( !theMap.IsEmpty() )
00502   {
00503     //======================================================================
00504     // 3) HAS initial vertex association
00505     //======================================================================
00506     switch ( theShape1.ShapeType() ) {
00507       // ----------------------------------------------------------------------
00508     case TopAbs_EDGE: { // TopAbs_EDGE
00509       // ----------------------------------------------------------------------
00510       if ( theMap.Extent() != 2 )
00511         RETURN_BAD_RESULT("Wrong map extent " << theMap.Extent() );
00512       TopoDS_Edge edge1 = TopoDS::Edge( theShape1 );
00513       TopoDS_Edge edge2 = TopoDS::Edge( theShape2 );
00514       if ( edge1.Orientation() >= TopAbs_INTERNAL ) edge1.Orientation( TopAbs_FORWARD );
00515       if ( edge2.Orientation() >= TopAbs_INTERNAL ) edge2.Orientation( TopAbs_FORWARD );
00516       TopoDS_Vertex VV1[2], VV2[2];
00517       TopExp::Vertices( edge1, VV1[0], VV1[1] );
00518       TopExp::Vertices( edge2, VV2[0], VV2[1] );
00519       int i1 = 0, i2 = 0;
00520       if ( theMap.IsBound( VV1[ i1 ] )) i1 = 1;
00521       if ( theMap.IsBound( VV2[ i2 ] )) i2 = 1;
00522       InsertAssociation( VV1[ i1 ], VV2[ i2 ], theMap );
00523       InsertAssociation( theShape1, theShape2, theMap );
00524       return true;
00525     }
00526       // ----------------------------------------------------------------------
00527     case TopAbs_FACE: { // TopAbs_FACE
00528       // ----------------------------------------------------------------------
00529       TopoDS_Face face1 = TopoDS::Face( theShape1 );
00530       TopoDS_Face face2 = TopoDS::Face( theShape2 );
00531       if ( face1.Orientation() >= TopAbs_INTERNAL ) face1.Orientation( TopAbs_FORWARD );
00532       if ( face2.Orientation() >= TopAbs_INTERNAL ) face2.Orientation( TopAbs_FORWARD );
00533 
00534       TopoDS_Vertex VV1[2], VV2[2];
00535       // find a not closed edge of face1 both vertices of which are associated
00536       int nbEdges = 0;
00537       TopExp_Explorer exp ( face1, TopAbs_EDGE );
00538       for ( ; VV2[ 1 ].IsNull() && exp.More(); exp.Next(), ++nbEdges ) {
00539         TopExp::Vertices( TopoDS::Edge( exp.Current() ), VV1[0], VV1[1] );
00540         if ( theMap.IsBound( VV1[0] ) ) {
00541           VV2[ 0 ] = TopoDS::Vertex( theMap( VV1[0] ));
00542           if ( theMap.IsBound( VV1[1] ) && !VV1[0].IsSame( VV1[1] ))
00543             VV2[ 1 ] = TopoDS::Vertex( theMap( VV1[1] ));
00544         }
00545       }
00546       if ( VV2[ 1 ].IsNull() ) { // 2 bound vertices not found
00547         if ( nbEdges > 1 ) {
00548           RETURN_BAD_RESULT("2 bound vertices not found" );
00549         } else {
00550           VV2[ 1 ] = VV2[ 0 ];
00551         }
00552       }
00553       list< TopoDS_Edge > edges1, edges2;
00554       int nbE = FindFaceAssociation( face1, VV1, face2, VV2, edges1, edges2 );
00555       if ( !nbE ) RETURN_BAD_RESULT("FindFaceAssociation() failed");
00556       FixAssocByPropagation( nbE, edges1, edges2, theMesh1, theMesh2 );
00557 
00558       list< TopoDS_Edge >::iterator eIt1 = edges1.begin();
00559       list< TopoDS_Edge >::iterator eIt2 = edges2.begin();
00560       for ( ; eIt1 != edges1.end(); ++eIt1, ++eIt2 )
00561       {
00562         InsertAssociation( *eIt1, *eIt2, theMap );
00563         VV1[0] = TopExp::FirstVertex( *eIt1, true );
00564         VV2[0] = TopExp::FirstVertex( *eIt2, true );
00565         InsertAssociation( VV1[0], VV2[0], theMap );
00566       }
00567       InsertAssociation( theShape1, theShape2, theMap );
00568       return true;
00569     }
00570       // ----------------------------------------------------------------------
00571     case TopAbs_SHELL: // TopAbs_SHELL, TopAbs_SOLID
00572     case TopAbs_SOLID: {
00573       // ----------------------------------------------------------------------
00574       TopoDS_Vertex VV1[2], VV2[2];
00575       // try to find a not closed edge of shape1 both vertices of which are associated
00576       TopoDS_Edge edge1;
00577       TopExp_Explorer exp ( theShape1, TopAbs_EDGE );
00578       for ( ; VV2[ 1 ].IsNull() && exp.More(); exp.Next() ) {
00579         edge1 = TopoDS::Edge( exp.Current() );
00580         if ( edge1.Orientation() >= TopAbs_INTERNAL ) edge1.Orientation( TopAbs_FORWARD );
00581         TopExp::Vertices( edge1 , VV1[0], VV1[1] );
00582         if ( theMap.IsBound( VV1[0] )) {
00583           VV2[ 0 ] = TopoDS::Vertex( theMap( VV1[0] ));
00584           if ( theMap.IsBound( VV1[1] ) && !VV1[0].IsSame( VV1[1] ))
00585             VV2[ 1 ] = TopoDS::Vertex( theMap( VV1[1] ));
00586         }
00587       }
00588       if ( VV2[ 1 ].IsNull() ) // 2 bound vertices not found
00589         RETURN_BAD_RESULT("2 bound vertices not found" );
00590       // get an edge2 of theShape2 corresponding to edge1
00591       TopoDS_Edge edge2 = GetEdgeByVertices( theMesh2, VV2[ 0 ], VV2[ 1 ]);
00592       if ( edge2.IsNull() )
00593         RETURN_BAD_RESULT("GetEdgeByVertices() failed");
00594 
00595       // build map of edge to faces if shapes are not sub-shapes of main ones
00596       bool isSubOfMain = false;
00597       if ( SMESHDS_SubMesh * sm = theMesh1->GetMeshDS()->MeshElements( theShape1 ))
00598         isSubOfMain = !sm->IsComplexSubmesh();
00599       else
00600         isSubOfMain = theMesh1->GetMeshDS()->ShapeToIndex( theShape1 );
00601       TAncestorMap e2f1, e2f2;
00602       const TAncestorMap& edgeToFace1 = isSubOfMain ? theMesh1->GetAncestorMap() : e2f1;
00603       const TAncestorMap& edgeToFace2 = isSubOfMain ? theMesh2->GetAncestorMap() : e2f2;
00604       if (!isSubOfMain) {
00605         TopExp::MapShapesAndAncestors( theShape1, TopAbs_EDGE, TopAbs_FACE, e2f1 );
00606         TopExp::MapShapesAndAncestors( theShape2, TopAbs_EDGE, TopAbs_FACE, e2f2 );
00607         if ( !edgeToFace1.Contains( edge1 ))
00608           RETURN_BAD_RESULT("edge1 does not belong to theShape1");
00609         if ( !edgeToFace2.Contains( edge2 ))
00610           RETURN_BAD_RESULT("edge2 does not belong to theShape2");
00611       }
00612       //
00613       // Look for 2 corresponing faces:
00614       //
00615       TopoDS_Shape F1, F2;
00616 
00617       // get a face sharing edge1 (F1)
00618       TopoDS_Shape FF2[2];
00619       TopTools_ListIteratorOfListOfShape ancestIt1( edgeToFace1.FindFromKey( edge1 ));
00620       for ( ; F1.IsNull() && ancestIt1.More(); ancestIt1.Next() )
00621         if ( ancestIt1.Value().ShapeType() == TopAbs_FACE )
00622           F1 = ancestIt1.Value().Oriented( TopAbs_FORWARD );
00623       if ( F1.IsNull() )
00624         RETURN_BAD_RESULT(" Face1 not found");
00625 
00626       // get 2 faces sharing edge2 (one of them is F2)
00627       TopTools_ListIteratorOfListOfShape ancestIt2( edgeToFace2.FindFromKey( edge2 ));
00628       for ( int i = 0; FF2[1].IsNull() && ancestIt2.More(); ancestIt2.Next() )
00629         if ( ancestIt2.Value().ShapeType() == TopAbs_FACE )
00630           FF2[ i++ ] = ancestIt2.Value().Oriented( TopAbs_FORWARD );
00631 
00632       // get oriented edge1 and edge2 from F1 and FF2[0]
00633       for ( exp.Init( F1, TopAbs_EDGE ); exp.More(); exp.Next() )
00634         if ( edge1.IsSame( exp.Current() )) {
00635           edge1 = TopoDS::Edge( exp.Current() );
00636           break;
00637         }
00638       for ( exp.Init( FF2[ 0 ], TopAbs_EDGE ); exp.More(); exp.Next() )
00639         if ( edge2.IsSame( exp.Current() )) {
00640           edge2 = TopoDS::Edge( exp.Current() );
00641           break;
00642         }
00643 
00644       // compare first vertices of edge1 and edge2
00645       TopExp::Vertices( edge1, VV1[0], VV1[1], true );
00646       TopExp::Vertices( edge2, VV2[0], VV2[1], true );
00647       F2 = FF2[ 0 ]; // (F2 !)
00648       if ( !VV1[ 0 ].IsSame( theMap( VV2[ 0 ], /*is2=*/true))) {
00649         edge2.Reverse();
00650         if ( FF2[ 1 ].IsNull() )
00651           F2.Reverse();
00652         else
00653           F2 = FF2[ 1 ];
00654       }
00655 
00656       TopTools_MapOfShape boundEdges;
00657 
00658       // association of face sub-shapes and neighbour faces
00659       list< pair < TopoDS_Face, TopoDS_Edge > > FE1, FE2;
00660       list< pair < TopoDS_Face, TopoDS_Edge > >::iterator fe1, fe2;
00661       FE1.push_back( make_pair( TopoDS::Face( F1 ), edge1 ));
00662       FE2.push_back( make_pair( TopoDS::Face( F2 ), edge2 ));
00663       for ( fe1 = FE1.begin(), fe2 = FE2.begin(); fe1 != FE1.end(); ++fe1, ++fe2 )
00664       {
00665         const TopoDS_Face& face1 = fe1->first;
00666         if ( theMap.IsBound( face1 ) ) continue;
00667         const TopoDS_Face& face2 = fe2->first;
00668         edge1 = fe1->second;
00669         edge2 = fe2->second;
00670         TopExp::Vertices( edge1, VV1[0], VV1[1], true );
00671         TopExp::Vertices( edge2, VV2[0], VV2[1], true );
00672         list< TopoDS_Edge > edges1, edges2;
00673         int nbE = FindFaceAssociation( face1, VV1, face2, VV2, edges1, edges2 );
00674         if ( !nbE ) RETURN_BAD_RESULT("FindFaceAssociation() failed");
00675         InsertAssociation( face1, face2, theMap ); // assoc faces
00676         MESSAGE("Assoc FACE " << theMesh1->GetMeshDS()->ShapeToIndex( face1 )<<
00677                 " to "        << theMesh2->GetMeshDS()->ShapeToIndex( face2 ));
00678         if ( nbE == 2 && (edge1.IsSame( edges1.front())) != (edge2.IsSame( edges2.front())))
00679         {
00680           Reverse( edges2, nbE );
00681         }
00682         list< TopoDS_Edge >::iterator eIt1 = edges1.begin();
00683         list< TopoDS_Edge >::iterator eIt2 = edges2.begin();
00684         for ( ; eIt1 != edges1.end(); ++eIt1, ++eIt2 )
00685         {
00686           if ( !boundEdges.Add( *eIt1 )) continue; // already associated
00687           InsertAssociation( *eIt1, *eIt2, theMap );  // assoc edges
00688           MESSAGE("Assoc edge " << theMesh1->GetMeshDS()->ShapeToIndex( *eIt1 )<<
00689                   " to "        << theMesh2->GetMeshDS()->ShapeToIndex( *eIt2 ));
00690           VV1[0] = TopExp::FirstVertex( *eIt1, true );
00691           VV2[0] = TopExp::FirstVertex( *eIt2, true );
00692           InsertAssociation( VV1[0], VV2[0], theMap ); // assoc vertices
00693           MESSAGE("Assoc vertex " << theMesh1->GetMeshDS()->ShapeToIndex( VV1[0] )<<
00694                   " to "          << theMesh2->GetMeshDS()->ShapeToIndex( VV2[0] ));
00695 
00696           // add adjacent faces to process
00697           TopoDS_Face nextFace1 = GetNextFace( edgeToFace1, *eIt1, face1 );
00698           TopoDS_Face nextFace2 = GetNextFace( edgeToFace2, *eIt2, face2 );
00699           if ( !nextFace1.IsNull() && !nextFace2.IsNull() ) {
00700             FE1.push_back( make_pair( nextFace1, *eIt1 ));
00701             FE2.push_back( make_pair( nextFace2, *eIt2 ));
00702           }
00703         }
00704       }
00705       InsertAssociation( theShape1, theShape2, theMap );
00706       return true;
00707     }
00708       // ----------------------------------------------------------------------
00709     case TopAbs_COMPOUND: { // GROUP
00710       // ----------------------------------------------------------------------
00711       // Maybe groups contain only one member
00712       TopoDS_Iterator it1( theShape1 ), it2( theShape2 );
00713       TopAbs_ShapeEnum memberType = it1.Value().ShapeType();
00714       int nbMembers = Count( theShape1, memberType, true );
00715       if ( nbMembers == 0 ) return true;
00716       if ( nbMembers == 1 ) {
00717         return FindSubShapeAssociation( it1.Value(), theMesh1, it2.Value(), theMesh2, theMap );
00718       }
00719       // Try to make shells of faces
00720       //
00721       BRep_Builder builder;
00722       TopoDS_Shell shell1, shell2;
00723       builder.MakeShell(shell1); builder.MakeShell(shell2);
00724       if ( memberType == TopAbs_FACE ) {
00725         // just add faces of groups to shells
00726         for (; it1.More(); it1.Next(), it2.Next() )
00727           builder.Add( shell1, it1.Value() ), builder.Add( shell2, it2.Value() );
00728       }
00729       else if ( memberType == TopAbs_EDGE ) {
00730         // Try to add faces sharing more than one edge of a group or
00731         // sharing all its vertices with the group
00732         TopTools_IndexedMapOfShape groupVertices[2];
00733         TopExp::MapShapes( theShape1, TopAbs_VERTEX, groupVertices[0]);
00734         TopExp::MapShapes( theShape2, TopAbs_VERTEX, groupVertices[1]);
00735         //
00736         TopTools_MapOfShape groupEdges[2], addedFaces[2];
00737         bool hasInitAssoc = (!theMap.IsEmpty()), initAssocOK = !hasInitAssoc;
00738         for (; it1.More(); it1.Next(), it2.Next() ) {
00739           groupEdges[0].Add( it1.Value() );
00740           groupEdges[1].Add( it2.Value() );
00741           if ( !initAssocOK ) {
00742             // for shell association there must be an edge with both vertices bound
00743             TopoDS_Vertex v1, v2;
00744             TopExp::Vertices( TopoDS::Edge( it1.Value().Oriented(TopAbs_FORWARD)), v1, v2 );
00745             initAssocOK = ( theMap.IsBound( v1 ) && theMap.IsBound( v2 ));
00746           }
00747         }
00748         for (int is2ndGroup = 0; initAssocOK && is2ndGroup < 2; ++is2ndGroup) {
00749           const TopoDS_Shape& group = is2ndGroup ? theShape2: theShape1;
00750           SMESH_Mesh*         mesh  = is2ndGroup ? theMesh2 : theMesh1;
00751           TopoDS_Shell&       shell = is2ndGroup ? shell2   : shell1;
00752           for ( TopoDS_Iterator it( group ); it.More(); it.Next() ) {
00753             const TopoDS_Edge& edge = TopoDS::Edge( it.Value() );
00754             TopoDS_Face face;
00755             for ( int iF = 0; iF < 2; ++iF ) { // loop on 2 faces sharing edge
00756               face = GetNextFace(mesh->GetAncestorMap(), edge, face);
00757               if ( !face.IsNull() ) {
00758                 int nbGroupEdges = 0;
00759                 for ( TopExp_Explorer f( face, TopAbs_EDGE ); f.More(); f.Next())
00760                   if ( groupEdges[ is2ndGroup ].Contains( f.Current() ))
00761                     if ( ++nbGroupEdges > 1 )
00762                       break;
00763                 bool add = (nbGroupEdges > 1 || Count( face, TopAbs_EDGE, true ) == 1 );
00764                 if ( !add ) {
00765                   add = true;
00766                   for ( TopExp_Explorer v( face, TopAbs_VERTEX ); add && v.More(); v.Next())
00767                     add = groupVertices[ is2ndGroup ].Contains( v.Current() );
00768                 }
00769                 if ( add && addedFaces[ is2ndGroup ].Add( face ))
00770                   builder.Add( shell, face );
00771               }
00772             }
00773           }
00774         }
00775       } else {
00776         RETURN_BAD_RESULT("Unexpected group type");
00777       }
00778       // Associate shells
00779       //
00780       int nbFaces1 = Count( shell1, TopAbs_FACE, 0 );
00781       int nbFaces2 = Count( shell2, TopAbs_FACE, 0 );
00782       if ( nbFaces1 != nbFaces2 )
00783         RETURN_BAD_RESULT("Different nb of faces found for shells");
00784       if ( nbFaces1 > 0 ) {
00785         bool ok = false;
00786         if ( nbFaces1 == 1 ) {
00787           TopoDS_Shape F1 = TopoDS_Iterator( shell1 ).Value();
00788           TopoDS_Shape F2 = TopoDS_Iterator( shell2 ).Value();
00789           ok = FindSubShapeAssociation( F1, theMesh1, F2, theMesh2, theMap );
00790         }
00791         else {
00792           ok = FindSubShapeAssociation(shell1, theMesh1, shell2, theMesh2, theMap );
00793         }
00794         // Check if all members are mapped 
00795         if ( ok ) {
00796           TopTools_MapOfShape boundMembers[2];
00797           TopoDS_Iterator mIt;
00798           for ( mIt.Initialize( theShape1 ); mIt.More(); mIt.Next())
00799             if ( theMap.IsBound( mIt.Value() )) {
00800               boundMembers[0].Add( mIt.Value() );
00801               boundMembers[1].Add( theMap( mIt.Value() ));
00802             }
00803           if ( boundMembers[0].Extent() != nbMembers ) {
00804             // make compounds of not bound members
00805             TopoDS_Compound comp[2];
00806             for ( int is2ndGroup = 0; is2ndGroup < 2; ++is2ndGroup ) {
00807               builder.MakeCompound( comp[is2ndGroup] );
00808               for ( mIt.Initialize( is2ndGroup ? theShape2:theShape1 ); mIt.More(); mIt.Next())
00809                 if ( ! boundMembers[ is2ndGroup ].Contains( mIt.Value() ))
00810                   builder.Add( comp[ is2ndGroup ], mIt.Value() );
00811             }
00812             // check if theMap contains initial association for the comp's
00813             bool hasInitialAssoc = false;
00814             if ( memberType == TopAbs_EDGE ) {
00815               for ( TopExp_Explorer v( comp[0], TopAbs_VERTEX ); v.More(); v.Next())
00816                 if ( theMap.IsBound( v.Current() )) {
00817                   hasInitialAssoc = true;
00818                   break;
00819                 }
00820             }
00821             if ( hasInitialAssoc == bool( !theMap.IsEmpty() ))
00822               ok = FindSubShapeAssociation( comp[0], theMesh1, comp[1], theMesh2, theMap );
00823             else {
00824               TShapeShapeMap tmpMap;
00825               ok = FindSubShapeAssociation( comp[0], theMesh1, comp[1], theMesh2, tmpMap );
00826               if ( ok ) {
00827                 TopTools_DataMapIteratorOfDataMapOfShapeShape mapIt( tmpMap._map1to2 );
00828                 for ( ; mapIt.More(); mapIt.Next() )
00829                   theMap.Bind( mapIt.Key(), mapIt.Value());
00830               }
00831             }
00832           }
00833         }
00834         return ok;
00835       }
00836       // Each edge of an edge group is shared by own faces
00837       // ------------------------------------------------------------------
00838       //
00839       // map vertices to edges sharing them, avoid doubling edges in lists
00840       TopTools_DataMapOfShapeListOfShape v2e[2];
00841       for (int isFirst = 0; isFirst < 2; ++isFirst ) {
00842         const TopoDS_Shape& group = isFirst ? theShape1 : theShape2;
00843         TopTools_DataMapOfShapeListOfShape& veMap = v2e[ isFirst ? 0 : 1 ];
00844         TopTools_MapOfShape addedEdges;
00845         for ( TopExp_Explorer e( group, TopAbs_EDGE ); e.More(); e.Next() ) {
00846           const TopoDS_Shape& edge = e.Current();
00847           if ( addedEdges.Add( edge )) {
00848             for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next()) {
00849               const TopoDS_Shape& vertex = v.Current();
00850               if ( !veMap.IsBound( vertex )) {
00851                 TopTools_ListOfShape l;
00852                 veMap.Bind( vertex, l );
00853               }
00854               veMap( vertex ).Append( edge );
00855             }
00856           }
00857         }   
00858       }
00859       while ( !v2e[0].IsEmpty() )
00860       {
00861         // find a bound vertex
00862         TopoDS_Vertex V[2];
00863         TopTools_DataMapIteratorOfDataMapOfShapeListOfShape v2eIt( v2e[0] );
00864         for ( ; v2eIt.More(); v2eIt.Next())
00865           if ( theMap.IsBound( v2eIt.Key() )) {
00866             V[0] = TopoDS::Vertex( v2eIt.Key() );
00867             V[1] = TopoDS::Vertex( theMap( V[0] ));
00868             break;
00869           }
00870         if ( V[0].IsNull() )
00871           RETURN_BAD_RESULT("No more bound vertices");
00872 
00873         while ( !V[0].IsNull() && v2e[0].IsBound( V[0] )) {
00874           TopTools_ListOfShape& edges0 = v2e[0]( V[0] );
00875           TopTools_ListOfShape& edges1 = v2e[1]( V[1] );
00876           int nbE0 = edges0.Extent(), nbE1 = edges1.Extent();
00877           if ( nbE0 != nbE1 )
00878             RETURN_BAD_RESULT("Different nb of edges: "<< nbE0 << " != " << nbE1);
00879 
00880           if ( nbE0 == 1 )
00881           {
00882             TopoDS_Edge e0 = TopoDS::Edge( edges0.First() );
00883             TopoDS_Edge e1 = TopoDS::Edge( edges1.First() );
00884             v2e[0].UnBind( V[0] );
00885             v2e[1].UnBind( V[1] );
00886             InsertAssociation( e0, e1, theMap );
00887             MESSAGE("Assoc edge " << theMesh1->GetMeshDS()->ShapeToIndex( e0 )<<
00888                     " to "        << theMesh2->GetMeshDS()->ShapeToIndex( e1 ));
00889             V[0] = GetNextVertex( e0, V[0] );
00890             V[1] = GetNextVertex( e1, V[1] );
00891             if ( !V[0].IsNull() ) {
00892               InsertAssociation( V[0], V[1], theMap );
00893               MESSAGE("Assoc vertex " << theMesh1->GetMeshDS()->ShapeToIndex( V[0] )<<
00894                       " to "          << theMesh2->GetMeshDS()->ShapeToIndex( V[1] ));
00895             }
00896           }
00897           else if ( nbE0 == 2 )
00898           {
00899             // one of edges must have both ends bound
00900             TopoDS_Vertex v0e0 = GetNextVertex( TopoDS::Edge( edges0.First() ), V[0] );
00901             TopoDS_Vertex v1e0 = GetNextVertex( TopoDS::Edge( edges0.Last() ),  V[0] );
00902             TopoDS_Vertex v0e1 = GetNextVertex( TopoDS::Edge( edges1.First() ), V[1] );
00903             TopoDS_Vertex v1e1 = GetNextVertex( TopoDS::Edge( edges1.Last() ),  V[1] );
00904             TopoDS_Shape e0b, e1b, e0n, e1n, v1b; // bound and not-bound
00905             TopoDS_Vertex v0n, v1n;
00906             if ( theMap.IsBound( v0e0 )) {
00907               v0n = v1e0; e0b = edges0.First(); e0n = edges0.Last(); v1b = theMap( v0e0 );
00908             } else if ( theMap.IsBound( v1e0 )) {
00909               v0n = v0e0; e0n = edges0.First(); e0b = edges0.Last(); v1b = theMap( v1e0 );
00910             } else {
00911               RETURN_BAD_RESULT("None of vertices bound");
00912             }
00913             if ( v1b.IsSame( v1e1 )) {
00914               v1n = v0e1; e1n = edges1.First(); e1b = edges1.Last();
00915             } else {
00916               v1n = v1e1; e1b = edges1.First(); e1n = edges1.Last();
00917             }
00918             InsertAssociation( e0b, e1b, theMap );
00919             InsertAssociation( e0n, e1n, theMap );
00920             InsertAssociation( v0n, v1n, theMap );
00921             MESSAGE("Assoc edge " << theMesh1->GetMeshDS()->ShapeToIndex( e0b )<<
00922                     " to "        << theMesh2->GetMeshDS()->ShapeToIndex( e1b ));
00923             MESSAGE("Assoc edge " << theMesh1->GetMeshDS()->ShapeToIndex( e0n )<<
00924                     " to "        << theMesh2->GetMeshDS()->ShapeToIndex( e1n ));
00925             MESSAGE("Assoc vertex " << theMesh1->GetMeshDS()->ShapeToIndex( v0n )<<
00926                     " to "          << theMesh2->GetMeshDS()->ShapeToIndex( v1n ));
00927             v2e[0].UnBind( V[0] );
00928             v2e[1].UnBind( V[1] );
00929             V[0] = v0n;
00930             V[1] = v1n;
00931           }
00932           else {
00933             RETURN_BAD_RESULT("Not implemented");
00934           }
00935         }
00936       } //while ( !v2e[0].IsEmpty() )
00937       return true;
00938     }
00939 
00940     default:
00941       RETURN_BAD_RESULT("Unexpected shape type");
00942 
00943     } // end switch by shape type
00944   } // end case of available initial vertex association
00945 
00946   //======================================================================
00947   // 4) NO INITIAL VERTEX ASSOCIATION
00948   //======================================================================
00949 
00950   switch ( theShape1.ShapeType() ) {
00951 
00952   case TopAbs_EDGE: {
00953     // ----------------------------------------------------------------------
00954     TopoDS_Edge edge1 = TopoDS::Edge( theShape1 );
00955     TopoDS_Edge edge2 = TopoDS::Edge( theShape2 );
00956     if ( IsPropagationPossible( theMesh1, theMesh2 ))
00957     {
00958       TopoDS_Edge prpEdge = GetPropagationEdge( theMesh1, edge2, edge1 ).second;
00959       if ( !prpEdge.IsNull() )
00960       {
00961         TopoDS_Vertex VV1[2], VV2[2];
00962         TopExp::Vertices( edge1,   VV1[0], VV1[1], true );
00963         TopExp::Vertices( prpEdge, VV2[0], VV2[1], true );
00964         InsertAssociation( VV1[ 0 ], VV2[ 0 ], theMap );
00965         InsertAssociation( VV1[ 1 ], VV2[ 1 ], theMap );
00966         if ( VV1[0].IsSame( VV1[1] ) || // one of edges is closed
00967              VV2[0].IsSame( VV2[1] ) )
00968         {
00969           InsertAssociation( edge1, prpEdge, theMap ); // insert with a proper orientation
00970         }
00971         InsertAssociation( theShape1, theShape2, theMap );
00972         return true; // done
00973       }
00974     }
00975     if ( SMESH_MesherHelper::IsClosedEdge( edge1 ) &&
00976          SMESH_MesherHelper::IsClosedEdge( edge2 ))
00977     {
00978       // TODO: find out a proper orientation (is it possible?)
00979       InsertAssociation( edge1, edge2, theMap ); // insert with a proper orientation
00980       InsertAssociation( TopExp::FirstVertex(edge1), TopExp::FirstVertex(edge2),
00981                          theMap );
00982       InsertAssociation( theShape1, theShape2, theMap );
00983       return true; // done
00984     }
00985     break; // try by vertex closeness
00986   }
00987 
00988   case TopAbs_FACE: {
00989     // ----------------------------------------------------------------------
00990     if ( IsPropagationPossible( theMesh1, theMesh2 )) // try by propagation in one mesh
00991     {
00992       TopoDS_Face face1 = TopoDS::Face(theShape1);
00993       TopoDS_Face face2 = TopoDS::Face(theShape2);
00994       if ( face1.Orientation() >= TopAbs_INTERNAL ) face1.Orientation( TopAbs_FORWARD );
00995       if ( face2.Orientation() >= TopAbs_INTERNAL ) face2.Orientation( TopAbs_FORWARD );
00996       TopoDS_Edge edge1, edge2;
00997       // get outer edge of theShape1
00998       edge1 = TopoDS::Edge( OuterShape( face1, TopAbs_EDGE ));
00999       // find out if any edge of face2 is a propagation edge of outer edge1
01000       map<int,TopoDS_Edge> propag_edges; // use map to find the closest propagation edge
01001       for ( TopExp_Explorer exp( face2, TopAbs_EDGE ); exp.More(); exp.Next() ) {
01002         edge2 = TopoDS::Edge( exp.Current() );
01003         pair<int,TopoDS_Edge> step_edge = GetPropagationEdge( theMesh1, edge2, edge1 );
01004         if ( !step_edge.second.IsNull() ) { // propagation found
01005           propag_edges.insert( step_edge );
01006           if ( step_edge.first == 1 ) break; // most close found
01007         }
01008       }
01009       if ( !propag_edges.empty() ) // propagation found
01010       {
01011         edge2 = propag_edges.begin()->second;
01012         TopoDS_Vertex VV1[2], VV2[2];
01013         TopExp::Vertices( edge1, VV1[0], VV1[1], true );
01014         TopExp::Vertices( edge2, VV2[0], VV2[1], true );
01015         list< TopoDS_Edge > edges1, edges2;
01016         int nbE = FindFaceAssociation( face1, VV1, face2, VV2, edges1, edges2 );
01017         if ( !nbE ) RETURN_BAD_RESULT("FindFaceAssociation() failed");
01018         // take care of proper association of propagated edges
01019         bool same1 = edge1.IsSame( edges1.front() );
01020         bool same2 = edge2.IsSame( edges2.front() );
01021         if ( same1 != same2 )
01022         {
01023           Reverse(edges2, nbE);
01024           if ( nbE != 2 ) // 2 degen edges of 4 (issue 0021144)
01025             edges2.splice( edges2.end(), edges2, edges2.begin());
01026         }
01027         // store association
01028         list< TopoDS_Edge >::iterator eIt1 = edges1.begin();
01029         list< TopoDS_Edge >::iterator eIt2 = edges2.begin();
01030         for ( ; eIt1 != edges1.end(); ++eIt1, ++eIt2 )
01031         {
01032           InsertAssociation( *eIt1, *eIt2, theMap );
01033           VV1[0] = TopExp::FirstVertex( *eIt1, true );
01034           VV2[0] = TopExp::FirstVertex( *eIt2, true );
01035           InsertAssociation( VV1[0], VV2[0], theMap );
01036         }
01037         InsertAssociation( theShape1, theShape2, theMap );
01038         return true;
01039       }
01040     }
01041     break; // try by vertex closeness
01042   }
01043   case TopAbs_COMPOUND: {
01044     // ----------------------------------------------------------------------
01045     if ( IsPropagationPossible( theMesh1, theMesh2 )) {
01046 
01047       // try to accosiate all using propagation
01048       if ( AssocGroupsByPropagation( theShape1, theShape2, *theMesh1, theMap ))
01049         return true;
01050 
01051       // find a boundary edge of theShape1
01052       TopoDS_Edge E = GetBoundaryEdge( theShape1, *theMesh1 );
01053       if ( E.IsNull() )
01054         break; // try by vertex closeness
01055 
01056       // find association for vertices of edge E
01057       TopoDS_Vertex VV1[2], VV2[2];
01058       for(TopExp_Explorer eexp(E, TopAbs_VERTEX); eexp.More(); eexp.Next()) {
01059         TopoDS_Vertex V1 = TopoDS::Vertex( eexp.Current() );
01060         // look for an edge ending in E whose one vertex is in theShape1
01061         // and the other, in theShape2
01062         const TopTools_ListOfShape& Ancestors = theMesh1->GetAncestors(V1);
01063         TopTools_ListIteratorOfListOfShape ita(Ancestors);
01064         for(; ita.More(); ita.Next()) {
01065           if( ita.Value().ShapeType() != TopAbs_EDGE ) continue;
01066           TopoDS_Edge edge = TopoDS::Edge(ita.Value());
01067           bool FromShape1 = false;
01068           for(TopExp_Explorer expe(theShape1, TopAbs_EDGE); expe.More(); expe.Next() ) {
01069             if(edge.IsSame(expe.Current())) {
01070               FromShape1 = true;
01071               break;
01072             }
01073           }
01074           if(!FromShape1) {
01075             // is it an edge between theShape1 and theShape2?
01076             TopExp_Explorer expv(edge, TopAbs_VERTEX);
01077             TopoDS_Vertex V2 = TopoDS::Vertex( expv.Current() );
01078             if(V2.IsSame(V1)) {
01079               expv.Next();
01080               V2 = TopoDS::Vertex( expv.Current() );
01081             }
01082             bool FromShape2 = false;
01083             for ( expv.Init( theShape2, TopAbs_VERTEX ); expv.More(); expv.Next()) {
01084               if ( V2.IsSame( expv.Current() )) {
01085                 FromShape2 = true;
01086                 break;
01087               }
01088             }
01089             if ( FromShape2 ) {
01090               if ( VV1[0].IsNull() )
01091                 VV1[0] = V1, VV2[0] = V2;
01092               else
01093                 VV1[1] = V1, VV2[1] = V2;
01094               break; // from loop on ancestors of V1
01095             }
01096           }
01097         }
01098       }
01099       if ( !VV1[1].IsNull() ) {
01100         InsertAssociation( VV1[0], VV2[0], theMap );
01101         InsertAssociation( VV1[1], VV2[1], theMap );
01102         return FindSubShapeAssociation( theShape1, theMesh1, theShape2, theMesh2, theMap);
01103       }
01104     }
01105     break; // try by vertex closeness
01106   }
01107   default:;
01108   }
01109 
01110   // 4.b) Find association by closeness of vertices
01111   // ----------------------------------------------
01112 
01113   TopTools_IndexedMapOfShape vMap1, vMap2;
01114   TopExp::MapShapes( theShape1, TopAbs_VERTEX, vMap1 );
01115   TopExp::MapShapes( theShape2, TopAbs_VERTEX, vMap2 );
01116   TopoDS_Vertex VV1[2], VV2[2];
01117 
01118   if ( vMap1.Extent() != vMap2.Extent() )
01119     RETURN_BAD_RESULT("Different nb of vertices");
01120 
01121   if ( vMap1.Extent() == 1 ) {
01122     InsertAssociation( vMap1(1), vMap2(1), theMap );
01123     if ( theShape1.ShapeType() == TopAbs_EDGE ) {
01124       InsertAssociation( theShape1, theShape2, theMap );
01125       return true;
01126     }
01127     return FindSubShapeAssociation( theShape1, theMesh1, theShape2, theMesh2, theMap);
01128   }
01129 
01130   // Try to associate by common vertices of an edge
01131   for ( int i = 1; i <= vMap1.Extent(); ++i )
01132   {
01133     const TopoDS_Shape& v1 = vMap1(i);
01134     if ( vMap2.Contains( v1 ))
01135     {
01136       // find an egde sharing v1 and sharing at the same time another common vertex
01137       PShapeIteratorPtr edgeIt = SMESH_MesherHelper::GetAncestors( v1, *theMesh1, TopAbs_EDGE);
01138       bool edgeFound = false;
01139       while ( edgeIt->more() && !edgeFound )
01140       {
01141         TopoDS_Edge edge = TopoDS::Edge( edgeIt->next()->Oriented(TopAbs_FORWARD));
01142         TopExp::Vertices(edge, VV1[0], VV1[1]);
01143         if ( !VV1[0].IsSame( VV1[1] ))
01144           edgeFound = ( vMap2.Contains( VV1[ v1.IsSame(VV1[0]) ? 1:0]));
01145       }
01146       if ( edgeFound )
01147       {
01148         InsertAssociation( VV1[0], VV1[0], theMap );
01149         InsertAssociation( VV1[1], VV1[1], theMap );
01150         if (FindSubShapeAssociation( theShape1, theMesh1, theShape2, theMesh2, theMap ))
01151           return true;
01152       }
01153     }
01154   }
01155 
01156   // Find transformation to make the shapes be of similar size at same location
01157 
01158   Bnd_Box box[2];
01159   for ( int i = 1; i <= vMap1.Extent(); ++i ) {
01160     box[ 0 ].Add( BRep_Tool::Pnt ( TopoDS::Vertex( vMap1( i ))));
01161     box[ 1 ].Add( BRep_Tool::Pnt ( TopoDS::Vertex( vMap2( i ))));
01162   }
01163 
01164   gp_Pnt gc[2]; // box center
01165   double x0,y0,z0, x1,y1,z1;
01166   box[0].Get( x0,y0,z0, x1,y1,z1 );
01167   gc[0] = 0.5 * ( gp_XYZ( x0,y0,z0 ) + gp_XYZ( x1,y1,z1 ));
01168   box[1].Get( x0,y0,z0, x1,y1,z1 );
01169   gc[1] = 0.5 * ( gp_XYZ( x0,y0,z0 ) + gp_XYZ( x1,y1,z1 ));
01170 
01171   // 1 -> 2
01172   gp_Vec vec01( gc[0], gc[1] );
01173   double scale = sqrt( box[1].SquareExtent() / box[0].SquareExtent() );
01174 
01175   // Find 2 closest vertices
01176 
01177   // get 2 linked vertices of shape 1 not belonging to an inner wire of a face
01178   TopoDS_Shape edge = getOuterEdge( theShape1, *theMesh1 );
01179   if ( edge.IsNull() || edge.ShapeType() != TopAbs_EDGE )
01180     RETURN_BAD_RESULT("Edge not found");
01181 
01182   TopExp::Vertices( TopoDS::Edge( edge.Oriented(TopAbs_FORWARD)), VV1[0], VV1[1]);
01183   if ( VV1[0].IsSame( VV1[1] ))
01184     RETURN_BAD_RESULT("Only closed edges");
01185 
01186   // find vertices closest to 2 linked vertices of shape 1
01187   for ( int i1 = 0; i1 < 2; ++i1 )
01188   {
01189     double dist2 = DBL_MAX;
01190     gp_Pnt p1 = BRep_Tool::Pnt( VV1[ i1 ]);
01191     p1.Translate( vec01 );
01192     p1.Scale( gc[1], scale );
01193     for ( int i2 = 1; i2 <= vMap2.Extent(); ++i2 )
01194     {
01195       TopoDS_Vertex V2 = TopoDS::Vertex( vMap2( i2 ));
01196       gp_Pnt p2 = BRep_Tool::Pnt ( V2 );
01197       double d2 = p1.SquareDistance( p2 );
01198       if ( d2 < dist2 && !V2.IsSame( VV2[ 0 ])) {
01199         VV2[ i1 ] = V2; dist2 = d2;
01200       }
01201     }
01202   }
01203 
01204   InsertAssociation( VV1[ 0 ], VV2[ 0 ], theMap );
01205   InsertAssociation( VV1[ 1 ], VV2[ 1 ], theMap );
01206   MESSAGE("Initial assoc VERT " << theMesh1->GetMeshDS()->ShapeToIndex( VV1[ 0 ] )<<
01207           " to "                << theMesh2->GetMeshDS()->ShapeToIndex( VV2[ 0 ] )<<
01208           "\nand         VERT " << theMesh1->GetMeshDS()->ShapeToIndex( VV1[ 1 ] )<<
01209           " to "                << theMesh2->GetMeshDS()->ShapeToIndex( VV2[ 1 ] ));
01210   if ( theShape1.ShapeType() == TopAbs_EDGE ) {
01211     InsertAssociation( theShape1, theShape2, theMap );
01212     return true;
01213   }
01214 
01215   return FindSubShapeAssociation( theShape1, theMesh1, theShape2, theMesh2, theMap );
01216 }
01217 
01218 //================================================================================
01229 //================================================================================
01230 
01231 int StdMeshers_ProjectionUtils::FindFaceAssociation(const TopoDS_Face&    face1,
01232                                                     TopoDS_Vertex         VV1[2],
01233                                                     const TopoDS_Face&    face2,
01234                                                     TopoDS_Vertex         VV2[2],
01235                                                     list< TopoDS_Edge > & edges1,
01236                                                     list< TopoDS_Edge > & edges2)
01237 {
01238   bool OK = false;
01239   list< int > nbEInW1, nbEInW2;
01240   int i_ok_wire_algo = -1;
01241   for ( int outer_wire_algo = 0; outer_wire_algo < 2 && !OK; ++outer_wire_algo )
01242   {
01243     edges1.clear();
01244     edges2.clear();
01245 
01246     if ( SMESH_Block::GetOrderedEdges( face1, VV1[0], edges1, nbEInW1, outer_wire_algo) !=
01247          SMESH_Block::GetOrderedEdges( face2, VV2[0], edges2, nbEInW2, outer_wire_algo) )
01248       CONT_BAD_RESULT("Different number of wires in faces ");
01249 
01250     if ( nbEInW1 != nbEInW2 )
01251       CONT_BAD_RESULT("Different number of edges in faces: " <<
01252                       nbEInW1.front() << " != " << nbEInW2.front());
01253 
01254     i_ok_wire_algo = outer_wire_algo;
01255 
01256     // Define if we need to reverse one of wires to make edges in lists match each other
01257 
01258     bool reverse = false;
01259 
01260     list< TopoDS_Edge >::iterator edgeIt;
01261     if ( !VV1[1].IsSame( TopExp::LastVertex( edges1.front(), true ))) {
01262       reverse = true;
01263       edgeIt = --edges1.end();
01264       // check if the second vertex belongs to the first or last edge in the wire
01265       if ( !VV1[1].IsSame( TopExp::FirstVertex( *edgeIt, true ))) {
01266         bool KO = true; // belongs to none
01267         if ( nbEInW1.size() > 1 ) { // several wires
01268           edgeIt = edges1.begin();
01269           std::advance( edgeIt, nbEInW1.front()-1 );
01270           KO = !VV1[1].IsSame( TopExp::FirstVertex( *edgeIt, true ));
01271         }
01272         if ( KO )
01273           CONT_BAD_RESULT("GetOrderedEdges() failed");
01274       }
01275     }
01276     edgeIt = --edges2.end();
01277     if ( !VV2[1].IsSame( TopExp::LastVertex( edges2.front(), true ))) {
01278       reverse = !reverse;
01279       // check if the second vertex belongs to the first or last edge in the wire
01280       if ( !VV2[1].IsSame( TopExp::FirstVertex( *edgeIt, true ))) {
01281         bool KO = true; // belongs to none
01282         if ( nbEInW2.size() > 1 ) { // several wires
01283           edgeIt = edges2.begin();
01284           std::advance( edgeIt, nbEInW2.front()-1 );
01285           KO = !VV2[1].IsSame( TopExp::FirstVertex( *edgeIt, true ));
01286         }
01287         if ( KO )
01288           CONT_BAD_RESULT("GetOrderedEdges() failed");
01289       }
01290     }
01291     if ( reverse )
01292     {
01293       Reverse( edges2 , nbEInW2.front());
01294       if (( VV1[1].IsSame( TopExp::LastVertex( edges1.front(), true ))) !=
01295           ( VV2[1].IsSame( TopExp::LastVertex( edges2.front(), true ))))
01296         CONT_BAD_RESULT("GetOrderedEdges() failed");
01297     }
01298     OK = true;
01299 
01300   } // loop algos getting an outer wire
01301   
01302   // Try to orient all (if !OK) or only internal wires (issue 0020996) by UV similarity
01303   if (( !OK || nbEInW1.size() > 1 ) && i_ok_wire_algo > -1 )
01304   {
01305     // Check that Vec(VV1[0],VV1[1]) in 2D on face1 is the same
01306     // as Vec(VV2[0],VV2[1]) on face2
01307     double vTol = BRep_Tool::Tolerance( VV1[0] );
01308     BRepAdaptor_Surface surface1( face1, false );
01309     double vTolUV =
01310       surface1.UResolution( vTol ) + surface1.VResolution( vTol ); // let's be tolerant
01311     gp_Pnt2d v0f1UV = BRep_Tool::Parameters( VV1[0], face1 );
01312     gp_Pnt2d v0f2UV = BRep_Tool::Parameters( VV2[0], face2 );
01313     gp_Pnt2d v1f1UV = BRep_Tool::Parameters( VV1[1], face1 );
01314     gp_Pnt2d v1f2UV = BRep_Tool::Parameters( VV2[1], face2 );
01315     gp_Vec2d v01f1Vec( v0f1UV, v1f1UV );
01316     gp_Vec2d v01f2Vec( v0f2UV, v1f2UV );
01317     if ( Abs( v01f1Vec.X()-v01f2Vec.X()) < vTolUV &&
01318          Abs( v01f1Vec.Y()-v01f2Vec.Y()) < vTolUV )
01319     {
01320       if ( !OK /*i_ok_wire_algo != 1*/ )
01321       {
01322         edges1.clear();
01323         edges2.clear();
01324         SMESH_Block::GetOrderedEdges( face1, VV1[0], edges1, nbEInW1, i_ok_wire_algo);
01325         SMESH_Block::GetOrderedEdges( face2, VV2[0], edges2, nbEInW2, i_ok_wire_algo);
01326       }
01327       gp_XY dUV = v0f2UV.XY() - v0f1UV.XY(); // UV shift between 2 faces
01328       // skip edges of the outer wire (if the outer wire is OK)
01329       list< int >::iterator nbEInW = nbEInW1.begin();
01330       list< TopoDS_Edge >::iterator edge1Beg = edges1.begin(), edge2Beg = edges2.begin();
01331       if ( OK )
01332       {
01333         for ( int i = 0; i < *nbEInW; ++i )
01334           ++edge1Beg, ++edge2Beg;
01335         ++nbEInW;
01336       }
01337       for ( ; nbEInW != nbEInW1.end(); ++nbEInW ) // loop on wires
01338       {
01339         // reach an end of edges of a current wire
01340         list< TopoDS_Edge >::iterator edge1End = edge1Beg, edge2End = edge2Beg;
01341         for ( int i = 0; i < *nbEInW; ++i )
01342           ++edge1End, ++edge2End;
01343         // rotate edges2 untill coincident with edges1 in 2D
01344         v0f1UV = BRep_Tool::Parameters( TopExp::FirstVertex(*edge1Beg,true), face1 );
01345         v1f1UV = BRep_Tool::Parameters( TopExp::LastVertex (*edge1Beg,true), face1 );
01346         v0f1UV.ChangeCoord() += dUV;
01347         v1f1UV.ChangeCoord() += dUV;
01348         int i = *nbEInW;
01349         while ( --i > 0 && !sameVertexUV( *edge2Beg, face2, 0, v0f1UV, vTolUV ))
01350           edges2.splice( edge2End, edges2, edge2Beg++ ); // move edge2Beg to place before edge2End
01351         if ( sameVertexUV( *edge2Beg, face2, 0, v0f1UV, vTolUV ))
01352         {
01353           if ( nbEInW == nbEInW1.begin() )
01354             OK = true; // OK is for the first wire
01355           // reverse edges2 if needed
01356           if ( !sameVertexUV( *edge2Beg, face2, 1, v1f1UV, vTolUV ))
01357           {
01358             Reverse( edges2 , *nbEInW, distance( edges2.begin(),edge2Beg ));
01359             // set correct edge2End
01360             edge2End = edges2.begin();
01361             std::advance( edge2End, std::accumulate( nbEInW1.begin(), nbEInW, *nbEInW));
01362           }
01363         }
01364         // prepare to the next wire loop
01365         edge1Beg = edge1End, edge2Beg = edge2End;
01366       }
01367     }
01368   }
01369 
01370   return OK ? nbEInW1.front() : 0;
01371 }
01372 
01373 //=======================================================================
01374 //function : InitVertexAssociation
01375 //purpose  : 
01376 //=======================================================================
01377 
01378 void StdMeshers_ProjectionUtils::InitVertexAssociation( const SMESH_Hypothesis* theHyp,
01379                                                         TShapeShapeMap &        theAssociationMap)
01380 {
01381   string hypName = theHyp->GetName();
01382   if ( hypName == "ProjectionSource1D" ) {
01383     const StdMeshers_ProjectionSource1D * hyp =
01384       static_cast<const StdMeshers_ProjectionSource1D*>( theHyp );
01385     if ( hyp->HasVertexAssociation() )
01386       InsertAssociation( hyp->GetTargetVertex(),hyp->GetSourceVertex(),theAssociationMap );
01387   }
01388   else if ( hypName == "ProjectionSource2D" ) {
01389     const StdMeshers_ProjectionSource2D * hyp =
01390       static_cast<const StdMeshers_ProjectionSource2D*>( theHyp );
01391     if ( hyp->HasVertexAssociation() ) {
01392       InsertAssociation( hyp->GetTargetVertex(1),hyp->GetSourceVertex(1),theAssociationMap);
01393       InsertAssociation( hyp->GetTargetVertex(2),hyp->GetSourceVertex(2),theAssociationMap);
01394     }
01395   }
01396   else if ( hypName == "ProjectionSource3D" ) {
01397     const StdMeshers_ProjectionSource3D * hyp =
01398       static_cast<const StdMeshers_ProjectionSource3D*>( theHyp );
01399     if ( hyp->HasVertexAssociation() ) {
01400       InsertAssociation( hyp->GetTargetVertex(1),hyp->GetSourceVertex(1),theAssociationMap);
01401       InsertAssociation( hyp->GetTargetVertex(2),hyp->GetSourceVertex(2),theAssociationMap);
01402     }
01403   }
01404 }
01405 
01406 //=======================================================================
01414 //=======================================================================
01415 
01416 bool StdMeshers_ProjectionUtils::InsertAssociation( const TopoDS_Shape& theShape1, // tgt
01417                                                     const TopoDS_Shape& theShape2, // src
01418                                                     TShapeShapeMap &    theAssociationMap)
01419 {
01420   if ( !theShape1.IsNull() && !theShape2.IsNull() ) {
01421     SHOW_SHAPE(theShape1,"Assoc ");
01422     SHOW_SHAPE(theShape2," to ");
01423     bool isNew = ( theAssociationMap.Bind( theShape1, theShape2 ));
01424     return isNew;
01425   }
01426   else {
01427     throw SALOME_Exception("StdMeshers_ProjectionUtils: attempt to associate NULL shape");
01428   }
01429   return false;
01430 }
01431 
01432 //=======================================================================
01440 //=======================================================================
01441 
01442 TopoDS_Edge StdMeshers_ProjectionUtils::GetEdgeByVertices( SMESH_Mesh*          theMesh,
01443                                                            const TopoDS_Vertex& theV1,
01444                                                            const TopoDS_Vertex& theV2)
01445 {
01446   if ( theMesh && !theV1.IsNull() && !theV2.IsNull() )
01447   {
01448     TopTools_ListIteratorOfListOfShape ancestorIt( theMesh->GetAncestors( theV1 ));
01449     for ( ; ancestorIt.More(); ancestorIt.Next() )
01450       if ( ancestorIt.Value().ShapeType() == TopAbs_EDGE )
01451         for ( TopExp_Explorer expV ( ancestorIt.Value(), TopAbs_VERTEX );
01452               expV.More();
01453               expV.Next() )
01454           if ( theV2.IsSame( expV.Current() ))
01455             return TopoDS::Edge( ancestorIt.Value() );
01456   }
01457   return TopoDS_Edge();
01458 }
01459 
01460 //================================================================================
01468 //================================================================================
01469 
01470 TopoDS_Face StdMeshers_ProjectionUtils::GetNextFace( const TAncestorMap& edgeToFaces,
01471                                                      const TopoDS_Edge&  edge,
01472                                                      const TopoDS_Face&  face)
01473 {
01474 //   if ( !edge.IsNull() && !face.IsNull() && edgeToFaces.Contains( edge ))
01475   if ( !edge.IsNull() && edgeToFaces.Contains( edge )) // PAL16202
01476   {
01477     TopTools_ListIteratorOfListOfShape ancestorIt( edgeToFaces.FindFromKey( edge ));
01478     for ( ; ancestorIt.More(); ancestorIt.Next() )
01479       if ( ancestorIt.Value().ShapeType() == TopAbs_FACE &&
01480            !face.IsSame( ancestorIt.Value() ))
01481         return TopoDS::Face( ancestorIt.Value() );
01482   }
01483   return TopoDS_Face();
01484 }
01485 
01486 //================================================================================
01490 //================================================================================
01491 
01492 TopoDS_Vertex StdMeshers_ProjectionUtils::GetNextVertex(const TopoDS_Edge&   edge,
01493                                                         const TopoDS_Vertex& vertex)
01494 {
01495   TopoDS_Vertex vF,vL;
01496   TopExp::Vertices(edge,vF,vL);
01497   if ( vF.IsSame( vL ))
01498     return TopoDS_Vertex();
01499   return vertex.IsSame( vF ) ? vL : vF; 
01500 }
01501 
01502 //================================================================================
01510 //================================================================================
01511 
01512 pair<int,TopoDS_Edge>
01513 StdMeshers_ProjectionUtils::GetPropagationEdge( SMESH_Mesh*        aMesh,
01514                                                 const TopoDS_Edge& theEdge,
01515                                                 const TopoDS_Edge& fromEdge)
01516 {
01517   TopTools_IndexedMapOfShape aChain;
01518   int step = 0;
01519 
01520   // List of edges, added to chain on the previous cycle pass
01521   TopTools_ListOfShape listPrevEdges;
01522   listPrevEdges.Append(fromEdge);
01523 
01524   // Collect all edges pass by pass
01525   while (listPrevEdges.Extent() > 0) {
01526     step++;
01527     // List of edges, added to chain on this cycle pass
01528     TopTools_ListOfShape listCurEdges;
01529 
01530     // Find the next portion of edges
01531     TopTools_ListIteratorOfListOfShape itE (listPrevEdges);
01532     for (; itE.More(); itE.Next()) {
01533       TopoDS_Shape anE = itE.Value();
01534 
01535       // Iterate on faces, having edge <anE>
01536       TopTools_ListIteratorOfListOfShape itA (aMesh->GetAncestors(anE));
01537       for (; itA.More(); itA.Next()) {
01538         TopoDS_Shape aW = itA.Value();
01539 
01540         // There are objects of different type among the ancestors of edge
01541         if (aW.ShapeType() == TopAbs_WIRE) {
01542           TopoDS_Shape anOppE;
01543 
01544           BRepTools_WireExplorer aWE (TopoDS::Wire(aW));
01545           Standard_Integer nb = 1, found = 0;
01546           TopTools_Array1OfShape anEdges (1,4);
01547           for (; aWE.More(); aWE.Next(), nb++) {
01548             if (nb > 4) {
01549               found = 0;
01550               break;
01551             }
01552             anEdges(nb) = aWE.Current();
01553             if (anEdges(nb).IsSame(anE)) found = nb;
01554           }
01555 
01556           if (nb == 5 && found > 0) {
01557             // Quadrangle face found, get an opposite edge
01558             Standard_Integer opp = found + 2;
01559             if (opp > 4) opp -= 4;
01560             anOppE = anEdges(opp);
01561 
01562             // add anOppE to aChain if ...
01563             if (!aChain.Contains(anOppE)) { // ... anOppE is not in aChain
01564               // Add found edge to the chain oriented so that to
01565               // have it co-directed with a forward MainEdge
01566               TopAbs_Orientation ori = anE.Orientation();
01567               if ( anEdges(opp).Orientation() == anEdges(found).Orientation() )
01568                 ori = TopAbs::Reverse( ori );
01569               anOppE.Orientation( ori );
01570               if ( anOppE.IsSame( theEdge ))
01571                 return make_pair( step, TopoDS::Edge( anOppE ));
01572               aChain.Add(anOppE);
01573               listCurEdges.Append(anOppE);
01574             }
01575           } // if (nb == 5 && found > 0)
01576         } // if (aF.ShapeType() == TopAbs_WIRE)
01577       } // for (; itF.More(); itF.Next())
01578     } // for (; itE.More(); itE.Next())
01579 
01580     listPrevEdges = listCurEdges;
01581   } // while (listPrevEdges.Extent() > 0)
01582 
01583   return make_pair( INT_MAX, TopoDS_Edge());
01584 }
01585 
01586 //================================================================================
01597 //================================================================================
01598 
01599 bool StdMeshers_ProjectionUtils::
01600 FindMatchingNodesOnFaces( const TopoDS_Face&     face1,
01601                           SMESH_Mesh*            mesh1,
01602                           const TopoDS_Face&     face2,
01603                           SMESH_Mesh*            mesh2,
01604                           const TShapeShapeMap & assocMap,
01605                           TNodeNodeMap &         node1To2Map)
01606 {
01607   SMESHDS_Mesh* meshDS1 = mesh1->GetMeshDS();
01608   SMESHDS_Mesh* meshDS2 = mesh2->GetMeshDS();
01609   
01610   SMESH_MesherHelper helper1( *mesh1 );
01611   SMESH_MesherHelper helper2( *mesh2 );
01612 
01613   // Get corresponding submeshes and roughly check match of meshes
01614 
01615   SMESHDS_SubMesh * SM2 = meshDS2->MeshElements( face2 );
01616   SMESHDS_SubMesh * SM1 = meshDS1->MeshElements( face1 );
01617   if ( !SM2 || !SM1 )
01618     RETURN_BAD_RESULT("Empty submeshes");
01619   if ( SM2->NbNodes()    != SM1->NbNodes() ||
01620        SM2->NbElements() != SM1->NbElements() )
01621     RETURN_BAD_RESULT("Different meshes on corresponding faces "
01622                       << meshDS1->ShapeToIndex( face1 ) << " and "
01623                       << meshDS2->ShapeToIndex( face2 ));
01624   if ( SM2->NbElements() == 0 )
01625     RETURN_BAD_RESULT("Empty submeshes");
01626 
01627   helper1.SetSubShape( face1 );
01628   helper2.SetSubShape( face2 );
01629   if ( helper1.HasSeam() != helper2.HasSeam() )
01630     RETURN_BAD_RESULT("Different faces' geometry");
01631 
01632   // Data to call SMESH_MeshEditor::FindMatchingNodes():
01633 
01634   // 1. Nodes of corresponding links:
01635 
01636   // get 2 matching edges, try to find not seam ones
01637   TopoDS_Edge edge1, edge2, seam1, seam2, anyEdge1, anyEdge2;
01638   TopExp_Explorer eE( OuterShape( face2, TopAbs_WIRE ), TopAbs_EDGE );
01639   do {
01640     // edge 2
01641     TopoDS_Edge e2 = TopoDS::Edge( eE.Current() );
01642     eE.Next();
01643     // edge 1
01644     if ( !assocMap.IsBound( e2, /*is2nd=*/true ))
01645       RETURN_BAD_RESULT("Association not found for edge " << meshDS2->ShapeToIndex( e2 ));
01646     TopoDS_Edge e1 = TopoDS::Edge( assocMap( e2, /*is2nd=*/true ));
01647     if ( !helper1.IsSubShape( e1, face1 ))
01648       RETURN_BAD_RESULT("Wrong association, edge " << meshDS1->ShapeToIndex( e1 ) <<
01649                         " isn't a sub-shape of face " << meshDS1->ShapeToIndex( face1 ));
01650     // check that there are nodes on edges
01651     SMESHDS_SubMesh * eSM1 = meshDS1->MeshElements( e1 );
01652     SMESHDS_SubMesh * eSM2 = meshDS2->MeshElements( e2 );
01653     bool nodesOnEdges = ( eSM1 && eSM2 && eSM1->NbNodes() && eSM2->NbNodes() );
01654     // check that the nodes on edges belong to faces
01655     // (as NETGEN ignores nodes on the degenerated geom edge)
01656     bool nodesOfFaces = false;
01657     if ( nodesOnEdges ) {
01658       const SMDS_MeshNode* n1 = eSM1->GetNodes()->next();
01659       const SMDS_MeshNode* n2 = eSM2->GetNodes()->next();
01660       nodesOfFaces = ( n1->GetInverseElementIterator(SMDSAbs_Face)->more() &&
01661                        n2->GetInverseElementIterator(SMDSAbs_Face)->more() );
01662     }
01663     if ( nodesOfFaces )
01664     {
01665       if ( helper2.IsRealSeam( e2 )) {
01666         seam1 = e1; seam2 = e2;
01667       }
01668       else {
01669         edge1 = e1; edge2 = e2;
01670       }
01671     }
01672     else {
01673       anyEdge1 = e1; anyEdge2 = e2;
01674     }
01675   } while ( edge2.IsNull() && eE.More() );
01676   //
01677   if ( edge2.IsNull() ) {
01678     edge1 = seam1; edge2 = seam2;
01679   }
01680   bool hasNodesOnEdge = (! edge2.IsNull() );
01681   if ( !hasNodesOnEdge ) {
01682     // 0020338 - nb segments == 1
01683     edge1 = anyEdge1; edge2 = anyEdge2;
01684   }
01685 
01686   // get 2 matching vertices
01687   TopoDS_Vertex V2 = TopExp::FirstVertex( TopoDS::Edge( edge2 ));
01688   if ( !assocMap.IsBound( V2, /*is2nd=*/true ))
01689     RETURN_BAD_RESULT("Association not found for vertex " << meshDS2->ShapeToIndex( V2 ));
01690   TopoDS_Vertex V1 = TopoDS::Vertex( assocMap( V2, /*is2nd=*/true ));
01691 
01692   // nodes on vertices
01693   const SMDS_MeshNode* vNode1 = SMESH_Algo::VertexNode( V1, meshDS1 );
01694   const SMDS_MeshNode* vNode2 = SMESH_Algo::VertexNode( V2, meshDS2 );
01695   if ( !vNode1 ) RETURN_BAD_RESULT("No node on vertex #" << meshDS1->ShapeToIndex( V1 ));
01696   if ( !vNode2 ) RETURN_BAD_RESULT("No node on vertex #" << meshDS2->ShapeToIndex( V2 ));
01697 
01698   // nodes on edges linked with nodes on vertices
01699   const SMDS_MeshNode* nullNode = 0;
01700   vector< const SMDS_MeshNode*> eNode1( 2, nullNode );
01701   vector< const SMDS_MeshNode*> eNode2( 2, nullNode );
01702   if ( hasNodesOnEdge )
01703   {
01704     int nbNodeToGet = 1;
01705     if ( helper1.IsClosedEdge( edge1 ) || helper2.IsClosedEdge( edge2 ) )
01706       nbNodeToGet = 2;
01707     for ( int is2 = 0; is2 < 2; ++is2 )
01708     {
01709       TopoDS_Edge &     edge  = is2 ? edge2 : edge1;
01710       SMESHDS_Mesh *    smDS  = is2 ? meshDS2 : meshDS1;
01711       SMESHDS_SubMesh* edgeSM = smDS->MeshElements( edge );
01712       // nodes linked with ones on vertices
01713       const SMDS_MeshNode*           vNode = is2 ? vNode2 : vNode1;
01714       vector< const SMDS_MeshNode*>& eNode = is2 ? eNode2 : eNode1;
01715       int nbGotNode = 0;
01716       SMDS_ElemIteratorPtr vElem = vNode->GetInverseElementIterator(SMDSAbs_Edge);
01717       while ( vElem->more() && nbGotNode != nbNodeToGet ) {
01718         const SMDS_MeshElement* elem = vElem->next();
01719         if ( edgeSM->Contains( elem ))
01720           eNode[ nbGotNode++ ] = 
01721             ( elem->GetNode(0) == vNode ) ? elem->GetNode(1) : elem->GetNode(0);
01722       }
01723       if ( nbGotNode > 1 ) // sort found nodes by param on edge
01724       {
01725         SMESH_MesherHelper* helper = is2 ? &helper2 : &helper1;
01726         double u0 = helper->GetNodeU( edge, eNode[ 0 ]);
01727         double u1 = helper->GetNodeU( edge, eNode[ 1 ]);
01728         if ( u0 > u1 ) std::swap( eNode[ 0 ], eNode[ 1 ]);
01729       }
01730       if ( nbGotNode == 0 )
01731         RETURN_BAD_RESULT("Found no nodes on edge " << smDS->ShapeToIndex( edge ) <<
01732                           " linked to " << vNode );
01733     }
01734   }
01735   else // 0020338 - nb segments == 1
01736   {
01737     // get 2 other matching vertices
01738     V2 = TopExp::LastVertex( TopoDS::Edge( edge2 ));
01739     if ( !assocMap.IsBound( V2, /*is2nd=*/true ))
01740       RETURN_BAD_RESULT("Association not found for vertex " << meshDS2->ShapeToIndex( V2 ));
01741     V1 = TopoDS::Vertex( assocMap( V2, /*is2nd=*/true ));
01742 
01743     // nodes on vertices
01744     eNode1[0] = SMESH_Algo::VertexNode( V1, meshDS1 );
01745     eNode2[0] = SMESH_Algo::VertexNode( V2, meshDS2 );
01746     if ( !eNode1[0] ) RETURN_BAD_RESULT("No node on vertex #" << meshDS1->ShapeToIndex( V1 ));
01747     if ( !eNode2[0] ) RETURN_BAD_RESULT("No node on vertex #" << meshDS2->ShapeToIndex( V2 ));
01748   }
01749 
01750   // 2. face sets
01751 
01752   set<const SMDS_MeshElement*> Elems1, Elems2;
01753   for ( int is2 = 0; is2 < 2; ++is2 )
01754   {
01755     set<const SMDS_MeshElement*> & elems = is2 ? Elems2 : Elems1;
01756     SMESHDS_SubMesh*                  sm = is2 ? SM2 : SM1;
01757     SMESH_MesherHelper*           helper = is2 ? &helper2 : &helper1;
01758     const TopoDS_Face &             face = is2 ? face2 : face1;
01759     SMDS_ElemIteratorPtr eIt = sm->GetElements();
01760 
01761     if ( !helper->IsRealSeam( is2 ? edge2 : edge1 ))
01762     {
01763       while ( eIt->more() ) elems.insert( eIt->next() );
01764     }
01765     else
01766     {
01767       // the only suitable edge is seam, i.e. it is a sphere.
01768       // FindMatchingNodes() will not know which way to go from any edge.
01769       // So we ignore all faces having nodes on edges or vertices except
01770       // one of faces sharing current start nodes
01771 
01772       // find a face to keep
01773       const SMDS_MeshElement* faceToKeep = 0;
01774       const SMDS_MeshNode* vNode = is2 ? vNode2 : vNode1;
01775       const SMDS_MeshNode* eNode = is2 ? eNode2[0] : eNode1[0];
01776       TIDSortedElemSet inSet, notInSet;
01777 
01778       const SMDS_MeshElement* f1 =
01779         SMESH_MeshEditor::FindFaceInSet( vNode, eNode, inSet, notInSet );
01780       if ( !f1 ) RETURN_BAD_RESULT("The first face on seam not found");
01781       notInSet.insert( f1 );
01782 
01783       const SMDS_MeshElement* f2 =
01784         SMESH_MeshEditor::FindFaceInSet( vNode, eNode, inSet, notInSet );
01785       if ( !f2 ) RETURN_BAD_RESULT("The second face on seam not found");
01786 
01787       // select a face with less UV of vNode
01788       const SMDS_MeshNode* notSeamNode[2] = {0, 0};
01789       for ( int iF = 0; iF < 2; ++iF ) {
01790         const SMDS_MeshElement* f = ( iF ? f2 : f1 );
01791         for ( int i = 0; !notSeamNode[ iF ] && i < f->NbNodes(); ++i ) {
01792           const SMDS_MeshNode* node = f->GetNode( i );
01793           if ( !helper->IsSeamShape( node->getshapeId() ))
01794             notSeamNode[ iF ] = node;
01795         }
01796       }
01797       gp_Pnt2d uv1 = helper->GetNodeUV( face, vNode, notSeamNode[0] );
01798       gp_Pnt2d uv2 = helper->GetNodeUV( face, vNode, notSeamNode[1] );
01799       if ( uv1.X() + uv1.Y() > uv2.X() + uv2.Y() )
01800         faceToKeep = f2;
01801       else
01802         faceToKeep = f1;
01803 
01804       // fill elem set
01805       elems.insert( faceToKeep );
01806       while ( eIt->more() ) {
01807         const SMDS_MeshElement* f = eIt->next();
01808         int nbNodes = f->NbNodes();
01809         if ( f->IsQuadratic() )
01810           nbNodes /= 2;
01811         bool onBnd = false;
01812         for ( int i = 0; !onBnd && i < nbNodes; ++i ) {
01813           const SMDS_MeshNode* node = f->GetNode( i );
01814           onBnd = ( node->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE);
01815         }
01816         if ( !onBnd )
01817           elems.insert( f );
01818       }
01819       // add also faces adjacent to faceToKeep
01820       int nbNodes = faceToKeep->NbNodes();
01821       if ( faceToKeep->IsQuadratic() ) nbNodes /= 2;
01822       notInSet.insert( f1 );
01823       notInSet.insert( f2 );
01824       for ( int i = 0; i < nbNodes; ++i ) {
01825         const SMDS_MeshNode* n1 = faceToKeep->GetNode( i );
01826         const SMDS_MeshNode* n2 = faceToKeep->GetNode(( i+1 ) % nbNodes );
01827         f1 = SMESH_MeshEditor::FindFaceInSet( n1, n2, inSet, notInSet );
01828         if ( f1 )
01829           elems.insert( f1 );
01830       }
01831     } // case on a sphere
01832   } // loop on 2 faces
01833 
01834   //  int quadFactor = (*Elems1.begin())->IsQuadratic() ? 2 : 1;
01835 
01836   node1To2Map.clear();
01837   int res = SMESH_MeshEditor::FindMatchingNodes( Elems1, Elems2,
01838                                                  vNode1, vNode2,
01839                                                  eNode1[0], eNode2[0],
01840                                                  node1To2Map);
01841   if ( res != SMESH_MeshEditor::SEW_OK )
01842     RETURN_BAD_RESULT("FindMatchingNodes() result " << res );
01843 
01844   // On a sphere, add matching nodes on the edge
01845 
01846   if ( helper1.IsRealSeam( edge1 ))
01847   {
01848     // sort nodes on edges by param on edge
01849     map< double, const SMDS_MeshNode* > u2nodesMaps[2];
01850     for ( int is2 = 0; is2 < 2; ++is2 )
01851     {
01852       TopoDS_Edge &     edge  = is2 ? edge2 : edge1;
01853       SMESHDS_Mesh *    smDS  = is2 ? meshDS2 : meshDS1;
01854       SMESHDS_SubMesh* edgeSM = smDS->MeshElements( edge );
01855       map< double, const SMDS_MeshNode* > & pos2nodes = u2nodesMaps[ is2 ];
01856 
01857       SMDS_NodeIteratorPtr nIt = edgeSM->GetNodes();
01858       while ( nIt->more() ) {
01859         const SMDS_MeshNode* node = nIt->next();
01860         const SMDS_EdgePosition* pos =
01861           static_cast<const SMDS_EdgePosition*>(node->GetPosition());
01862         pos2nodes.insert( make_pair( pos->GetUParameter(), node ));
01863       }
01864       if ( pos2nodes.size() != edgeSM->NbNodes() )
01865         RETURN_BAD_RESULT("Equal params of nodes on edge "
01866                           << smDS->ShapeToIndex( edge ) << " of face " << is2 );
01867     }
01868     if ( u2nodesMaps[0].size() != u2nodesMaps[1].size() )
01869       RETURN_BAD_RESULT("Different nb of new nodes on edges or wrong params");
01870 
01871     // compare edge orientation
01872     double u1 = helper1.GetNodeU( edge1, vNode1 );
01873     double u2 = helper2.GetNodeU( edge2, vNode2 );
01874     bool isFirst1 = ( u1 < u2nodesMaps[0].begin()->first );
01875     bool isFirst2 = ( u2 < u2nodesMaps[1].begin()->first );
01876     bool reverse ( isFirst1 != isFirst2 );
01877 
01878     // associate matching nodes
01879     map< double, const SMDS_MeshNode* >::iterator u_Node1, u_Node2, end1;
01880     map< double, const SMDS_MeshNode* >::reverse_iterator uR_Node2;
01881     u_Node1 = u2nodesMaps[0].begin();
01882     u_Node2 = u2nodesMaps[1].begin();
01883     uR_Node2 = u2nodesMaps[1].rbegin();
01884     end1 = u2nodesMaps[0].end();
01885     for ( ; u_Node1 != end1; ++u_Node1 ) {
01886       const SMDS_MeshNode* n1 = u_Node1->second;
01887       const SMDS_MeshNode* n2 = ( reverse ? (uR_Node2++)->second : (u_Node2++)->second );
01888       node1To2Map.insert( make_pair( n1, n2 ));
01889     }
01890 
01891     // associate matching nodes on the last vertices
01892     V2 = TopExp::LastVertex( TopoDS::Edge( edge2 ));
01893     if ( !assocMap.IsBound( V2, /*is2nd=*/true ))
01894       RETURN_BAD_RESULT("Association not found for vertex " << meshDS2->ShapeToIndex( V2 ));
01895     V1 = TopoDS::Vertex( assocMap( V2, /*is2nd=*/true ));
01896     vNode1 = SMESH_Algo::VertexNode( V1, meshDS1 );
01897     vNode2 = SMESH_Algo::VertexNode( V2, meshDS2 );
01898     if ( !vNode1 ) RETURN_BAD_RESULT("No node on vertex #" << meshDS1->ShapeToIndex( V1 ));
01899     if ( !vNode2 ) RETURN_BAD_RESULT("No node on vertex #" << meshDS2->ShapeToIndex( V2 ));
01900     node1To2Map.insert( make_pair( vNode1, vNode2 ));
01901   }
01902 
01903 // don't know why this condition is usually true :(
01904 //   if ( node1To2Map.size() * quadFactor < SM1->NbNodes() )
01905 //     MESSAGE("FindMatchingNodes() found too few node pairs starting from nodes ("
01906 //             << vNode1->GetID() << " - " << eNode1[0]->GetID() << ") ("
01907 //             << vNode2->GetID() << " - " << eNode2[0]->GetID() << "):"
01908 //             << node1To2Map.size() * quadFactor << " < " << SM1->NbNodes());
01909   
01910   return true;
01911 }
01912 
01913 //================================================================================
01920 //================================================================================
01921 
01922 TopoDS_Shape StdMeshers_ProjectionUtils::OuterShape( const TopoDS_Face& face,
01923                                                      TopAbs_ShapeEnum   type)
01924 {
01925   TopExp_Explorer exp( BRepTools::OuterWire( face ), type );
01926   if ( exp.More() )
01927     return exp.Current();
01928   return TopoDS_Shape();
01929 }
01930 
01931 //================================================================================
01938 //================================================================================
01939 
01940 bool StdMeshers_ProjectionUtils::MakeComputed(SMESH_subMesh * sm, const int iterationNb)
01941 {
01942   if ( iterationNb > 10 )
01943     RETURN_BAD_RESULT("Infinite recursive projection");
01944   if ( !sm )
01945     RETURN_BAD_RESULT("NULL submesh");
01946   if ( sm->IsMeshComputed() )
01947     return true;
01948 
01949   SMESH_Mesh* mesh = sm->GetFather();
01950   SMESH_Gen* gen   = mesh->GetGen();
01951   SMESH_Algo* algo = sm->GetAlgo();
01952   if ( !algo )
01953   {
01954     if ( sm->GetSubShape().ShapeType() != TopAbs_COMPOUND )
01955       RETURN_BAD_RESULT("No algo assigned to submesh " << sm->GetId());
01956     // group
01957     bool computed = true;
01958     for ( TopoDS_Iterator grMember( sm->GetSubShape() ); grMember.More(); grMember.Next())
01959       if ( SMESH_subMesh* grSub = mesh->GetSubMesh( grMember.Value() ))
01960         if ( !MakeComputed( grSub, iterationNb + 1 ))
01961           computed = false;
01962     return computed;
01963   }
01964 
01965   string algoType = algo->GetName();
01966   if ( algoType.substr(0, 11) != "Projection_")
01967     return gen->Compute( *mesh, sm->GetSubShape() );
01968 
01969   // try to compute source mesh
01970 
01971   const list <const SMESHDS_Hypothesis *> & hyps =
01972     algo->GetUsedHypothesis( *mesh, sm->GetSubShape() );
01973 
01974   TopoDS_Shape srcShape;
01975   SMESH_Mesh* srcMesh = 0;
01976   list <const SMESHDS_Hypothesis*>::const_iterator hIt = hyps.begin();
01977   for ( ; srcShape.IsNull() && hIt != hyps.end(); ++hIt ) {
01978     string hypName = (*hIt)->GetName();
01979     if ( hypName == "ProjectionSource1D" ) {
01980       const StdMeshers_ProjectionSource1D * hyp =
01981         static_cast<const StdMeshers_ProjectionSource1D*>( *hIt );
01982       srcShape = hyp->GetSourceEdge();
01983       srcMesh = hyp->GetSourceMesh();
01984     }
01985     else if ( hypName == "ProjectionSource2D" ) {
01986       const StdMeshers_ProjectionSource2D * hyp =
01987         static_cast<const StdMeshers_ProjectionSource2D*>( *hIt );
01988       srcShape = hyp->GetSourceFace();
01989       srcMesh = hyp->GetSourceMesh();
01990     }
01991     else if ( hypName == "ProjectionSource3D" ) {
01992       const StdMeshers_ProjectionSource3D * hyp =
01993         static_cast<const StdMeshers_ProjectionSource3D*>( *hIt );
01994       srcShape = hyp->GetSource3DShape();
01995       srcMesh = hyp->GetSourceMesh();
01996     }
01997   }
01998   if ( srcShape.IsNull() ) // no projection source defined
01999     return gen->Compute( *mesh, sm->GetSubShape() );
02000 
02001   if ( srcShape.IsSame( sm->GetSubShape() ))
02002     RETURN_BAD_RESULT("Projection from self");
02003     
02004   if ( !srcMesh )
02005     srcMesh = mesh;
02006 
02007   if ( MakeComputed( srcMesh->GetSubMesh( srcShape ), iterationNb + 1 ))
02008     return gen->Compute( *mesh, sm->GetSubShape() );
02009 
02010   return false;
02011 }
02012 
02013 //================================================================================
02020 //================================================================================
02021 
02022 int StdMeshers_ProjectionUtils::Count(const TopoDS_Shape&    shape,
02023                                       const TopAbs_ShapeEnum type,
02024                                       const bool             ignoreSame)
02025 {
02026   if ( ignoreSame ) {
02027     TopTools_IndexedMapOfShape map;
02028     TopExp::MapShapes( shape, type, map );
02029     return map.Extent();
02030   }
02031   else {
02032     int nb = 0;
02033     for ( TopExp_Explorer exp( shape, type ); exp.More(); exp.Next() )
02034       ++nb;
02035     return nb;
02036   }
02037 }
02038 
02039 //================================================================================
02043 //================================================================================
02044 
02045 TopoDS_Edge StdMeshers_ProjectionUtils::GetBoundaryEdge(const TopoDS_Shape& edgeContainer,
02046                                                         const SMESH_Mesh&   mesh)
02047 {
02048   TopTools_IndexedMapOfShape facesOfEdgeContainer, facesNearEdge;
02049   TopExp::MapShapes( edgeContainer, TopAbs_FACE, facesOfEdgeContainer );
02050 
02051   if ( !facesOfEdgeContainer.IsEmpty() ) 
02052     for ( TopExp_Explorer exp(edgeContainer, TopAbs_EDGE); exp.More(); exp.Next() )
02053     {
02054       const TopoDS_Edge& edge = TopoDS::Edge( exp.Current() );
02055       facesNearEdge.Clear();
02056       PShapeIteratorPtr faceIt = SMESH_MesherHelper::GetAncestors( edge, mesh, TopAbs_FACE );
02057       while ( const TopoDS_Shape* face = faceIt->next() )
02058         if ( facesOfEdgeContainer.Contains( *face ))
02059           if ( facesNearEdge.Add( *face ) && facesNearEdge.Extent() > 1 )
02060             break;
02061       if ( facesNearEdge.Extent() == 1 )
02062         return edge;
02063     }
02064 
02065   return TopoDS_Edge();
02066 }
02067 
02068 
02069 namespace { // Definition of event listeners
02070 
02071   SMESH_subMeshEventListener* GetSrcSubMeshListener();
02072 
02073   //================================================================================
02078   //================================================================================
02079 
02080   struct HypModifWaiter: SMESH_subMeshEventListener
02081   {
02082     HypModifWaiter():SMESH_subMeshEventListener(false,// won't be deleted by submesh
02083                                                 "StdMeshers_ProjectionUtils::HypModifWaiter") {}
02084     void ProcessEvent(const int event, const int eventType, SMESH_subMesh* subMesh,
02085                       EventListenerData*, const SMESH_Hypothesis*)
02086     {
02087       if ( event     == SMESH_subMesh::MODIF_HYP &&
02088            eventType == SMESH_subMesh::ALGO_EVENT)
02089       {
02090         // delete current source listener
02091         subMesh->DeleteEventListener( GetSrcSubMeshListener() );
02092         // let algo set a new one
02093         if ( SMESH_Algo* algo = subMesh->GetAlgo() )
02094           algo->SetEventListener( subMesh );
02095       }
02096     }
02097   };
02098   //================================================================================
02102   //================================================================================
02103 
02104   SMESH_subMeshEventListener* GetHypModifWaiter() {
02105     static HypModifWaiter aHypModifWaiter;
02106     return &aHypModifWaiter;
02107   }
02108   //================================================================================
02112   //================================================================================
02113 
02114   SMESH_subMeshEventListener* GetSrcSubMeshListener() {
02115     static SMESH_subMeshEventListener srcListener(false, // won't be deleted by submesh
02116                                                   "StdMeshers_ProjectionUtils::SrcSubMeshListener");
02117     return &srcListener;
02118   }
02119 }
02120 
02121 //================================================================================
02128 //================================================================================
02129 
02130 void StdMeshers_ProjectionUtils::SetEventListener(SMESH_subMesh* subMesh,
02131                                                   TopoDS_Shape   srcShape,
02132                                                   SMESH_Mesh*    srcMesh)
02133 {
02134   // Set listener that resets an event listener on source submesh when
02135   // "ProjectionSource*D" hypothesis is modified since source shape can be changed
02136   subMesh->SetEventListener( GetHypModifWaiter(),0,subMesh);
02137 
02138   // Set an event listener to submesh of the source shape
02139   if ( !srcShape.IsNull() )
02140   {
02141     if ( !srcMesh )
02142       srcMesh = subMesh->GetFather();
02143 
02144     SMESH_subMesh* srcShapeSM = srcMesh->GetSubMesh( srcShape );
02145 
02146     if ( srcShapeSM != subMesh ) {
02147       if ( srcShapeSM->GetSubMeshDS() &&
02148            srcShapeSM->GetSubMeshDS()->IsComplexSubmesh() )
02149       {  // source shape is a group
02150         TopExp_Explorer it(srcShapeSM->GetSubShape(), // explore the group into sub-shapes...
02151                            subMesh->GetSubShape().ShapeType()); // ...of target shape type
02152         for (; it.More(); it.Next())
02153         {
02154           SMESH_subMesh* srcSM = srcMesh->GetSubMesh( it.Current() );
02155           if ( srcSM != subMesh )
02156           {
02157             SMESH_subMeshEventListenerData* data =
02158               srcSM->GetEventListenerData(GetSrcSubMeshListener());
02159             if ( data )
02160               data->mySubMeshes.push_back( subMesh );
02161             else
02162               data = SMESH_subMeshEventListenerData::MakeData( subMesh );
02163             subMesh->SetEventListener ( GetSrcSubMeshListener(), data, srcSM );
02164           }
02165         }
02166       }
02167       else
02168       {
02169         subMesh->SetEventListener( GetSrcSubMeshListener(),
02170                                    SMESH_subMeshEventListenerData::MakeData( subMesh ),
02171                                    srcShapeSM );
02172       }
02173     }
02174   }
02175 }