Salome HOME
f531ba16ec9675970c6ef53363f0a3e692b8972d
[modules/geom.git] / src / GEOMAlgo / GEOMAlgo_GetInPlaceAPI.cxx
1 // Copyright (C) 2007-2014  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License, or (at your option) any later version.
10 //
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
14 // Lesser General Public License for more details.
15 //
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
19 //
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22 // File:     GEOMAlgo_GetInPlaceAPI.cxx
23 // Created:
24 // Author:   Sergey KHROMOV
25
26
27 #include <GEOMAlgo_GetInPlaceAPI.hxx>
28 #include <GEOMAlgo_GetInPlace.hxx>
29 #include <GEOM_Function.hxx>
30 #include <GEOM_Object.hxx>
31 #include <GEOMUtils.hxx>
32
33 #include <Bnd_Box.hxx>
34 #include <BRepAdaptor_Surface.hxx>
35 #include <BRepBndLib.hxx>
36 #include <BRepBuilderAPI_MakeVertex.hxx>
37 #include <BRepExtrema_DistShapeShape.hxx>
38 #include <BRepGProp.hxx>
39 #include <BRep_Tool.hxx>
40 #include <Geom2d_Curve.hxx>
41 #include <GProp_GProps.hxx>
42 #include <gp_Pnt.hxx>
43 #include <Precision.hxx>
44 #include <TDataStd_IntegerArray.hxx>
45 #include <TopExp.hxx>
46 #include <TopExp_Explorer.hxx>
47 #include <TopoDS.hxx>
48 #include <TopoDS_Vertex.hxx>
49 #include <TopTools_MapOfShape.hxx>
50
51
52 //=======================================================================
53 //function : GetInPlace
54 //purpose  : 
55 //=======================================================================
56 Standard_Boolean GEOMAlgo_GetInPlaceAPI::GetInPlace
57                       (const TopoDS_Shape        &theWhere,
58                        const TopoDS_Shape        &theWhat,
59                              GEOMAlgo_GetInPlace &theGIP)
60 {
61   if (theWhere.IsNull() || theWhat.IsNull()) {
62     return Standard_False;
63   }
64
65   // Compute confusion tolerance.
66   Standard_Real    aTolConf = Precision::Confusion();
67   Standard_Integer i;
68
69   for (i = 0; i < 2; ++i) {
70     TopExp_Explorer anExp(i == 0 ? theWhere : theWhat, TopAbs_VERTEX);
71
72     for (; anExp.More(); anExp.Next()) {
73       const TopoDS_Vertex aVtx = TopoDS::Vertex(anExp.Current());
74       const Standard_Real aTolVtx = BRep_Tool::Tolerance(aVtx);
75
76       if (aTolVtx > aTolConf) {
77         aTolConf = aTolVtx;
78       }
79     }
80   }
81
82   // Compute mass tolerance.
83   Bnd_Box       aBoundingBox;
84   Standard_Real aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
85   Standard_Real aMassTol;
86
87   BRepBndLib::Add(theWhere, aBoundingBox);
88   BRepBndLib::Add(theWhat,  aBoundingBox);
89   aBoundingBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
90   aMassTol = Max(aXmax - aXmin, aYmax - aYmin);
91   aMassTol = Max(aMassTol, aZmax - aZmin);
92   aMassTol *= aTolConf;
93
94   // Searching for the sub-shapes inside the ShapeWhere shape
95   theGIP.SetTolerance(aTolConf);
96   theGIP.SetTolMass(aMassTol);
97   theGIP.SetTolCG(aTolConf);
98
99   theGIP.SetArgument(theWhat);
100   theGIP.SetShapeWhere(theWhere);
101
102   theGIP.Perform();
103
104   int iErr = theGIP.ErrorStatus();
105
106   if (iErr) {
107     return Standard_False;
108   }
109
110   return Standard_True;
111 }
112
113 //=======================================================================
114 //function : GetInPlaceOld
115 //purpose  : 
116 //=======================================================================
117 Standard_Integer GEOMAlgo_GetInPlaceAPI::GetInPlaceOld
118             (const TopoDS_Shape         &theWhere,
119              const TopoDS_Shape         &theWhat,
120                    TopTools_ListOfShape &theShapesInPlace)
121 {
122   theShapesInPlace.Clear();
123
124   if (theWhere.IsNull() || theWhat.IsNull()) {
125     // Error: aWhere and aWhat TopoDS_Shape are Null.
126     return 1;
127   }
128
129   TopoDS_Shape     aPntShape;
130   TopoDS_Vertex    aVertex;
131   bool             isFound = false;
132   TopAbs_ShapeEnum iType   = TopAbs_SOLID;
133   //Standard_Real    aWhat_Mass = 0., aWhere_Mass = 0.;
134   Standard_Real    tab_aWhat[4],    tab_aWhere[4];
135   Standard_Real    dl_l = 1e-3;
136   Standard_Real    min_l, Tol_0D, Tol_1D, Tol_2D, Tol_3D, Tol_Mass;
137   Standard_Real    aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
138   Bnd_Box          BoundingBox;
139   gp_Pnt           aPnt, aPnt_aWhat, tab_Pnt[2];
140   GProp_GProps     aProps;
141
142   iType = GEOMUtils::GetTypeOfSimplePart(theWhat);
143   if (iType == TopAbs_SHAPE) {
144     // Error: An attempt to extract a shape of not supported type.
145     return 2;
146   }
147
148   TopExp_Explorer Exp_aWhat  ( theWhat,  iType );
149   TopExp_Explorer Exp_aWhere ( theWhere, iType );
150   TopExp_Explorer Exp_Edge   ( theWhere, TopAbs_EDGE );
151
152   // Find the shortest edge in theShapeWhere shape
153   BRepBndLib::Add(theWhere, BoundingBox);
154   BoundingBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
155   min_l = fabs(aXmax - aXmin);
156   if( min_l < fabs(aYmax - aYmin) ) min_l = fabs(aYmax - aYmin);
157   if( min_l < fabs(aZmax - aZmin) ) min_l = fabs(aZmax - aZmin);
158   min_l /= dl_l;
159   // Mantis issue 0020908 BEGIN
160   if (!Exp_Edge.More()) {
161     min_l = Precision::Confusion();
162   }
163   // Mantis issue 0020908 END
164   for ( Standard_Integer nbEdge = 0; Exp_Edge.More(); Exp_Edge.Next(), nbEdge++ ) {
165     TopExp_Explorer Exp_Vertex( Exp_Edge.Current(), TopAbs_VERTEX);
166     for ( Standard_Integer nbVertex = 0; Exp_Vertex.More(); Exp_Vertex.Next(), nbVertex++ ) {
167       aPnt = BRep_Tool::Pnt( TopoDS::Vertex( Exp_Vertex.Current() ) );
168       tab_Pnt[nbVertex] = aPnt;
169     }
170     if ( ! tab_Pnt[0].IsEqual(tab_Pnt[1], dl_l) ) {
171       BRepGProp::LinearProperties(Exp_Edge.Current(), aProps);
172       if ( aProps.Mass() < min_l ) min_l = aProps.Mass();
173     }
174   }
175
176   // Compute tolerances
177   Tol_0D = dl_l;
178   Tol_1D = dl_l * min_l;
179   Tol_2D = dl_l * ( min_l * min_l) * ( 2. + dl_l);
180   Tol_3D = dl_l * ( min_l * min_l * min_l ) * ( 3. + (3 * dl_l) + (dl_l * dl_l) );
181
182   if (Tol_0D < Precision::Confusion()) Tol_0D = Precision::Confusion();
183   if (Tol_1D < Precision::Confusion()) Tol_1D = Precision::Confusion();
184   if (Tol_2D < Precision::Confusion()) Tol_2D = Precision::Confusion();
185   if (Tol_3D < Precision::Confusion()) Tol_3D = Precision::Confusion();
186
187   Tol_Mass = Tol_3D;
188   if ( iType == TopAbs_VERTEX )    Tol_Mass = Tol_0D;
189   else if ( iType == TopAbs_EDGE ) Tol_Mass = Tol_1D;
190   else if ( iType == TopAbs_FACE ) Tol_Mass = Tol_2D;
191
192   // Searching for the sub-shapes inside the ShapeWhere shape
193   TopTools_MapOfShape map_aWhere;
194   for ( Exp_aWhere.ReInit(); Exp_aWhere.More(); Exp_aWhere.Next() ) {
195     if (!map_aWhere.Add(Exp_aWhere.Current()))
196       continue; // skip repeated shape to avoid mass addition
197     GetShapeProperties( Exp_aWhere.Current(), tab_aWhere, aPnt );
198     for ( Exp_aWhat.ReInit(); Exp_aWhat.More(); Exp_aWhat.Next() ) {
199       GetShapeProperties( Exp_aWhat.Current(), tab_aWhat, aPnt_aWhat );
200       if ( fabs(tab_aWhat[3] - tab_aWhere[3]) <= Tol_Mass && aPnt_aWhat.Distance(aPnt) <= Tol_1D )
201         isFound = true;
202       else {
203         if ( (tab_aWhat[3] - tab_aWhere[3]) > Tol_Mass ) {
204           aPntShape = BRepBuilderAPI_MakeVertex( aPnt ).Shape();
205           aVertex   = TopoDS::Vertex( aPntShape );
206           BRepExtrema_DistShapeShape aWhereDistance ( aVertex, Exp_aWhere.Current() );
207           BRepExtrema_DistShapeShape aWhatDistance  ( aVertex, Exp_aWhat.Current() );
208           if ( aWhereDistance.IsDone() && aWhatDistance.IsDone() &&
209                fabs(aWhereDistance.Value() - aWhatDistance.Value()) <= Tol_1D )
210           {
211             // 0020162: "EDF 961 GEOM : Getinplace is getting additionnal orthogonal faces"
212             // aVertex must be projected to the same point on Where and on What
213             gp_Pnt pOnWhat  = aWhatDistance.PointOnShape2(1);
214             gp_Pnt pOnWhere = aWhereDistance.PointOnShape2(1);
215             isFound = ( pOnWhat.Distance(pOnWhere) <= Tol_1D );
216             if ( isFound && iType == TopAbs_FACE )
217             {
218               // check normals at pOnWhat and pOnWhere
219               const double angleTol = M_PI/180.;
220               gp_Vec normToWhat  = GetNormal( TopoDS::Face(Exp_aWhat.Current()), aWhatDistance);
221               gp_Vec normToWhere = GetNormal( TopoDS::Face(Exp_aWhere.Current()), aWhereDistance);
222               if ( normToWhat * normToWhere < 0 )
223                 normToWhat.Reverse();
224               isFound = ( normToWhat.Angle( normToWhere ) < angleTol );
225             }
226           }
227         }
228       }
229       if ( isFound ) {
230         theShapesInPlace.Append(Exp_aWhere.Current());
231         //aWhere_Mass += tab_aWhere[3];
232         isFound = false;
233         break;
234       }
235     }
236   }
237
238   if (theShapesInPlace.Extent() == 0) {
239     // Not found any Results
240     return 3;
241   }
242
243   return 0;
244 }
245
246 //=======================================================================
247 //function : GetNormal
248 //purpose  : 
249 //=======================================================================
250 gp_Vec GEOMAlgo_GetInPlaceAPI::GetNormal
251                          (const TopoDS_Face                &theFace,
252                           const BRepExtrema_DistShapeShape &theExtrema)
253 {
254   gp_Vec defaultNorm(1,0,0); // to have same normals on different faces
255   try {
256     // get UV at extrema point
257     Standard_Real u,v, f,l;
258     switch ( theExtrema.SupportTypeShape2(1) ) {
259     case BRepExtrema_IsInFace: {
260       theExtrema.ParOnFaceS2(1, u, v );
261       break;
262     }
263     case BRepExtrema_IsOnEdge: {
264       TopoDS_Edge edge = TopoDS::Edge( theExtrema.SupportOnShape2(1));
265       Handle(Geom2d_Curve) pcurve =
266         BRep_Tool::CurveOnSurface(edge, theFace, f,l);
267
268       theExtrema.ParOnEdgeS2( 1, u );
269       gp_Pnt2d uv = pcurve->Value( u );
270       u = uv.Coord(1);
271       v = uv.Coord(2);
272       break;
273     }
274     case BRepExtrema_IsVertex: return defaultNorm;
275     }
276     // get derivatives
277     BRepAdaptor_Surface surface( theFace, false );
278     gp_Vec du, dv; gp_Pnt p;
279     surface.D1( u, v, p, du, dv );
280
281     return du ^ dv;
282
283   } catch (Standard_Failure ) {
284   }
285   return defaultNorm;
286 }
287
288 //=======================================================================
289 //function : GetShapeProperties
290 //purpose  : 
291 //=======================================================================
292 void GEOMAlgo_GetInPlaceAPI::GetShapeProperties(const TopoDS_Shape  &theShape,
293                                                       Standard_Real  theTab[],
294                                                       gp_Pnt        &theVertex)
295 {
296   GProp_GProps  aProps;
297   gp_Pnt        aCenterMass;
298   Standard_Real aShapeSize;
299
300   if    (theShape.ShapeType() == TopAbs_VERTEX) {
301     aCenterMass = BRep_Tool::Pnt(TopoDS::Vertex(theShape));
302   } else if (theShape.ShapeType() == TopAbs_EDGE) {
303     BRepGProp::LinearProperties(theShape,  aProps);
304   } else if (theShape.ShapeType() == TopAbs_FACE) {
305     BRepGProp::SurfaceProperties(theShape, aProps);
306   } else {
307     BRepGProp::VolumeProperties(theShape,  aProps);
308   }
309
310   if (theShape.ShapeType() == TopAbs_VERTEX) {
311     aShapeSize = 1;
312   } else {
313     aCenterMass = aProps.CentreOfMass();
314     aShapeSize  = aProps.Mass();
315   }
316
317   theVertex = aCenterMass;
318   theTab[0] = theVertex.X();
319   theTab[1] = theVertex.Y();
320   theTab[2] = theVertex.Z();
321   theTab[3] = aShapeSize;
322 }
323
324 //=======================================================================
325 //function : GetInPlaceByHistory
326 //purpose  : 
327 //=======================================================================
328 Standard_Boolean GEOMAlgo_GetInPlaceAPI::GetInPlaceByHistory
329                       (const Handle(GEOM_Function)      &theWhereFunction,
330                        const TopTools_IndexedMapOfShape &theWhereIndices,
331                        const TopoDS_Shape               &theWhat,
332                              TopTools_ListOfShape       &theShapesInPlace)
333 {
334   if (theWhereFunction.IsNull() || theWhat.IsNull())
335     return Standard_False;
336
337   if (theWhereIndices.Contains(theWhat)) {
338     // entity was not changed by the operation
339     theShapesInPlace.Append(theWhat);
340
341     return Standard_True;
342   }
343
344   // try to find in history
345   TDF_Label aHistoryLabel = theWhereFunction->GetHistoryEntry(Standard_False);
346
347   // search in history for all argument shapes
348   Standard_Boolean isFound = Standard_False;
349   Standard_Boolean isGood = Standard_False;
350
351   TDF_LabelSequence aLabelSeq;
352   theWhereFunction->GetDependency(aLabelSeq);
353   Standard_Integer nbArg = aLabelSeq.Length();
354
355   for (Standard_Integer iarg = 1; iarg <= nbArg && !isFound; iarg++) {
356
357     TDF_Label anArgumentRefLabel = aLabelSeq.Value(iarg);
358
359     Handle(GEOM_Object) anArgumentObject = GEOM_Object::GetReferencedObject(anArgumentRefLabel);
360     TopoDS_Shape anArgumentShape = anArgumentObject->GetValue();
361
362     TopTools_IndexedMapOfShape anArgumentIndices;
363     TopExp::MapShapes(anArgumentShape, anArgumentIndices);
364
365     if (anArgumentIndices.Contains(theWhat)) {
366       isFound = Standard_True;
367       Standard_Integer aWhatIndex = anArgumentIndices.FindIndex(theWhat);
368
369       // Find corresponding label in history
370       TDF_Label anArgumentHistoryLabel =
371         theWhereFunction->GetArgumentHistoryEntry(anArgumentRefLabel, Standard_False);
372       if (anArgumentHistoryLabel.IsNull()) {
373         // Lost History of operation argument. Possibly, all its entities was removed.
374         isGood = Standard_True;
375       }
376       else {
377         TDF_Label aWhatHistoryLabel = anArgumentHistoryLabel.FindChild(aWhatIndex, Standard_False);
378
379         if (aWhatHistoryLabel.IsNull()) {
380           // Removed entity ? Compound ? Compsolid ? Shell ? Wire
381           isGood = Standard_False;
382         } else {
383           Handle(TDataStd_IntegerArray) anIntegerArray;
384           if (!aWhatHistoryLabel.FindAttribute(TDataStd_IntegerArray::GetID(), anIntegerArray)) {
385             //Error: Empty modifications history for the sought shape.
386             isGood = Standard_False;
387           }
388           else {
389             isGood = Standard_True;
390             Standard_Integer imod, aModifLen = anIntegerArray->Array()->Length();
391             for (imod = 1; imod <= aModifLen; imod++) {
392               const Standard_Integer anIndex =
393                 anIntegerArray->Array()->Value(imod);
394
395               theShapesInPlace.Append(theWhereIndices.FindKey(anIndex));
396             }
397           }
398         }
399       }
400     }
401   }
402
403   isFound = isGood;
404
405   if (!isFound) {
406     // try compound/compsolid/shell/wire element by element
407     Standard_Boolean isFoundAny = Standard_False;
408     TopTools_MapOfShape mapShape;
409
410     if (theWhat.ShapeType() == TopAbs_COMPOUND ||
411         theWhat.ShapeType() == TopAbs_COMPSOLID) {
412       // recursive processing of compound/compsolid
413       TopoDS_Iterator anIt (theWhat, Standard_True, Standard_True);
414       for (; anIt.More(); anIt.Next()) {
415         if (mapShape.Add(anIt.Value())) {
416           TopoDS_Shape curWhat = anIt.Value();
417           isFoundAny = GetInPlaceByHistory(theWhereFunction, theWhereIndices, curWhat, theShapesInPlace);
418           if (isFoundAny) isFound = Standard_True;
419         }
420       }
421     }
422     else if (theWhat.ShapeType() == TopAbs_SHELL) {
423       // try to replace a shell by its faces images
424       TopExp_Explorer anExp (theWhat, TopAbs_FACE);
425       for (; anExp.More(); anExp.Next()) {
426         if (mapShape.Add(anExp.Current())) {
427           TopoDS_Shape curWhat = anExp.Current();
428           isFoundAny = GetInPlaceByHistory(theWhereFunction, theWhereIndices, curWhat, theShapesInPlace);
429           if (isFoundAny) isFound = Standard_True;
430         }
431       }
432     }
433     else if (theWhat.ShapeType() == TopAbs_WIRE) {
434       // try to replace a wire by its edges images
435       TopExp_Explorer anExp (theWhat, TopAbs_EDGE);
436       for (; anExp.More(); anExp.Next()) {
437         if (mapShape.Add(anExp.Current())) {
438           TopoDS_Shape curWhat = anExp.Current();
439           isFoundAny = GetInPlaceByHistory(theWhereFunction, theWhereIndices, curWhat, theShapesInPlace);
440           if (isFoundAny) isFound = Standard_True;
441         }
442       }
443     }
444     else {
445       // Removed entity
446     }
447   }
448
449   return isFound;
450 }