0024428: Implementation of LGPL license
[occt.git] / src / BOPAlgo / BOPAlgo_PaveFiller_6.cxx
1 // Created by: Peter KURNEV
2 // Copyright (c) 2010-2014 OPEN CASCADE SAS
3 // Copyright (c) 2007-2010 CEA/DEN, EDF R&D, OPEN CASCADE
4 // Copyright (c) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN, CEDRAT,
5 //                         EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 //
7 // This file is part of Open CASCADE Technology software library.
8 //
9 // This library is free software; you can redistribute it and / or modify it
10 // under the terms of the GNU Lesser General Public version 2.1 as published
11 // by the Free Software Foundation, with special exception defined in the file
12 // OCCT_LGPL_EXCEPTION.txt. Consult the file LICENSE_LGPL_21.txt included in OCCT
13 // distribution for complete text of the license and disclaimer of any warranty.
14 //
15 // Alternatively, this file may be used under the terms of Open CASCADE
16 // commercial license or contractual agreement.
17
18 #include <BOPAlgo_PaveFiller.ixx>
19
20 #include <Precision.hxx>
21 #include <NCollection_IncAllocator.hxx>
22 #include <Bnd_Box.hxx>
23
24 #include <Geom_Curve.hxx>
25 #include <Geom2d_Curve.hxx>
26
27 #include <GeomAPI_ProjectPointOnCurve.hxx>
28 #include <GeomAPI_ProjectPointOnSurf.hxx>
29
30 #include <TopoDS_Edge.hxx>
31 #include <TopoDS_Face.hxx>
32 #include <TopoDS_Vertex.hxx>
33 #include <TopoDS_Compound.hxx>
34
35 #include <TopExp_Explorer.hxx>
36
37 #include <BRep_Builder.hxx>
38 #include <BRep_Tool.hxx>
39
40 #include <BRepBndLib.hxx>
41 #include <BRepTools.hxx>
42
43 #include <BRepAdaptor_Curve.hxx>
44 #include <BRepAdaptor_Surface.hxx>
45
46 #include <IntTools_FaceFace.hxx>
47 #include <IntTools_SequenceOfCurves.hxx>
48 #include <IntTools_SequenceOfPntOn2Faces.hxx>
49 #include <IntTools_Curve.hxx>
50 #include <IntTools_PntOn2Faces.hxx>
51 #include <IntTools_Tools.hxx>
52
53 #include <IntSurf_ListOfPntOn2S.hxx>
54 #include <IntSurf_PntOn2S.hxx>
55
56 #include <BOPTools_AlgoTools.hxx>
57 #include <BOPTools_AlgoTools3D.hxx>
58
59 #include <BOPCol_MapOfInteger.hxx>
60 #include <BOPCol_ListOfShape.hxx>
61 #include <BOPCol_DataMapOfShapeInteger.hxx>
62 #include <BOPCol_ListOfInteger.hxx>
63 #include <BOPCol_IndexedMapOfInteger.hxx>
64 #include <BOPCol_DataMapOfIntegerReal.hxx>
65
66 #include <BOPInt_Context.hxx>
67 #include <BOPInt_Tools.hxx>
68
69 #include <BOPDS_Interf.hxx>
70 #include <BOPDS_Iterator.hxx>
71 #include <BOPDS_Curve.hxx>
72 #include <BOPDS_Point.hxx>
73 #include <BOPDS_FaceInfo.hxx>
74 #include <BOPDS_Curve.hxx>
75 #include <BOPDS_MapOfPaveBlock.hxx>
76 #include <BOPDS_PaveBlock.hxx>
77 #include <BOPDS_VectorOfCurve.hxx>
78 #include <BOPDS_VectorOfPoint.hxx>
79 #include <BOPDS_ShapeInfo.hxx>
80 #include <BOPDS_PaveBlock.hxx>
81 #include <BOPDS_ListOfPave.hxx>
82 #include <BOPDS_ListOfPaveBlock.hxx>
83 #include <BOPDS_CoupleOfPaveBlocks.hxx>
84 #include <BOPDS_FaceInfo.hxx>
85 #include <BOPDS_CommonBlock.hxx>
86
87 #include <BOPAlgo_Tools.hxx>
88 #include <BRepBuilderAPI_MakeVertex.hxx>
89 #include <TopExp.hxx>
90 #include <BOPInt_ShrunkRange.hxx>
91 #include <BOPDS_DataMapOfPaveBlockListOfPaveBlock.hxx>
92
93 static void ToleranceFF(const BRepAdaptor_Surface& aBAS1,
94                         const BRepAdaptor_Surface& aBAS2,
95                         Standard_Real& aTolFF);
96
97 //=======================================================================
98 //function : PerformFF
99 //purpose  : 
100 //=======================================================================
101 void BOPAlgo_PaveFiller::PerformFF()
102 {
103   Standard_Integer iSize;
104   Standard_Boolean bValid;
105   //
106   myErrorStatus=0;
107   //
108   myIterator->Initialize(TopAbs_FACE, TopAbs_FACE);
109   iSize=myIterator->ExpectedLength();
110   if (!iSize) {
111     return; 
112   }
113   //
114   Standard_Boolean bJustAdd, bApp, bCompC2D1, bCompC2D2, bIsDone;
115   Standard_Boolean bToSplit, bTangentFaces;
116   Standard_Integer nF1, nF2, aNbCurves, aNbPoints, iX, i, iP, iC, aNbLP;
117   Standard_Real aApproxTol, aTolR3D, aTolR2D, aTolFF;
118   BRepAdaptor_Surface aBAS1, aBAS2;
119   BOPCol_MapOfInteger aMI;
120   //
121   BOPDS_VectorOfInterfFF& aFFs=myDS->InterfFF();
122   aFFs.SetStartSize(iSize);
123   aFFs.SetIncrement(iSize);
124   aFFs.Init();
125   //
126   bApp=mySectionAttribute.Approximation();
127   bCompC2D1=mySectionAttribute.PCurveOnS1();
128   bCompC2D2=mySectionAttribute.PCurveOnS2();
129   aApproxTol=1.e-7;
130   bToSplit = Standard_False;
131   //
132   for (; myIterator->More(); myIterator->Next()) {
133     myIterator->Value(nF1, nF2, bJustAdd);
134     if(bJustAdd) {
135       continue;
136     }
137     //
138     const TopoDS_Face& aF1=(*(TopoDS_Face *)(&myDS->Shape(nF1)));
139     const TopoDS_Face& aF2=(*(TopoDS_Face *)(&myDS->Shape(nF2)));
140     //
141     aBAS1.Initialize(aF1, Standard_False);
142     aBAS2.Initialize(aF2, Standard_False);
143     //
144     if (aBAS1.GetType() == GeomAbs_Plane && 
145         aBAS2.GetType() == GeomAbs_Plane) {
146       Standard_Boolean bToIntersect;
147       //
148       if (aMI.Add(nF1)) {
149         myDS->UpdateFaceInfoOn(nF1);
150         myDS->UpdateFaceInfoIn(nF1);
151       }
152       if (aMI.Add(nF2)) {
153         myDS->UpdateFaceInfoOn(nF2);
154         myDS->UpdateFaceInfoIn(nF2);
155       }
156       //
157       bToIntersect = CheckPlanes(nF1, nF2);
158       if (!bToIntersect) {
159         myDS->AddInterf(nF1, nF2);
160         iX=aFFs.Append()-1;
161         BOPDS_InterfFF& aFF=aFFs(iX);
162         aFF.SetIndices(nF1, nF2);
163         aFF.Init(0, 0);
164         continue;
165       }
166     }
167     //
168     IntTools_FaceFace aFaceFace;
169     //
170     IntSurf_ListOfPntOn2S aListOfPnts;
171     GetEFPnts(nF1, nF2, aListOfPnts);
172     aNbLP = aListOfPnts.Extent();
173     if (aNbLP) {
174       aFaceFace.SetList(aListOfPnts);
175     }
176
177     aFaceFace.SetParameters(bApp, bCompC2D1, bCompC2D2, aApproxTol);
178     //
179     aFaceFace.Perform(aF1, aF2);
180     //
181     bIsDone=aFaceFace.IsDone();
182     if (bIsDone) {
183       aTolR3D=aFaceFace.TolReached3d();
184       aTolR2D=aFaceFace.TolReached2d();
185       bTangentFaces=aFaceFace.TangentFaces();
186       //
187       ToleranceFF(aBAS1, aBAS2, aTolFF);
188       //
189       if (aTolR3D < aTolFF){
190         aTolR3D=aTolFF;
191       }
192       if (aTolR2D < 1.e-7){
193         aTolR2D=1.e-7;
194       }
195       //
196       aFaceFace.PrepareLines3D(bToSplit);
197       //
198       const IntTools_SequenceOfCurves& aCvsX=aFaceFace.Lines();
199       const IntTools_SequenceOfPntOn2Faces& aPntsX=aFaceFace.Points();
200       //
201       aNbCurves=aCvsX.Length();
202       aNbPoints=aPntsX.Length();
203       //
204       myDS->AddInterf(nF1, nF2);
205       //
206       iX=aFFs.Append()-1;
207       BOPDS_InterfFF& aFF=aFFs(iX);
208       aFF.SetIndices(nF1, nF2);
209       //
210       aFF.SetTolR3D(aTolR3D);
211       aFF.SetTolR2D(aTolR2D);
212       aFF.SetTangentFaces(bTangentFaces);
213       //
214       // Curves, Points 
215       aFF.Init(aNbCurves, aNbPoints);
216       //
217       // Curves
218       BOPDS_VectorOfCurve& aVNC=aFF.ChangeCurves();
219       for (i=1; i<=aNbCurves; ++i) {
220         Bnd_Box aBox;
221         //
222         const IntTools_Curve& aIC=aCvsX(i);
223         const Handle(Geom_Curve)& aC3D= aIC.Curve();
224         bValid=BOPInt_Tools::CheckCurve(aC3D, aTolR3D, aBox);
225         if (bValid) {
226           iC=aVNC.Append()-1;
227           BOPDS_Curve& aNC=aVNC(iC);
228           aNC.SetCurve(aIC);
229           aNC.SetBox(aBox);
230         }
231       }
232       //
233       // Points
234       BOPDS_VectorOfPoint& aVNP=aFF.ChangePoints();
235       for (i=1; i<=aNbPoints; ++i) {
236         const IntTools_PntOn2Faces& aPi=aPntsX(i);
237         const gp_Pnt& aP=aPi.P1().Pnt();
238         //
239         iP=aVNP.Append()-1;
240         BOPDS_Point& aNP=aVNP(iP);
241         aNP.SetPnt(aP);
242       }
243     //}// if (aNbCs || aNbPs)
244     }// if (bIsDone) {
245     else {// 904/L1
246       iX=aFFs.Append()-1;
247       BOPDS_InterfFF& aFF=aFFs(iX);
248       aFF.SetIndices(nF1, nF2);
249       aNbCurves=0;
250       aNbPoints=0;
251       aFF.Init(aNbCurves, aNbPoints);
252     }
253   }// for (; myIterator->More(); myIterator->Next()) {
254 }
255 //=======================================================================
256 //function : MakeBlocks
257 //purpose  : 
258 //=======================================================================
259 void BOPAlgo_PaveFiller::MakeBlocks()
260 {
261   Standard_Integer aNbFF;
262   //
263   myErrorStatus=0;
264   //
265   BOPDS_VectorOfInterfFF& aFFs=myDS->InterfFF();
266   aNbFF=aFFs.Extent();
267   if (!aNbFF) {
268     return;
269   }
270   //
271   Standard_Boolean bExist, bValid2D;
272   Standard_Integer i, nF1, nF2, aNbC, aNbP, j;
273   Standard_Integer nV1, nV2;
274   Standard_Real aTolR3D, aTolR2D, aT1, aT2, aTol;
275   Handle(NCollection_IncAllocator) aAllocator;
276   BOPDS_ListIteratorOfListOfPaveBlock aItLPB;
277   TopoDS_Edge aES;
278   Handle(BOPDS_PaveBlock) aPBOut;
279   //
280   //-----------------------------------------------------scope f
281   aAllocator=new NCollection_IncAllocator();
282   //
283   BOPCol_ListOfInteger aLSE(aAllocator);
284   BOPCol_MapOfInteger aMVOnIn(100, aAllocator), aMF(100, aAllocator),
285                       aMVStick(100,aAllocator), aMVEF(100, aAllocator),
286                       aMVB(100, aAllocator), aMI(100, aAllocator);
287   BOPDS_MapOfPaveBlock aMPBOnIn(100, aAllocator),
288                        aMPBAdd(100, aAllocator);
289   BOPDS_ListOfPaveBlock aLPB(aAllocator);
290   BOPDS_IndexedDataMapOfShapeCoupleOfPaveBlocks aMSCPB(100, aAllocator); 
291   BOPCol_DataMapOfShapeInteger aMVI(100, aAllocator);
292   BOPDS_DataMapOfPaveBlockListOfPaveBlock aDMExEdges(100, aAllocator);
293   BOPCol_DataMapOfIntegerReal aMVTol(100, aAllocator);
294   BOPCol_DataMapIteratorOfDataMapOfIntegerReal aItMV;
295   BOPCol_DataMapOfIntegerInteger aDMI(100, aAllocator);
296   //
297   for (i=0; i<aNbFF; ++i) {
298     BOPDS_InterfFF& aFF=aFFs(i);
299     aFF.Indices(nF1, nF2);
300     //
301     BOPDS_VectorOfPoint& aVP=aFF.ChangePoints();
302     aNbP=aVP.Extent();
303     BOPDS_VectorOfCurve& aVC=aFF.ChangeCurves();
304     aNbC=aVC.Extent();
305     if (!aNbP && !aNbC) {
306       continue;
307     }
308     //
309     const TopoDS_Face& aF1=(*(TopoDS_Face *)(&myDS->Shape(nF1)));
310     const TopoDS_Face& aF2=(*(TopoDS_Face *)(&myDS->Shape(nF2)));
311     //
312     aTolR3D=aFF.TolR3D();
313     aTolR2D=aFF.TolR2D();
314     //
315     // Update face info
316     if (aMF.Add(nF1)) {
317       myDS->UpdateFaceInfoOn(nF1);
318     }
319     if (aMF.Add(nF2)) {
320       myDS->UpdateFaceInfoOn(nF2);
321     }
322    
323     BOPDS_FaceInfo& aFI1=myDS->ChangeFaceInfo(nF1);
324     BOPDS_FaceInfo& aFI2=myDS->ChangeFaceInfo(nF2);
325     //
326     aMVOnIn.Clear();
327     aMPBOnIn.Clear();
328     aMVB.Clear();
329     aMVTol.Clear();
330     //
331     myDS->VerticesOnIn(nF1, nF2, aMVOnIn, aMPBOnIn);
332     myDS->SharedEdges(nF1, nF2, aLSE, aAllocator);
333     
334     // 1. Treat Points
335     for (j=0; j<aNbP; ++j) {
336       TopoDS_Vertex aV;
337       BOPDS_CoupleOfPaveBlocks aCPB;
338       //
339       BOPDS_Point& aNP=aVP.ChangeValue(j);
340       const gp_Pnt& aP=aNP.Pnt();
341       //
342       bExist=IsExistingVertex(aP, aTolR3D, aMVOnIn);
343       if (!bExist) {
344         BOPTools_AlgoTools::MakeNewVertex(aP, aTolR3D, aV);
345         //
346         aCPB.SetIndexInterf(i);
347         aCPB.SetIndex(j);
348         aMSCPB.Add(aV, aCPB);
349       }
350     }
351     //
352     // 2. Treat Curves
353     aMVStick.Clear();
354     aMVEF.Clear();
355     GetStickVertices(nF1, nF2, aMVStick, aMVEF, aMI);
356     //
357     for (j=0; j<aNbC; ++j) {
358       BOPDS_Curve& aNC=aVC.ChangeValue(j);
359       const IntTools_Curve& aIC=aNC.Curve();
360       // DEBt
361       aNC.InitPaveBlock1();
362       //
363       PutPavesOnCurve(aMVOnIn, aTolR3D, aNC, nF1, nF2, aMI, aMVEF, aMVTol);
364       //
365       PutStickPavesOnCurve(aF1, aF2, aMI, aNC, aMVStick, aMVTol);
366       //904/F7
367       if (aNbC == 1) {
368         PutEFPavesOnCurve(aNC, aMI, aMVEF, aMVTol);
369       }
370       //
371       if (aIC.HasBounds()) {
372         PutBoundPaveOnCurve(aF1, aF2, aTolR3D, aNC, aMVOnIn, aMVB);
373       }
374     }//for (j=0; j<aNbC; ++j) {
375     //
376     // Put closing pave if needed
377     for (j=0; j<aNbC; ++j) {
378       BOPDS_Curve& aNC=aVC.ChangeValue(j);
379       PutClosingPaveOnCurve (aNC);
380     }
381     //
382     // 3. Make section edges
383     for (j=0; j<aNbC; ++j) {
384       BOPDS_Curve& aNC=aVC.ChangeValue(j);
385       const IntTools_Curve& aIC=aNC.Curve();
386       //
387       BOPDS_ListOfPaveBlock& aLPBC=aNC.ChangePaveBlocks();
388       Handle(BOPDS_PaveBlock)& aPB1=aNC.ChangePaveBlock1();
389       //
390       aLPB.Clear();
391       aPB1->Update(aLPB, Standard_False);
392       //
393       aItLPB.Initialize(aLPB);
394       for (; aItLPB.More(); aItLPB.Next()) {
395         Handle(BOPDS_PaveBlock)& aPB=aItLPB.ChangeValue();
396         aPB->Indices(nV1, nV2);
397         aPB->Range  (aT1, aT2);
398         //
399         if (fabs(aT1 - aT2) < Precision::PConfusion()) {
400           continue;
401         }
402         //
403         bValid2D=myContext->IsValidBlockForFaces(aT1, aT2, aIC, aF1, aF2, aTolR3D);
404         if (!bValid2D) {
405           continue;
406         }
407         //
408         bExist=IsExistingPaveBlock(aPB, aNC, aTolR3D, aLSE);
409         if (bExist) {
410           continue;
411         }
412         //
413         bExist=IsExistingPaveBlock(aPB, aNC, aTolR3D, aMPBOnIn, aPBOut);
414         if (bExist) {
415           if (aMPBAdd.Add(aPBOut)) {
416             Standard_Boolean bInBothFaces = Standard_True;
417             if (!myDS->IsCommonBlock(aPBOut)) {
418               Standard_Integer nE;
419               Standard_Real aTolE;
420               //
421               nE = aPBOut->Edge();
422               const TopoDS_Edge& aE = *(TopoDS_Edge*)&myDS->Shape(nE);
423               aTolE = BRep_Tool::Tolerance(aE);
424               if (aTolR3D > aTolE) {
425                 myDS->UpdateEdgeTolerance(nE, aTolR3D);
426               }
427               bInBothFaces = Standard_False;
428             } else {
429               bInBothFaces = (aFI1.PaveBlocksOn().Contains(aPBOut) ||
430                               aFI1.PaveBlocksIn().Contains(aPBOut))&&
431                              (aFI2.PaveBlocksOn().Contains(aPBOut) ||
432                               aFI2.PaveBlocksIn().Contains(aPBOut));
433             }
434             if (!bInBothFaces) {
435               PreparePostTreatFF(i, aPBOut, aMSCPB, aMVI, aVC);
436             }
437           }
438           continue;
439         }
440         //
441         // Make Edge
442         const TopoDS_Vertex& aV1=(*(TopoDS_Vertex *)(&myDS->Shape(nV1)));
443         const TopoDS_Vertex& aV2=(*(TopoDS_Vertex *)(&myDS->Shape(nV2)));
444         //
445         BOPTools_AlgoTools::MakeEdge (aIC, aV1, aT1, aV2, aT2, aTolR3D, aES);
446         BOPTools_AlgoTools::MakePCurve(aES, aF1, aF2, aIC, 
447                                        mySectionAttribute.PCurveOnS1(),
448                                        mySectionAttribute.PCurveOnS2());
449         //
450         if (BOPTools_AlgoTools::IsMicroEdge(aES, myContext)) {
451           continue;
452         }
453         //
454         // Append the Pave Block to the Curve j
455         aLPBC.Append(aPB);
456         //
457         // Keep info for post treatment 
458         BOPDS_CoupleOfPaveBlocks aCPB;
459         aCPB.SetIndexInterf(i);
460         aCPB.SetIndex(j);
461         aCPB.SetPaveBlock1(aPB);
462         //
463         aMSCPB.Add(aES, aCPB);
464         aMVI.Bind(aV1, nV1);
465         aMVI.Bind(aV2, nV2);
466         //
467         aMVTol.UnBind(nV1);
468         aMVTol.UnBind(nV2);
469       }
470       //
471       aLPBC.RemoveFirst();
472     }//for (j=0; j<aNbC; ++j) {
473     //back to previous tolerance values for unused vertices
474     aItMV.Initialize(aMVTol);
475     for (; aItMV.More(); aItMV.Next()) {
476       nV1 = aItMV.Key();
477       aTol = aItMV.Value();
478       //
479       const TopoDS_Vertex& aV = *(TopoDS_Vertex*)&myDS->Shape(nV1);
480       const Handle(BRep_TVertex)& TV = *((Handle(BRep_TVertex)*)&aV.TShape());
481       TV->Tolerance(aTol);
482     }
483     //
484     ProcessExistingPaveBlocks(i, aMPBOnIn, aMSCPB, aMVI, aMVB, aMPBAdd);
485   }//for (i=0; i<aNbFF; ++i) {
486   // 
487   // post treatment
488   myErrorStatus=PostTreatFF(aMSCPB, aMVI, aDMExEdges, aDMI, aAllocator);
489   if (myErrorStatus) {
490     return;
491   }
492   //
493   // update face info
494   UpdateFaceInfo(aDMExEdges);
495   //Update all pave blocks
496   UpdatePaveBlocks(aDMI);
497   //-----------------------------------------------------scope t
498   aMF.Clear();
499   aMVStick.Clear();
500   aMPBOnIn.Clear();
501   aMVOnIn.Clear();
502   aDMExEdges.Clear();
503   aMI.Clear();
504   aDMI.Clear();
505   aAllocator.Nullify();
506 }
507
508 //=======================================================================
509 //function : PostTreatFF
510 //purpose  : 
511 //=======================================================================
512   Standard_Integer BOPAlgo_PaveFiller::PostTreatFF
513     (BOPDS_IndexedDataMapOfShapeCoupleOfPaveBlocks& theMSCPB,
514      BOPCol_DataMapOfShapeInteger& aMVI,
515      BOPDS_DataMapOfPaveBlockListOfPaveBlock& aDMExEdges,
516      BOPCol_DataMapOfIntegerInteger& aDMI,
517      Handle(NCollection_BaseAllocator)& theAllocator)
518 {
519   Standard_Integer iRet, aNbS;
520   //
521   iRet=0;
522   aNbS=theMSCPB.Extent();
523   if (!aNbS) {
524     return iRet;
525   }
526   //
527   Standard_Boolean bHasPaveBlocks, bOld;
528   Standard_Integer iErr, nSx, nVSD, iX, iP, iC, j, nV, iV = 0, iE, k;
529   Standard_Integer jx, aNbLPBx;
530   Standard_Real aT;
531   TopAbs_ShapeEnum aType;
532   TopoDS_Shape aV, aE;
533   BOPCol_ListIteratorOfListOfShape aItLS;
534   BOPDS_ListIteratorOfListOfPaveBlock aItLPB;
535   BOPDS_PDS aPDS;
536   Handle(BOPDS_PaveBlock) aPB1;
537   BOPDS_Pave aPave[2], aPave1[2];
538   BOPDS_ShapeInfo aSI;
539   //
540   BOPCol_ListOfShape aLS(theAllocator);
541   BOPAlgo_PaveFiller aPF(theAllocator);
542   //
543   BOPDS_VectorOfInterfFF& aFFs=myDS->InterfFF();
544   //
545   // <-DEB f
546   //
547   // 0
548   if (aNbS==1) {
549     const TopoDS_Shape& aS=theMSCPB.FindKey(1);
550     const BOPDS_CoupleOfPaveBlocks &aCPB=theMSCPB.FindFromIndex(1);
551     
552     //
553     aType=aS.ShapeType();
554     if (aType==TopAbs_VERTEX) {
555       aSI.SetShapeType(aType);
556       aSI.SetShape(aS);
557       iV=myDS->Append(aSI);
558       //
559       iX=aCPB.IndexInterf();
560       iP=aCPB.Index();
561       BOPDS_InterfFF& aFF=aFFs(iX); 
562       BOPDS_VectorOfPoint& aVNP=aFF.ChangePoints();
563       BOPDS_Point& aNP=aVNP(iP);
564       aNP.SetIndex(iV);
565     }
566     else if (aType==TopAbs_EDGE) {
567       aPB1=aCPB.PaveBlock1();
568       //
569       if (aPB1->HasEdge()) {
570         BOPDS_ListOfPaveBlock aLPBx;
571         aLPBx.Append(aPB1);
572         aDMExEdges.Bind(aPB1, aLPBx);
573       } else {
574         aSI.SetShapeType(aType);
575         aSI.SetShape(aS);
576         iE=myDS->Append(aSI);
577         //
578         aPB1->SetEdge(iE);
579       }
580     }
581     return iRet;
582   }
583   //
584   // 1 prepare arguments
585   for (k=1; k<=aNbS; ++k) {
586     const TopoDS_Shape& aS=theMSCPB.FindKey(k);
587     aLS.Append(aS);
588   }
589   //
590   // 2 Fuse shapes
591   aPF.SetArguments(aLS);
592   aPF.Perform();
593   iErr=aPF.ErrorStatus();
594   if (iErr) {
595     iRet=1;
596     return iRet;
597   }
598   aPDS=aPF.PDS();
599   //
600   aItLS.Initialize(aLS);
601   for (; aItLS.More(); aItLS.Next()) {
602     const TopoDS_Shape& aSx=aItLS.Value();
603     nSx=aPDS->Index(aSx);
604     const BOPDS_ShapeInfo& aSIx=aPDS->ShapeInfo(nSx);
605     //
606     aType=aSIx.ShapeType();
607     //
608     if (aType==TopAbs_VERTEX) {
609       if (aPDS->HasShapeSD(nSx, nVSD)) {
610         aV=aPDS->Shape(nVSD);
611       }
612       else {
613         aV=aSx;
614       }
615       // index of new vertex in theDS -> iV
616       if (!aMVI.IsBound(aV)) {
617         aSI.SetShapeType(aType);
618         aSI.SetShape(aV);
619         iV=myDS->Append(aSI);
620         //
621         aMVI.Bind(aV, iV);
622       }
623       else {
624         iV=aMVI.Find(aV);
625       }
626       // update FF interference
627       const BOPDS_CoupleOfPaveBlocks &aCPB=theMSCPB.FindFromKey(aSx);
628       iX=aCPB.IndexInterf();
629       iP=aCPB.Index();
630       BOPDS_InterfFF& aFF=aFFs(iX);
631       BOPDS_VectorOfPoint& aVNP=aFF.ChangePoints();
632       BOPDS_Point& aNP=aVNP(iP);
633       aNP.SetIndex(iV);
634     }//if (aType==TopAbs_VERTEX) {
635     //
636     else if (aType==TopAbs_EDGE) {
637       bHasPaveBlocks=aPDS->HasPaveBlocks(nSx);
638       const BOPDS_CoupleOfPaveBlocks &aCPB=theMSCPB.FindFromKey(aSx);
639       iX=aCPB.IndexInterf();
640       iC=aCPB.Index();
641       aPB1=aCPB.PaveBlock1();
642       //
643       bOld = aPB1->HasEdge();
644       if (bOld) {
645         BOPDS_ListOfPaveBlock aLPBx;
646         aDMExEdges.Bind(aPB1, aLPBx);
647       }
648       //
649       if (!bHasPaveBlocks) {
650         if (bOld) {
651           aDMExEdges.ChangeFind(aPB1).Append(aPB1);
652         } else {
653           aSI.SetShapeType(aType);
654           aSI.SetShape(aSx);
655           iE=myDS->Append(aSI);
656           //
657           aPB1->SetEdge(iE);
658         }
659       }
660       else {
661         BOPDS_InterfFF& aFF=aFFs(iX);
662         BOPDS_VectorOfCurve& aVNC=aFF.ChangeCurves();
663         BOPDS_Curve& aNC=aVNC(iC);
664         BOPDS_ListOfPaveBlock& aLPBC=aNC.ChangePaveBlocks();
665         //
666         const BOPDS_ListOfPaveBlock& aLPBx=aPDS->PaveBlocks(nSx);
667         aNbLPBx=aLPBx.Extent();
668         //
669         if (bOld && !aNbLPBx) {
670           aDMExEdges.ChangeFind(aPB1).Append(aPB1);
671           continue;
672         }
673         //
674         if (!bOld) {
675           aItLPB.Initialize(aLPBC);
676           for (; aItLPB.More(); aItLPB.Next()) {
677             const Handle(BOPDS_PaveBlock)& aPBC=aItLPB.Value();
678             if (aPBC==aPB1) {
679               aLPBC.Remove(aItLPB);
680               break;
681             }
682           } 
683         }
684         //
685         if (!aNbLPBx) {
686           aE=aSx;
687           //
688           if (!aMVI.IsBound(aE)) {
689             aSI.SetShapeType(aType);
690             aSI.SetShape(aE);
691             iE=myDS->Append(aSI);
692             aMVI.Bind(aE, iE);
693           }
694           else {
695             iE=aMVI.Find(aE);
696           }
697           // append new PaveBlock to aLPBC
698           aPB1->SetEdge(iE);
699           aLPBC.Append(aPB1);
700         } // if (!aNbLPBx) {
701         //
702         else {
703           aItLPB.Initialize(aLPBx);
704           if (bOld) {
705             aPave1[0] = aPB1->Pave1();
706             aPave1[1] = aPB1->Pave2();
707           }
708           for (; aItLPB.More(); aItLPB.Next()) {
709             const Handle(BOPDS_PaveBlock)& aPBx=aItLPB.Value();
710             const Handle(BOPDS_PaveBlock) aPBRx=aPDS->RealPaveBlock(aPBx);
711             //
712             // update vertices of paves
713             aPave[0]=aPBx->Pave1();
714             aPave[1]=aPBx->Pave2();
715             for (j=0; j<2; ++j) {
716               jx = 0;
717               if (bOld) {
718                 aT = aPave[j].Parameter();
719                 if (aT == aPave1[0].Parameter()) {
720                   jx = 1;
721                 } else if (aT == aPave1[1].Parameter()) {
722                   jx = 2;
723                 }
724                 //
725                 if (jx) {
726                   iV = aPave1[jx-1].Index();
727                 }
728               } 
729               if (!jx) {
730                 nV=aPave[j].Index();
731                 aV=aPDS->Shape(nV);
732                 //
733                 if (!aMVI.IsBound(aV)) {// index of new vertex in theDS -> iV
734                   aSI.SetShapeType(TopAbs_VERTEX);
735                   aSI.SetShape(aV);
736                   iV=myDS->Append(aSI);
737                   aMVI.Bind(aV, iV);
738                 }
739                 else {
740                   iV=aMVI.Find(aV);
741                 }
742               }
743               const BOPDS_Pave& aP1 = !j ? aPB1->Pave1() : aPB1->Pave2();
744               if (aP1.Parameter() == aPave[j].Parameter() && aP1.Index() != iV) {
745                 aDMI.Bind(aP1.Index(), iV);
746               }
747               //
748               aPave[j].SetIndex(iV);
749             }
750             //
751             // add edge
752             aE=aPDS->Shape(aPBRx->Edge());
753             //
754             if (!aMVI.IsBound(aE)) {
755               aSI.SetShapeType(aType);
756               aSI.SetShape(aE);
757               iE=myDS->Append(aSI);
758               aMVI.Bind(aE, iE);
759             }
760             else {
761               iE=aMVI.Find(aE);
762             }
763             // append new PaveBlock to aLPBC
764             Handle(BOPDS_PaveBlock) aPBC=new BOPDS_PaveBlock();
765             //
766             aPBC->SetPave1(aPave[0]);
767             aPBC->SetPave2(aPave[1]);
768             aPBC->SetEdge(iE);
769             if (bOld) {
770               aPBC->SetOriginalEdge(aPB1->OriginalEdge());
771               aDMExEdges.ChangeFind(aPB1).Append(aPBC);
772             }
773             else {
774               aLPBC.Append(aPBC);
775             }
776           }
777         }
778       }
779     }//else if (aType==TopAbs_EDGE)
780   }//for (; aItLS.More(); aItLS.Next()) {
781   return iRet;
782 }
783
784 //=======================================================================
785 //function : UpdateFaceInfo
786 //purpose  : 
787 //=======================================================================
788 void BOPAlgo_PaveFiller::UpdateFaceInfo(BOPDS_DataMapOfPaveBlockListOfPaveBlock& theDME)
789 {
790   Standard_Integer i, j, nV1, nF1, nF2, 
791                    aNbFF, aNbC, aNbP, aNbS, aNbPBIn;
792   BOPDS_IndexedMapOfPaveBlock aMPBCopy;
793   BOPDS_ListIteratorOfListOfPaveBlock aItLPB;
794   //
795   BOPDS_VectorOfInterfFF& aFFs=myDS->InterfFF();
796   aNbFF=aFFs.Extent();
797   //
798   //1. Sections (curves, points);
799   for (i=0; i<aNbFF; ++i) {
800     BOPDS_InterfFF& aFF=aFFs(i);
801     aFF.Indices(nF1, nF2);
802     //
803     BOPDS_FaceInfo& aFI1=myDS->ChangeFaceInfo(nF1);
804     BOPDS_FaceInfo& aFI2=myDS->ChangeFaceInfo(nF2);
805     //
806     BOPDS_VectorOfCurve& aVNC=aFF.ChangeCurves();
807     aNbC=aVNC.Extent();
808     for (j=0; j<aNbC; ++j) {
809       BOPDS_Curve& aNC=aVNC(j);
810       BOPDS_ListOfPaveBlock& aLPBC=aNC.ChangePaveBlocks();
811       aItLPB.Initialize(aLPBC);
812       //
813       if (aItLPB.More() && theDME.IsBound(aLPBC.First())) {
814         const Handle(BOPDS_PaveBlock)& aPB=aLPBC.First();
815         BOPDS_ListOfPaveBlock& aLPB = theDME.ChangeFind(aPB);
816         UpdateExistingPaveBlocks(aPB, aLPB, nF1, nF2);
817         aLPBC.Clear();
818         continue;
819       }
820       //
821       for(; aItLPB.More(); aItLPB.Next()) {
822         const Handle(BOPDS_PaveBlock)& aPB=aItLPB.Value();
823         aFI1.ChangePaveBlocksSc().Add(aPB);
824         aFI2.ChangePaveBlocksSc().Add(aPB);
825       }
826     }
827     // VerticesSc
828     const BOPDS_VectorOfPoint& aVNP=aFF.Points();
829     aNbP=aVNP.Extent();
830     for (j=0; j<aNbP; ++j) {
831       const BOPDS_Point& aNP=aVNP(j);
832       nV1=aNP.Index();
833       aFI1.ChangeVerticesSc().Add(nV1);
834       aFI2.ChangeVerticesSc().Add(nV1);
835     }
836   }
837   //
838   //2. PaveBlocksIn
839   if (theDME.IsEmpty()) {
840     return;
841   }
842   //
843   aNbS=myDS->NbSourceShapes();
844   for (i=0; i<aNbS; ++i) {
845     const BOPDS_ShapeInfo& aSI=myDS->ShapeInfo(i);
846     if (aSI.ShapeType()!=TopAbs_FACE) {
847       continue;
848     }
849     if(!myDS->HasFaceInfo(i)) {
850       continue;
851     }
852     BOPDS_FaceInfo& aFI=myDS->ChangeFaceInfo(i);
853     //
854     BOPDS_IndexedMapOfPaveBlock& aMPBIn=aFI.ChangePaveBlocksIn();
855     aMPBCopy.Assign(aMPBIn);
856     aMPBIn.Clear();
857     //
858     aNbPBIn=aMPBCopy.Extent();
859     for (j=1; j<=aNbPBIn; ++j) {
860       const Handle(BOPDS_PaveBlock)& aPB = aMPBCopy(j);
861       if (theDME.IsBound(aPB)) {
862         const BOPDS_ListOfPaveBlock& aLPB = theDME.Find(aPB);
863         aItLPB.Initialize(aLPB);
864         for (; aItLPB.More(); aItLPB.Next()) {
865           const Handle(BOPDS_PaveBlock)& aPB1 = aItLPB.Value();
866           aMPBIn.Add(aPB1);
867         }
868       } else {
869         aMPBIn.Add(aPB);
870       }
871     }//for (j=1; j<=aNbPBIn; ++j) {
872   }//for (i=0; i<aNbS; ++i) {
873 }
874 //=======================================================================
875 //function : IsExistingVertex
876 //purpose  : 
877 //=======================================================================
878   Standard_Boolean BOPAlgo_PaveFiller::IsExistingVertex
879     (const gp_Pnt& aP,
880      const Standard_Real theTolR3D,
881      const BOPCol_MapOfInteger& aMVOnIn)const
882 {
883   Standard_Boolean bRet;
884   Standard_Integer nV, iFlag;
885   Standard_Real aTolV;
886   gp_Pnt aPV;
887   Bnd_Box aBoxP;
888   BOPCol_MapIteratorOfMapOfInteger aIt;
889   //
890   bRet=Standard_True;
891   //
892   aBoxP.Add(aP);
893   aBoxP.Enlarge(theTolR3D);
894   //
895   aIt.Initialize(aMVOnIn);
896   for (; aIt.More(); aIt.Next()) {
897     Bnd_Box aBoxV;
898     //
899     nV=aIt.Value();
900     const TopoDS_Vertex& aV=(*(TopoDS_Vertex *)(&myDS->Shape(nV)));
901     aPV=BRep_Tool::Pnt(aV);
902     aTolV=BRep_Tool::Tolerance(aV);
903     aBoxV.Add(aP);
904     aBoxV.Enlarge(aTolV);
905     //
906     if (!aBoxP.IsOut(aBoxV)) {
907       iFlag=BOPTools_AlgoTools::ComputeVV(aV, aP, theTolR3D);
908       if (!iFlag) {
909         return bRet;
910       }
911     }
912   }
913   return !bRet;
914 }
915 //=======================================================================
916 //function : IsExistingPaveBlock
917 //purpose  : 
918 //=======================================================================
919   Standard_Boolean BOPAlgo_PaveFiller::IsExistingPaveBlock
920     (const Handle(BOPDS_PaveBlock)& thePB,
921      const BOPDS_Curve& theNC,
922      const Standard_Real theTolR3D,
923      const BOPCol_ListOfInteger& theLSE)
924 {
925   Standard_Boolean bRet=Standard_True;
926   //
927   if (theLSE.IsEmpty()) {
928     return !bRet;
929   } 
930   //
931   Standard_Real aT1, aT2, aTm, aTx, aTol;
932   Standard_Integer nE, iFlag;
933   gp_Pnt aPm;
934   Bnd_Box aBoxPm;
935   BOPCol_ListIteratorOfListOfInteger aItLI;
936   //
937   thePB->Range(aT1, aT2);
938   aTm=IntTools_Tools::IntermediatePoint (aT1, aT2);
939   theNC.Curve().D0(aTm, aPm);
940   aBoxPm.Add(aPm);
941   aBoxPm.Enlarge(theTolR3D);
942   //
943   aItLI.Initialize(theLSE);
944   for (; aItLI.More(); aItLI.Next()) {
945     nE=aItLI.Value();
946     const BOPDS_ShapeInfo& aSIE=myDS->ChangeShapeInfo(nE);
947     const Bnd_Box& aBoxE=aSIE.Box();
948     if (!aBoxE.IsOut(aBoxPm)) {
949       const TopoDS_Edge& aE=(*(TopoDS_Edge *)(&aSIE.Shape()));
950       aTol = BRep_Tool::Tolerance(aE);
951       aTol = aTol > theTolR3D ? aTol : theTolR3D;
952       iFlag=myContext->ComputePE(aPm, aTol, aE, aTx);
953       if (!iFlag) {
954         return bRet;
955       }
956     }
957   }
958   return !bRet;
959 }
960
961 //=======================================================================
962 //function : IsExistingPaveBlock
963 //purpose  : 
964 //=======================================================================
965   Standard_Boolean BOPAlgo_PaveFiller::IsExistingPaveBlock
966     (const Handle(BOPDS_PaveBlock)& thePB,
967      const BOPDS_Curve& theNC,
968      const Standard_Real theTolR3D,
969      const BOPDS_MapOfPaveBlock& theMPBOnIn,
970      Handle(BOPDS_PaveBlock&) aPBOut)
971 {
972   Standard_Boolean bRet;
973   Standard_Real aT1, aT2, aTm, aTx;
974   Standard_Integer nSp, iFlag1, iFlag2, nV11, nV12, nV21, nV22;
975   gp_Pnt aP1, aPm, aP2;
976   Bnd_Box aBoxP1, aBoxPm, aBoxP2;
977   BOPDS_MapIteratorOfMapOfPaveBlock aIt;
978   //
979   bRet=Standard_False;
980   const IntTools_Curve& aIC=theNC.Curve();
981   //
982   thePB->Range(aT1, aT2);
983   thePB->Indices(nV11, nV12);
984   //first point
985   aIC.D0(aT1, aP1);
986   aBoxP1.Add(aP1);
987   aBoxP1.Enlarge(theTolR3D);
988   //intermediate point
989   aTm=IntTools_Tools::IntermediatePoint (aT1, aT2);
990   aIC.D0(aTm, aPm);
991   aBoxPm.Add(aPm);
992   aBoxPm.Enlarge(theTolR3D);
993   //last point
994   aIC.D0(aT2, aP2);
995   aBoxP2.Add(aP2);
996   aBoxP2.Enlarge(theTolR3D);
997   //
998   aIt.Initialize(theMPBOnIn);
999   for (; aIt.More(); aIt.Next()) {
1000     const Handle(BOPDS_PaveBlock)& aPB=aIt.Value();
1001     aPB->Indices(nV21, nV22);
1002     nSp=aPB->Edge();
1003     const BOPDS_ShapeInfo& aSISp=myDS->ChangeShapeInfo(nSp);
1004     const TopoDS_Edge& aSp=(*(TopoDS_Edge *)(&aSISp.Shape()));
1005     const Bnd_Box& aBoxSp=aSISp.Box();
1006     //
1007     iFlag1 = (nV11 == nV21 || nV11 == nV22) ? 2 : (!aBoxSp.IsOut(aBoxP1) ? 1 : 0);
1008     iFlag2 = (nV12 == nV21 || nV12 == nV22) ? 2 : (!aBoxSp.IsOut(aBoxP2) ? 1 : 0);
1009     if (iFlag1 && iFlag2) {
1010       if (aBoxSp.IsOut(aBoxPm) || myContext->ComputePE(aPm, theTolR3D, aSp, aTx)) {
1011         continue;
1012       }
1013       //
1014       if (iFlag1 == 1) {
1015         iFlag1 = !myContext->ComputePE(aP1, theTolR3D, aSp, aTx);
1016       }
1017       //
1018       if (iFlag2 == 1) {
1019         iFlag2 = !myContext->ComputePE(aP2, theTolR3D, aSp, aTx);
1020       }
1021       //
1022       if (iFlag1 && iFlag2) {
1023         aPBOut = aPB;
1024         bRet=Standard_True;
1025         break;
1026       }
1027     }
1028   }
1029   return bRet;
1030 }
1031
1032 //=======================================================================
1033 //function : PutBoundPaveOnCurve
1034 //purpose  : 
1035 //=======================================================================
1036   void BOPAlgo_PaveFiller::PutBoundPaveOnCurve(const TopoDS_Face& aF1,
1037                                                const TopoDS_Face& aF2,
1038                                                const Standard_Real aTolR3D,
1039                                                BOPDS_Curve& aNC,
1040                                                BOPCol_MapOfInteger& aMVOnIn,
1041                                                BOPCol_MapOfInteger& aMVB)
1042 {
1043   Standard_Boolean bVF;
1044   Standard_Integer nV, iFlag, nVn, j, aNbEP;
1045   Standard_Real aT[2], aTmin, aTmax, aTV, aTol, aTolVnew;
1046   gp_Pnt aP[2];
1047   TopoDS_Vertex aVn;
1048   BOPDS_ListIteratorOfListOfPave aItLP;
1049   BOPDS_Pave aPn, aPMM[2];
1050   //
1051   aTolVnew = Precision::Confusion();
1052   //
1053   const IntTools_Curve& aIC=aNC.Curve();
1054   aIC.Bounds(aT[0], aT[1], aP[0], aP[1]);
1055   //
1056   Handle(BOPDS_PaveBlock)& aPB=aNC.ChangePaveBlock1();
1057   const BOPDS_ListOfPave& aLP=aPB->ExtPaves();
1058   //
1059   aNbEP=aLP.Extent();
1060   if (aNbEP) {
1061     aTmin=1.e10;
1062     aTmax=-aTmin;
1063     //
1064     aItLP.Initialize(aLP);
1065     for (; aItLP.More(); aItLP.Next()) {
1066       const BOPDS_Pave& aPv=aItLP.Value();
1067       aPv.Contents(nV, aTV);
1068       if (aTV<aTmin) {
1069         aPMM[0]=aPv;
1070         aTmin=aTV;
1071       }
1072       if (aTV>aTmax) {
1073         aPMM[1]=aPv;
1074         aTmax=aTV;
1075       }
1076     }
1077   }
1078   //
1079   for (j=0; j<2; ++j) {
1080     //if curve is closed, process only one bound
1081     if (j && aP[1].IsEqual(aP[0], aTolVnew)) {
1082       continue;
1083     }
1084     //
1085     iFlag=1;
1086     //
1087     if (aNbEP) {
1088       Bnd_Box aBoxP;
1089       //
1090       aBoxP.Set(aP[j]);
1091       aTol = aTolR3D+Precision::Confusion();
1092       aBoxP.Enlarge(aTol);
1093       const BOPDS_Pave& aPV=aPMM[j];
1094       nV=aPV.Index();
1095       const BOPDS_ShapeInfo& aSIV=myDS->ShapeInfo(nV);
1096       const TopoDS_Vertex& aV=(*(TopoDS_Vertex *)(&aSIV.Shape()));
1097       const Bnd_Box& aBoxV=aSIV.Box();
1098       if (!aBoxP.IsOut(aBoxV)){
1099         iFlag=BOPTools_AlgoTools::ComputeVV(aV, aP[j], aTol);
1100       }
1101     }
1102     if (iFlag) {
1103       // 900/L5
1104       bVF=myContext->IsValidPointForFaces (aP[j], aF1, aF2, aTolR3D);
1105       if (!bVF) {
1106         continue;
1107       }
1108       //
1109       BOPDS_ShapeInfo aSIVn;
1110       //
1111       BOPTools_AlgoTools::MakeNewVertex(aP[j], aTolR3D, aVn);
1112       aSIVn.SetShapeType(TopAbs_VERTEX);
1113       aSIVn.SetShape(aVn);
1114       //
1115       nVn=myDS->Append(aSIVn);
1116       //
1117       aPn.SetIndex(nVn);
1118       aPn.SetParameter(aT[j]);
1119       aPB->AppendExtPave(aPn);
1120       //
1121       aVn=(*(TopoDS_Vertex *)(&myDS->Shape(nVn)));
1122       BOPTools_AlgoTools::UpdateVertex (aIC, aT[j], aVn);
1123       //
1124       aTolVnew = BRep_Tool::Tolerance(aVn);
1125       //
1126       BOPDS_ShapeInfo& aSIDS=myDS->ChangeShapeInfo(nVn);
1127       Bnd_Box& aBoxDS=aSIDS.ChangeBox();
1128       BRepBndLib::Add(aVn, aBoxDS);
1129       aMVOnIn.Add(nVn);
1130       aMVB.Add(nVn);
1131     }
1132   }
1133 }
1134
1135 //=======================================================================
1136 //function : PutPavesOnCurve
1137 //purpose  : 
1138 //=======================================================================
1139   void BOPAlgo_PaveFiller::PutPavesOnCurve(const BOPCol_MapOfInteger& aMVOnIn,
1140                                            const Standard_Real aTolR3D,
1141                                            BOPDS_Curve& aNC,
1142                                            const Standard_Integer nF1,
1143                                            const Standard_Integer nF2,
1144                                            const BOPCol_MapOfInteger& aMI,
1145                                            const BOPCol_MapOfInteger& aMVEF,
1146                                            BOPCol_DataMapOfIntegerReal& aMVTol)
1147 {
1148   Standard_Boolean bInBothFaces;
1149   Standard_Integer nV;
1150   BOPCol_MapIteratorOfMapOfInteger aIt;
1151   //
1152   const Bnd_Box& aBoxC=aNC.Box();
1153   //
1154   //Put EF vertices first
1155   aIt.Initialize(aMVEF);
1156   for (; aIt.More(); aIt.Next()) {
1157     nV=aIt.Value();
1158     PutPaveOnCurve(nV, aTolR3D, aNC, aMI, aMVTol, 2);
1159   }
1160   //Put all other vertices
1161   aIt.Initialize(aMVOnIn);
1162   for (; aIt.More(); aIt.Next()) {
1163     nV=aIt.Value();
1164     if (aMVEF.Contains(nV)) {
1165       continue;
1166     }
1167     //
1168     const BOPDS_ShapeInfo& aSIV=myDS->ShapeInfo(nV);
1169     const Bnd_Box& aBoxV=aSIV.Box();
1170     //
1171     if (aBoxC.IsOut(aBoxV)){
1172       continue;
1173     }
1174     if (!myDS->IsNewShape(nV)) {
1175       const BOPDS_FaceInfo& aFI1 = myDS->FaceInfo(nF1);
1176       const BOPDS_FaceInfo& aFI2 = myDS->FaceInfo(nF2);
1177       //
1178       bInBothFaces = (aFI1.VerticesOn().Contains(nV) ||
1179                       aFI1.VerticesIn().Contains(nV))&&
1180                      (aFI2.VerticesOn().Contains(nV) ||
1181                       aFI2.VerticesIn().Contains(nV));
1182       if (!bInBothFaces) {
1183         continue;
1184       }
1185     }
1186     //
1187     PutPaveOnCurve(nV, aTolR3D, aNC, aMI, aMVTol, 1);
1188   }
1189 }
1190
1191 //=======================================================================
1192 //function : ExtendedTolerance
1193 //purpose  : 
1194 //=======================================================================
1195 Standard_Boolean BOPAlgo_PaveFiller::ExtendedTolerance(const Standard_Integer nV,
1196                                                        const BOPCol_MapOfInteger& aMI,
1197                                                        Standard_Real& aTolVExt,
1198                                                        const Standard_Integer aType)
1199 {
1200   Standard_Boolean bFound = Standard_False;
1201   if (!(myDS->IsNewShape(nV))) {
1202     return bFound;
1203   }
1204   //
1205   Standard_Integer i, k, aNbLines, aNbInt;
1206   Standard_Real aT11, aT12, aD1, aD2, aD;
1207   TopoDS_Vertex aV;
1208   gp_Pnt aPV, aP11, aP12;
1209   //
1210   k = 0;
1211   aNbInt = 2;
1212   if (aType == 1) {
1213     aNbInt = 1;
1214   } else if (aType == 2) {
1215     k = 1;
1216   }
1217   //
1218   aV = (*(TopoDS_Vertex *)(&myDS->Shape(nV)));
1219   aPV=BRep_Tool::Pnt(aV);
1220   //
1221   BOPDS_VectorOfInterfEE& aEEs=myDS->InterfEE();
1222   BOPDS_VectorOfInterfEF& aEFs=myDS->InterfEF();
1223   //
1224   for (; k<aNbInt; ++k) {
1225     aNbLines = !k ? aEEs.Extent() : aEFs.Extent();
1226     for (i = 0; i < aNbLines; ++i) {
1227       BOPDS_Interf *aInt = !k ? (BOPDS_Interf*) (&aEEs(i)) :
1228         (BOPDS_Interf*) (&aEFs(i));
1229       if (aInt->IndexNew() == nV) {
1230         if (aMI.Contains(aInt->Index1()) && aMI.Contains(aInt->Index2())) {
1231           const IntTools_CommonPrt& aComPrt = !k ? aEEs(i).CommonPart() :
1232             aEFs(i).CommonPart();
1233           //
1234           const TopoDS_Edge& aE1=aComPrt.Edge1();
1235           aComPrt.Range1(aT11, aT12);
1236           BOPTools_AlgoTools::PointOnEdge(aE1, aT11, aP11);
1237           BOPTools_AlgoTools::PointOnEdge(aE1, aT12, aP12);
1238           aD1=aPV.Distance(aP11);
1239           aD2=aPV.Distance(aP12);
1240           aD=(aD1>aD2)? aD1 : aD2;
1241           if (aD>aTolVExt) {
1242             aTolVExt=aD;
1243           }
1244           return !bFound;
1245         }//if (aMI.Contains(aEF.Index1()) && aMI.Contains(aEF.Index2())) {
1246       }//if (aInt->IndexNew() == nV) {
1247     }//for (i = 0; i < aNbLines; ++i) {
1248   }//for (k=0; k<2; ++k) {
1249   return bFound;
1250 }
1251
1252 //=======================================================================
1253 //function : GetEFPnts
1254 //purpose  : 
1255 //=======================================================================
1256   void BOPAlgo_PaveFiller::GetEFPnts(const Standard_Integer nF1,
1257                                      const Standard_Integer nF2,
1258                                      IntSurf_ListOfPntOn2S& aListOfPnts)
1259 {
1260   Standard_Integer nE, nF, nFOpposite, aNbEFs, i;
1261   Standard_Real U1, U2, V1, V2, f, l;
1262   BOPCol_MapOfInteger aMI;
1263   //
1264   //collect indexes of all shapes from nF1 and nF2.
1265   GetFullShapeMap(nF1, aMI);
1266   GetFullShapeMap(nF2, aMI);
1267   //
1268   BOPDS_VectorOfInterfEF& aEFs=myDS->InterfEF();
1269   aNbEFs = aEFs.Extent();
1270   //
1271   for(i = 0; i < aNbEFs; ++i) {
1272     const BOPDS_InterfEF& aEF = aEFs(i);
1273     if (aEF.HasIndexNew()) {
1274       aEF.Indices(nE, nFOpposite);
1275       if(aMI.Contains(nE) && aMI.Contains(nFOpposite)) {
1276         const IntTools_CommonPrt& aCP = aEF.CommonPart();
1277         Standard_Real aPar = aCP.VertexParameter1();
1278         const TopoDS_Edge& aE = (*(TopoDS_Edge*)(&myDS->Shape(nE)));
1279         const TopoDS_Face& aFOpposite = (*(TopoDS_Face*)(&myDS->Shape(nFOpposite)));
1280         //
1281         const Handle(Geom_Curve)& aCurve = BRep_Tool::Curve(aE, f, l);
1282         //
1283         nF = (nFOpposite == nF1) ? nF2 : nF1;
1284         const TopoDS_Face& aF = (*(TopoDS_Face*)(&myDS->Shape(nF)));
1285         Handle(Geom2d_Curve) aPCurve = BRep_Tool::CurveOnSurface(aE, aF, f, l);
1286         //
1287         GeomAPI_ProjectPointOnSurf& aProj = myContext->ProjPS(aFOpposite);
1288         //
1289         gp_Pnt aPoint;
1290         aCurve->D0(aPar, aPoint);
1291         IntSurf_PntOn2S aPnt;
1292         if(!aPCurve.IsNull()) {
1293           gp_Pnt2d aP2d = aPCurve->Value(aPar);
1294           aProj.Perform(aPoint);
1295           if(aProj.IsDone()) {
1296             aProj.LowerDistanceParameters(U1,V1);
1297             if (nF == nF1) {
1298               aPnt.SetValue(aP2d.X(),aP2d.Y(),U1,V1);
1299             } else {
1300               aPnt.SetValue(U1,V1,aP2d.X(),aP2d.Y());
1301             }
1302             aListOfPnts.Append(aPnt);
1303           }
1304         }
1305         else {
1306           GeomAPI_ProjectPointOnSurf& aProj1 = myContext->ProjPS(aF);
1307           aProj1.Perform(aPoint);
1308           aProj.Perform(aPoint);
1309           if(aProj1.IsDone() && aProj.IsDone()){
1310             aProj1.LowerDistanceParameters(U1,V1);
1311             aProj.LowerDistanceParameters(U2,V2);
1312             if (nF == nF1) {
1313               aPnt.SetValue(U1,V1,U2,V2);
1314             } else {
1315               aPnt.SetValue(U2,V2,U1,V1);
1316             }
1317             aListOfPnts.Append(aPnt);
1318           }
1319         }
1320       }
1321     }
1322   }
1323 }
1324
1325 //=======================================================================
1326 //function : ProcessUnUsedVertices
1327 //purpose  : 
1328 //=======================================================================
1329   void BOPAlgo_PaveFiller::PutEFPavesOnCurve(BOPDS_Curve& aNC,
1330                                              const BOPCol_MapOfInteger& aMI,
1331                                              const BOPCol_MapOfInteger& aMVEF,
1332                                              BOPCol_DataMapOfIntegerReal& aMVTol)
1333 {
1334   if (!aMVEF.Extent()) {
1335     return;
1336   }
1337   //
1338   const IntTools_Curve& aIC=aNC.Curve();
1339   GeomAbs_CurveType aTypeC;
1340   aTypeC=aIC.Type();
1341   if (!(aTypeC==GeomAbs_BezierCurve || aTypeC==GeomAbs_BSplineCurve)) {
1342     return;
1343   }
1344   //
1345   Standard_Integer nV;
1346   BOPCol_MapOfInteger aMV;
1347   //
1348   aMV.Assign(aMVEF);
1349   RemoveUsedVertices(aNC, aMV);
1350   if (!aMV.Extent()) {
1351     return;
1352   }
1353   //
1354   Standard_Real aDist;
1355   BOPDS_Pave aPave;
1356   //
1357   const Handle(Geom_Curve)& aC3D=aIC.Curve();
1358   GeomAPI_ProjectPointOnCurve& aProjPT = myContext->ProjPT(aC3D);
1359   //
1360   BOPCol_MapIteratorOfMapOfInteger aItMI;
1361   aItMI.Initialize(aMV);
1362   for (; aItMI.More(); aItMI.Next()) {
1363     nV = aItMI.Value();
1364     const TopoDS_Vertex& aV = (*(TopoDS_Vertex *)(&myDS->Shape(nV)));
1365     gp_Pnt aPV = BRep_Tool::Pnt(aV);
1366     aProjPT.Perform(aPV);
1367     Standard_Integer aNbPoints = aProjPT.NbPoints();
1368     if (aNbPoints) {
1369       aDist = aProjPT.LowerDistance();
1370       PutPaveOnCurve(nV, aDist, aNC, aMI, aMVTol);
1371     }
1372   }
1373 }
1374
1375 //=======================================================================
1376 //function : ProcessUnUsedVertices
1377 //purpose  : 
1378 //=======================================================================
1379   void BOPAlgo_PaveFiller::PutStickPavesOnCurve(const TopoDS_Face& aF1,
1380                                                 const TopoDS_Face& aF2,
1381                                                 const BOPCol_MapOfInteger& aMI,
1382                                                 BOPDS_Curve& aNC,
1383                                                 const BOPCol_MapOfInteger& aMVStick,
1384                                                 BOPCol_DataMapOfIntegerReal& aMVTol)
1385 {
1386   BOPCol_MapOfInteger aMV;
1387   aMV.Assign(aMVStick);
1388   RemoveUsedVertices(aNC, aMV);
1389   //
1390   if (!aMV.Extent()) {
1391     return;
1392   }
1393   //
1394   GeomAbs_SurfaceType aType1, aType2;
1395   Handle(Geom_Surface) aS1=BRep_Tool::Surface(aF1);
1396   Handle(Geom_Surface) aS2=BRep_Tool::Surface(aF2);
1397   GeomAdaptor_Surface aGAS1(aS1);
1398   GeomAdaptor_Surface aGAS2(aS2);
1399   //
1400   aType1=aGAS1.GetType();
1401   aType2=aGAS2.GetType();
1402   //
1403   if (aType1==GeomAbs_Torus  || aType2==GeomAbs_Torus) {
1404     GeomAbs_CurveType aTypeC;
1405     //
1406     const IntTools_Curve& aIC=aNC.Curve();
1407     aTypeC=aIC.Type();
1408     if (aTypeC==GeomAbs_BezierCurve || aTypeC==GeomAbs_BSplineCurve) {
1409       Handle(Geom2d_Curve) aC2D[2];
1410       //
1411       aC2D[0]=aIC.FirstCurve2d();
1412       aC2D[1]=aIC.SecondCurve2d();
1413       if (!aC2D[0].IsNull() && !aC2D[1].IsNull()) {
1414         Standard_Integer nV, m, n;
1415         Standard_Real aTC[2], aD, aD2, u, v, aDT2, aScPr, aDScPr;
1416         gp_Pnt aPC[2], aPV;
1417         gp_Dir aDN[2];
1418         gp_Pnt2d aP2D;
1419         BOPCol_MapIteratorOfMapOfInteger aItMI, aItMI1;
1420         //
1421         aDT2=2e-7;     // the rich criteria
1422         aDScPr=5.e-9;  // the creasing criteria
1423         aIC.Bounds(aTC[0], aTC[1], aPC[0], aPC[1]);
1424         //
1425         aItMI.Initialize(aMV);
1426         for (; aItMI.More(); aItMI.Next()) {
1427           nV = aItMI.Value();
1428           const TopoDS_Vertex& aV=*((TopoDS_Vertex*)&myDS->Shape(nV));
1429           aPV=BRep_Tool::Pnt(aV);
1430           //
1431           for (m=0; m<2; ++m) {
1432             aD2=aPC[m].SquareDistance(aPV);
1433             if (aD2>aDT2) {// no rich
1434               continue; 
1435             }
1436             //
1437             for (n=0; n<2; ++n) {
1438               Handle(Geom_Surface)& aS=(!n)? aS1 : aS2;
1439               aC2D[n]->D0(aTC[m], aP2D);
1440               aP2D.Coord(u, v);
1441               BOPTools_AlgoTools3D::GetNormalToSurface(aS, u, v, aDN[n]);
1442             }
1443             // 
1444             aScPr=aDN[0]*aDN[1];
1445             if (aScPr<0.) {
1446               aScPr=-aScPr;
1447             }
1448             aScPr=1.-aScPr;
1449             //
1450             if (aScPr>aDScPr) {
1451               continue;
1452             }
1453             //
1454             // The intersection curve aIC is vanishing curve (the crease)
1455             aD=sqrt(aD2);
1456             //
1457             PutPaveOnCurve(nV, aD, aNC, aMI, aMVTol);
1458           }
1459         }//for (jVU=1; jVU=aNbVU; ++jVU) {
1460       }
1461     }//if (aTypeC==GeomAbs_BezierCurve || aTypeC==GeomAbs_BSplineCurve) {
1462   }//if(aType1==GeomAbs_Torus  || aType2==GeomAbs_Torus) {
1463 }
1464
1465 //=======================================================================
1466 //function : GetStickVertices
1467 //purpose  : 
1468 //=======================================================================
1469   void BOPAlgo_PaveFiller::GetStickVertices(const Standard_Integer nF1,
1470                                             const Standard_Integer nF2,
1471                                             BOPCol_MapOfInteger& aMVStick,
1472                                             BOPCol_MapOfInteger& aMVEF,
1473                                             BOPCol_MapOfInteger& aMI)
1474 {
1475   Standard_Integer nS1, nS2, nVNew, aTypeInt, i;
1476   //
1477   BOPDS_VectorOfInterfVV& aVVs=myDS->InterfVV();
1478   BOPDS_VectorOfInterfVE& aVEs=myDS->InterfVE();
1479   BOPDS_VectorOfInterfEE& aEEs=myDS->InterfEE();
1480   BOPDS_VectorOfInterfVF& aVFs=myDS->InterfVF();
1481   BOPDS_VectorOfInterfEF& aEFs=myDS->InterfEF();
1482   //
1483   Standard_Integer aNbLines[5] = {aVVs.Extent(), aVEs.Extent(), aEEs.Extent(),
1484                                   aVFs.Extent(), aEFs.Extent()};
1485   //collect indices of all shapes from nF1 and nF2.
1486   aMI.Clear();
1487   GetFullShapeMap(nF1, aMI);
1488   GetFullShapeMap(nF2, aMI);
1489   //
1490   //collect VV, VE, EE, VF interferences
1491   for (aTypeInt = 0; aTypeInt < 4; ++aTypeInt) {
1492     for (i = 0; i < aNbLines[aTypeInt]; ++i) {
1493       BOPDS_Interf* aInt = (aTypeInt==0) ? (BOPDS_Interf*)(&aVVs(i)) : 
1494         ((aTypeInt==1) ? (BOPDS_Interf*)(&aVEs(i)) :
1495          ((aTypeInt==2) ? (BOPDS_Interf*)(&aEEs(i)) : (BOPDS_Interf*)(&aVFs(i))));
1496       if (aInt->HasIndexNew()) {
1497         aInt->Indices(nS1, nS2);
1498         if(aMI.Contains(nS1) && aMI.Contains(nS2)) {
1499           nVNew = aInt->IndexNew();
1500           aMVStick.Add(nVNew);
1501         }
1502       }
1503     }
1504   }
1505   //collect EF interferences
1506   for (i = 0; i < aNbLines[4]; ++i) {
1507     const BOPDS_InterfEF& aInt = aEFs(i);
1508     if (aInt.HasIndexNew()) {
1509       aInt.Indices(nS1, nS2);
1510       if(aMI.Contains(nS1) && aMI.Contains(nS2)) {
1511         nVNew = aInt.IndexNew();
1512         aMVStick.Add(nVNew);
1513         aMVEF.Add(nVNew);
1514       }
1515     }
1516   }
1517 }
1518
1519 //=======================================================================
1520 // function: GetFullShapeMap
1521 // purpose: 
1522 //=======================================================================
1523 void BOPAlgo_PaveFiller::GetFullShapeMap(const Standard_Integer nF,
1524                                          BOPCol_MapOfInteger& aMI)
1525 {
1526   BOPCol_ListIteratorOfListOfInteger aIt;
1527   Standard_Integer nS;
1528   //
1529   const BOPDS_ShapeInfo& aSI = myDS->ShapeInfo(nF);
1530   const BOPCol_ListOfInteger& aLI = aSI.SubShapes();
1531   //
1532   aMI.Add(nF);
1533   aIt.Initialize(aLI);
1534   for (; aIt.More(); aIt.Next()) {
1535     nS = aIt.Value();
1536     aMI.Add(nS);
1537   }
1538 }
1539
1540 //=======================================================================
1541 // function: RemoveUsedVertices
1542 // purpose: 
1543 //=======================================================================
1544 void BOPAlgo_PaveFiller::RemoveUsedVertices(BOPDS_Curve& aNC,
1545                                             BOPCol_MapOfInteger& aMV)
1546 {
1547   if (!aMV.Extent()) {
1548     return;
1549   }
1550   Standard_Integer nV;
1551   BOPDS_Pave aPave;
1552   BOPDS_ListIteratorOfListOfPave aItLP;
1553   //
1554   Handle(BOPDS_PaveBlock)& aPB=aNC.ChangePaveBlock1();
1555   const BOPDS_ListOfPave& aLP = aPB->ExtPaves();
1556   aItLP.Initialize(aLP);
1557   for (;aItLP.More();aItLP.Next()) {
1558     aPave = aItLP.Value();
1559     nV = aPave.Index();
1560     aMV.Remove(nV);
1561   }
1562 }
1563
1564 //=======================================================================
1565 //function : PutPaveOnCurve
1566 //purpose  : 
1567 //=======================================================================
1568   void BOPAlgo_PaveFiller::PutPaveOnCurve(const Standard_Integer nV,
1569                                           const Standard_Real aTolR3D,
1570                                           BOPDS_Curve& aNC,
1571                                           const BOPCol_MapOfInteger& aMI,
1572                                           BOPCol_DataMapOfIntegerReal& aMVTol,
1573                                           const Standard_Integer iCheckExtend)
1574 {
1575   Standard_Boolean bIsVertexOnLine;
1576   Standard_Real aT, aTol, aTolNew;
1577   BOPDS_Pave aPave;
1578   //
1579   const TopoDS_Vertex aV = (*(TopoDS_Vertex *)(&myDS->Shape(nV)));
1580   Handle(BOPDS_PaveBlock)& aPB=aNC.ChangePaveBlock1();
1581   const IntTools_Curve& aIC = aNC.Curve();
1582   //
1583   bIsVertexOnLine=myContext->IsVertexOnLine(aV, aIC, aTolR3D, aT);
1584   if (!bIsVertexOnLine && iCheckExtend) {
1585     aTol = BRep_Tool::Tolerance(aV);
1586     //
1587     ExtendedTolerance(nV, aMI, aTol, iCheckExtend);
1588     bIsVertexOnLine=myContext->IsVertexOnLine(aV, aTol, aIC, aTolR3D, aT);
1589   }
1590   //
1591   if (bIsVertexOnLine) {
1592     aPave.SetIndex(nV);
1593     aPave.SetParameter(aT);
1594     //
1595     aPB->AppendExtPave(aPave);
1596     //
1597     aTol = BRep_Tool::Tolerance(aV);
1598     //
1599     BOPTools_AlgoTools::UpdateVertex (aIC, aT, aV);
1600     //
1601     if (!aMVTol.IsBound(nV)) {
1602       aTolNew = BRep_Tool::Tolerance(aV);
1603       if (aTolNew > aTol) {
1604         aMVTol.Bind(nV, aTol);
1605       }
1606     }
1607     // 
1608     BOPDS_ShapeInfo& aSIDS=myDS->ChangeShapeInfo(nV);
1609     Bnd_Box& aBoxDS=aSIDS.ChangeBox();
1610     BRepBndLib::Add(aV, aBoxDS);
1611   }
1612 }
1613
1614 //=======================================================================
1615 //function : ProcessOldPaveBlocks
1616 //purpose  : 
1617 //=======================================================================
1618   void BOPAlgo_PaveFiller::ProcessExistingPaveBlocks
1619     (const Standard_Integer theInt,
1620      const BOPDS_MapOfPaveBlock& aMPBOnIn,
1621      BOPDS_IndexedDataMapOfShapeCoupleOfPaveBlocks& aMSCPB,
1622      BOPCol_DataMapOfShapeInteger& aMVI,
1623      const BOPCol_MapOfInteger& aMVB,
1624      BOPDS_MapOfPaveBlock& aMPB)
1625 {
1626   Standard_Integer nV, nE, iFlag;
1627   Standard_Real aT;
1628   BOPCol_MapIteratorOfMapOfInteger aItB;
1629   BOPDS_MapIteratorOfMapOfPaveBlock aItPB;
1630   //
1631   BOPDS_VectorOfInterfFF& aFFs=myDS->InterfFF();
1632   BOPDS_InterfFF& aFF = aFFs(theInt);
1633   BOPDS_VectorOfCurve& aVC=aFF.ChangeCurves();
1634   //  
1635   aItB.Initialize(aMVB);
1636   for (; aItB.More(); aItB.Next()) {
1637     nV = aItB.Value();
1638     const BOPDS_ShapeInfo& aSIV=myDS->ShapeInfo(nV);
1639     const Bnd_Box& aBoxV=aSIV.Box();
1640     const TopoDS_Vertex& aV = *(TopoDS_Vertex*)&aSIV.Shape();
1641     if (!aMVI.IsBound(aV)) {
1642       continue;
1643     }
1644     //
1645     aItPB.Initialize(aMPBOnIn);
1646     for (; aItPB.More(); aItPB.Next()) {
1647       const Handle(BOPDS_PaveBlock)& aPB = aItPB.Value();
1648       if (aPB->Pave1().Index() == nV || aPB->Pave2().Index() == nV) {
1649         continue;
1650       }
1651       //
1652       if (aMPB.Contains(aPB)) {
1653         continue;
1654       }
1655       nE=aPB->Edge();
1656       const BOPDS_ShapeInfo& aSIE=myDS->ShapeInfo(nE);
1657       const Bnd_Box& aBoxE=aSIE.Box();
1658       //
1659       if (!aBoxV.IsOut(aBoxE)) {
1660         const TopoDS_Edge& aE = *(TopoDS_Edge*)&aSIE.Shape();
1661         //
1662         iFlag=myContext->ComputeVE (aV, aE, aT);
1663         if (!iFlag) {
1664           aMPB.Add(aPB);
1665           //
1666           PreparePostTreatFF(theInt, aPB, aMSCPB, aMVI, aVC);
1667         }
1668       }
1669     }
1670   }
1671 }
1672
1673 //=======================================================================
1674 //function : UpdateExistingPaveBlocks
1675 //purpose  : 
1676 //=======================================================================
1677   void BOPAlgo_PaveFiller::UpdateExistingPaveBlocks
1678     (const Handle(BOPDS_PaveBlock)& aPBf,
1679      BOPDS_ListOfPaveBlock& aLPB,
1680      const Standard_Integer nF1,
1681      const Standard_Integer nF2) 
1682 {
1683   Standard_Integer nE;
1684   Standard_Boolean bCB;
1685   Handle(BOPDS_PaveBlock) aPB, aPB1, aPB2, aPB2n;
1686   Handle(BOPDS_CommonBlock) aCB;
1687   BOPDS_ListIteratorOfListOfPaveBlock aIt, aIt1, aIt2;
1688   BOPDS_IndexedMapOfPaveBlock aMPB;
1689   //
1690   //remove micro edges from aLPB
1691   aIt.Initialize(aLPB);
1692   for (; aIt.More();) {
1693     aPB = aIt.Value();
1694     const TopoDS_Edge& aE = *(TopoDS_Edge*)&myDS->Shape(aPB->Edge());
1695     if (BOPTools_AlgoTools::IsMicroEdge(aE, myContext)) {
1696       aLPB.Remove(aIt);
1697       continue;
1698     }
1699     aIt.Next();
1700   }
1701   //
1702   if (!aLPB.Extent()) {
1703     return;
1704   }
1705   //update face info
1706   myDS->UpdateFaceInfoOn(nF1);
1707   //
1708   myDS->UpdateFaceInfoOn(nF2);
1709   //
1710   BOPDS_FaceInfo& aFI1 = myDS->ChangeFaceInfo(nF1);
1711   BOPDS_FaceInfo& aFI2 = myDS->ChangeFaceInfo(nF2);
1712   //
1713   BOPDS_IndexedMapOfPaveBlock& aMPBOn1 = aFI1.ChangePaveBlocksOn();
1714   BOPDS_IndexedMapOfPaveBlock& aMPBIn1 = aFI1.ChangePaveBlocksIn();
1715   BOPDS_IndexedMapOfPaveBlock& aMPBOn2 = aFI2.ChangePaveBlocksOn();
1716   BOPDS_IndexedMapOfPaveBlock& aMPBIn2 = aFI2.ChangePaveBlocksIn();
1717   //
1718   // remove old pave blocks
1719   const Handle(BOPDS_CommonBlock)& aCB1 = myDS->CommonBlock(aPBf);
1720   bCB = !aCB1.IsNull();
1721   BOPDS_ListOfPaveBlock aLPB1;
1722   //
1723   if (bCB) {
1724     aLPB1.Assign(aCB1->PaveBlocks());
1725   } else {
1726     aLPB1.Append(aPBf);
1727   }
1728   aIt1.Initialize(aLPB1);
1729   for (; aIt1.More(); aIt1.Next()) {
1730     aPB1 = aIt1.Value();
1731     nE = aPB1->OriginalEdge();
1732     //
1733     BOPDS_ListOfPaveBlock& aLPB2 = myDS->ChangePaveBlocks(nE);
1734     aIt2.Initialize(aLPB2);
1735     for (; aIt2.More(); aIt2.Next()) {
1736       aPB2 = aIt2.Value();
1737       if (aPB1 == aPB2) {
1738         aLPB2.Remove(aIt2);
1739         break;
1740       }
1741     }
1742   }
1743   //
1744   if (bCB) {
1745     //create new pave blocks
1746     const BOPCol_ListOfInteger& aFaces = aCB1->Faces();
1747     aIt.Initialize(aLPB);
1748     for (; aIt.More(); aIt.Next()) {
1749       Handle(BOPDS_PaveBlock)& aPB = aIt.ChangeValue();
1750       //
1751       aCB = new BOPDS_CommonBlock;
1752       aIt1.Initialize(aLPB1);
1753       for (; aIt1.More(); aIt1.Next()) {
1754         aPB2 = aIt1.Value();
1755         nE = aPB2->OriginalEdge();
1756         //
1757         aPB2n = new BOPDS_PaveBlock;
1758         aPB2n->SetPave1(aPB->Pave1());
1759         aPB2n->SetPave2(aPB->Pave2());
1760         aPB2n->SetEdge(aPB->Edge());
1761         aPB2n->SetOriginalEdge(nE);
1762         aCB->AddPaveBlock(aPB2n);
1763         myDS->SetCommonBlock(aPB2n, aCB);
1764         myDS->ChangePaveBlocks(nE).Append(aPB2n);
1765       }
1766       aCB->AddFaces(aFaces);
1767       myDS->SortPaveBlocks(aCB);
1768       //
1769       aPB=aCB->PaveBlocks().First();
1770     }
1771   } 
1772   //
1773   aIt.Initialize(aLPB);
1774   for (; aIt.More(); aIt.Next()) {
1775     Handle(BOPDS_PaveBlock)& aPB = aIt.ChangeValue();
1776     nE = aPB->OriginalEdge();
1777     //
1778     Standard_Integer nF = (aMPBOn1.Contains(aPBf) || 
1779                            aMPBIn1.Contains(aPBf)) ? nF2 : nF1;
1780     const TopoDS_Face& aF = *(TopoDS_Face*)&myDS->Shape(nF);
1781     IntTools_Range aShrR(aPB->Pave1().Parameter(), aPB->Pave2().Parameter());
1782     const TopoDS_Edge& aE = *(TopoDS_Edge*)&myDS->Shape(aPB->Edge());
1783     //
1784     Standard_Boolean bCom = BOPTools_AlgoTools::IsBlockInOnFace(aShrR, aF, aE, myContext);
1785     if (bCom) {
1786       if (bCB) {
1787         aCB = myDS->CommonBlock(aPB);
1788         aCB->AddFace(nF);
1789       } else {
1790         aCB = new BOPDS_CommonBlock;
1791         aCB->AddPaveBlock(aPB);
1792         aCB->AddFace(nF1);
1793         aCB->AddFace(nF2);
1794         //
1795         myDS->SetCommonBlock(aPB, aCB);
1796       }
1797       aMPB.Add(aPB);
1798     }
1799     if (!bCB) {
1800       myDS->ChangePaveBlocks(nE).Append(aPB);
1801     }
1802   }
1803   //
1804   Standard_Integer i, aNbPB;
1805   Standard_Boolean bIn1, bIn2;
1806   //
1807   bIn1 = aMPBOn1.Contains(aPBf) || aMPBIn1.Contains(aPBf);
1808   bIn2 = aMPBOn2.Contains(aPBf) || aMPBIn2.Contains(aPBf);
1809   //
1810   aNbPB=aMPB.Extent();
1811   for (i=1; i<=aNbPB; ++i) {
1812     aPB = aMPB(i);
1813     if (!bIn1) {
1814       aMPBIn1.Add(aPB);
1815     }
1816     //
1817     if (!bIn2) {
1818       aMPBIn2.Add(aPB);
1819     }
1820   }
1821 }
1822
1823 //=======================================================================
1824 // function: PutClosingPaveOnCurve
1825 // purpose:
1826 //=======================================================================
1827   void BOPAlgo_PaveFiller::PutClosingPaveOnCurve(BOPDS_Curve& aNC)
1828 {
1829   Standard_Boolean bIsClosed, bHasBounds, bAdded;
1830   Standard_Integer nVC, j;
1831   Standard_Real aT[2], aTC, dT, aTx;
1832   gp_Pnt aP[2] ; 
1833   BOPDS_Pave aPVx;
1834   BOPDS_ListIteratorOfListOfPave aItLP;
1835   //
1836   const IntTools_Curve& aIC=aNC.Curve();
1837   const Handle(Geom_Curve)& aC3D=aIC.Curve();
1838   if(aC3D.IsNull()) {
1839     return;
1840   }
1841   //
1842   bIsClosed=IntTools_Tools::IsClosed(aC3D);
1843   if (!bIsClosed) {
1844     return;
1845   }
1846   //
1847   bHasBounds=aIC.HasBounds ();
1848   if (!bHasBounds){
1849     return;
1850   }
1851   // 
1852   bAdded=Standard_False;
1853   dT=Precision::PConfusion();
1854   aIC.Bounds (aT[0], aT[1], aP[0], aP[1]);
1855   //
1856   Handle(BOPDS_PaveBlock)& aPB=aNC.ChangePaveBlock1();
1857   BOPDS_ListOfPave& aLP=aPB->ChangeExtPaves();
1858   //
1859   aItLP.Initialize(aLP);
1860   for (; aItLP.More() && !bAdded; aItLP.Next()) {
1861     const BOPDS_Pave& aPC=aItLP.Value();
1862     nVC=aPC.Index();
1863     aTC=aPC.Parameter();
1864     //
1865     for (j=0; j<2; ++j) {
1866       if (fabs(aTC-aT[j]) < dT) {
1867         aTx=(!j) ? aT[1] : aT[0];
1868         aPVx.SetIndex(nVC);
1869         aPVx.SetParameter(aTx);
1870         aLP.Append(aPVx);
1871         //
1872         bAdded=Standard_True;
1873         break;
1874       }
1875     }
1876   }
1877 }
1878
1879 //=======================================================================
1880 //function : PreparePostTreatFF
1881 //purpose  : 
1882 //=======================================================================
1883   void BOPAlgo_PaveFiller::PreparePostTreatFF
1884     (const Standard_Integer aInt,
1885      const Handle(BOPDS_PaveBlock)& aPB,
1886      BOPDS_IndexedDataMapOfShapeCoupleOfPaveBlocks& aMSCPB,
1887      BOPCol_DataMapOfShapeInteger& aMVI,
1888      BOPDS_VectorOfCurve& aVC) 
1889 {
1890   Standard_Integer nV1, nV2;
1891   //
1892   Standard_Integer iC=aVC.Append()-1;
1893   BOPDS_ListOfPaveBlock& aLPBC = aVC(iC).ChangePaveBlocks();
1894   aLPBC.Append(aPB);
1895   //
1896   aPB->Indices(nV1, nV2);
1897   const TopoDS_Vertex& aV1=(*(TopoDS_Vertex *)(&myDS->Shape(nV1)));
1898   const TopoDS_Vertex& aV2=(*(TopoDS_Vertex *)(&myDS->Shape(nV2)));
1899   const TopoDS_Edge& aE = *(TopoDS_Edge*)&myDS->Shape(aPB->Edge());
1900   // Keep info for post treatment 
1901   BOPDS_CoupleOfPaveBlocks aCPB;
1902   aCPB.SetIndexInterf(aInt);
1903   aCPB.SetIndex(iC);
1904   aCPB.SetPaveBlock1(aPB);
1905   //
1906   aMSCPB.Add(aE, aCPB);
1907   aMVI.Bind(aV1, nV1);
1908   aMVI.Bind(aV2, nV2);
1909 }
1910
1911 //=======================================================================
1912 //function : CheckPlanes
1913 //purpose  : 
1914 //=======================================================================
1915 Standard_Boolean 
1916   BOPAlgo_PaveFiller::CheckPlanes(const Standard_Integer nF1,
1917                                   const Standard_Integer nF2)const
1918 {
1919   Standard_Boolean bToIntersect;
1920   Standard_Integer i, nV2, iCnt;
1921   BOPCol_MapIteratorOfMapOfInteger aIt;
1922   //
1923   bToIntersect=Standard_False;
1924   //
1925   const BOPDS_FaceInfo& aFI1=myDS->ChangeFaceInfo(nF1);
1926   const BOPDS_FaceInfo& aFI2=myDS->ChangeFaceInfo(nF2);
1927   //
1928   const BOPCol_MapOfInteger& aMVIn1=aFI1.VerticesIn();
1929   const BOPCol_MapOfInteger& aMVOn1=aFI1.VerticesOn();
1930   //
1931   iCnt=0;
1932   for (i=0; (i<2 && !bToIntersect); ++i) {
1933     const BOPCol_MapOfInteger& aMV2=(!i) ? aFI2.VerticesIn() 
1934       : aFI2.VerticesOn();
1935     //
1936     aIt.Initialize(aMV2);
1937     for (; aIt.More(); aIt.Next()) {
1938       nV2=aIt.Value();
1939       if (aMVIn1.Contains(nV2) || aMVOn1.Contains(nV2)) {
1940         ++iCnt;
1941         if (iCnt>1) {
1942           bToIntersect=!bToIntersect;
1943           break;
1944         }
1945       }
1946     }
1947   }
1948   //
1949   return bToIntersect;
1950 }
1951 //=======================================================================
1952 //function : UpdatePaveBlocks
1953 //purpose  : 
1954 //=======================================================================
1955 void BOPAlgo_PaveFiller::UpdatePaveBlocks(const BOPCol_DataMapOfIntegerInteger& aDMI)
1956 {
1957   if (aDMI.IsEmpty()) {
1958     return;
1959   }
1960   //
1961   Standard_Integer nSp, aNbPBP, nV[2], i, j;
1962   Standard_Real aT[2];
1963   Standard_Boolean bCB, bRebuild;
1964   BOPDS_ListIteratorOfListOfPaveBlock aItPB;
1965   BOPDS_MapOfPaveBlock aMPB;
1966   //
1967   BOPDS_VectorOfListOfPaveBlock& aPBP=myDS->ChangePaveBlocksPool();
1968   aNbPBP = aPBP.Extent();
1969   for (i=0; i<aNbPBP; ++i) {
1970     BOPDS_ListOfPaveBlock& aLPB=aPBP(i);
1971     //
1972     aItPB.Initialize(aLPB);
1973     for (; aItPB.More(); aItPB.Next()) {
1974       Handle(BOPDS_PaveBlock) aPB=aItPB.Value();
1975       const Handle(BOPDS_CommonBlock)& aCB=myDS->CommonBlock(aPB);
1976       bCB = !aCB.IsNull();
1977       if (bCB) {
1978         aPB=aCB->PaveBlock1();
1979       }
1980       //
1981       if (aMPB.Add(aPB)) {
1982         bRebuild = Standard_False;
1983         aPB->Indices(nV[0], nV[1]);
1984         aPB->Range(aT[0], aT[1]);
1985         //
1986         for (j = 0; j < 2; ++j) {
1987           if (aDMI.IsBound(nV[j])) {
1988             BOPDS_Pave aPave;
1989             //
1990             nV[j] = aDMI.Find(nV[j]);
1991             aPave.SetIndex(nV[j]);
1992             aPave.SetParameter(aT[j]);
1993             //
1994             bRebuild = Standard_True;
1995             if (!j) {
1996               aPB->SetPave1(aPave);
1997             } else {
1998               aPB->SetPave2(aPave);
1999             }
2000           }
2001         }
2002         //
2003         if (bRebuild) {
2004           nSp = SplitEdge(aPB->Edge(), nV[0], aT[0], nV[1], aT[1]);
2005           if (bCB) {
2006             aCB->SetEdge(nSp);
2007           }
2008           else {
2009             aPB->SetEdge(nSp);
2010           }
2011         }// if (bRebuild) {
2012       }// if (aMPB.Add(aPB)) {
2013     }// for (; aItPB.More(); aItPB.Next()) {
2014   }// for (i=0; i<aNbPBP; ++i) {
2015   aMPB.Clear();
2016 }
2017 //=======================================================================
2018 //function : ToleranceFF
2019 //purpose  : Computes the TolFF according to the tolerance value and 
2020 //           types of the faces.
2021 //=======================================================================
2022   void ToleranceFF(const BRepAdaptor_Surface& aBAS1,
2023                    const BRepAdaptor_Surface& aBAS2,
2024                    Standard_Real& aTolFF)
2025 {
2026   Standard_Real aTol1, aTol2;
2027   Standard_Boolean isAna1, isAna2;
2028   //
2029   aTol1 = aBAS1.Tolerance();
2030   aTol2 = aBAS2.Tolerance();
2031   aTolFF = Max(aTol1, aTol2);
2032   //
2033   isAna1 = (aBAS1.GetType() == GeomAbs_Plane ||
2034             aBAS1.GetType() == GeomAbs_Cylinder ||
2035             aBAS1.GetType() == GeomAbs_Cone ||
2036             aBAS1.GetType() == GeomAbs_Sphere ||
2037             aBAS1.GetType() == GeomAbs_Torus);
2038   //
2039   isAna2 = (aBAS2.GetType() == GeomAbs_Plane ||
2040             aBAS2.GetType() == GeomAbs_Cylinder ||
2041             aBAS2.GetType() == GeomAbs_Cone ||
2042             aBAS2.GetType() == GeomAbs_Sphere ||
2043             aBAS2.GetType() == GeomAbs_Torus);
2044   //
2045   if (!isAna1 || !isAna2) {
2046     aTolFF =  Max(aTolFF, 5.e-6);
2047   }
2048 }