Salome HOME
9d21bb9ee910f961df9c5d12f7fd55d9207357c8
[modules/geom.git] / src / GEOMAlgo / GEOMAlgo_GetInPlaceAPI.cxx
1 // Copyright (C) 2007-2023  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_Object.hxx>
30 #include <GEOMUtils.hxx>
31
32 #include <Bnd_Box.hxx>
33 #include <BRepAdaptor_Surface.hxx>
34 #include <BRepBndLib.hxx>
35 #include <BRepBuilderAPI_MakeVertex.hxx>
36 #include <BRepExtrema_DistShapeShape.hxx>
37 #include <BRepGProp.hxx>
38 #include <BRep_Tool.hxx>
39 #include <Geom2d_Curve.hxx>
40 #include <GProp_GProps.hxx>
41 #include <gp_Pnt.hxx>
42 #include <Precision.hxx>
43 #include <TDataStd_IntegerArray.hxx>
44 #include <TopExp.hxx>
45 #include <TopExp_Explorer.hxx>
46 #include <TopoDS.hxx>
47 #include <TopoDS_Vertex.hxx>
48 #include <TColStd_MapOfInteger.hxx>
49
50 //=======================================================================
51 //function : GetInPlace
52 //purpose  : 
53 //=======================================================================
54 Standard_Boolean GEOMAlgo_GetInPlaceAPI::GetInPlace
55                       (const TopoDS_Shape        &theWhere,
56                        const TopoDS_Shape        &theWhat,
57                              GEOMAlgo_GetInPlace &theGIP)
58 {
59   if (theWhere.IsNull() || theWhat.IsNull()) {
60     return Standard_False;
61   }
62
63   // Compute confusion tolerance.
64   Standard_Real    aTolConf = Precision::Confusion();
65   Standard_Integer i;
66
67   for (i = 0; i < 2; ++i) {
68     TopExp_Explorer anExp(i == 0 ? theWhere : theWhat, TopAbs_VERTEX);
69
70     for (; anExp.More(); anExp.Next()) {
71       const TopoDS_Vertex aVtx = TopoDS::Vertex(anExp.Current());
72       const Standard_Real aTolVtx = BRep_Tool::Tolerance(aVtx);
73
74       if (aTolVtx > aTolConf) {
75         aTolConf = aTolVtx;
76       }
77     }
78   }
79
80   // Compute mass tolerance.
81   Bnd_Box       aBoundingBox;
82   Standard_Real aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
83   Standard_Real aMassTol;
84
85   BRepBndLib::Add(theWhere, aBoundingBox);
86   BRepBndLib::Add(theWhat,  aBoundingBox);
87   aBoundingBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
88   aMassTol = Max(aXmax - aXmin, aYmax - aYmin);
89   aMassTol = Max(aMassTol, aZmax - aZmin);
90   aMassTol *= aTolConf;
91
92   // Searching for the sub-shapes inside the ShapeWhere shape
93   theGIP.SetTolerance(aTolConf);
94   theGIP.SetTolMass(aMassTol);
95   theGIP.SetTolCG(aTolConf);
96
97   theGIP.SetArgument(theWhat);
98   theGIP.SetShapeWhere(theWhere);
99
100   theGIP.Perform();
101
102   int iErr = theGIP.ErrorStatus();
103
104   if (iErr) {
105     return Standard_False;
106   }
107
108   return Standard_True;
109 }
110
111 //=======================================================================
112 //function : GetInPlaceOld
113 //purpose  : 
114 //=======================================================================
115 Standard_Integer GEOMAlgo_GetInPlaceAPI::GetInPlaceOld
116             (const TopoDS_Shape         &theWhere,
117              const TopoDS_Shape         &theWhat,
118                    TopTools_ListOfShape &theShapesInPlace)
119 {
120   theShapesInPlace.Clear();
121
122   if (theWhere.IsNull() || theWhat.IsNull()) {
123     // Error: aWhere and aWhat TopoDS_Shape are Null.
124     return 1;
125   }
126
127   // Check shape type.
128   TopAbs_ShapeEnum iType = GEOMUtils::GetTypeOfSimplePart(theWhat);
129
130   if (iType == TopAbs_SHAPE) {
131     // Error: An attempt to extract a shape of not supported type.
132     return 2;
133   }
134
135   // Compute confusion tolerance.
136   Standard_Real    aTolConf = Precision::Confusion();
137   Standard_Integer i;
138
139   for (i = 0; i < 2; ++i) {
140     TopExp_Explorer anExp(i == 0 ? theWhere : theWhat, TopAbs_VERTEX);
141
142     for (; anExp.More(); anExp.Next()) {
143       const TopoDS_Vertex aVtx = TopoDS::Vertex(anExp.Current());
144       const Standard_Real aTolVtx = BRep_Tool::Tolerance(aVtx);
145
146       if (aTolVtx > aTolConf) {
147         aTolConf = aTolVtx;
148       }
149     }
150   }
151
152   // Compute mass tolerance.
153   Bnd_Box       aBoundingBox;
154   Standard_Real aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
155   Standard_Real aMassTol;
156
157   BRepBndLib::Add(theWhere, aBoundingBox);
158   BRepBndLib::Add(theWhat,  aBoundingBox);
159   aBoundingBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
160   aMassTol = Max(aXmax - aXmin, aYmax - aYmin);
161   aMassTol = Max(aMassTol, aZmax - aZmin);
162   aMassTol *= aTolConf;
163
164   // Compute the result.
165   TopExp_Explorer     Exp_aWhat  (theWhat,  iType);
166   TopExp_Explorer     Exp_aWhere (theWhere, iType);
167   Standard_Real       tab_aWhat[4], tab_aWhere[4];
168   gp_Pnt              aPnt, aPnt_aWhat;
169   TopoDS_Shape        aPntShape;
170   TopoDS_Vertex       aVertex;
171   bool                isFound = false;
172   TopTools_MapOfShape map_aWhere;
173
174   for (; Exp_aWhere.More(); Exp_aWhere.Next()) {
175     if (!map_aWhere.Add(Exp_aWhere.Current()))
176       continue; // skip repeated shape to avoid mass addition
177     GetShapeProperties( Exp_aWhere.Current(), tab_aWhere, aPnt );
178     for ( Exp_aWhat.ReInit(); Exp_aWhat.More(); Exp_aWhat.Next() ) {
179       GetShapeProperties( Exp_aWhat.Current(), tab_aWhat, aPnt_aWhat );
180       if (fabs(tab_aWhat[3] - tab_aWhere[3]) <= aMassTol && aPnt_aWhat.Distance(aPnt) <= aTolConf)
181         isFound = true;
182       else {
183         if (tab_aWhat[3] > tab_aWhere[3]) {
184           aPntShape = BRepBuilderAPI_MakeVertex( aPnt ).Shape();
185           aVertex   = TopoDS::Vertex( aPntShape );
186           BRepExtrema_DistShapeShape aWhereDistance ( aVertex, Exp_aWhere.Current() );
187           BRepExtrema_DistShapeShape aWhatDistance  ( aVertex, Exp_aWhat.Current() );
188           if (aWhereDistance.IsDone() && aWhatDistance.IsDone() &&
189               fabs(aWhereDistance.Value() - aWhatDistance.Value()) <= aTolConf)
190           {
191             // 0020162: "EDF 961 GEOM : Getinplace is getting additional orthogonal faces"
192             // aVertex must be projected to the same point on Where and on What
193             gp_Pnt pOnWhat  = aWhatDistance.PointOnShape2(1);
194             gp_Pnt pOnWhere = aWhereDistance.PointOnShape2(1);
195             isFound = (pOnWhat.Distance(pOnWhere) <= aTolConf);
196             if ( isFound && iType == TopAbs_FACE )
197             {
198               // check normals at pOnWhat and pOnWhere
199               const double angleTol = M_PI/180.;
200               gp_Vec normToWhat  = GetNormal( TopoDS::Face(Exp_aWhat.Current()), aWhatDistance);
201               gp_Vec normToWhere = GetNormal( TopoDS::Face(Exp_aWhere.Current()), aWhereDistance);
202               if ( normToWhat * normToWhere < 0 )
203                 normToWhat.Reverse();
204               isFound = ( normToWhat.Angle( normToWhere ) < angleTol );
205             }
206           }
207         }
208       }
209       if ( isFound ) {
210         theShapesInPlace.Append(Exp_aWhere.Current());
211         //aWhere_Mass += tab_aWhere[3];
212         isFound = false;
213         break;
214       }
215     }
216   }
217
218   if (theShapesInPlace.Extent() == 0) {
219     // Not found any Results
220     return 3;
221   }
222
223   return 0;
224 }
225
226 //=======================================================================
227 //function : GetNormal
228 //purpose  : 
229 //=======================================================================
230 gp_Vec GEOMAlgo_GetInPlaceAPI::GetNormal
231                          (const TopoDS_Face                &theFace,
232                           const BRepExtrema_DistShapeShape &theExtrema)
233 {
234   gp_Vec defaultNorm(1,0,0); // to have same normals on different faces
235   try {
236     // get UV at extrema point
237     Standard_Real u,v, f,l;
238     switch ( theExtrema.SupportTypeShape2(1) ) {
239     case BRepExtrema_IsInFace: {
240       theExtrema.ParOnFaceS2(1, u, v );
241       break;
242     }
243     case BRepExtrema_IsOnEdge: {
244       TopoDS_Edge edge = TopoDS::Edge( theExtrema.SupportOnShape2(1));
245       Handle(Geom2d_Curve) pcurve =
246         BRep_Tool::CurveOnSurface(edge, theFace, f,l);
247
248       theExtrema.ParOnEdgeS2( 1, u );
249       gp_Pnt2d uv = pcurve->Value( u );
250       u = uv.Coord(1);
251       v = uv.Coord(2);
252       break;
253     }
254     case BRepExtrema_IsVertex: return defaultNorm;
255     }
256     // get derivatives
257     BRepAdaptor_Surface surface( theFace, false );
258     gp_Vec du, dv; gp_Pnt p;
259     surface.D1( u, v, p, du, dv );
260
261     return du ^ dv;
262
263   } catch (Standard_Failure&) {
264   }
265   return defaultNorm;
266 }
267
268 //=======================================================================
269 //function : GetShapeProperties
270 //purpose  : 
271 //=======================================================================
272 void GEOMAlgo_GetInPlaceAPI::GetShapeProperties(const TopoDS_Shape  &theShape,
273                                                       Standard_Real  theTab[],
274                                                       gp_Pnt        &theVertex)
275 {
276   GProp_GProps  aProps;
277   gp_Pnt        aCenterMass;
278   Standard_Real aShapeSize;
279
280   if    (theShape.ShapeType() == TopAbs_VERTEX) {
281     aCenterMass = BRep_Tool::Pnt(TopoDS::Vertex(theShape));
282   } else if (theShape.ShapeType() == TopAbs_EDGE) {
283     BRepGProp::LinearProperties(theShape,  aProps);
284   } else if (theShape.ShapeType() == TopAbs_FACE) {
285     BRepGProp::SurfaceProperties(theShape, aProps);
286   } else {
287     BRepGProp::VolumeProperties(theShape,  aProps);
288   }
289
290   if (theShape.ShapeType() == TopAbs_VERTEX) {
291     aShapeSize = 1;
292   } else {
293     aCenterMass = aProps.CentreOfMass();
294     aShapeSize  = aProps.Mass();
295   }
296
297   theVertex = aCenterMass;
298   theTab[0] = theVertex.X();
299   theTab[1] = theVertex.Y();
300   theTab[2] = theVertex.Z();
301   theTab[3] = aShapeSize;
302 }
303
304 //=======================================================================
305 //function : GetInPlaceByHistory
306 //purpose  : 
307 //=======================================================================
308 Standard_Boolean GEOMAlgo_GetInPlaceAPI::GetInPlaceByHistory
309                       (const Handle(GEOM_Function)      &theWhereFunction,
310                        const TopTools_IndexedMapOfShape &theWhereIndices,
311                        const TopoDS_Shape               &theWhat,
312                              TopTools_ListOfShape       &theShapesInPlace)
313 {
314   if (theWhereFunction.IsNull() || theWhat.IsNull())
315     return Standard_False;
316
317   if (theWhereIndices.Contains(theWhat)) {
318     // entity was not changed by the operation
319     theShapesInPlace.Append(theWhat);
320
321     return Standard_True;
322   }
323
324   // try to find in history
325   //TDF_Label aHistoryLabel = theWhereFunction->GetHistoryEntry(Standard_False);
326
327   // search in history for all argument shapes
328   Standard_Boolean isFound = Standard_False;
329   Standard_Boolean isGood = Standard_False;
330
331   TDF_LabelSequence aLabelSeq;
332   theWhereFunction->GetDependency(aLabelSeq);
333   Standard_Integer nbArg = aLabelSeq.Length();
334
335   for (Standard_Integer iarg = 1; iarg <= nbArg && !isFound; iarg++) {
336
337     TDF_Label anArgumentRefLabel = aLabelSeq.Value(iarg);
338
339     Handle(GEOM_Object) anArgumentObject = GEOM_Object::GetReferencedObject(anArgumentRefLabel);
340     TopoDS_Shape anArgumentShape = anArgumentObject->GetValue();
341
342     TopTools_IndexedMapOfShape anArgumentIndices;
343     TopExp::MapShapes(anArgumentShape, anArgumentIndices);
344
345     if (anArgumentIndices.Contains(theWhat)) {
346       isFound = Standard_True;
347       Standard_Integer aWhatIndex = anArgumentIndices.FindIndex(theWhat);
348
349       // Find corresponding label in history
350       TDF_Label anArgumentHistoryLabel =
351         theWhereFunction->GetArgumentHistoryEntry(anArgumentRefLabel, Standard_False);
352       if (anArgumentHistoryLabel.IsNull()) {
353         // Lost History of operation argument. Possibly, all its entities was removed.
354         isGood = Standard_True;
355       }
356       else {
357         TDF_Label aWhatHistoryLabel = anArgumentHistoryLabel.FindChild(aWhatIndex, Standard_False);
358
359         if (aWhatHistoryLabel.IsNull()) {
360           // Removed entity ? Compound ? Compsolid ? Shell ? Wire
361           isGood = Standard_False;
362         } else {
363           Handle(TDataStd_IntegerArray) anIntegerArray;
364           if (!aWhatHistoryLabel.FindAttribute(TDataStd_IntegerArray::GetID(), anIntegerArray)) {
365             //Error: Empty modifications history for the sought shape.
366             isGood = Standard_False;
367           }
368           else {
369             isGood = Standard_True;
370             Standard_Integer imod, aModifLen = anIntegerArray->Array()->Length();
371             for (imod = 1; imod <= aModifLen; imod++) {
372               const Standard_Integer anIndex =
373                 anIntegerArray->Array()->Value(imod);
374
375               theShapesInPlace.Append(theWhereIndices.FindKey(anIndex));
376             }
377           }
378         }
379       }
380     }
381   }
382
383   isFound = isGood;
384
385   if (!isFound) {
386     // try compound/compsolid/shell/wire element by element
387     Standard_Boolean isFoundAny = Standard_False;
388     TopTools_MapOfShape mapShape;
389
390     if (theWhat.ShapeType() == TopAbs_COMPOUND ||
391         theWhat.ShapeType() == TopAbs_COMPSOLID) {
392       // recursive processing of compound/compsolid
393       TopoDS_Iterator anIt (theWhat, Standard_True, Standard_True);
394       for (; anIt.More(); anIt.Next()) {
395         if (mapShape.Add(anIt.Value())) {
396           TopoDS_Shape curWhat = anIt.Value();
397           isFoundAny = GetInPlaceByHistory(theWhereFunction, theWhereIndices, curWhat, theShapesInPlace);
398           if (isFoundAny) isFound = Standard_True;
399         }
400       }
401     }
402     else if (theWhat.ShapeType() == TopAbs_SHELL) {
403       // try to replace a shell by its faces images
404       TopExp_Explorer anExp (theWhat, TopAbs_FACE);
405       for (; anExp.More(); anExp.Next()) {
406         if (mapShape.Add(anExp.Current())) {
407           TopoDS_Shape curWhat = anExp.Current();
408           isFoundAny = GetInPlaceByHistory(theWhereFunction, theWhereIndices, curWhat, theShapesInPlace);
409           if (isFoundAny) isFound = Standard_True;
410         }
411       }
412     }
413     else if (theWhat.ShapeType() == TopAbs_WIRE) {
414       // try to replace a wire by its edges images
415       TopExp_Explorer anExp (theWhat, TopAbs_EDGE);
416       for (; anExp.More(); anExp.Next()) {
417         if (mapShape.Add(anExp.Current())) {
418           TopoDS_Shape curWhat = anExp.Current();
419           isFoundAny = GetInPlaceByHistory(theWhereFunction, theWhereIndices, curWhat, theShapesInPlace);
420           if (isFoundAny) isFound = Standard_True;
421         }
422       }
423     }
424     else {
425       // Removed entity
426     }
427   }
428
429   return isFound;
430 }
431
432 //=======================================================================
433 //function : GetInPlaceByHistory
434 //purpose  : 
435 //=======================================================================
436 Standard_Boolean
437 GEOMAlgo_GetInPlaceAPI::GetInPlaceMap (const Handle(GEOM_Function)       & theWhereFunction,
438                                        const TopoDS_Shape                & theWhat,
439                                        std::vector< std::vector< int > > & theResVec)
440 {
441   //theResVec.clear();
442
443   if (theWhereFunction.IsNull() || theWhat.IsNull() || theWhereFunction->GetValue().IsNull() )
444     return Standard_False;
445   TopoDS_Shape theWhere = theWhereFunction->GetValue();
446
447   TopTools_IndexedMapOfShape whereIndices, whatIndices;
448   TopExp::MapShapes( theWhere, whereIndices );
449   TopExp::MapShapes( theWhat,  whatIndices );
450
451   theResVec.resize( whatIndices.Extent() + 1 );
452
453   // first, try by history
454
455   Standard_Boolean isFound = Standard_False;
456
457   TDF_LabelSequence aLabelSeq;
458   theWhereFunction->GetDependency(aLabelSeq);
459   Standard_Integer nbArg = aLabelSeq.Length();
460
461   for (Standard_Integer iarg = 1; iarg <= nbArg && !isFound; iarg++)
462   {
463     TDF_Label anArgumentRefLabel = aLabelSeq.Value(iarg);
464
465     Handle(GEOM_Object) anArgumentObject = GEOM_Object::GetReferencedObject(anArgumentRefLabel);
466     TopoDS_Shape anArgumentShape = anArgumentObject->GetValue();
467
468     TopTools_IndexedMapOfShape anArgumentIndices;
469     TopExp::MapShapes(anArgumentShape, anArgumentIndices);
470
471     if (( isFound = anArgumentIndices.Contains(theWhat)))
472     {
473       // Find corresponding label in history
474       TDF_Label anArgumentHistoryLabel =
475         theWhereFunction->GetArgumentHistoryEntry(anArgumentRefLabel, Standard_False);
476       if ( anArgumentHistoryLabel.IsNull())
477       {
478         // Lost History of operation argument. Possibly, theWhat was removed from theWhere
479         isFound = false;
480         break;
481       }
482       else
483       {
484         Standard_Integer aWhatIndex = anArgumentIndices.FindIndex(theWhat);
485         for ( int i = 0, iSubWhat = aWhatIndex; i < whatIndices.Extent(); ++i, ++iSubWhat )
486         {
487           TDF_Label aHistoryLabel = anArgumentHistoryLabel.FindChild( iSubWhat, Standard_False );
488           if ( !aHistoryLabel.IsNull() )
489           {
490             Handle(TDataStd_IntegerArray) anIntegerArray;
491             if (aHistoryLabel.FindAttribute(TDataStd_IntegerArray::GetID(), anIntegerArray))
492             {
493               Standard_Integer imod, aModifLen = anIntegerArray->Array()->Length();
494               for ( imod = 1; imod <= aModifLen; imod++ )
495               {
496                 Standard_Integer     whereIndex = anIntegerArray->Array()->Value( imod );
497                 const TopoDS_Shape& argSubShape = anArgumentIndices( iSubWhat );
498                 Standard_Integer      whatIndex = whatIndices.FindIndex( argSubShape );
499                 theResVec[ whatIndex ].push_back( whereIndex );
500               }
501             }
502           }
503         }
504       }
505     }
506   }
507
508   if ( !isFound ) // use GetInPlace()
509   {
510     GEOMAlgo_GetInPlace gip;
511     if ( ! GetInPlace( theWhere, theWhat, gip ))
512       return false;
513
514     const TopTools_DataMapOfShapeListOfShape& img = gip.Images();
515     TopTools_DataMapIteratorOfDataMapOfShapeListOfShape imgIt( img );
516     for ( ; imgIt.More(); imgIt.Next() )
517     {
518       const TopoDS_Shape&              whatSub = imgIt.Key();
519       const TopTools_ListOfShape& whereSubList = imgIt.Value();
520       int whatID = whatIndices.FindIndex( whatSub );
521       if ( whatID == 0 ) continue;
522       TopTools_ListIteratorOfListOfShape whereIt( whereSubList );
523       for ( ; whereIt.More(); whereIt.Next() )
524       {
525         int whereID = whereIndices.FindIndex( whereIt.Value() );
526         if ( whereID && whatSub.ShapeType() == whereIt.Value().ShapeType() )
527           theResVec[ whatID ].push_back( whereID );
528       }
529       isFound = true;
530     }
531   }
532
533   if ( isFound )
534   {
535     // check that all sub-shapes are found and restore missing history of wires, shells etc.
536     for ( int iSubWhat = 1; iSubWhat <= whatIndices.Extent(); ++iSubWhat )
537     {
538       if ( !theResVec[ iSubWhat ].empty() )
539         continue;
540       const TopoDS_Shape& whatSubShape = whatIndices( iSubWhat );
541       switch ( whatSubShape.ShapeType() )
542       {
543       case TopAbs_COMPOUND:
544       case TopAbs_COMPSOLID:
545         continue; // not possible?
546       case TopAbs_SHELL:
547       case TopAbs_WIRE:
548       {
549         // find a corresponding sub-shape in theWhere
550         TColStd_MapOfInteger whereIDs;
551         TopAbs_ShapeEnum subType =
552           whatSubShape.ShapeType() == TopAbs_WIRE ? TopAbs_EDGE : TopAbs_FACE;
553         for ( TopExp_Explorer subIt( whatSubShape, subType ); subIt.More(); subIt.Next() )
554         {
555           int whatID = whatIndices.FindIndex( subIt.Current() );
556           std::vector< int > & whereIDsVec = theResVec[ whatID ];
557           for ( size_t i = 0; i < whereIDsVec.size(); ++i )
558             whereIDs.Add( whereIDsVec[i] );
559         }
560         // look for a Where sub-shape including all whereIDs
561         TopExp_Explorer whereIt( theWhere, whatSubShape.ShapeType() );
562         for ( ; whereIt.More(); whereIt.Next() )
563         {
564           bool isAllIn = true;
565           for ( TopoDS_Iterator it( whereIt.Current() ); it.More() && isAllIn; it.Next() )
566           {
567             int whereID = whereIndices.FindIndex( it.Value() );
568             isAllIn = whereIDs.Contains( whereID );
569           }
570           if ( isAllIn )
571           {
572             theResVec[ iSubWhat ].push_back( whereIndices.FindIndex( whereIt.Current() ));
573             break;
574           }
575         }
576       }
577       default:; // removed sub-shape
578       }
579     }
580   }
581
582   return isFound;
583 }