Version: 6.3.1

src/StdMeshers/StdMeshers_ProjectionUtils.cxx

Go to the documentation of this file.
00001 // Copyright (C) 2007-2011  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;
00311       TopExp_Explorer edgeExp1( theGroup1, TopAbs_EDGE );
00312       for ( ; bndEdge.IsNull() && edgeExp1.More(); edgeExp1.Next())
00313         if ( HERE::IsBoundaryEdge( TopoDS::Edge( edgeExp1.Current()), theGroup1, theMesh ))
00314           bndEdge = edgeExp1.Current();
00315       if ( bndEdge.IsNull() )
00316         return false;
00317 
00318       list< TopoDS_Shape > edges(1, bndEdge);
00319       list< TopoDS_Shape >::iterator edge1 = edges.begin();
00320       for ( ; edge1 != edges.end(); ++edge1 )
00321       {
00322         // there must be one or zero not associated faces between ancestors of edge
00323         // belonging to theGroup1
00324         TopoDS_Shape face1;
00325         ancestIt.Initialize( theMesh.GetAncestors( *edge1 ) );
00326         for ( ; ancestIt.More() && face1.IsNull(); ancestIt.Next() ) {
00327           if ( ancestIt.Value().ShapeType() == TopAbs_FACE &&
00328                !theMap.IsBound( ancestIt.Value() ) &&
00329                subshapes1.Contains( ancestIt.Value() ))
00330             face1 = ancestIt.Value();
00331 
00332           // add edges of face1 to start searching for adjacent faces from
00333           for ( TopExp_Explorer e(face1, TopAbs_EDGE); e.More(); e.Next())
00334             if ( !edge1->IsSame( e.Current() ))
00335               edges.push_back( e.Current() );
00336         }
00337         if ( !face1.IsNull() ) {
00338           // find the corresponding face of theGroup2
00339           TopoDS_Shape edge2 = theMap( *edge1 );
00340           TopoDS_Shape face2;
00341           ancestIt.Initialize( theMesh.GetAncestors( edge2 ) );
00342           for ( ; ancestIt.More() && face2.IsNull(); ancestIt.Next() ) {
00343             if ( ancestIt.Value().ShapeType() == TopAbs_FACE &&
00344                  !theMap.IsBound( ancestIt.Value() ) &&
00345                  subshapes2.Contains( ancestIt.Value() ))
00346               face2 = ancestIt.Value();
00347           }
00348           if ( face2.IsNull() )
00349             return false;
00350 
00351           HERE::InsertAssociation( face1, face2, theMap);
00352         }
00353       }
00354     }
00355     return true;
00356   }
00357 
00358   //================================================================================
00363   //================================================================================
00364 
00365   bool sameVertexUV( const TopoDS_Edge& edge,
00366                      const TopoDS_Face& face,
00367                      const int&         vIndex,
00368                      const gp_Pnt2d&    uv,
00369                      const double&      tol2d )
00370   {
00371     TopoDS_Vertex VV[2];
00372     TopExp::Vertices( edge, VV[0], VV[1], true);
00373     gp_Pnt2d v1UV = BRep_Tool::Parameters( VV[vIndex], face);
00374     double dist2d = v1UV.Distance( uv );
00375     return dist2d < tol2d;
00376   }
00377 
00378 } // namespace
00379 
00380 //=======================================================================
00391 //=======================================================================
00392 
00393 bool StdMeshers_ProjectionUtils::FindSubShapeAssociation(const TopoDS_Shape& theShape1,
00394                                                          SMESH_Mesh*         theMesh1,
00395                                                          const TopoDS_Shape& theShape2,
00396                                                          SMESH_Mesh*         theMesh2,
00397                                                          TShapeShapeMap &    theMap)
00398 {
00399   // Structure of this long function is following
00400   // 1) Group->group projection: theShape1 is a group member,
00401   //    theShape2 is a group. We find a group theShape1 is in and recall self.
00402   // 2) Accosiate same shapes with different location (partners).
00403   // 3) If vertex association is given, perform accosiation according to shape type:
00404   //       switch ( ShapeType ) {
00405   //         case TopAbs_EDGE:
00406   //         case ...:
00407   //       }
00408   // 4) else try to accosiate in different ways:
00409   //       a) accosiate shapes by propagation and other simple cases
00410   //            switch ( ShapeType ) {
00411   //            case TopAbs_EDGE:
00412   //            case ...:
00413   //            }
00414   //       b) find association of a couple of vertices and recall self.
00415   //
00416 
00417   theMeshDS[0] = theMesh1->GetMeshDS(); // debug
00418   theMeshDS[1] = theMesh2->GetMeshDS();
00419 
00420   // =================================================================================
00421   // 1) Is it the case of associating a group member -> another group? (PAL16202, 16203)
00422   // =================================================================================
00423   if ( theShape1.ShapeType() != theShape2.ShapeType() ) {
00424     TopoDS_Shape group1, group2;
00425     if ( theShape1.ShapeType() == TopAbs_COMPOUND ) {
00426       group1 = theShape1;
00427       group2 = FindGroupContaining( theShape2, theMesh2, group1 );
00428     }
00429     else if ( theShape2.ShapeType() == TopAbs_COMPOUND ) {
00430       group2 = theShape2;
00431       group1 = FindGroupContaining( theShape1, theMesh1, group2 );
00432     }
00433     if ( group1.IsNull() || group2.IsNull() )
00434       RETURN_BAD_RESULT("Different shape types");
00435     // Associate compounds
00436     return FindSubShapeAssociation(group1, theMesh1, group2, theMesh2, theMap );
00437   }
00438 
00439   bool bidirect = ( !theShape1.IsSame( theShape2 ));
00440 
00441   // ============
00442   // 2) Is partner?
00443   // ============
00444   bool partner = theShape1.IsPartner( theShape2 );
00445   TopTools_DataMapIteratorOfDataMapOfShapeShape vvIt( theMap );
00446   for ( ; partner && vvIt.More(); vvIt.Next() )
00447     partner = vvIt.Key().IsPartner( vvIt.Value() );
00448 
00449   if ( partner ) // Same shape with different location
00450   {
00451     // recursively associate all subshapes of theShape1 and theShape2
00452     typedef list< pair< TopoDS_Shape, TopoDS_Shape > > TShapePairsList;
00453     TShapePairsList shapesQueue( 1, make_pair( theShape1, theShape2 ));
00454     TShapePairsList::iterator s1_s2 = shapesQueue.begin();
00455     for ( ; s1_s2 != shapesQueue.end(); ++s1_s2 )
00456     {
00457       InsertAssociation( s1_s2->first, s1_s2->second, theMap, bidirect);
00458       TopoDS_Iterator s1It( s1_s2->first), s2It( s1_s2->second );
00459       for ( ; s1It.More(); s1It.Next(), s2It.Next() )
00460         shapesQueue.push_back( make_pair( s1It.Value(), s2It.Value() ));
00461     }
00462     return true;
00463   }
00464 
00465   if ( !theMap.IsEmpty() )
00466   {
00467     //======================================================================
00468     // 3) HAS initial vertex association
00469     //======================================================================
00470     switch ( theShape1.ShapeType() ) {
00471       // ----------------------------------------------------------------------
00472     case TopAbs_EDGE: { // TopAbs_EDGE
00473       // ----------------------------------------------------------------------
00474       if ( theMap.Extent() != 2 )
00475         RETURN_BAD_RESULT("Wrong map extent " << theMap.Extent() );
00476       TopoDS_Edge edge1 = TopoDS::Edge( theShape1 );
00477       TopoDS_Edge edge2 = TopoDS::Edge( theShape2 );
00478       if ( edge1.Orientation() >= TopAbs_INTERNAL ) edge1.Orientation( TopAbs_FORWARD );
00479       if ( edge2.Orientation() >= TopAbs_INTERNAL ) edge2.Orientation( TopAbs_FORWARD );
00480       TopoDS_Vertex VV1[2], VV2[2];
00481       TopExp::Vertices( edge1, VV1[0], VV1[1] );
00482       TopExp::Vertices( edge2, VV2[0], VV2[1] );
00483       int i1 = 0, i2 = 0;
00484       if ( theMap.IsBound( VV1[ i1 ] )) i1 = 1;
00485       if ( theMap.IsBound( VV2[ i2 ] )) i2 = 1;
00486       InsertAssociation( VV1[ i1 ], VV2[ i2 ], theMap, bidirect);
00487       InsertAssociation( theShape1, theShape2, theMap, bidirect );
00488       return true;
00489     }
00490       // ----------------------------------------------------------------------
00491     case TopAbs_FACE: { // TopAbs_FACE
00492       // ----------------------------------------------------------------------
00493       TopoDS_Face face1 = TopoDS::Face( theShape1 );
00494       TopoDS_Face face2 = TopoDS::Face( theShape2 );
00495       if ( face1.Orientation() >= TopAbs_INTERNAL ) face1.Orientation( TopAbs_FORWARD );
00496       if ( face2.Orientation() >= TopAbs_INTERNAL ) face2.Orientation( TopAbs_FORWARD );
00497 
00498       TopoDS_Vertex VV1[2], VV2[2];
00499       // find a not closed edge of face1 both vertices of which are associated
00500       int nbEdges = 0;
00501       TopExp_Explorer exp ( face1, TopAbs_EDGE );
00502       for ( ; VV2[ 1 ].IsNull() && exp.More(); exp.Next(), ++nbEdges ) {
00503         TopExp::Vertices( TopoDS::Edge( exp.Current() ), VV1[0], VV1[1] );
00504         if ( theMap.IsBound( VV1[0] ) ) {
00505           VV2[ 0 ] = TopoDS::Vertex( theMap( VV1[0] ));
00506           if ( theMap.IsBound( VV1[1] ) && !VV1[0].IsSame( VV1[1] ))
00507             VV2[ 1 ] = TopoDS::Vertex( theMap( VV1[1] ));
00508         }
00509       }
00510       if ( VV2[ 1 ].IsNull() ) { // 2 bound vertices not found
00511         if ( nbEdges > 1 ) {
00512           RETURN_BAD_RESULT("2 bound vertices not found" );
00513         } else {
00514           VV2[ 1 ] = VV2[ 0 ];
00515         }
00516       }
00517       list< TopoDS_Edge > edges1, edges2;
00518       int nbE = FindFaceAssociation( face1, VV1, face2, VV2, edges1, edges2 );
00519       if ( !nbE ) RETURN_BAD_RESULT("FindFaceAssociation() failed");
00520       FixAssocByPropagation( nbE, edges1, edges2, theMesh1, theMesh2 );
00521 
00522       list< TopoDS_Edge >::iterator eIt1 = edges1.begin();
00523       list< TopoDS_Edge >::iterator eIt2 = edges2.begin();
00524       for ( ; eIt1 != edges1.end(); ++eIt1, ++eIt2 )
00525       {
00526         InsertAssociation( *eIt1, *eIt2, theMap, bidirect);
00527         VV1[0] = TopExp::FirstVertex( *eIt1, true );
00528         VV2[0] = TopExp::FirstVertex( *eIt2, true );
00529         InsertAssociation( VV1[0], VV2[0], theMap, bidirect);
00530       }
00531       InsertAssociation( theShape1, theShape2, theMap, bidirect );
00532       return true;
00533     }
00534       // ----------------------------------------------------------------------
00535     case TopAbs_SHELL: // TopAbs_SHELL, TopAbs_SOLID
00536     case TopAbs_SOLID: {
00537       // ----------------------------------------------------------------------
00538       TopoDS_Vertex VV1[2], VV2[2];
00539       // try to find a not closed edge of shape1 both vertices of which are associated
00540       TopoDS_Edge edge1;
00541       TopExp_Explorer exp ( theShape1, TopAbs_EDGE );
00542       for ( ; VV2[ 1 ].IsNull() && exp.More(); exp.Next() ) {
00543         edge1 = TopoDS::Edge( exp.Current() );
00544         if ( edge1.Orientation() >= TopAbs_INTERNAL ) edge1.Orientation( TopAbs_FORWARD );
00545         TopExp::Vertices( edge1 , VV1[0], VV1[1] );
00546         if ( theMap.IsBound( VV1[0] )) {
00547           VV2[ 0 ] = TopoDS::Vertex( theMap( VV1[0] ));
00548           if ( theMap.IsBound( VV1[1] ) && !VV1[0].IsSame( VV1[1] ))
00549             VV2[ 1 ] = TopoDS::Vertex( theMap( VV1[1] ));
00550         }
00551       }
00552       if ( VV2[ 1 ].IsNull() ) // 2 bound vertices not found
00553         RETURN_BAD_RESULT("2 bound vertices not found" );
00554       // get an edge2 of theShape2 corresponding to edge1
00555       TopoDS_Edge edge2 = GetEdgeByVertices( theMesh2, VV2[ 0 ], VV2[ 1 ]);
00556       if ( edge2.IsNull() )
00557         RETURN_BAD_RESULT("GetEdgeByVertices() failed");
00558 
00559       // build map of edge to faces if shapes are not subshapes of main ones
00560       bool isSubOfMain = false;
00561       if ( SMESHDS_SubMesh * sm = theMesh1->GetMeshDS()->MeshElements( theShape1 ))
00562         isSubOfMain = !sm->IsComplexSubmesh();
00563       else
00564         isSubOfMain = theMesh1->GetMeshDS()->ShapeToIndex( theShape1 );
00565       TAncestorMap e2f1, e2f2;
00566       const TAncestorMap& edgeToFace1 = isSubOfMain ? theMesh1->GetAncestorMap() : e2f1;
00567       const TAncestorMap& edgeToFace2 = isSubOfMain ? theMesh2->GetAncestorMap() : e2f2;
00568       if (!isSubOfMain) {
00569         TopExp::MapShapesAndAncestors( theShape1, TopAbs_EDGE, TopAbs_FACE, e2f1 );
00570         TopExp::MapShapesAndAncestors( theShape2, TopAbs_EDGE, TopAbs_FACE, e2f2 );
00571         if ( !edgeToFace1.Contains( edge1 ))
00572           RETURN_BAD_RESULT("edge1 does not belong to theShape1");
00573         if ( !edgeToFace2.Contains( edge2 ))
00574           RETURN_BAD_RESULT("edge2 does not belong to theShape2");
00575       }
00576       //
00577       // Look for 2 corresponing faces:
00578       //
00579       TopoDS_Shape F1, F2;
00580 
00581       // get a face sharing edge1 (F1)
00582       TopoDS_Shape FF2[2];
00583       TopTools_ListIteratorOfListOfShape ancestIt1( edgeToFace1.FindFromKey( edge1 ));
00584       for ( ; F1.IsNull() && ancestIt1.More(); ancestIt1.Next() )
00585         if ( ancestIt1.Value().ShapeType() == TopAbs_FACE )
00586           F1 = ancestIt1.Value().Oriented( TopAbs_FORWARD );
00587       if ( F1.IsNull() )
00588         RETURN_BAD_RESULT(" Face1 not found");
00589 
00590       // get 2 faces sharing edge2 (one of them is F2)
00591       TopTools_ListIteratorOfListOfShape ancestIt2( edgeToFace2.FindFromKey( edge2 ));
00592       for ( int i = 0; FF2[1].IsNull() && ancestIt2.More(); ancestIt2.Next() )
00593         if ( ancestIt2.Value().ShapeType() == TopAbs_FACE )
00594           FF2[ i++ ] = ancestIt2.Value().Oriented( TopAbs_FORWARD );
00595 
00596       // get oriented edge1 and edge2 from F1 and FF2[0]
00597       for ( exp.Init( F1, TopAbs_EDGE ); exp.More(); exp.Next() )
00598         if ( edge1.IsSame( exp.Current() )) {
00599           edge1 = TopoDS::Edge( exp.Current() );
00600           break;
00601         }
00602       for ( exp.Init( FF2[ 0 ], TopAbs_EDGE ); exp.More(); exp.Next() )
00603         if ( edge2.IsSame( exp.Current() )) {
00604           edge2 = TopoDS::Edge( exp.Current() );
00605           break;
00606         }
00607 
00608       // compare first vertices of edge1 and edge2
00609       TopExp::Vertices( edge1, VV1[0], VV1[1], true );
00610       TopExp::Vertices( edge2, VV2[0], VV2[1], true );
00611       F2 = FF2[ 0 ]; // (F2 !)
00612       if ( !VV1[ 0 ].IsSame( theMap( VV2[ 0 ]))) {
00613         edge2.Reverse();
00614         if ( FF2[ 1 ].IsNull() )
00615           F2.Reverse();
00616         else
00617           F2 = FF2[ 1 ];
00618       }
00619 
00620       TopTools_MapOfShape boundEdges;
00621 
00622       // association of face subshapes and neighbour faces
00623       list< pair < TopoDS_Face, TopoDS_Edge > > FE1, FE2;
00624       list< pair < TopoDS_Face, TopoDS_Edge > >::iterator fe1, fe2;
00625       FE1.push_back( make_pair( TopoDS::Face( F1 ), edge1 ));
00626       FE2.push_back( make_pair( TopoDS::Face( F2 ), edge2 ));
00627       for ( fe1 = FE1.begin(), fe2 = FE2.begin(); fe1 != FE1.end(); ++fe1, ++fe2 )
00628       {
00629         const TopoDS_Face& face1 = fe1->first;
00630         if ( theMap.IsBound( face1 ) ) continue;
00631         const TopoDS_Face& face2 = fe2->first;
00632         edge1 = fe1->second;
00633         edge2 = fe2->second;
00634         TopExp::Vertices( edge1, VV1[0], VV1[1], true );
00635         TopExp::Vertices( edge2, VV2[0], VV2[1], true );
00636         list< TopoDS_Edge > edges1, edges2;
00637         int nbE = FindFaceAssociation( face1, VV1, face2, VV2, edges1, edges2 );
00638         if ( !nbE ) RETURN_BAD_RESULT("FindFaceAssociation() failed");
00639         InsertAssociation( face1, face2, theMap, bidirect); // assoc faces
00640         MESSAGE("Assoc FACE " << theMesh1->GetMeshDS()->ShapeToIndex( face1 )<<
00641                 " to "        << theMesh2->GetMeshDS()->ShapeToIndex( face2 ));
00642         if ( nbE == 2 && (edge1.IsSame( edges1.front())) != (edge2.IsSame( edges2.front())))
00643         {
00644           Reverse( edges2, nbE );
00645         }
00646         list< TopoDS_Edge >::iterator eIt1 = edges1.begin();
00647         list< TopoDS_Edge >::iterator eIt2 = edges2.begin();
00648         for ( ; eIt1 != edges1.end(); ++eIt1, ++eIt2 )
00649         {
00650           if ( !boundEdges.Add( *eIt1 )) continue; // already associated
00651           InsertAssociation( *eIt1, *eIt2, theMap, bidirect);  // assoc edges
00652           MESSAGE("Assoc edge " << theMesh1->GetMeshDS()->ShapeToIndex( *eIt1 )<<
00653                   " to "        << theMesh2->GetMeshDS()->ShapeToIndex( *eIt2 ));
00654           VV1[0] = TopExp::FirstVertex( *eIt1, true );
00655           VV2[0] = TopExp::FirstVertex( *eIt2, true );
00656           InsertAssociation( VV1[0], VV2[0], theMap, bidirect); // assoc vertices
00657           MESSAGE("Assoc vertex " << theMesh1->GetMeshDS()->ShapeToIndex( VV1[0] )<<
00658                   " to "          << theMesh2->GetMeshDS()->ShapeToIndex( VV2[0] ));
00659 
00660           // add adjacent faces to process
00661           TopoDS_Face nextFace1 = GetNextFace( edgeToFace1, *eIt1, face1 );
00662           TopoDS_Face nextFace2 = GetNextFace( edgeToFace2, *eIt2, face2 );
00663           if ( !nextFace1.IsNull() && !nextFace2.IsNull() ) {
00664             FE1.push_back( make_pair( nextFace1, *eIt1 ));
00665             FE2.push_back( make_pair( nextFace2, *eIt2 ));
00666           }
00667         }
00668       }
00669       InsertAssociation( theShape1, theShape2, theMap, bidirect );
00670       return true;
00671     }
00672       // ----------------------------------------------------------------------
00673     case TopAbs_COMPOUND: { // GROUP
00674       // ----------------------------------------------------------------------
00675       // Maybe groups contain only one member
00676       TopoDS_Iterator it1( theShape1 ), it2( theShape2 );
00677       TopAbs_ShapeEnum memberType = it1.Value().ShapeType();
00678       int nbMembers = Count( theShape1, memberType, true );
00679       if ( nbMembers == 0 ) return true;
00680       if ( nbMembers == 1 ) {
00681         return FindSubShapeAssociation( it1.Value(), theMesh1, it2.Value(), theMesh2, theMap );
00682       }
00683       // Try to make shells of faces
00684       //
00685       BRep_Builder builder;
00686       TopoDS_Shell shell1, shell2;
00687       builder.MakeShell(shell1); builder.MakeShell(shell2);
00688       if ( memberType == TopAbs_FACE ) {
00689         // just add faces of groups to shells
00690         for (; it1.More(); it1.Next(), it2.Next() )
00691           builder.Add( shell1, it1.Value() ), builder.Add( shell2, it2.Value() );
00692       }
00693       else if ( memberType == TopAbs_EDGE ) {
00694         // Try to add faces sharing more than one edge of a group or
00695         // sharing all its vertices with the group
00696         TopTools_IndexedMapOfShape groupVertices[2];
00697         TopExp::MapShapes( theShape1, TopAbs_VERTEX, groupVertices[0]);
00698         TopExp::MapShapes( theShape2, TopAbs_VERTEX, groupVertices[1]);
00699         //
00700         TopTools_MapOfShape groupEdges[2], addedFaces[2];
00701         bool hasInitAssoc = (!theMap.IsEmpty()), initAssocOK = !hasInitAssoc;
00702         for (; it1.More(); it1.Next(), it2.Next() ) {
00703           groupEdges[0].Add( it1.Value() );
00704           groupEdges[1].Add( it2.Value() );
00705           if ( !initAssocOK ) {
00706             // for shell association there must be an edge with both vertices bound
00707             TopoDS_Vertex v1, v2;
00708             TopExp::Vertices( TopoDS::Edge( it1.Value().Oriented(TopAbs_FORWARD)), v1, v2 );
00709             initAssocOK = ( theMap.IsBound( v1 ) && theMap.IsBound( v2 ));
00710           }
00711         }
00712         for (int is2ndGroup = 0; initAssocOK && is2ndGroup < 2; ++is2ndGroup) {
00713           const TopoDS_Shape& group = is2ndGroup ? theShape2: theShape1;
00714           SMESH_Mesh*         mesh  = is2ndGroup ? theMesh2 : theMesh1;
00715           TopoDS_Shell&       shell = is2ndGroup ? shell2   : shell1;
00716           for ( TopoDS_Iterator it( group ); it.More(); it.Next() ) {
00717             const TopoDS_Edge& edge = TopoDS::Edge( it.Value() );
00718             TopoDS_Face face;
00719             for ( int iF = 0; iF < 2; ++iF ) { // loop on 2 faces sharing edge
00720               face = GetNextFace(mesh->GetAncestorMap(), edge, face);
00721               if ( !face.IsNull() ) {
00722                 int nbGroupEdges = 0;
00723                 for ( TopExp_Explorer f( face, TopAbs_EDGE ); f.More(); f.Next())
00724                   if ( groupEdges[ is2ndGroup ].Contains( f.Current() ))
00725                     if ( ++nbGroupEdges > 1 )
00726                       break;
00727                 bool add = (nbGroupEdges > 1 || Count( face, TopAbs_EDGE, true ) == 1 );
00728                 if ( !add ) {
00729                   add = true;
00730                   for ( TopExp_Explorer v( face, TopAbs_VERTEX ); add && v.More(); v.Next())
00731                     add = groupVertices[ is2ndGroup ].Contains( v.Current() );
00732                 }
00733                 if ( add && addedFaces[ is2ndGroup ].Add( face ))
00734                   builder.Add( shell, face );
00735               }
00736             }
00737           }
00738         }
00739       } else {
00740         RETURN_BAD_RESULT("Unexpected group type");
00741       }
00742       // Associate shells
00743       //
00744       int nbFaces1 = Count( shell1, TopAbs_FACE, 0 );
00745       int nbFaces2 = Count( shell2, TopAbs_FACE, 0 );
00746       if ( nbFaces1 != nbFaces2 )
00747         RETURN_BAD_RESULT("Different nb of faces found for shells");
00748       if ( nbFaces1 > 0 ) {
00749         bool ok = false;
00750         if ( nbFaces1 == 1 ) {
00751           TopoDS_Shape F1 = TopoDS_Iterator( shell1 ).Value();
00752           TopoDS_Shape F2 = TopoDS_Iterator( shell2 ).Value();
00753           ok = FindSubShapeAssociation( F1, theMesh1, F2, theMesh2, theMap );
00754         }
00755         else {
00756           ok = FindSubShapeAssociation(shell1, theMesh1, shell2, theMesh2, theMap );
00757         }
00758         // Check if all members are mapped 
00759         if ( ok ) {
00760           TopTools_MapOfShape boundMembers[2];
00761           TopoDS_Iterator mIt;
00762           for ( mIt.Initialize( theShape1 ); mIt.More(); mIt.Next())
00763             if ( theMap.IsBound( mIt.Value() )) {
00764               boundMembers[0].Add( mIt.Value() );
00765               boundMembers[1].Add( theMap( mIt.Value() ));
00766             }
00767           if ( boundMembers[0].Extent() != nbMembers ) {
00768             // make compounds of not bound members
00769             TopoDS_Compound comp[2];
00770             for ( int is2ndGroup = 0; is2ndGroup < 2; ++is2ndGroup ) {
00771               builder.MakeCompound( comp[is2ndGroup] );
00772               for ( mIt.Initialize( is2ndGroup ? theShape2:theShape1 ); mIt.More(); mIt.Next())
00773                 if ( ! boundMembers[ is2ndGroup ].Contains( mIt.Value() ))
00774                   builder.Add( comp[ is2ndGroup ], mIt.Value() );
00775             }
00776             // check if theMap contains initial association for the comp's
00777             bool hasInitialAssoc = false;
00778             if ( memberType == TopAbs_EDGE ) {
00779               for ( TopExp_Explorer v( comp[0], TopAbs_VERTEX ); v.More(); v.Next())
00780                 if ( theMap.IsBound( v.Current() )) {
00781                   hasInitialAssoc = true;
00782                   break;
00783                 }
00784             }
00785             if ( hasInitialAssoc == bool( !theMap.IsEmpty() ))
00786               ok = FindSubShapeAssociation( comp[0], theMesh1, comp[1], theMesh2, theMap );
00787             else {
00788               TShapeShapeMap tmpMap;
00789               ok = FindSubShapeAssociation( comp[0], theMesh1, comp[1], theMesh2, tmpMap );
00790               if ( ok ) {
00791                 TopTools_DataMapIteratorOfDataMapOfShapeShape mapIt( tmpMap );
00792                 for ( ; mapIt.More(); mapIt.Next() )
00793                   theMap.Bind( mapIt.Key(), mapIt.Value());
00794               }
00795             }
00796           }
00797         }
00798         return ok;
00799       }
00800       // Each edge of an edge group is shared by own faces
00801       // ------------------------------------------------------------------
00802       //
00803       // map vertices to edges sharing them, avoid doubling edges in lists
00804       TopTools_DataMapOfShapeListOfShape v2e[2];
00805       for (int isFirst = 0; isFirst < 2; ++isFirst ) {
00806         const TopoDS_Shape& group = isFirst ? theShape1 : theShape2;
00807         TopTools_DataMapOfShapeListOfShape& veMap = v2e[ isFirst ? 0 : 1 ];
00808         TopTools_MapOfShape addedEdges;
00809         for ( TopExp_Explorer e( group, TopAbs_EDGE ); e.More(); e.Next() ) {
00810           const TopoDS_Shape& edge = e.Current();
00811           if ( addedEdges.Add( edge )) {
00812             for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next()) {
00813               const TopoDS_Shape& vertex = v.Current();
00814               if ( !veMap.IsBound( vertex )) {
00815                 TopTools_ListOfShape l;
00816                 veMap.Bind( vertex, l );
00817               }
00818               veMap( vertex ).Append( edge );
00819             }
00820           }
00821         }   
00822       }
00823       while ( !v2e[0].IsEmpty() )
00824       {
00825         // find a bound vertex
00826         TopoDS_Vertex V[2];
00827         TopTools_DataMapIteratorOfDataMapOfShapeListOfShape v2eIt( v2e[0] );
00828         for ( ; v2eIt.More(); v2eIt.Next())
00829           if ( theMap.IsBound( v2eIt.Key() )) {
00830             V[0] = TopoDS::Vertex( v2eIt.Key() );
00831             V[1] = TopoDS::Vertex( theMap( V[0] ));
00832             break;
00833           }
00834         if ( V[0].IsNull() )
00835           RETURN_BAD_RESULT("No more bound vertices");
00836 
00837         while ( !V[0].IsNull() && v2e[0].IsBound( V[0] )) {
00838           TopTools_ListOfShape& edges0 = v2e[0]( V[0] );
00839           TopTools_ListOfShape& edges1 = v2e[1]( V[1] );
00840           int nbE0 = edges0.Extent(), nbE1 = edges1.Extent();
00841           if ( nbE0 != nbE1 )
00842             RETURN_BAD_RESULT("Different nb of edges: "<< nbE0 << " != " << nbE1);
00843 
00844           if ( nbE0 == 1 )
00845           {
00846             TopoDS_Edge e0 = TopoDS::Edge( edges0.First() );
00847             TopoDS_Edge e1 = TopoDS::Edge( edges1.First() );
00848             v2e[0].UnBind( V[0] );
00849             v2e[1].UnBind( V[1] );
00850             InsertAssociation( e0, e1, theMap, bidirect );
00851             MESSAGE("Assoc edge " << theMesh1->GetMeshDS()->ShapeToIndex( e0 )<<
00852                     " to "        << theMesh2->GetMeshDS()->ShapeToIndex( e1 ));
00853             V[0] = GetNextVertex( e0, V[0] );
00854             V[1] = GetNextVertex( e1, V[1] );
00855             if ( !V[0].IsNull() ) {
00856               InsertAssociation( V[0], V[1], theMap, bidirect );
00857               MESSAGE("Assoc vertex " << theMesh1->GetMeshDS()->ShapeToIndex( V[0] )<<
00858                       " to "          << theMesh2->GetMeshDS()->ShapeToIndex( V[1] ));
00859             }
00860           }
00861           else if ( nbE0 == 2 )
00862           {
00863             // one of edges must have both ends bound
00864             TopoDS_Vertex v0e0 = GetNextVertex( TopoDS::Edge( edges0.First() ), V[0] );
00865             TopoDS_Vertex v1e0 = GetNextVertex( TopoDS::Edge( edges0.Last() ),  V[0] );
00866             TopoDS_Vertex v0e1 = GetNextVertex( TopoDS::Edge( edges1.First() ), V[1] );
00867             TopoDS_Vertex v1e1 = GetNextVertex( TopoDS::Edge( edges1.Last() ),  V[1] );
00868             TopoDS_Shape e0b, e1b, e0n, e1n, v1b; // bound and not-bound
00869             TopoDS_Vertex v0n, v1n;
00870             if ( theMap.IsBound( v0e0 )) {
00871               v0n = v1e0; e0b = edges0.First(); e0n = edges0.Last(); v1b = theMap( v0e0 );
00872             } else if ( theMap.IsBound( v1e0 )) {
00873               v0n = v0e0; e0n = edges0.First(); e0b = edges0.Last(); v1b = theMap( v1e0 );
00874             } else {
00875               RETURN_BAD_RESULT("None of vertices bound");
00876             }
00877             if ( v1b.IsSame( v1e1 )) {
00878               v1n = v0e1; e1n = edges1.First(); e1b = edges1.Last();
00879             } else {
00880               v1n = v1e1; e1b = edges1.First(); e1n = edges1.Last();
00881             }
00882             InsertAssociation( e0b, e1b, theMap, bidirect );
00883             InsertAssociation( e0n, e1n, theMap, bidirect );
00884             InsertAssociation( v0n, v1n, theMap, bidirect );
00885             MESSAGE("Assoc edge " << theMesh1->GetMeshDS()->ShapeToIndex( e0b )<<
00886                     " to "        << theMesh2->GetMeshDS()->ShapeToIndex( e1b ));
00887             MESSAGE("Assoc edge " << theMesh1->GetMeshDS()->ShapeToIndex( e0n )<<
00888                     " to "        << theMesh2->GetMeshDS()->ShapeToIndex( e1n ));
00889             MESSAGE("Assoc vertex " << theMesh1->GetMeshDS()->ShapeToIndex( v0n )<<
00890                     " to "          << theMesh2->GetMeshDS()->ShapeToIndex( v1n ));
00891             v2e[0].UnBind( V[0] );
00892             v2e[1].UnBind( V[1] );
00893             V[0] = v0n;
00894             V[1] = v1n;
00895           }
00896           else {
00897             RETURN_BAD_RESULT("Not implemented");
00898           }
00899         }
00900       } //while ( !v2e[0].IsEmpty() )
00901       return true;
00902     }
00903 
00904     default:
00905       RETURN_BAD_RESULT("Unexpected shape type");
00906 
00907     } // end switch by shape type
00908   } // end case of available initial vertex association
00909 
00910   //======================================================================
00911   // 4) NO INITIAL VERTEX ASSOCIATION
00912   //======================================================================
00913 
00914   switch ( theShape1.ShapeType() ) {
00915 
00916   case TopAbs_EDGE: {
00917     // ----------------------------------------------------------------------
00918     TopoDS_Edge edge1 = TopoDS::Edge( theShape1 );
00919     TopoDS_Edge edge2 = TopoDS::Edge( theShape2 );
00920     if ( IsPropagationPossible( theMesh1, theMesh2 ))
00921     {
00922       TopoDS_Edge prpEdge = GetPropagationEdge( theMesh1, edge2, edge1 ).second;
00923       if ( !prpEdge.IsNull() )
00924       {
00925         TopoDS_Vertex VV1[2], VV2[2];
00926         TopExp::Vertices( edge1,   VV1[0], VV1[1], true );
00927         TopExp::Vertices( prpEdge, VV2[0], VV2[1], true );
00928         InsertAssociation( VV1[ 0 ], VV2[ 0 ], theMap, bidirect);
00929         InsertAssociation( VV1[ 1 ], VV2[ 1 ], theMap, bidirect);
00930         if ( VV1[0].IsSame( VV1[1] ) || // one of edges is closed
00931              VV2[0].IsSame( VV2[1] ) )
00932         {
00933           InsertAssociation( edge1, prpEdge, theMap, bidirect); // insert with a proper orientation
00934         }
00935         InsertAssociation( theShape1, theShape2, theMap, bidirect );
00936         return true; // done
00937       }
00938     }
00939     if ( SMESH_MesherHelper::IsClosedEdge( edge1 ) &&
00940          SMESH_MesherHelper::IsClosedEdge( edge2 ))
00941     {
00942       // TODO: find out a proper orientation (is it possible?)
00943       InsertAssociation( edge1, edge2, theMap, bidirect); // insert with a proper orientation
00944       InsertAssociation( TopExp::FirstVertex(edge1), TopExp::FirstVertex(edge2),
00945                          theMap, bidirect);
00946       InsertAssociation( theShape1, theShape2, theMap, bidirect );
00947       return true; // done
00948     }
00949     break; // try by vertex closeness
00950   }
00951 
00952   case TopAbs_FACE: {
00953     // ----------------------------------------------------------------------
00954     if ( IsPropagationPossible( theMesh1, theMesh2 )) // try by propagation in one mesh
00955     {
00956       TopoDS_Face face1 = TopoDS::Face(theShape1);
00957       TopoDS_Face face2 = TopoDS::Face(theShape2);
00958       if ( face1.Orientation() >= TopAbs_INTERNAL ) face1.Orientation( TopAbs_FORWARD );
00959       if ( face2.Orientation() >= TopAbs_INTERNAL ) face2.Orientation( TopAbs_FORWARD );
00960       TopoDS_Edge edge1, edge2;
00961       // get outer edge of theShape1
00962       edge1 = TopoDS::Edge( OuterShape( face1, TopAbs_EDGE ));
00963       // find out if any edge of face2 is a propagation edge of outer edge1
00964       map<int,TopoDS_Edge> propag_edges; // use map to find the closest propagation edge
00965       for ( TopExp_Explorer exp( face2, TopAbs_EDGE ); exp.More(); exp.Next() ) {
00966         edge2 = TopoDS::Edge( exp.Current() );
00967         pair<int,TopoDS_Edge> step_edge = GetPropagationEdge( theMesh1, edge2, edge1 );
00968         if ( !step_edge.second.IsNull() ) { // propagation found
00969           propag_edges.insert( step_edge );
00970           if ( step_edge.first == 1 ) break; // most close found
00971         }
00972       }
00973       if ( !propag_edges.empty() ) // propagation found
00974       {
00975         edge2 = propag_edges.begin()->second;
00976         TopoDS_Vertex VV1[2], VV2[2];
00977         TopExp::Vertices( edge1, VV1[0], VV1[1], true );
00978         TopExp::Vertices( edge2, VV2[0], VV2[1], true );
00979         list< TopoDS_Edge > edges1, edges2;
00980         int nbE = FindFaceAssociation( face1, VV1, face2, VV2, edges1, edges2 );
00981         if ( !nbE ) RETURN_BAD_RESULT("FindFaceAssociation() failed");
00982         // take care of proper association of propagated edges
00983         bool same1 = edge1.IsSame( edges1.front() );
00984         bool same2 = edge2.IsSame( edges2.front() );
00985         if ( same1 != same2 )
00986         {
00987           Reverse(edges2, nbE);
00988           if ( nbE != 2 ) // 2 degen edges of 4 (issue 0021144)
00989             edges2.splice( edges2.end(), edges2, edges2.begin());
00990         }
00991         // store association
00992         list< TopoDS_Edge >::iterator eIt1 = edges1.begin();
00993         list< TopoDS_Edge >::iterator eIt2 = edges2.begin();
00994         for ( ; eIt1 != edges1.end(); ++eIt1, ++eIt2 )
00995         {
00996           InsertAssociation( *eIt1, *eIt2, theMap, bidirect);
00997           VV1[0] = TopExp::FirstVertex( *eIt1, true );
00998           VV2[0] = TopExp::FirstVertex( *eIt2, true );
00999           InsertAssociation( VV1[0], VV2[0], theMap, bidirect);
01000         }
01001         InsertAssociation( theShape1, theShape2, theMap, bidirect );
01002         return true;
01003       }
01004     }
01005     break; // try by vertex closeness
01006   }
01007   case TopAbs_COMPOUND: {
01008     // ----------------------------------------------------------------------
01009     if ( IsPropagationPossible( theMesh1, theMesh2 )) {
01010 
01011       // try to accosiate all using propagation
01012       if ( AssocGroupsByPropagation( theShape1, theShape2, *theMesh1, theMap ))
01013         return true;
01014 
01015       // find a boundary edge for theShape1
01016       TopoDS_Edge E;
01017       for(TopExp_Explorer exp(theShape1, TopAbs_EDGE); exp.More(); exp.Next() ) {
01018         E = TopoDS::Edge( exp.Current() );
01019         if ( IsBoundaryEdge( E, theShape1, *theMesh1 ))
01020           break;
01021         else
01022           E.Nullify();
01023       }
01024       if ( E.IsNull() )
01025         break; // try by vertex closeness
01026 
01027       // find association for vertices of edge E
01028       TopoDS_Vertex VV1[2], VV2[2];
01029       for(TopExp_Explorer eexp(E, TopAbs_VERTEX); eexp.More(); eexp.Next()) {
01030         TopoDS_Vertex V1 = TopoDS::Vertex( eexp.Current() );
01031         // look for an edge ending in E whose one vertex is in theShape1
01032         // and the other, in theShape2
01033         const TopTools_ListOfShape& Ancestors = theMesh1->GetAncestors(V1);
01034         TopTools_ListIteratorOfListOfShape ita(Ancestors);
01035         for(; ita.More(); ita.Next()) {
01036           if( ita.Value().ShapeType() != TopAbs_EDGE ) continue;
01037           TopoDS_Edge edge = TopoDS::Edge(ita.Value());
01038           bool FromShape1 = false;
01039           for(TopExp_Explorer expe(theShape1, TopAbs_EDGE); expe.More(); expe.Next() ) {
01040             if(edge.IsSame(expe.Current())) {
01041               FromShape1 = true;
01042               break;
01043             }
01044           }
01045           if(!FromShape1) {
01046             // is it an edge between theShape1 and theShape2?
01047             TopExp_Explorer expv(edge, TopAbs_VERTEX);
01048             TopoDS_Vertex V2 = TopoDS::Vertex( expv.Current() );
01049             if(V2.IsSame(V1)) {
01050               expv.Next();
01051               V2 = TopoDS::Vertex( expv.Current() );
01052             }
01053             bool FromShape2 = false;
01054             for ( expv.Init( theShape2, TopAbs_VERTEX ); expv.More(); expv.Next()) {
01055               if ( V2.IsSame( expv.Current() )) {
01056                 FromShape2 = true;
01057                 break;
01058               }
01059             }
01060             if ( FromShape2 ) {
01061               if ( VV1[0].IsNull() )
01062                 VV1[0] = V1, VV2[0] = V2;
01063               else
01064                 VV1[1] = V1, VV2[1] = V2;
01065               break; // from loop on ancestors of V1
01066             }
01067           }
01068         }
01069       }
01070       if ( !VV1[1].IsNull() ) {
01071         InsertAssociation( VV1[0], VV2[0], theMap, bidirect);
01072         InsertAssociation( VV1[1], VV2[1], theMap, bidirect);
01073         return FindSubShapeAssociation( theShape1, theMesh1, theShape2, theMesh2, theMap);
01074       }
01075     }
01076     break; // try by vertex closeness
01077   }
01078   default:;
01079   }
01080 
01081   // Find association by closeness of vertices
01082   // ------------------------------------------
01083 
01084   TopTools_IndexedMapOfShape vMap1, vMap2;
01085   TopExp::MapShapes( theShape1, TopAbs_VERTEX, vMap1 );
01086   TopExp::MapShapes( theShape2, TopAbs_VERTEX, vMap2 );
01087   TopoDS_Vertex VV1[2], VV2[2];
01088 
01089   if ( vMap1.Extent() != vMap2.Extent() )
01090     RETURN_BAD_RESULT("Different nb of vertices");
01091 
01092   if ( vMap1.Extent() == 1 ) {
01093     InsertAssociation( vMap1(1), vMap2(1), theMap, bidirect);
01094     if ( theShape1.ShapeType() == TopAbs_EDGE ) {
01095       InsertAssociation( theShape1, theShape2, theMap, bidirect );
01096       return true;
01097     }
01098     return FindSubShapeAssociation( theShape1, theMesh1, theShape2, theMesh2, theMap);
01099   }
01100 
01101   // Try to associate by common vertices of an edge
01102   for ( int i = 1; i <= vMap1.Extent(); ++i )
01103   {
01104     const TopoDS_Shape& v1 = vMap1(i);
01105     if ( vMap2.Contains( v1 ))
01106     {
01107       // find an egde sharing v1 and sharing at the same time another common vertex
01108       PShapeIteratorPtr edgeIt = SMESH_MesherHelper::GetAncestors( v1, *theMesh1, TopAbs_EDGE);
01109       bool edgeFound = false;
01110       while ( edgeIt->more() && !edgeFound )
01111       {
01112         TopoDS_Edge edge = TopoDS::Edge( edgeIt->next()->Oriented(TopAbs_FORWARD));
01113         TopExp::Vertices(edge, VV1[0], VV1[1]);
01114         if ( !VV1[0].IsSame( VV1[1] ))
01115           edgeFound = ( vMap2.Contains( VV1[ v1.IsSame(VV1[0]) ? 1:0]));
01116       }
01117       if ( edgeFound )
01118       {
01119         InsertAssociation( VV1[0], VV1[0], theMap, bidirect );
01120         InsertAssociation( VV1[1], VV1[1], theMap, bidirect );
01121         if (FindSubShapeAssociation( theShape1, theMesh1, theShape2, theMesh2, theMap ))
01122           return true;
01123       }
01124     }
01125   }
01126 
01127   // Find transformation to make the shapes be of similar size at same location
01128 
01129   Bnd_Box box[2];
01130   for ( int i = 1; i <= vMap1.Extent(); ++i ) {
01131     box[ 0 ].Add( BRep_Tool::Pnt ( TopoDS::Vertex( vMap1( i ))));
01132     box[ 1 ].Add( BRep_Tool::Pnt ( TopoDS::Vertex( vMap2( i ))));
01133   }
01134 
01135   gp_Pnt gc[2]; // box center
01136   double x0,y0,z0, x1,y1,z1;
01137   box[0].Get( x0,y0,z0, x1,y1,z1 );
01138   gc[0] = 0.5 * ( gp_XYZ( x0,y0,z0 ) + gp_XYZ( x1,y1,z1 ));
01139   box[1].Get( x0,y0,z0, x1,y1,z1 );
01140   gc[1] = 0.5 * ( gp_XYZ( x0,y0,z0 ) + gp_XYZ( x1,y1,z1 ));
01141 
01142   // 1 -> 2
01143   gp_Vec vec01( gc[0], gc[1] );
01144   double scale = sqrt( box[1].SquareExtent() / box[0].SquareExtent() );
01145 
01146   // Find 2 closest vertices
01147 
01148   // get 2 linked vertices of shape 1 not belonging to an inner wire of a face
01149   TopoDS_Shape edge = theShape1;
01150   TopExp_Explorer expF( theShape1, TopAbs_FACE ), expE;
01151   if ( expF.More() ) {
01152     for ( ; expF.More(); expF.Next() ) {
01153       edge.Nullify();
01154       TopoDS_Shape wire = OuterShape( TopoDS::Face( expF.Current() ), TopAbs_WIRE );
01155       for ( expE.Init( wire, TopAbs_EDGE ); edge.IsNull() && expE.More(); expE.Next() )
01156         if ( !SMESH_MesherHelper::IsClosedEdge( TopoDS::Edge( expE.Current() )))
01157           edge = expE.Current();
01158       if ( !edge.IsNull() )
01159         break;
01160     }
01161   } else if (edge.ShapeType() != TopAbs_EDGE) { // no faces
01162     edge.Nullify();
01163     for ( expE.Init( theShape1, TopAbs_EDGE ); edge.IsNull() && expE.More(); expE.Next() )
01164       if ( !SMESH_MesherHelper::IsClosedEdge( TopoDS::Edge( expE.Current() )))
01165         edge = expE.Current();
01166   }
01167   if ( edge.IsNull() || edge.ShapeType() != TopAbs_EDGE )
01168     RETURN_BAD_RESULT("Edge not found");
01169 
01170   TopExp::Vertices( TopoDS::Edge( edge.Oriented(TopAbs_FORWARD)), VV1[0], VV1[1]);
01171   if ( VV1[0].IsSame( VV1[1] ))
01172     RETURN_BAD_RESULT("Only closed edges");
01173 
01174   // find vertices closest to 2 linked vertices of shape 1
01175   for ( int i1 = 0; i1 < 2; ++i1 )
01176   {
01177     double dist2 = DBL_MAX;
01178     gp_Pnt p1 = BRep_Tool::Pnt( VV1[ i1 ]);
01179     p1.Translate( vec01 );
01180     p1.Scale( gc[1], scale );
01181     for ( int i2 = 1; i2 <= vMap2.Extent(); ++i2 )
01182     {
01183       TopoDS_Vertex V2 = TopoDS::Vertex( vMap2( i2 ));
01184       gp_Pnt p2 = BRep_Tool::Pnt ( V2 );
01185       double d2 = p1.SquareDistance( p2 );
01186       if ( d2 < dist2 && !V2.IsSame( VV2[ 0 ])) {
01187         VV2[ i1 ] = V2; dist2 = d2;
01188       }
01189     }
01190   }
01191 
01192   InsertAssociation( VV1[ 0 ], VV2[ 0 ], theMap, bidirect);
01193   InsertAssociation( VV1[ 1 ], VV2[ 1 ], theMap, bidirect);
01194   MESSAGE("Initial assoc VERT " << theMesh1->GetMeshDS()->ShapeToIndex( VV1[ 0 ] )<<
01195           " to "                << theMesh2->GetMeshDS()->ShapeToIndex( VV2[ 0 ] )<<
01196           "\nand         VERT " << theMesh1->GetMeshDS()->ShapeToIndex( VV1[ 1 ] )<<
01197           " to "                << theMesh2->GetMeshDS()->ShapeToIndex( VV2[ 1 ] ));
01198   if ( theShape1.ShapeType() == TopAbs_EDGE ) {
01199     InsertAssociation( theShape1, theShape2, theMap, bidirect );
01200     return true;
01201   }
01202 
01203   return FindSubShapeAssociation( theShape1, theMesh1, theShape2, theMesh2, theMap );
01204 }
01205 
01206 //================================================================================
01217 //================================================================================
01218 
01219 int StdMeshers_ProjectionUtils::FindFaceAssociation(const TopoDS_Face&    face1,
01220                                                     TopoDS_Vertex         VV1[2],
01221                                                     const TopoDS_Face&    face2,
01222                                                     TopoDS_Vertex         VV2[2],
01223                                                     list< TopoDS_Edge > & edges1,
01224                                                     list< TopoDS_Edge > & edges2)
01225 {
01226   bool OK = false;
01227   list< int > nbEInW1, nbEInW2;
01228   int i_ok_wire_algo = -1;
01229   for ( int outer_wire_algo = 0; outer_wire_algo < 2 && !OK; ++outer_wire_algo )
01230   {
01231     edges1.clear();
01232     edges2.clear();
01233 
01234     if ( SMESH_Block::GetOrderedEdges( face1, VV1[0], edges1, nbEInW1, outer_wire_algo) !=
01235          SMESH_Block::GetOrderedEdges( face2, VV2[0], edges2, nbEInW2, outer_wire_algo) )
01236       CONT_BAD_RESULT("Different number of wires in faces ");
01237 
01238     if ( nbEInW1 != nbEInW2 )
01239       CONT_BAD_RESULT("Different number of edges in faces: " <<
01240                       nbEInW1.front() << " != " << nbEInW2.front());
01241 
01242     i_ok_wire_algo = outer_wire_algo;
01243 
01244     // Define if we need to reverse one of wires to make edges in lists match each other
01245 
01246     bool reverse = false;
01247 
01248     list< TopoDS_Edge >::iterator edgeIt;
01249     if ( !VV1[1].IsSame( TopExp::LastVertex( edges1.front(), true ))) {
01250       reverse = true;
01251       edgeIt = --edges1.end();
01252       // check if the second vertex belongs to the first or last edge in the wire
01253       if ( !VV1[1].IsSame( TopExp::FirstVertex( *edgeIt, true ))) {
01254         bool KO = true; // belongs to none
01255         if ( nbEInW1.size() > 1 ) { // several wires
01256           edgeIt = edges1.begin();
01257           std::advance( edgeIt, nbEInW1.front()-1 );
01258           KO = !VV1[1].IsSame( TopExp::FirstVertex( *edgeIt, true ));
01259         }
01260         if ( KO )
01261           CONT_BAD_RESULT("GetOrderedEdges() failed");
01262       }
01263     }
01264     edgeIt = --edges2.end();
01265     if ( !VV2[1].IsSame( TopExp::LastVertex( edges2.front(), true ))) {
01266       reverse = !reverse;
01267       // check if the second vertex belongs to the first or last edge in the wire
01268       if ( !VV2[1].IsSame( TopExp::FirstVertex( *edgeIt, true ))) {
01269         bool KO = true; // belongs to none
01270         if ( nbEInW2.size() > 1 ) { // several wires
01271           edgeIt = edges2.begin();
01272           std::advance( edgeIt, nbEInW2.front()-1 );
01273           KO = !VV2[1].IsSame( TopExp::FirstVertex( *edgeIt, true ));
01274         }
01275         if ( KO )
01276           CONT_BAD_RESULT("GetOrderedEdges() failed");
01277       }
01278     }
01279     if ( reverse )
01280     {
01281       Reverse( edges2 , nbEInW2.front());
01282       if (( VV1[1].IsSame( TopExp::LastVertex( edges1.front(), true ))) !=
01283           ( VV2[1].IsSame( TopExp::LastVertex( edges2.front(), true ))))
01284         CONT_BAD_RESULT("GetOrderedEdges() failed");
01285     }
01286     OK = true;
01287 
01288   } // loop algos getting an outer wire
01289   
01290   // Try to orient all (if !OK) or only internal wires (issue 0020996) by UV similarity
01291   if (( !OK || nbEInW1.size() > 1 ) && i_ok_wire_algo > -1 )
01292   {
01293     // Check that Vec(VV1[0],VV1[1]) in 2D on face1 is the same
01294     // as Vec(VV2[0],VV2[1]) on face2
01295     double vTol = BRep_Tool::Tolerance( VV1[0] );
01296     BRepAdaptor_Surface surface1( face1, false );
01297     double vTolUV =
01298       surface1.UResolution( vTol ) + surface1.VResolution( vTol ); // let's be tolerant
01299     gp_Pnt2d v0f1UV = BRep_Tool::Parameters( VV1[0], face1 );
01300     gp_Pnt2d v0f2UV = BRep_Tool::Parameters( VV2[0], face2 );
01301     gp_Pnt2d v1f1UV = BRep_Tool::Parameters( VV1[1], face1 );
01302     gp_Pnt2d v1f2UV = BRep_Tool::Parameters( VV2[1], face2 );
01303     gp_Vec2d v01f1Vec( v0f1UV, v1f1UV );
01304     gp_Vec2d v01f2Vec( v0f2UV, v1f2UV );
01305     if ( Abs( v01f1Vec.X()-v01f2Vec.X()) < vTolUV &&
01306          Abs( v01f1Vec.Y()-v01f2Vec.Y()) < vTolUV )
01307     {
01308       if ( !OK /*i_ok_wire_algo != 1*/ )
01309       {
01310         edges1.clear();
01311         edges2.clear();
01312         SMESH_Block::GetOrderedEdges( face1, VV1[0], edges1, nbEInW1, i_ok_wire_algo);
01313         SMESH_Block::GetOrderedEdges( face2, VV2[0], edges2, nbEInW2, i_ok_wire_algo);
01314       }
01315       gp_XY dUV = v0f2UV.XY() - v0f1UV.XY(); // UV shift between 2 faces
01316       // skip edges of the outer wire (if the outer wire is OK)
01317       list< int >::iterator nbEInW = nbEInW1.begin();
01318       list< TopoDS_Edge >::iterator edge1Beg = edges1.begin(), edge2Beg = edges2.begin();
01319       if ( OK )
01320       {
01321         for ( int i = 0; i < *nbEInW; ++i )
01322           ++edge1Beg, ++edge2Beg;
01323         ++nbEInW;
01324       }
01325       for ( ; nbEInW != nbEInW1.end(); ++nbEInW ) // loop on wires
01326       {
01327         // reach an end of edges of a current wire
01328         list< TopoDS_Edge >::iterator edge1End = edge1Beg, edge2End = edge2Beg;
01329         for ( int i = 0; i < *nbEInW; ++i )
01330           ++edge1End, ++edge2End;
01331         // rotate edges2 untill coincident with edges1 in 2D
01332         v0f1UV = BRep_Tool::Parameters( TopExp::FirstVertex(*edge1Beg,true), face1 );
01333         v1f1UV = BRep_Tool::Parameters( TopExp::LastVertex (*edge1Beg,true), face1 );
01334         v0f1UV.ChangeCoord() += dUV;
01335         v1f1UV.ChangeCoord() += dUV;
01336         int i = *nbEInW;
01337         while ( --i > 0 && !sameVertexUV( *edge2Beg, face2, 0, v0f1UV, vTolUV ))
01338           edges2.splice( edge2End, edges2, edge2Beg++ ); // move edge2Beg to place before edge2End
01339         if ( sameVertexUV( *edge2Beg, face2, 0, v0f1UV, vTolUV ))
01340         {
01341           if ( nbEInW == nbEInW1.begin() )
01342             OK = true; // OK is for the first wire
01343           // reverse edges2 if needed
01344           if ( !sameVertexUV( *edge2Beg, face2, 1, v1f1UV, vTolUV ))
01345           {
01346             Reverse( edges2 , *nbEInW, distance( edges2.begin(),edge2Beg ));
01347             // set correct edge2End
01348             edge2End = edges2.begin();
01349             std::advance( edge2End, std::accumulate( nbEInW1.begin(), nbEInW, *nbEInW));
01350           }
01351         }
01352         // prepare to the next wire loop
01353         edge1Beg = edge1End, edge2Beg = edge2End;
01354       }
01355     }
01356   }
01357 
01358   return OK ? nbEInW1.front() : 0;
01359 }
01360 
01361 //=======================================================================
01362 //function : InitVertexAssociation
01363 //purpose  : 
01364 //=======================================================================
01365 
01366 void StdMeshers_ProjectionUtils::InitVertexAssociation( const SMESH_Hypothesis* theHyp,
01367                                                         TShapeShapeMap &        theAssociationMap,
01368                                                         const TopoDS_Shape&     theTargetShape)
01369 {
01370   string hypName = theHyp->GetName();
01371   if ( hypName == "ProjectionSource1D" ) {
01372     const StdMeshers_ProjectionSource1D * hyp =
01373       static_cast<const StdMeshers_ProjectionSource1D*>( theHyp );
01374     if ( hyp->HasVertexAssociation() )
01375       InsertAssociation( hyp->GetSourceVertex(),hyp->GetTargetVertex(),theAssociationMap);
01376   }
01377   else if ( hypName == "ProjectionSource2D" ) {
01378     const StdMeshers_ProjectionSource2D * hyp =
01379       static_cast<const StdMeshers_ProjectionSource2D*>( theHyp );
01380     if ( hyp->HasVertexAssociation() ) {
01381       InsertAssociation( hyp->GetSourceVertex(1),hyp->GetTargetVertex(1),theAssociationMap);
01382       InsertAssociation( hyp->GetSourceVertex(2),hyp->GetTargetVertex(2),theAssociationMap);
01383     }
01384   }
01385   else if ( hypName == "ProjectionSource3D" ) {
01386     const StdMeshers_ProjectionSource3D * hyp =
01387       static_cast<const StdMeshers_ProjectionSource3D*>( theHyp );
01388     if ( hyp->HasVertexAssociation() ) {
01389       InsertAssociation( hyp->GetSourceVertex(1),hyp->GetTargetVertex(1),theAssociationMap);
01390       InsertAssociation( hyp->GetSourceVertex(2),hyp->GetTargetVertex(2),theAssociationMap);
01391     }
01392   }
01393 }
01394 
01395 //=======================================================================
01403 //=======================================================================
01404 
01405 bool StdMeshers_ProjectionUtils::InsertAssociation( const TopoDS_Shape& theShape1,
01406                                                     const TopoDS_Shape& theShape2,
01407                                                     TShapeShapeMap &    theAssociationMap,
01408                                                     const bool          theBidirectional)
01409 {
01410   if ( !theShape1.IsNull() && !theShape2.IsNull() ) {
01411     SHOW_SHAPE(theShape1,"Assoc ");
01412     SHOW_SHAPE(theShape2," to ");
01413     bool isNew = ( theAssociationMap.Bind( theShape1, theShape2 ));
01414     if ( theBidirectional )
01415       theAssociationMap.Bind( theShape2, theShape1 );
01416     return isNew;
01417   }
01418   else {
01419     throw SALOME_Exception("StdMeshers_ProjectionUtils: attempt to associate NULL shape");
01420   }
01421   return false;
01422 }
01423 
01424 //=======================================================================
01432 //=======================================================================
01433 
01434 TopoDS_Edge StdMeshers_ProjectionUtils::GetEdgeByVertices( SMESH_Mesh*          theMesh,
01435                                                            const TopoDS_Vertex& theV1,
01436                                                            const TopoDS_Vertex& theV2)
01437 {
01438   if ( theMesh && !theV1.IsNull() && !theV2.IsNull() )
01439   {
01440     TopTools_ListIteratorOfListOfShape ancestorIt( theMesh->GetAncestors( theV1 ));
01441     for ( ; ancestorIt.More(); ancestorIt.Next() )
01442       if ( ancestorIt.Value().ShapeType() == TopAbs_EDGE )
01443         for ( TopExp_Explorer expV ( ancestorIt.Value(), TopAbs_VERTEX );
01444               expV.More();
01445               expV.Next() )
01446           if ( theV2.IsSame( expV.Current() ))
01447             return TopoDS::Edge( ancestorIt.Value() );
01448   }
01449   return TopoDS_Edge();
01450 }
01451 
01452 //================================================================================
01460 //================================================================================
01461 
01462 TopoDS_Face StdMeshers_ProjectionUtils::GetNextFace( const TAncestorMap& edgeToFaces,
01463                                                      const TopoDS_Edge&  edge,
01464                                                      const TopoDS_Face&  face)
01465 {
01466 //   if ( !edge.IsNull() && !face.IsNull() && edgeToFaces.Contains( edge ))
01467   if ( !edge.IsNull() && edgeToFaces.Contains( edge )) // PAL16202
01468   {
01469     TopTools_ListIteratorOfListOfShape ancestorIt( edgeToFaces.FindFromKey( edge ));
01470     for ( ; ancestorIt.More(); ancestorIt.Next() )
01471       if ( ancestorIt.Value().ShapeType() == TopAbs_FACE &&
01472            !face.IsSame( ancestorIt.Value() ))
01473         return TopoDS::Face( ancestorIt.Value() );
01474   }
01475   return TopoDS_Face();
01476 }
01477 
01478 //================================================================================
01482 //================================================================================
01483 
01484 TopoDS_Vertex StdMeshers_ProjectionUtils::GetNextVertex(const TopoDS_Edge&   edge,
01485                                                         const TopoDS_Vertex& vertex)
01486 {
01487   TopoDS_Vertex vF,vL;
01488   TopExp::Vertices(edge,vF,vL);
01489   if ( vF.IsSame( vL ))
01490     return TopoDS_Vertex();
01491   return vertex.IsSame( vF ) ? vL : vF; 
01492 }
01493 
01494 //================================================================================
01502 //================================================================================
01503 
01504 pair<int,TopoDS_Edge>
01505 StdMeshers_ProjectionUtils::GetPropagationEdge( SMESH_Mesh*        aMesh,
01506                                                 const TopoDS_Edge& theEdge,
01507                                                 const TopoDS_Edge& fromEdge)
01508 {
01509   TopTools_IndexedMapOfShape aChain;
01510   int step = 0;
01511 
01512   // List of edges, added to chain on the previous cycle pass
01513   TopTools_ListOfShape listPrevEdges;
01514   listPrevEdges.Append(fromEdge);
01515 
01516   // Collect all edges pass by pass
01517   while (listPrevEdges.Extent() > 0) {
01518     step++;
01519     // List of edges, added to chain on this cycle pass
01520     TopTools_ListOfShape listCurEdges;
01521 
01522     // Find the next portion of edges
01523     TopTools_ListIteratorOfListOfShape itE (listPrevEdges);
01524     for (; itE.More(); itE.Next()) {
01525       TopoDS_Shape anE = itE.Value();
01526 
01527       // Iterate on faces, having edge <anE>
01528       TopTools_ListIteratorOfListOfShape itA (aMesh->GetAncestors(anE));
01529       for (; itA.More(); itA.Next()) {
01530         TopoDS_Shape aW = itA.Value();
01531 
01532         // There are objects of different type among the ancestors of edge
01533         if (aW.ShapeType() == TopAbs_WIRE) {
01534           TopoDS_Shape anOppE;
01535 
01536           BRepTools_WireExplorer aWE (TopoDS::Wire(aW));
01537           Standard_Integer nb = 1, found = 0;
01538           TopTools_Array1OfShape anEdges (1,4);
01539           for (; aWE.More(); aWE.Next(), nb++) {
01540             if (nb > 4) {
01541               found = 0;
01542               break;
01543             }
01544             anEdges(nb) = aWE.Current();
01545             if (anEdges(nb).IsSame(anE)) found = nb;
01546           }
01547 
01548           if (nb == 5 && found > 0) {
01549             // Quadrangle face found, get an opposite edge
01550             Standard_Integer opp = found + 2;
01551             if (opp > 4) opp -= 4;
01552             anOppE = anEdges(opp);
01553 
01554             // add anOppE to aChain if ...
01555             if (!aChain.Contains(anOppE)) { // ... anOppE is not in aChain
01556               // Add found edge to the chain oriented so that to
01557               // have it co-directed with a forward MainEdge
01558               TopAbs_Orientation ori = anE.Orientation();
01559               if ( anEdges(opp).Orientation() == anEdges(found).Orientation() )
01560                 ori = TopAbs::Reverse( ori );
01561               anOppE.Orientation( ori );
01562               if ( anOppE.IsSame( theEdge ))
01563                 return make_pair( step, TopoDS::Edge( anOppE ));
01564               aChain.Add(anOppE);
01565               listCurEdges.Append(anOppE);
01566             }
01567           } // if (nb == 5 && found > 0)
01568         } // if (aF.ShapeType() == TopAbs_WIRE)
01569       } // for (; itF.More(); itF.Next())
01570     } // for (; itE.More(); itE.Next())
01571 
01572     listPrevEdges = listCurEdges;
01573   } // while (listPrevEdges.Extent() > 0)
01574 
01575   return make_pair( INT_MAX, TopoDS_Edge());
01576 }
01577 
01578 //================================================================================
01589 //================================================================================
01590 
01591 bool StdMeshers_ProjectionUtils::
01592 FindMatchingNodesOnFaces( const TopoDS_Face&     face1,
01593                           SMESH_Mesh*            mesh1,
01594                           const TopoDS_Face&     face2,
01595                           SMESH_Mesh*            mesh2,
01596                           const TShapeShapeMap & assocMap,
01597                           TNodeNodeMap &         node1To2Map)
01598 {
01599   SMESHDS_Mesh* meshDS1 = mesh1->GetMeshDS();
01600   SMESHDS_Mesh* meshDS2 = mesh2->GetMeshDS();
01601   
01602   SMESH_MesherHelper helper1( *mesh1 );
01603   SMESH_MesherHelper helper2( *mesh2 );
01604 
01605   // Get corresponding submeshes and roughly check match of meshes
01606 
01607   SMESHDS_SubMesh * SM2 = meshDS2->MeshElements( face2 );
01608   SMESHDS_SubMesh * SM1 = meshDS1->MeshElements( face1 );
01609   if ( !SM2 || !SM1 )
01610     RETURN_BAD_RESULT("Empty submeshes");
01611   if ( SM2->NbNodes()    != SM1->NbNodes() ||
01612        SM2->NbElements() != SM1->NbElements() )
01613     RETURN_BAD_RESULT("Different meshes on corresponding faces "
01614                       << meshDS1->ShapeToIndex( face1 ) << " and "
01615                       << meshDS2->ShapeToIndex( face2 ));
01616   if ( SM2->NbElements() == 0 )
01617     RETURN_BAD_RESULT("Empty submeshes");
01618 
01619   helper1.SetSubShape( face1 );
01620   helper2.SetSubShape( face2 );
01621   if ( helper1.HasSeam() != helper2.HasSeam() )
01622     RETURN_BAD_RESULT("Different faces' geometry");
01623 
01624   // Data to call SMESH_MeshEditor::FindMatchingNodes():
01625 
01626   // 1. Nodes of corresponding links:
01627 
01628   // get 2 matching edges, try to find not seam ones
01629   TopoDS_Edge edge1, edge2, seam1, seam2, anyEdge1, anyEdge2;
01630   TopExp_Explorer eE( OuterShape( face2, TopAbs_WIRE ), TopAbs_EDGE );
01631   do {
01632     // edge 2
01633     TopoDS_Edge e2 = TopoDS::Edge( eE.Current() );
01634     eE.Next();
01635     // edge 1
01636     if ( !assocMap.IsBound( e2 ))
01637       RETURN_BAD_RESULT("Association not found for edge " << meshDS2->ShapeToIndex( e2 ));
01638     TopoDS_Edge e1 = TopoDS::Edge( assocMap( e2 ));
01639     if ( !helper1.IsSubShape( e1, face1 ))
01640       RETURN_BAD_RESULT("Wrong association, edge " << meshDS1->ShapeToIndex( e1 ) <<
01641                         " isn't a subshape of face " << meshDS1->ShapeToIndex( face1 ));
01642     // check that there are nodes on edges
01643     SMESHDS_SubMesh * eSM1 = meshDS1->MeshElements( e1 );
01644     SMESHDS_SubMesh * eSM2 = meshDS2->MeshElements( e2 );
01645     bool nodesOnEdges = ( eSM1 && eSM2 && eSM1->NbNodes() && eSM2->NbNodes() );
01646     // check that the nodes on edges belong to faces
01647     // (as NETGEN ignores nodes on the degenerated geom edge)
01648     bool nodesOfFaces = false;
01649     if ( nodesOnEdges ) {
01650       const SMDS_MeshNode* n1 = eSM1->GetNodes()->next();
01651       const SMDS_MeshNode* n2 = eSM2->GetNodes()->next();
01652       nodesOfFaces = ( n1->GetInverseElementIterator(SMDSAbs_Face)->more() &&
01653                        n2->GetInverseElementIterator(SMDSAbs_Face)->more() );
01654     }
01655     if ( nodesOfFaces )
01656     {
01657       if ( helper2.IsRealSeam( e2 )) {
01658         seam1 = e1; seam2 = e2;
01659       }
01660       else {
01661         edge1 = e1; edge2 = e2;
01662       }
01663     }
01664     else {
01665       anyEdge1 = e1; anyEdge2 = e2;
01666     }
01667   } while ( edge2.IsNull() && eE.More() );
01668   //
01669   if ( edge2.IsNull() ) {
01670     edge1 = seam1; edge2 = seam2;
01671   }
01672   bool hasNodesOnEdge = (! edge2.IsNull() );
01673   if ( !hasNodesOnEdge ) {
01674     // 0020338 - nb segments == 1
01675     edge1 = anyEdge1; edge2 = anyEdge2;
01676   }
01677 
01678   // get 2 matching vertices
01679   TopoDS_Vertex V2 = TopExp::FirstVertex( TopoDS::Edge( edge2 ));
01680   if ( !assocMap.IsBound( V2 ))
01681     RETURN_BAD_RESULT("Association not found for vertex " << meshDS2->ShapeToIndex( V2 ));
01682   TopoDS_Vertex V1 = TopoDS::Vertex( assocMap( V2 ));
01683 
01684   // nodes on vertices
01685   const SMDS_MeshNode* vNode1 = SMESH_Algo::VertexNode( V1, meshDS1 );
01686   const SMDS_MeshNode* vNode2 = SMESH_Algo::VertexNode( V2, meshDS2 );
01687   if ( !vNode1 ) RETURN_BAD_RESULT("No node on vertex #" << meshDS1->ShapeToIndex( V1 ));
01688   if ( !vNode2 ) RETURN_BAD_RESULT("No node on vertex #" << meshDS2->ShapeToIndex( V2 ));
01689 
01690   // nodes on edges linked with nodes on vertices
01691   const SMDS_MeshNode* nullNode = 0;
01692   vector< const SMDS_MeshNode*> eNode1( 2, nullNode );
01693   vector< const SMDS_MeshNode*> eNode2( 2, nullNode );
01694   if ( hasNodesOnEdge )
01695   {
01696     int nbNodeToGet = 1;
01697     if ( helper1.IsClosedEdge( edge1 ) || helper2.IsClosedEdge( edge2 ) )
01698       nbNodeToGet = 2;
01699     for ( int is2 = 0; is2 < 2; ++is2 )
01700     {
01701       TopoDS_Edge &     edge  = is2 ? edge2 : edge1;
01702       SMESHDS_Mesh *    smDS  = is2 ? meshDS2 : meshDS1;
01703       SMESHDS_SubMesh* edgeSM = smDS->MeshElements( edge );
01704       // nodes linked with ones on vertices
01705       const SMDS_MeshNode*           vNode = is2 ? vNode2 : vNode1;
01706       vector< const SMDS_MeshNode*>& eNode = is2 ? eNode2 : eNode1;
01707       int nbGotNode = 0;
01708       SMDS_ElemIteratorPtr vElem = vNode->GetInverseElementIterator(SMDSAbs_Edge);
01709       while ( vElem->more() && nbGotNode != nbNodeToGet ) {
01710         const SMDS_MeshElement* elem = vElem->next();
01711         if ( edgeSM->Contains( elem ))
01712           eNode[ nbGotNode++ ] = 
01713             ( elem->GetNode(0) == vNode ) ? elem->GetNode(1) : elem->GetNode(0);
01714       }
01715       if ( nbGotNode > 1 ) // sort found nodes by param on edge
01716       {
01717         SMESH_MesherHelper* helper = is2 ? &helper2 : &helper1;
01718         double u0 = helper->GetNodeU( edge, eNode[ 0 ]);
01719         double u1 = helper->GetNodeU( edge, eNode[ 1 ]);
01720         if ( u0 > u1 ) std::swap( eNode[ 0 ], eNode[ 1 ]);
01721       }
01722       if ( nbGotNode == 0 )
01723         RETURN_BAD_RESULT("Found no nodes on edge " << smDS->ShapeToIndex( edge ) <<
01724                           " linked to " << vNode );
01725     }
01726   }
01727   else // 0020338 - nb segments == 1
01728   {
01729     // get 2 other matching vertices
01730     V2 = TopExp::LastVertex( TopoDS::Edge( edge2 ));
01731     if ( !assocMap.IsBound( V2 ))
01732       RETURN_BAD_RESULT("Association not found for vertex " << meshDS2->ShapeToIndex( V2 ));
01733     V1 = TopoDS::Vertex( assocMap( V2 ));
01734 
01735     // nodes on vertices
01736     eNode1[0] = SMESH_Algo::VertexNode( V1, meshDS1 );
01737     eNode2[0] = SMESH_Algo::VertexNode( V2, meshDS2 );
01738     if ( !eNode1[0] ) RETURN_BAD_RESULT("No node on vertex #" << meshDS1->ShapeToIndex( V1 ));
01739     if ( !eNode2[0] ) RETURN_BAD_RESULT("No node on vertex #" << meshDS2->ShapeToIndex( V2 ));
01740   }
01741 
01742   // 2. face sets
01743 
01744   set<const SMDS_MeshElement*> Elems1, Elems2;
01745   for ( int is2 = 0; is2 < 2; ++is2 )
01746   {
01747     set<const SMDS_MeshElement*> & elems = is2 ? Elems2 : Elems1;
01748     SMESHDS_SubMesh*                  sm = is2 ? SM2 : SM1;
01749     SMESH_MesherHelper*           helper = is2 ? &helper2 : &helper1;
01750     const TopoDS_Face &             face = is2 ? face2 : face1;
01751     SMDS_ElemIteratorPtr eIt = sm->GetElements();
01752 
01753     if ( !helper->IsRealSeam( is2 ? edge2 : edge1 ))
01754     {
01755       while ( eIt->more() ) elems.insert( eIt->next() );
01756     }
01757     else
01758     {
01759       // the only suitable edge is seam, i.e. it is a sphere.
01760       // FindMatchingNodes() will not know which way to go from any edge.
01761       // So we ignore all faces having nodes on edges or vertices except
01762       // one of faces sharing current start nodes
01763 
01764       // find a face to keep
01765       const SMDS_MeshElement* faceToKeep = 0;
01766       const SMDS_MeshNode* vNode = is2 ? vNode2 : vNode1;
01767       const SMDS_MeshNode* eNode = is2 ? eNode2[0] : eNode1[0];
01768       TIDSortedElemSet inSet, notInSet;
01769 
01770       const SMDS_MeshElement* f1 =
01771         SMESH_MeshEditor::FindFaceInSet( vNode, eNode, inSet, notInSet );
01772       if ( !f1 ) RETURN_BAD_RESULT("The first face on seam not found");
01773       notInSet.insert( f1 );
01774 
01775       const SMDS_MeshElement* f2 =
01776         SMESH_MeshEditor::FindFaceInSet( vNode, eNode, inSet, notInSet );
01777       if ( !f2 ) RETURN_BAD_RESULT("The second face on seam not found");
01778 
01779       // select a face with less UV of vNode
01780       const SMDS_MeshNode* notSeamNode[2] = {0, 0};
01781       for ( int iF = 0; iF < 2; ++iF ) {
01782         const SMDS_MeshElement* f = ( iF ? f2 : f1 );
01783         for ( int i = 0; !notSeamNode[ iF ] && i < f->NbNodes(); ++i ) {
01784           const SMDS_MeshNode* node = f->GetNode( i );
01785           if ( !helper->IsSeamShape( node->getshapeId() ))
01786             notSeamNode[ iF ] = node;
01787         }
01788       }
01789       gp_Pnt2d uv1 = helper->GetNodeUV( face, vNode, notSeamNode[0] );
01790       gp_Pnt2d uv2 = helper->GetNodeUV( face, vNode, notSeamNode[1] );
01791       if ( uv1.X() + uv1.Y() > uv2.X() + uv2.Y() )
01792         faceToKeep = f2;
01793       else
01794         faceToKeep = f1;
01795 
01796       // fill elem set
01797       elems.insert( faceToKeep );
01798       while ( eIt->more() ) {
01799         const SMDS_MeshElement* f = eIt->next();
01800         int nbNodes = f->NbNodes();
01801         if ( f->IsQuadratic() )
01802           nbNodes /= 2;
01803         bool onBnd = false;
01804         for ( int i = 0; !onBnd && i < nbNodes; ++i ) {
01805           const SMDS_MeshNode* node = f->GetNode( i );
01806           onBnd = ( node->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE);
01807         }
01808         if ( !onBnd )
01809           elems.insert( f );
01810       }
01811       // add also faces adjacent to faceToKeep
01812       int nbNodes = faceToKeep->NbNodes();
01813       if ( faceToKeep->IsQuadratic() ) nbNodes /= 2;
01814       notInSet.insert( f1 );
01815       notInSet.insert( f2 );
01816       for ( int i = 0; i < nbNodes; ++i ) {
01817         const SMDS_MeshNode* n1 = faceToKeep->GetNode( i );
01818         const SMDS_MeshNode* n2 = faceToKeep->GetNode( i+1 % nbNodes );
01819         f1 = SMESH_MeshEditor::FindFaceInSet( n1, n2, inSet, notInSet );
01820         if ( f1 )
01821           elems.insert( f1 );
01822       }
01823     } // case on a sphere
01824   } // loop on 2 faces
01825 
01826   //  int quadFactor = (*Elems1.begin())->IsQuadratic() ? 2 : 1;
01827 
01828   node1To2Map.clear();
01829   int res = SMESH_MeshEditor::FindMatchingNodes( Elems1, Elems2,
01830                                                  vNode1, vNode2,
01831                                                  eNode1[0], eNode2[0],
01832                                                  node1To2Map);
01833   if ( res != SMESH_MeshEditor::SEW_OK )
01834     RETURN_BAD_RESULT("FindMatchingNodes() result " << res );
01835 
01836   // On a sphere, add matching nodes on the edge
01837 
01838   if ( helper1.IsRealSeam( edge1 ))
01839   {
01840     // sort nodes on edges by param on edge
01841     map< double, const SMDS_MeshNode* > u2nodesMaps[2];
01842     for ( int is2 = 0; is2 < 2; ++is2 )
01843     {
01844       TopoDS_Edge &     edge  = is2 ? edge2 : edge1;
01845       SMESHDS_Mesh *    smDS  = is2 ? meshDS2 : meshDS1;
01846       SMESHDS_SubMesh* edgeSM = smDS->MeshElements( edge );
01847       map< double, const SMDS_MeshNode* > & pos2nodes = u2nodesMaps[ is2 ];
01848 
01849       SMDS_NodeIteratorPtr nIt = edgeSM->GetNodes();
01850       while ( nIt->more() ) {
01851         const SMDS_MeshNode* node = nIt->next();
01852         const SMDS_EdgePosition* pos =
01853           static_cast<const SMDS_EdgePosition*>(node->GetPosition());
01854         pos2nodes.insert( make_pair( pos->GetUParameter(), node ));
01855       }
01856       if ( pos2nodes.size() != edgeSM->NbNodes() )
01857         RETURN_BAD_RESULT("Equal params of nodes on edge "
01858                           << smDS->ShapeToIndex( edge ) << " of face " << is2 );
01859     }
01860     if ( u2nodesMaps[0].size() != u2nodesMaps[1].size() )
01861       RETURN_BAD_RESULT("Different nb of new nodes on edges or wrong params");
01862 
01863     // compare edge orientation
01864     double u1 = helper1.GetNodeU( edge1, vNode1 );
01865     double u2 = helper2.GetNodeU( edge2, vNode2 );
01866     bool isFirst1 = ( u1 < u2nodesMaps[0].begin()->first );
01867     bool isFirst2 = ( u2 < u2nodesMaps[1].begin()->first );
01868     bool reverse ( isFirst1 != isFirst2 );
01869 
01870     // associate matching nodes
01871     map< double, const SMDS_MeshNode* >::iterator u_Node1, u_Node2, end1;
01872     map< double, const SMDS_MeshNode* >::reverse_iterator uR_Node2;
01873     u_Node1 = u2nodesMaps[0].begin();
01874     u_Node2 = u2nodesMaps[1].begin();
01875     uR_Node2 = u2nodesMaps[1].rbegin();
01876     end1 = u2nodesMaps[0].end();
01877     for ( ; u_Node1 != end1; ++u_Node1 ) {
01878       const SMDS_MeshNode* n1 = u_Node1->second;
01879       const SMDS_MeshNode* n2 = ( reverse ? (uR_Node2++)->second : (u_Node2++)->second );
01880       node1To2Map.insert( make_pair( n1, n2 ));
01881     }
01882 
01883     // associate matching nodes on the last vertices
01884     V2 = TopExp::LastVertex( TopoDS::Edge( edge2 ));
01885     if ( !assocMap.IsBound( V2 ))
01886       RETURN_BAD_RESULT("Association not found for vertex " << meshDS2->ShapeToIndex( V2 ));
01887     V1 = TopoDS::Vertex( assocMap( V2 ));
01888     vNode1 = SMESH_Algo::VertexNode( V1, meshDS1 );
01889     vNode2 = SMESH_Algo::VertexNode( V2, meshDS2 );
01890     if ( !vNode1 ) RETURN_BAD_RESULT("No node on vertex #" << meshDS1->ShapeToIndex( V1 ));
01891     if ( !vNode2 ) RETURN_BAD_RESULT("No node on vertex #" << meshDS2->ShapeToIndex( V2 ));
01892     node1To2Map.insert( make_pair( vNode1, vNode2 ));
01893   }
01894 
01895 // don't know why this condition is usually true :(
01896 //   if ( node1To2Map.size() * quadFactor < SM1->NbNodes() )
01897 //     MESSAGE("FindMatchingNodes() found too few node pairs starting from nodes ("
01898 //             << vNode1->GetID() << " - " << eNode1[0]->GetID() << ") ("
01899 //             << vNode2->GetID() << " - " << eNode2[0]->GetID() << "):"
01900 //             << node1To2Map.size() * quadFactor << " < " << SM1->NbNodes());
01901   
01902   return true;
01903 }
01904 
01905 //================================================================================
01912 //================================================================================
01913 
01914 TopoDS_Shape StdMeshers_ProjectionUtils::OuterShape( const TopoDS_Face& face,
01915                                                      TopAbs_ShapeEnum   type)
01916 {
01917   TopExp_Explorer exp( BRepTools::OuterWire( face ), type );
01918   if ( exp.More() )
01919     return exp.Current();
01920   return TopoDS_Shape();
01921 }
01922 
01923 //================================================================================
01930 //================================================================================
01931 
01932 bool StdMeshers_ProjectionUtils::MakeComputed(SMESH_subMesh * sm, const int iterationNb)
01933 {
01934   if ( iterationNb > 10 )
01935     RETURN_BAD_RESULT("Infinite recursive projection");
01936   if ( !sm )
01937     RETURN_BAD_RESULT("NULL submesh");
01938   if ( sm->IsMeshComputed() )
01939     return true;
01940 
01941   SMESH_Mesh* mesh = sm->GetFather();
01942   SMESH_Gen* gen   = mesh->GetGen();
01943   SMESH_Algo* algo = gen->GetAlgo( *mesh, sm->GetSubShape() );
01944   if ( !algo )
01945   {
01946     if ( sm->GetSubShape().ShapeType() != TopAbs_COMPOUND )
01947       RETURN_BAD_RESULT("No algo assigned to submesh " << sm->GetId());
01948     // group
01949     bool computed = true;
01950     for ( TopoDS_Iterator grMember( sm->GetSubShape() ); grMember.More(); grMember.Next())
01951       if ( SMESH_subMesh* grSub = mesh->GetSubMesh( grMember.Value() ))
01952         if ( !MakeComputed( grSub, iterationNb + 1 ))
01953           computed = false;
01954     return computed;
01955   }
01956 
01957   string algoType = algo->GetName();
01958   if ( algoType.substr(0, 11) != "Projection_")
01959     return gen->Compute( *mesh, sm->GetSubShape() );
01960 
01961   // try to compute source mesh
01962 
01963   const list <const SMESHDS_Hypothesis *> & hyps =
01964     algo->GetUsedHypothesis( *mesh, sm->GetSubShape() );
01965 
01966   TopoDS_Shape srcShape;
01967   SMESH_Mesh* srcMesh = 0;
01968   list <const SMESHDS_Hypothesis*>::const_iterator hIt = hyps.begin();
01969   for ( ; srcShape.IsNull() && hIt != hyps.end(); ++hIt ) {
01970     string hypName = (*hIt)->GetName();
01971     if ( hypName == "ProjectionSource1D" ) {
01972       const StdMeshers_ProjectionSource1D * hyp =
01973         static_cast<const StdMeshers_ProjectionSource1D*>( *hIt );
01974       srcShape = hyp->GetSourceEdge();
01975       srcMesh = hyp->GetSourceMesh();
01976     }
01977     else if ( hypName == "ProjectionSource2D" ) {
01978       const StdMeshers_ProjectionSource2D * hyp =
01979         static_cast<const StdMeshers_ProjectionSource2D*>( *hIt );
01980       srcShape = hyp->GetSourceFace();
01981       srcMesh = hyp->GetSourceMesh();
01982     }
01983     else if ( hypName == "ProjectionSource3D" ) {
01984       const StdMeshers_ProjectionSource3D * hyp =
01985         static_cast<const StdMeshers_ProjectionSource3D*>( *hIt );
01986       srcShape = hyp->GetSource3DShape();
01987       srcMesh = hyp->GetSourceMesh();
01988     }
01989   }
01990   if ( srcShape.IsNull() ) // no projection source defined
01991     return gen->Compute( *mesh, sm->GetSubShape() );
01992 
01993   if ( srcShape.IsSame( sm->GetSubShape() ))
01994     RETURN_BAD_RESULT("Projection from self");
01995     
01996   if ( !srcMesh )
01997     srcMesh = mesh;
01998 
01999   if ( MakeComputed( srcMesh->GetSubMesh( srcShape ), iterationNb + 1 ))
02000     return gen->Compute( *mesh, sm->GetSubShape() );
02001 
02002   return false;
02003 }
02004 
02005 //================================================================================
02012 //================================================================================
02013 
02014 int StdMeshers_ProjectionUtils::Count(const TopoDS_Shape&    shape,
02015                                       const TopAbs_ShapeEnum type,
02016                                       const bool             ignoreSame)
02017 {
02018   if ( ignoreSame ) {
02019     TopTools_IndexedMapOfShape map;
02020     TopExp::MapShapes( shape, type, map );
02021     return map.Extent();
02022   }
02023   else {
02024     int nb = 0;
02025     for ( TopExp_Explorer exp( shape, type ); exp.More(); exp.Next() )
02026       ++nb;
02027     return nb;
02028   }
02029 }
02030 
02031 //================================================================================
02035 //================================================================================
02036 
02037 bool StdMeshers_ProjectionUtils::IsBoundaryEdge(const TopoDS_Edge&  edge,
02038                                                 const TopoDS_Shape& edgeContainer,
02039                                                 SMESH_Mesh&         mesh)
02040 {
02041   TopTools_IndexedMapOfShape facesOfEdgeContainer, facesNearEdge;
02042   TopExp::MapShapes( edgeContainer, TopAbs_FACE, facesOfEdgeContainer );
02043 
02044   const TopTools_ListOfShape& EAncestors = mesh.GetAncestors(edge);
02045   TopTools_ListIteratorOfListOfShape itea(EAncestors);
02046   for(; itea.More(); itea.Next()) {
02047     if( itea.Value().ShapeType() == TopAbs_FACE &&
02048         facesOfEdgeContainer.Contains( itea.Value() ))
02049     {
02050       facesNearEdge.Add( itea.Value() );
02051       if ( facesNearEdge.Extent() > 1 )
02052         return false;
02053     }
02054   }
02055   return ( facesNearEdge.Extent() == 1 );
02056 }
02057 
02058 
02059 namespace {
02060 
02061   SMESH_subMeshEventListener* GetSrcSubMeshListener();
02062 
02063   //================================================================================
02068   //================================================================================
02069 
02070   struct HypModifWaiter: SMESH_subMeshEventListener
02071   {
02072     HypModifWaiter():SMESH_subMeshEventListener(0){} // won't be deleted by submesh
02073 
02074     void ProcessEvent(const int event, const int eventType, SMESH_subMesh* subMesh,
02075                       EventListenerData*, const SMESH_Hypothesis*)
02076     {
02077       if ( event     == SMESH_subMesh::MODIF_HYP &&
02078            eventType == SMESH_subMesh::ALGO_EVENT)
02079       {
02080         // delete current source listener
02081         subMesh->DeleteEventListener( GetSrcSubMeshListener() );
02082         // let algo set a new one
02083         SMESH_Gen* gen = subMesh->GetFather()->GetGen();
02084         if ( SMESH_Algo* algo = gen->GetAlgo( *subMesh->GetFather(),
02085                                               subMesh->GetSubShape() ))
02086           algo->SetEventListener( subMesh );
02087       }
02088     }
02089   };
02090   //================================================================================
02094   //================================================================================
02095 
02096   SMESH_subMeshEventListener* GetHypModifWaiter() {
02097     static HypModifWaiter aHypModifWaiter;
02098     return &aHypModifWaiter;
02099   }
02100   //================================================================================
02104   //================================================================================
02105 
02106   SMESH_subMeshEventListener* GetSrcSubMeshListener() {
02107     static SMESH_subMeshEventListener srcListener(0); // won't be deleted by submesh
02108     return &srcListener;
02109   }
02110 }
02111 
02112 //================================================================================
02119 //================================================================================
02120 
02121 void StdMeshers_ProjectionUtils::SetEventListener(SMESH_subMesh* subMesh,
02122                                                   TopoDS_Shape   srcShape,
02123                                                   SMESH_Mesh*    srcMesh)
02124 {
02125   // Set listener that resets an event listener on source submesh when
02126   // "ProjectionSource*D" hypothesis is modified since source shape can be changed
02127   subMesh->SetEventListener( GetHypModifWaiter(),0,subMesh);
02128 
02129   // Set an event listener to submesh of the source shape
02130   if ( !srcShape.IsNull() )
02131   {
02132     if ( !srcMesh )
02133       srcMesh = subMesh->GetFather();
02134 
02135     SMESH_subMesh* srcShapeSM = srcMesh->GetSubMesh( srcShape );
02136 
02137     if ( srcShapeSM != subMesh ) {
02138       if ( srcShapeSM->GetSubMeshDS() &&
02139            srcShapeSM->GetSubMeshDS()->IsComplexSubmesh() )
02140       {  // source shape is a group
02141         TopExp_Explorer it(srcShapeSM->GetSubShape(), // explore the group into subshapes...
02142                            subMesh->GetSubShape().ShapeType()); // ...of target shape type
02143         for (; it.More(); it.Next())
02144         {
02145           SMESH_subMesh* srcSM = srcMesh->GetSubMesh( it.Current() );
02146           if ( srcSM != subMesh )
02147           {
02148             SMESH_subMeshEventListenerData* data =
02149               srcSM->GetEventListenerData(GetSrcSubMeshListener());
02150             if ( data )
02151               data->mySubMeshes.push_back( subMesh );
02152             else
02153               data = SMESH_subMeshEventListenerData::MakeData( subMesh );
02154             subMesh->SetEventListener ( GetSrcSubMeshListener(), data, srcSM );
02155           }
02156         }
02157       }
02158       else
02159       {
02160         subMesh->SetEventListener( GetSrcSubMeshListener(),
02161                                    SMESH_subMeshEventListenerData::MakeData( subMesh ),
02162                                    srcShapeSM );
02163       }
02164     }
02165   }
02166 }
Copyright © 2007-2011 CEA/DEN, EDF R&D, OPEN CASCADE
Copyright © 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN, CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS