Salome HOME
correct integration mistake
[modules/geom.git] / src / GEOMImpl / GEOMImpl_PipeDriver.cxx
1 //  Copyright (C) 2007-2008  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.
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 #include <Standard_Stream.hxx>
23
24 #include <GEOMImpl_PipeDriver.hxx>
25
26 #include <GEOMImpl_IShapesOperations.hxx>
27 #include <GEOMImpl_IPipeDiffSect.hxx>
28 #include <GEOMImpl_IPipeShellSect.hxx>
29 #include <GEOMImpl_IPipeBiNormal.hxx>
30 #include <GEOMImpl_IPipe.hxx>
31 #include <GEOMImpl_GlueDriver.hxx>
32 #include <GEOMImpl_Types.hxx>
33 #include <GEOM_Function.hxx>
34
35 //#include <GEOMAlgo_GlueAnalyser.hxx>
36
37 #include <ShapeAnalysis_FreeBounds.hxx>
38 #include <ShapeAnalysis_Edge.hxx>
39 #include <ShapeFix_Face.hxx>
40 #include <ShapeFix_Shell.hxx>
41 #include <ShapeFix_Shape.hxx>
42 #include <ShapeFix_ShapeTolerance.hxx>
43
44 #include <BRep_Tool.hxx>
45 #include <BRep_Builder.hxx>
46 #include <BRepBuilderAPI_MakeWire.hxx>
47 #include <BRepBuilderAPI_Sewing.hxx>
48 #include <BRepCheck_Analyzer.hxx>
49 #include <BRepOffsetAPI_MakePipe.hxx>
50 #include <BRepOffsetAPI_MakePipeShell.hxx>
51 #include <GProp_GProps.hxx>
52 #include <BRepGProp.hxx>
53 #include <BRepBuilderAPI_MakeFace.hxx>
54 #include <BRepBuilderAPI_Copy.hxx>
55
56 #include <TopAbs.hxx>
57 #include <TopExp.hxx>
58 #include <TopExp_Explorer.hxx>
59 #include <TopoDS.hxx>
60 #include <TopoDS_Wire.hxx>
61 #include <TopoDS_Edge.hxx>
62 #include <TopoDS_Shape.hxx>
63 #include <TopoDS_Solid.hxx>
64 #include <TopoDS_Shell.hxx>
65 #include <TopoDS_Face.hxx>
66 #include <TopoDS_Compound.hxx>
67 #include <TopTools_SequenceOfShape.hxx>
68 #include <TopTools_IndexedDataMapOfShapeShape.hxx>
69 #include <TopTools_IndexedDataMapOfShapeListOfShape.hxx>
70 #include <TopTools_ListIteratorOfListOfShape.hxx>
71
72 #include <GeomAPI_ProjectPointOnCurve.hxx>
73 #include <GeomAPI_Interpolate.hxx>
74 #include <Geom_TrimmedCurve.hxx>
75 #include <Geom_Plane.hxx>
76 #include <Geom_RectangularTrimmedSurface.hxx>
77 #include <Geom_BezierSurface.hxx>
78 #include <Geom_Line.hxx>
79 #include <Geom_Conic.hxx>
80 #include <Geom_BSplineCurve.hxx>
81 #include <Geom_BSplineSurface.hxx>
82 #include <GeomFill_BSplineCurves.hxx>
83 #include <GeomConvert_ApproxCurve.hxx>
84 #include <GeomConvert.hxx>
85
86 #include <TColgp_SequenceOfPnt.hxx>
87 #include <TColgp_HArray1OfPnt.hxx>
88 #include <TColgp_Array2OfPnt.hxx>
89 #include <TColStd_HSequenceOfTransient.hxx>
90
91 #include <Precision.hxx>
92 #include <Standard_NullObject.hxx>
93 #include <Standard_TypeMismatch.hxx>
94 #include <Standard_ConstructionError.hxx>
95
96 #include "utilities.h"
97
98 //=======================================================================
99 //function : GetID
100 //purpose  :
101 //=======================================================================
102 const Standard_GUID& GEOMImpl_PipeDriver::GetID()
103 {
104   static Standard_GUID aPipeDriver("FF1BBB19-5D14-4df2-980B-3A668264EA16");
105   return aPipeDriver;
106 }
107
108 //=======================================================================
109 //function : GEOMImpl_PipeDriver
110 //purpose  :
111 //=======================================================================
112 GEOMImpl_PipeDriver::GEOMImpl_PipeDriver()
113 {
114 }
115
116 //=======================================================================
117 //function : FillForOtherEdges
118 //purpose  : auxilary for CreatePipeForShellSections()
119 //=======================================================================
120 static bool FillForOtherEdges(const TopoDS_Shape& F1,
121                               const TopoDS_Shape& E1,
122                               const TopoDS_Shape& V1,
123                               TopTools_IndexedDataMapOfShapeShape& FF)
124 {
125   //cout<<"FillForOtherEdges"<<endl;
126   // find other pairs for vertexes and edges
127   // creating map of vertex edges for both faces
128   TopTools_IndexedDataMapOfShapeListOfShape aMapVertEdge1;
129   TopExp::MapShapesAndAncestors(F1, TopAbs_VERTEX, TopAbs_EDGE, aMapVertEdge1);
130   if (!FF.Contains(F1))
131     MESSAGE("    FillForOtherEdges: map FF not contains key F1");
132   if (!FF.Contains(E1))
133     MESSAGE("    FillForOtherEdges: map FF not contains key E1");
134   if (!FF.Contains(V1))
135     MESSAGE("    FillForOtherEdges: map FF not contains key V1");
136   const TopoDS_Shape& F2 = FF.FindFromKey(F1);
137   const TopoDS_Shape& E2 = FF.FindFromKey(E1);
138   const TopoDS_Shape& V2 = FF.FindFromKey(V1);
139   TopTools_IndexedDataMapOfShapeListOfShape aMapVertEdge2;
140   TopExp::MapShapesAndAncestors(F2, TopAbs_VERTEX, TopAbs_EDGE, aMapVertEdge2);
141
142   TopoDS_Edge ES1 = TopoDS::Edge(E1);
143   TopoDS_Edge ES2 = TopoDS::Edge(E2);
144   TopoDS_Shape VS1 = V1;
145   TopoDS_Shape VS2 = V2;
146
147   ShapeAnalysis_Edge sae;
148   while(1) {
149     if(!aMapVertEdge1.Contains(VS1))
150       MESSAGE ("    FillForOtherEdges: map aMapVertEdge1 not contains key VS1");
151     const TopTools_ListOfShape& aList1 = aMapVertEdge1.FindFromKey(VS1);
152     //TopoDS_Shape E1next;
153     TopTools_ListIteratorOfListOfShape anIter1(aList1);
154     if(anIter1.Value().IsSame(ES1)) {
155       anIter1.Next();
156     }
157     //E1next = anIter1.Value();
158     if(!aMapVertEdge2.Contains(VS2))
159       MESSAGE ("    FillForOtherEdges: map aMapVertEdge2 not contains key VS2");
160     const TopTools_ListOfShape& aList2 = aMapVertEdge2.FindFromKey(VS2);
161     //TopoDS_Shape E2next;
162     TopTools_ListIteratorOfListOfShape anIter2(aList2);
163     if(anIter2.Value().IsSame(ES2)) {
164       anIter2.Next();
165     }
166     //E2next = anIter2.Value();
167     //ES1 = TopoDS::Edge(E1next);
168     //ES2 = TopoDS::Edge(E2next);
169     ES1 = TopoDS::Edge(anIter1.Value());
170     ES2 = TopoDS::Edge(anIter2.Value());
171     if(!FF.Contains(ES1)) {
172       FF.Add(ES1,ES2);
173     }
174     if(VS1.IsSame(sae.FirstVertex(ES1)))
175       VS1 = sae.LastVertex(ES1);
176     else
177       VS1 = sae.FirstVertex(ES1);
178     if(VS2.IsSame(sae.FirstVertex(ES2)))
179       VS2 = sae.LastVertex(ES2);
180     else
181       VS2 = sae.FirstVertex(ES2);
182     if(VS1.IsSame(V1))
183       break;
184     if(!FF.Contains(VS1)) {
185       FF.Add(VS1,VS2);
186     }
187   }
188
189   return true;
190 }
191
192 //=======================================================================
193 //function : FillCorrespondingEdges
194 //purpose  : auxilary for CreatePipeForShellSections()
195 //=======================================================================
196 static bool FillCorrespondingEdges(const TopoDS_Shape& FS1,
197                                    const TopoDS_Shape& FS2,
198                                    const TopoDS_Vertex& aLoc1,
199                                    const TopoDS_Vertex& aLoc2,
200                                    const TopoDS_Wire& aWirePath,
201                                    TopTools_IndexedDataMapOfShapeShape& FF)
202 {
203   //cout<<"FillCorrespondingEdges"<<endl;
204   // find corresponding edges
205   TopExp_Explorer expw1(FS1,TopAbs_WIRE);
206   TopoDS_Wire aWire1 = TopoDS::Wire(expw1.Current());
207   //exp = TopExp_Explorer(FS2,TopAbs_WIRE);
208   TopExp_Explorer expw2(FS2,TopAbs_WIRE);
209   TopoDS_Wire aWire2 = TopoDS::Wire(expw2.Current());
210   BRepOffsetAPI_MakePipeShell aBuilder(aWirePath);
211   aBuilder.Add(aWire1, aLoc1);
212   aBuilder.Add(aWire2, aLoc2);
213   if(!aBuilder.IsReady()) {
214     return false;
215   }
216   aBuilder.Build();
217   TopoDS_Shape aShape = aBuilder.Shape();
218   /*
219   TopoDS_Compound C;
220   BRep_Builder B;
221   B.MakeCompound(C);
222   B.Add(C,aShape);
223   B.Add(C,FS1);
224   B.Add(C,FS2);
225   BRepTools::Write(C,"/dn02/users_Linux/skl/work/Bugs/14857/comp.brep");
226   */
227   ShapeAnalysis_Edge sae;
228   double tol = Max( BRep_Tool::Tolerance(TopoDS::Face(FS1)),
229                     BRep_Tool::Tolerance(TopoDS::Face(FS2)) );
230   TopTools_MapOfShape Vs1,Vs2;
231   TopExp_Explorer exp;
232   exp.Init( FS1, TopAbs_EDGE );
233   TopoDS_Edge E1 = TopoDS::Edge(exp.Current());
234   TopoDS_Vertex V11 = sae.FirstVertex(E1);
235   TopoDS_Vertex V21 = sae.LastVertex(E1);
236   gp_Pnt P11 = BRep_Tool::Pnt(V11);
237   gp_Pnt P21 = BRep_Tool::Pnt(V21);
238   //cout<<"P11("<<P11.X()<<","<<P11.Y()<<","<<P11.Z()<<")"<<endl;
239   //cout<<"P21("<<P21.X()<<","<<P21.Y()<<","<<P21.Z()<<")"<<endl;
240   // find corresponding vertexes from created shape
241   TopoDS_Vertex VN11,VN21;
242   for( exp.Init( aShape, TopAbs_VERTEX ); exp.More(); exp.Next() ) {
243     TopoDS_Vertex V = TopoDS::Vertex(exp.Current());
244     gp_Pnt P = BRep_Tool::Pnt(V);
245     if(P.Distance(P11)<tol) {
246       VN11 = V;
247     }
248     if(P.Distance(P21)<tol) {
249       VN21 = V;
250     }
251   }
252   // find edge contains VN11 and VN21 and corresponding vertexes
253   TopoDS_Vertex VN12,VN22;
254   for( exp.Init( aShape, TopAbs_FACE ); exp.More(); exp.Next() ) {
255     TopoDS_Shape F = exp.Current();
256     TopExp_Explorer expe;
257     bool IsFind = false;
258     for( expe.Init( F, TopAbs_EDGE ); expe.More(); expe.Next() ) {
259       TopoDS_Edge E = TopoDS::Edge(expe.Current());
260       TopoDS_Vertex VF = sae.FirstVertex(E);
261       TopoDS_Vertex VL = sae.LastVertex(E);
262       if( (VF.IsSame(VN11) && VL.IsSame(VN21)) || (VF.IsSame(VN21) && VL.IsSame(VN11)) ) {
263         IsFind = true;
264         break;
265       }
266     }
267     if(IsFind) {
268       for( expe.Init( F, TopAbs_EDGE ); expe.More(); expe.Next() ) {
269         TopoDS_Edge E = TopoDS::Edge(expe.Current());
270         TopoDS_Vertex VF = sae.FirstVertex(E);
271         TopoDS_Vertex VL = sae.LastVertex(E);
272         if( VF.IsSame(VN11) && !VL.IsSame(VN21) )
273           VN12 = VL;
274         if( VL.IsSame(VN11) && !VF.IsSame(VN21) )
275           VN12 = VF;
276         if( VF.IsSame(VN21) && !VL.IsSame(VN11) )
277           VN22 = VL;
278         if( VL.IsSame(VN21) && !VF.IsSame(VN11) )
279           VN22 = VF;
280       }
281       break;
282     }
283   }
284   // find vertexes from FS2 corresponded to VN12 and VN22
285   // and find edge from FS2 contains V12 and V22,
286   // this edge will be corresponded to edge E1
287   TopoDS_Vertex V12,V22;
288   gp_Pnt PN12 = BRep_Tool::Pnt(VN12);
289   gp_Pnt PN22 = BRep_Tool::Pnt(VN22);
290   //cout<<"PN12("<<PN12.X()<<","<<PN12.Y()<<","<<PN12.Z()<<")"<<endl;
291   //cout<<"PN22("<<PN22.X()<<","<<PN22.Y()<<","<<PN22.Z()<<")"<<endl;
292   TopoDS_Edge E2;
293   TopExp_Explorer expe;
294   for( expe.Init( FS2, TopAbs_EDGE ); expe.More(); expe.Next() ) {
295     TopoDS_Edge E = TopoDS::Edge(expe.Current());
296     TopoDS_Vertex VF = sae.FirstVertex(E);
297     TopoDS_Vertex VL = sae.LastVertex(E);
298     gp_Pnt PF = BRep_Tool::Pnt(VF);
299     gp_Pnt PL = BRep_Tool::Pnt(VL);
300     if( PF.Distance(PN12)<tol && PL.Distance(PN22)<tol ) {
301       V12 = VF;
302       V22 = VL;
303       E2 = E;
304       break;
305     }
306     if( PF.Distance(PN22)<tol && PL.Distance(PN12)<tol ) {
307       V12 = VL;
308       V22 = VF;
309       E2 = E;
310       break;
311     }
312   }
313   FF.Add(V11,V12);
314   FF.Add(V21,V22);
315   FF.Add(E1,E2);
316
317   // find other pairs for vertexes and edges
318   // creating map of vertex edges for both faces
319   return FillForOtherEdges(FS1,E1,V21,FF);
320
321   //return true;
322 }
323
324 //=======================================================================
325 //function : FillCorrespondingEdges
326 //purpose  : auxilary for CreatePipeShellsWithoutPath()
327 //=======================================================================
328 static bool FillCorrespondingEdges(const TopoDS_Shape& FS1,
329                                    const TopoDS_Shape& FS2,
330                                    const TopoDS_Vertex& aLoc1,
331                                    const TopoDS_Vertex& aLoc2,
332                                    TopTools_IndexedDataMapOfShapeShape& FF)
333 {
334   //cout<<"FillCorrespondingEdges"<<endl;
335
336   gp_Pnt P1 = BRep_Tool::Pnt(aLoc1);
337   gp_Pnt P2 = BRep_Tool::Pnt(aLoc2);
338   gp_Vec aDir(P1,P2);
339
340   ShapeAnalysis_Edge sae;
341   double tol = Max( BRep_Tool::Tolerance(TopoDS::Face(FS1)),
342                     BRep_Tool::Tolerance(TopoDS::Face(FS2)) );
343   TopTools_MapOfShape Vs1,Vs2;
344
345   TopoDS_Vertex V11=aLoc1, V12=aLoc2, V21, V22;
346   TopoDS_Edge E1,E2;
347
348   TopExp_Explorer exp1;
349   for( exp1.Init(FS1,TopAbs_EDGE); exp1.More(); exp1.Next() ) {
350     E1 = TopoDS::Edge(exp1.Current());
351     TopoDS_Vertex V1 = sae.FirstVertex(E1);
352     TopoDS_Vertex V2 = sae.LastVertex(E1);
353     gp_Pnt Ptmp1 = BRep_Tool::Pnt(V1);
354     gp_Pnt Ptmp2 = BRep_Tool::Pnt(V2);
355     //cout<<"P11("<<P11.X()<<","<<P11.Y()<<","<<P11.Z()<<")"<<endl;
356     //cout<<"P21("<<P21.X()<<","<<P21.Y()<<","<<P21.Z()<<")"<<endl;
357     if(P1.Distance(Ptmp1)<tol) {
358       V21 = V2;
359       break;
360     }
361     if(P1.Distance(Ptmp2)<tol) {
362       V21 = V1;
363       break;
364     }
365   }
366
367   TopoDS_Edge E21,E22;
368   TopoDS_Vertex VE21,VE22;
369   int nbe=0;
370   for( exp1.Init(FS2,TopAbs_EDGE); exp1.More() && nbe<2; exp1.Next() ) {
371     TopoDS_Edge E = TopoDS::Edge(exp1.Current());
372     TopoDS_Vertex V1 = sae.FirstVertex(E);
373     TopoDS_Vertex V2 = sae.LastVertex(E);
374     gp_Pnt Ptmp1 = BRep_Tool::Pnt(V1);
375     gp_Pnt Ptmp2 = BRep_Tool::Pnt(V2);
376     if(P2.Distance(Ptmp1)<tol) {
377       if(nbe==0) {
378         E21 = E;
379         VE21 = V2;
380         nbe++;
381       }
382       else if(nbe==1) {
383         E22 = E;
384         VE22 = V2;
385         nbe++;
386       }
387     }
388     if(P2.Distance(Ptmp2)<tol) {
389       if(nbe==0) {
390         E21 = E;
391         VE21 = V1;
392         nbe++;
393       }
394       else if(nbe==1) {
395         E22 = E;
396         VE22 = V1;
397         nbe++;
398       }
399     }
400   }
401
402   gp_Pnt PV21 = BRep_Tool::Pnt(V21);
403   gp_Pnt PE21 = BRep_Tool::Pnt(VE21);
404   gp_Pnt PE22 = BRep_Tool::Pnt(VE22);
405   gp_Vec aDir1(PV21,PE21);
406   gp_Vec aDir2(PV21,PE22);
407   double ang1 = aDir.Angle(aDir1);
408   double ang2 = aDir.Angle(aDir2);
409   if(fabs(ang1)<fabs(ang2)) {
410     E2 = E21;
411     V22 = VE21;
412   }
413   else {
414     E2 = E22;
415     V22 = VE22;
416   }
417
418   FF.Add(V11,V12);
419   FF.Add(V21,V22);
420   FF.Add(E1,E2);
421
422   // find other pairs for vertexes and edges
423   return FillForOtherEdges(FS1,E1,V21,FF);
424 }
425
426 //=======================================================================
427 //function : FindNextPairOfFaces
428 //purpose  : auxilary for CreatePipeForShellSections()
429 //=======================================================================
430 static void FindNextPairOfFaces(const TopoDS_Shape& aCurFace,
431                                 TopTools_IndexedDataMapOfShapeListOfShape& aMapEdgeFaces1,
432                                 TopTools_IndexedDataMapOfShapeListOfShape& aMapEdgeFaces2,
433                                 TopTools_IndexedDataMapOfShapeShape& FF,
434                                 GEOMImpl_IPipe* aCI)
435 {
436   //cout<<"FindNextPairOfFaces"<<endl;
437   TopExp_Explorer anExp;
438   for ( anExp.Init( aCurFace, TopAbs_EDGE ); anExp.More(); anExp.Next() ) {
439     TopoDS_Shape E1 = anExp.Current();
440     if(!FF.Contains(E1)) {
441       if(aCI) delete aCI;
442       Standard_ConstructionError::Raise("FindNextPairOfFaces: Can not find edge in map");
443     }
444     if(!FF.Contains(E1))
445       MESSAGE ("    FindNextPairOfFaces: map FF not contains key E1");
446     const TopoDS_Shape& E2 = FF.FindFromKey(E1);
447     TopExp_Explorer anExpV;
448     anExpV.Init( E1, TopAbs_VERTEX );
449     TopoDS_Shape V1 = anExpV.Current();
450     if(!FF.Contains(V1)) {
451       if(aCI) delete aCI;
452       Standard_ConstructionError::Raise("FindNextPairOfFaces: Can not find vertex in map");
453     }
454
455     if(!aMapEdgeFaces1.Contains(E1))
456       MESSAGE ("    FindNextPairOfFaces: map aMapEdgeFaces1 not contains key E1");
457     const TopTools_ListOfShape& aList1 = aMapEdgeFaces1.FindFromKey(E1);
458     if(aList1.Extent()<2)
459       continue;
460     TopTools_ListIteratorOfListOfShape anIter(aList1);
461     if(anIter.Value().IsEqual(aCurFace)) {
462       anIter.Next();
463     }
464     TopoDS_Shape F1other = anIter.Value();
465     if(FF.Contains(F1other))
466       continue;
467
468     if(!FF.Contains(aCurFace))
469       MESSAGE ("    FindNextPairOfFaces: map FF not contains key aCurFace");
470     const TopoDS_Shape& F2 = FF.FindFromKey(aCurFace);
471     if(!aMapEdgeFaces2.Contains(E2))
472       MESSAGE ("    FindNextPairOfFaces: map aMapEdgeFaces2 not contains key E2");
473     const TopTools_ListOfShape& aList2 = aMapEdgeFaces2.FindFromKey(E2);
474     if(aList2.Extent()<2) {
475       if(aCI) delete aCI;
476       Standard_ConstructionError::Raise("FindNextPairOfFaces: Can not find corresponding face");
477     }
478     TopTools_ListIteratorOfListOfShape anIter2(aList2);
479     if(anIter2.Value().IsEqual(F2)) {
480       anIter2.Next();
481     }
482     TopoDS_Shape F2other = anIter2.Value();
483     FF.Add(F1other,F2other);
484
485     // add pairs of edges to FF
486     bool stat =  FillForOtherEdges(F1other,E1,V1,FF);
487     if( !stat ) {
488       if(aCI) delete aCI;
489       Standard_ConstructionError::Raise("FindNextPairOfFaces: Can not mapping other egdes");
490     }
491
492     FindNextPairOfFaces(F1other, aMapEdgeFaces1, aMapEdgeFaces2, FF, aCI);
493   }
494 }
495
496 //=======================================================================
497 //function : FindFirstPairFaces
498 //purpose  : auxilary for Execute()
499 //=======================================================================
500 static void FindFirstPairFaces(const TopoDS_Shape& S1, const TopoDS_Shape& S2,
501                                TopoDS_Vertex& V1, TopoDS_Vertex& V2,
502                                TopoDS_Shape& FS1, TopoDS_Shape& FS2)
503 {
504   //cout<<"FindFirstPairFaces"<<endl;
505
506   // check if vertexes are subshapes of sections
507   gp_Pnt P1 = BRep_Tool::Pnt(V1);
508   gp_Pnt P2 = BRep_Tool::Pnt(V2);
509   TopoDS_Vertex V1new,V2new;
510   TopExp_Explorer exp;
511   double mindist = 1.e10;
512   for( exp.Init( S1, TopAbs_VERTEX ); exp.More(); exp.Next() ) {
513     TopoDS_Vertex V = TopoDS::Vertex(exp.Current());
514     gp_Pnt P = BRep_Tool::Pnt(V);
515     double dist = P1.Distance(P);
516     if(dist<mindist) {
517       mindist = dist;
518       V1new = V;
519     }
520   }
521   mindist = 1.e10;
522   for( exp.Init( S2, TopAbs_VERTEX ); exp.More(); exp.Next() ) {
523     TopoDS_Vertex V = TopoDS::Vertex(exp.Current());
524     gp_Pnt P = BRep_Tool::Pnt(V);
525     double dist = P2.Distance(P);
526     if(dist<mindist) {
527       mindist = dist;
528       V2new = V;
529     }
530   }
531
532   //gp_Pnt P1new = BRep_Tool::Pnt(V1new);
533   //gp_Pnt P2new = BRep_Tool::Pnt(V2new);
534   //cout<<"  P1("<<P1.X()<<","<<P1.Y()<<","<<P1.Z()<<")"<<endl;
535   //cout<<"  P2("<<P2.X()<<","<<P2.Y()<<","<<P2.Z()<<")"<<endl;
536   //cout<<"  P1new("<<P1new.X()<<","<<P1new.Y()<<","<<P1new.Z()<<")"<<endl;
537   //cout<<"  P2new("<<P2new.X()<<","<<P2new.Y()<<","<<P2new.Z()<<")"<<endl;
538
539   // replace vertexes if it is needed
540   if(!V1.IsSame(V1new)) {
541     V1 = V1new;
542     P1 = BRep_Tool::Pnt(V1);
543     MESSAGE ("  replace V1");
544   }
545   else
546     MESSAGE ("  not replace V1");
547   if(!V2.IsSame(V2new)) {
548     V2 = V2new;
549     P2 = BRep_Tool::Pnt(V2);
550     MESSAGE ("  replace V2");
551   }
552   else
553     MESSAGE ("  not replace V2");
554
555   TopTools_IndexedDataMapOfShapeListOfShape aMapVertFaces1;
556   TopExp::MapShapesAndAncestors(S1, TopAbs_VERTEX, TopAbs_FACE, aMapVertFaces1);
557   TopTools_IndexedDataMapOfShapeListOfShape aMapVertFaces2;
558   TopExp::MapShapesAndAncestors(S2, TopAbs_VERTEX, TopAbs_FACE, aMapVertFaces2);
559
560   if(!aMapVertFaces1.Contains(V1))
561     MESSAGE ("    FindFirstPairFaces: map aMapVertFaces1 not contains key V1");
562   const TopTools_ListOfShape& aList1 = aMapVertFaces1.FindFromKey(V1);
563   TopTools_ListIteratorOfListOfShape anIter1(aList1);
564   FS1 = anIter1.Value();
565   // find middle point
566   double x1=0., y1=0., z1=0.;
567   int nbv1=0;
568   for( exp.Init( FS1, TopAbs_VERTEX ); exp.More(); exp.Next() ) {
569     TopoDS_Vertex V = TopoDS::Vertex(exp.Current());
570     gp_Pnt P = BRep_Tool::Pnt(V);
571     x1 += P.X();
572     y1 += P.Y();
573     z1 += P.Z();
574     nbv1++;
575   }
576   gp_Pnt PM1(x1/nbv1, y1/nbv1, z1/nbv1);
577
578   TColgp_SequenceOfPnt Ps;
579   TopTools_SequenceOfShape Fs;
580   if(!aMapVertFaces2.Contains(V2))
581     MESSAGE ("    FindFirstPairFaces: map aMapVertFaces2 not contains key V2");
582   const TopTools_ListOfShape& aList2 = aMapVertFaces2.FindFromKey(V2);
583   TopTools_ListIteratorOfListOfShape anIter2(aList2);
584   for(; anIter2.More(); anIter2.Next()) {
585     TopoDS_Shape F = anIter2.Value();
586     double x2=0., y2=0., z2=0.;
587     int nbv2=0;
588     for( exp.Init( F, TopAbs_VERTEX ); exp.More(); exp.Next() ) {
589       TopoDS_Vertex V = TopoDS::Vertex(exp.Current());
590       gp_Pnt P = BRep_Tool::Pnt(V);
591       x2 += P.X();
592       y2 += P.Y();
593       z2 += P.Z();
594       nbv2++;
595     }
596     gp_Pnt PM(x2/nbv1, y2/nbv1, z2/nbv1);
597     Fs.Append(F);
598     Ps.Append(PM);
599   }
600
601   gp_Vec aDir(P1,P2);
602   int i=1;
603   double MinAng = PI;
604   int numface = 0;
605   for(; i<=Fs.Length(); i++) {
606     gp_Vec tmpDir(PM1,Ps(i));
607     double ang = fabs(aDir.Angle(tmpDir));
608     if(ang<MinAng) {
609       MinAng = ang;
610       numface = i;
611     }
612   }
613   FS2 = Fs(numface);
614 }
615
616 //=======================================================================
617 //function : CreatePipeForShellSections
618 //purpose  : auxilary for Execute()
619 //=======================================================================
620 static TopoDS_Shape CreatePipeForShellSections(const TopoDS_Wire& aWirePath,
621                                                GEOMImpl_IPipe* aCI)
622 {
623   //cout<<"CreatePipeForShellSections"<<endl;
624   //TopoDS_Shape res;
625   int i,j;
626   BRep_Builder B;
627
628   GEOMImpl_IPipeShellSect* aCIDS = (GEOMImpl_IPipeShellSect*)aCI;
629   Handle(TColStd_HSequenceOfTransient) aBasesObjs = aCIDS->GetBases();
630   Handle(TColStd_HSequenceOfTransient) aSubBasesObjs = aCIDS->GetSubBases();
631   Handle(TColStd_HSequenceOfTransient) aLocObjs = aCIDS->GetLocations();
632   Standard_Boolean aWithContact = (aCIDS->GetWithContactMode());
633   Standard_Boolean aWithCorrect = (aCIDS->GetWithCorrectionMode());
634
635   Standard_Integer nbBases = aBasesObjs->Length(),
636     nbSubBases = (aSubBasesObjs.IsNull() ? 0 :aSubBasesObjs->Length()),
637     nbLocs = (aLocObjs.IsNull() ? 0 :aLocObjs->Length());
638
639   if( nbLocs != nbBases) {
640     if(aCI) delete aCI;
641     Standard_ConstructionError::Raise("Number of sections is not equal to number of locations ");
642   }
643   if( nbSubBases && nbSubBases != nbBases) {
644     if(aCI) delete aCI;
645     Standard_ConstructionError::Raise("Number of sections is not equal to number of subsections ");
646   }
647
648   //BRepOffsetAPI_MakePipeShell aBuilder(aWirePath);
649
650   TopTools_SequenceOfShape VLocs;
651   for(i=1; i<=nbBases; i++) {
652     Handle(Standard_Transient) anItemLoc = aLocObjs->Value(i);
653     if(anItemLoc.IsNull())
654       continue;
655     Handle(GEOM_Function) aRefLoc = Handle(GEOM_Function)::DownCast(anItemLoc);
656     TopoDS_Shape aShapeLoc = aRefLoc->GetValue();
657     if(aShapeLoc.IsNull() || aShapeLoc.ShapeType() != TopAbs_VERTEX)
658       continue;
659     VLocs.Append(aShapeLoc);
660   }
661   nbLocs = VLocs.Length();
662   if( nbLocs != nbBases) {
663     if(aCI) delete aCI;
664     Standard_ConstructionError::Raise("One of location shapes is not a vertex");
665   }
666   // split wire path by location points
667   TColgp_SequenceOfPnt PLocs;
668   for(i=1; i<=nbLocs; i++) {
669     TopoDS_Vertex V = TopoDS::Vertex(VLocs.Value(i));
670     PLocs.Append(BRep_Tool::Pnt(V));
671   }
672
673   TopTools_SequenceOfShape Edges;
674   TopTools_SequenceOfShape Wires;
675   ShapeAnalysis_Edge sae;
676
677   if(nbLocs==2) {
678     TopExp_Explorer anExp;
679     for ( anExp.Init( aWirePath, TopAbs_EDGE ); anExp.More(); anExp.Next() ) {
680       Edges.Append(anExp.Current());
681     }
682     Standard_Integer Num1 = 0;
683     Standard_Integer Num2 = 0;
684     for(i=1; i<=Edges.Length(); i++) {
685       TopoDS_Edge E = TopoDS::Edge(Edges.Value(i));
686       double tol = BRep_Tool::Tolerance(E);
687       TopoDS_Vertex V1 = sae.FirstVertex(E);
688       TopoDS_Vertex V2 = sae.LastVertex(E);
689       gp_Pnt P1 = BRep_Tool::Pnt(V1);
690       gp_Pnt P2 = BRep_Tool::Pnt(V2);
691       if( P1.Distance(PLocs.First()) < tol ) {
692         Num1 = i;
693       }
694       if( P2.Distance(PLocs.Last()) < tol ) {
695         Num2 = i;
696       }
697     }
698     if( Num1>0 && Num2>0 ) {
699       TopoDS_Wire W;
700       B.MakeWire(W);
701       for(i=Num1; i<=Num2; i++) {
702         B.Add(W,Edges.Value(i));
703       }
704       Wires.Append(W);
705     }
706     else {
707       Wires.Append(aWirePath);
708     }
709   }
710   else {
711     TopExp_Explorer anExp;
712     for ( anExp.Init( aWirePath, TopAbs_EDGE ); anExp.More(); anExp.Next() ) {
713       Edges.Append(anExp.Current());
714     }
715     TopoDS_Edge edge = TopoDS::Edge(Edges.First());
716     double tol = BRep_Tool::Tolerance(edge);
717     TopoDS_Vertex VF = sae.FirstVertex(edge);
718     gp_Pnt PF = BRep_Tool::Pnt(VF);
719     //cout<<"PF("<<PF.X()<<","<<PF.Y()<<","<<PF.Z()<<")"<<endl;
720     if( PF.Distance(PLocs.First()) > tol ) {
721       if(aCI) delete aCI;
722       Standard_ConstructionError::Raise
723         ("First location shapes is not coincided with first vertex of aWirePath");
724     }
725     VLocs.ChangeValue(1) = VF;
726     edge = TopoDS::Edge(Edges.Last());
727     tol = BRep_Tool::Tolerance(edge);
728     TopoDS_Vertex VL = sae.LastVertex(edge);
729     gp_Pnt PL = BRep_Tool::Pnt(VL);
730     if( PL.Distance(PLocs.Last()) > tol ) {
731       if(aCI) delete aCI;
732       Standard_ConstructionError::Raise
733         ("Last location shapes is not coincided with last vertex of aWirePath");
734     }
735     VLocs.ChangeValue(nbLocs) = VL;
736     int jcurr = 2;
737     TopTools_SequenceOfShape tmpEdges;
738     for(i=1; i<=Edges.Length() && jcurr<nbLocs; i++) {
739       TopoDS_Edge E = TopoDS::Edge(Edges.Value(i));
740       tol = BRep_Tool::Tolerance(E);
741       TopoDS_Vertex V1 = sae.FirstVertex(E);
742       TopoDS_Vertex V2 = sae.LastVertex(E);
743       gp_Pnt P1 = BRep_Tool::Pnt(V1);
744       gp_Pnt P2 = BRep_Tool::Pnt(V2);
745       if( P2.Distance(PLocs.Value(jcurr)) < tol ) {
746         // make wire from current edge and add created
747         // wire to Wires
748         TopoDS_Wire W;
749         B.MakeWire(W);
750         for(j=1; j<=tmpEdges.Length(); j++)
751           B.Add(W,tmpEdges.Value(j));
752         B.Add(W,E);
753         Wires.Append(W);
754         VLocs.ChangeValue(jcurr) = V2;
755         jcurr++;
756         tmpEdges.Clear();
757       }
758       else {
759         // find distance between E and aLocs(jcurr)
760         double fp,lp;
761         Handle(Geom_Curve) C = BRep_Tool::Curve(E,fp,lp);
762         GeomAPI_ProjectPointOnCurve PPC (PLocs.Value(jcurr),C);
763         if( PPC.NbPoints()>0 &&
764             PLocs.Value(jcurr).Distance(PPC.Point(1)) < tol ) {
765           double param = PPC.Parameter(1);
766           gp_Pnt PC1;
767           C->D0(param,PC1);
768           // split current edge
769           Handle(Geom_TrimmedCurve) tc1 = new Geom_TrimmedCurve(C,fp,param);
770           Handle(Geom_TrimmedCurve) tc2 = new Geom_TrimmedCurve(C,param,lp);
771           TopoDS_Edge E1,E2;
772           gp_Pnt Pfp;
773           C->D0(fp,Pfp);
774           if(Pfp.Distance(P1)<tol) {
775             B.MakeEdge(E1,tc1,tol);
776             B.Add(E1,V1);
777             TopoDS_Shape tmpV = VLocs.Value(jcurr).Oriented(TopAbs_REVERSED);
778             B.Add(E1,TopoDS::Vertex(tmpV));
779             tmpEdges.Append(E1);
780             B.MakeEdge(E2,tc2,tol);
781             tmpV = VLocs.Value(jcurr).Oriented(TopAbs_FORWARD);
782             B.Add(E2,TopoDS::Vertex(tmpV));
783             B.Add(E2,V2);
784           }
785           else {
786             B.MakeEdge(E1,tc2,tol);
787             TopoDS_Shape tmpV = VLocs.Value(jcurr).Oriented(TopAbs_FORWARD);
788             B.Add(E1,TopoDS::Vertex(tmpV));
789             B.Add(E1,V1);
790             E1.Reverse();
791             tmpEdges.Append(E1);
792             B.MakeEdge(E2,tc1,tol);
793             B.Add(E2,V2);
794             tmpV = VLocs.Value(jcurr).Oriented(TopAbs_REVERSED);
795             B.Add(E2,TopoDS::Vertex(tmpV));
796             E2.Reverse();
797           }
798           // create wire from tmpEdges
799           TopoDS_Wire W;
800           B.MakeWire(W);
801           for(j=1; j<=tmpEdges.Length(); j++)
802             B.Add(W,tmpEdges.Value(j));
803           Wires.Append(W);
804           jcurr++;
805           tmpEdges.Clear();
806           Edges.Remove(i);
807           Edges.InsertAfter(i-1,E1);
808           Edges.InsertAfter(i,E2);
809         }
810         else {
811           tmpEdges.Append(E);
812         }
813       }
814     }
815     // create wire from other edges
816     TopoDS_Wire W;
817     B.MakeWire(W);
818     for(; i<=Edges.Length(); i++)
819       B.Add(W,Edges.Value(i));
820     Wires.Append(W);
821     //cout<<"Wires.Length()="<<Wires.Length()<<endl;
822   }
823
824   if( Wires.Length() != nbLocs-1 ) {
825     if(aCI) delete aCI;
826     Standard_ConstructionError::Raise
827       ("One of location shapes is not lied on the path");
828   }
829
830   //TopTools_SequenceOfShape aSeqBases;
831   //TopTools_SequenceOfShape aSeqSubBases;
832   //TopTools_SequenceOfShape aSeqFaces;
833   TopoDS_Compound aComp;
834   B.MakeCompound(aComp);
835   for (i = 1; i < nbBases; i++) {
836     TopoDS_Wire WPath = TopoDS::Wire(Wires.Value(i));
837     // 1 section
838     Handle(Standard_Transient) anItem1 = aBasesObjs->Value(i);
839     if(anItem1.IsNull())
840       continue;
841     Handle(GEOM_Function) aRefBase1 = Handle(GEOM_Function)::DownCast(anItem1);
842     if(aRefBase1.IsNull())
843       continue;
844     TopoDS_Shape aShBase1 = aRefBase1->GetValue();
845     if(aShBase1.IsNull())
846       continue;
847     TopAbs_ShapeEnum aType1 = aShBase1.ShapeType();
848     // 2 section
849     Handle(Standard_Transient) anItem2 = aBasesObjs->Value(i+1);
850     if(anItem2.IsNull())
851       continue;
852     Handle(GEOM_Function) aRefBase2 = Handle(GEOM_Function)::DownCast(anItem2);
853     if(aRefBase2.IsNull())
854       continue;
855     TopoDS_Shape aShBase2 = aRefBase2->GetValue();
856     if(aShBase2.IsNull())
857       continue;
858     TopAbs_ShapeEnum aType2 = aShBase2.ShapeType();
859
860     //BRepTools::Write(aShBase1,"/dn02/users_Linux/skl/work/Bugs/14857/base1.brep");
861
862     bool OkSec = ( aType1==TopAbs_SHELL || aType1==TopAbs_FACE ) &&
863                  ( aType2==TopAbs_SHELL || aType2==TopAbs_FACE );
864     if( !OkSec ) {
865       if(aCI) delete aCI;
866       Standard_ConstructionError::Raise("One of section shapes has invalid type");
867     }
868
869     bool CreateFewSolids = false;
870     // compare sections
871     TopExp_Explorer anExp;
872     Standard_Integer nbf1 = 0;
873     for ( anExp.Init( aShBase1, TopAbs_FACE ); anExp.More(); anExp.Next() ) {
874       nbf1++;
875     }
876     Standard_Integer nbf2 = 0;
877     for ( anExp.Init( aShBase2, TopAbs_FACE ); anExp.More(); anExp.Next() ) {
878       nbf2++;
879     }
880     if(nbf1==nbf2) {
881       CreateFewSolids = true;
882     }
883
884     /*
885     // check orientation of sections
886     bool NeedReverse = false;
887     {
888       // first section
889       anExp.Init( aShBase1, TopAbs_FACE );
890       TopoDS_Shape aFace = anExp.Current();
891       TColgp_SequenceOfPnt aPnts;
892       double xc=0, yc=0, zc=0;
893       for ( anExp.Init( aFace, TopAbs_VERTEX ); anExp.More(); anExp.Next() ) {
894         TopoDS_Vertex V = TopoDS::Vertex(anExp.Current());
895         aPnts.Append(BRep_Tool::Pnt(V));
896         xc += aPnts.Last().X();
897         yc += aPnts.Last().Y();
898         zc += aPnts.Last().Z();
899       }
900       gp_Pnt PC( xc/aPnts.Length(), yc/aPnts.Length(), zc/aPnts.Length() );
901       gp_Vec V1(PC,aPnts.Value(1));
902       gp_Vec V2(PC,aPnts.Value(2));
903       gp_Vec VN = V1.Crossed(V2);
904       for(int ip=2; ip<aPnts.Length(); ip++) {
905         V1 = gp_Vec(PC,aPnts.Value(ip));
906         V2 = gp_Vec(PC,aPnts.Value(ip+1));
907         VN.Add(V1.Crossed(V2));
908       }
909       gp_Vec PathNorm;
910       gp_Pnt PLoc = BRep_Tool::Pnt(TopoDS::Vertex(VLocs(i)));
911       TopExp_Explorer WE;
912       for ( WE.Init( WPath, TopAbs_EDGE ); WE.More(); WE.Next() ) {
913         TopoDS_Edge edge = TopoDS::Edge(WE.Current());
914         double tol = BRep_Tool::Tolerance(edge);
915         TopoDS_Vertex VF = sae.FirstVertex(edge);
916         gp_Pnt PF = BRep_Tool::Pnt(VF);
917         if( PF.Distance(PLoc) < tol ) {
918           double fp,lp;
919           Handle(Geom_Curve) C = BRep_Tool::Curve(edge,fp,lp);
920           gp_Pnt P1,P2;
921           C->D0(fp,P1);
922           if( P1.Distance(PLoc) < tol ) {
923             C->D0(fp+(lp-fp)/100,P2);
924           }
925           else {
926             C->D0(lp,P1);
927             C->D0(lp+(fp-lp)/100,P2);
928           }
929           PathNorm = gp_Vec(P1,P2);
930           break;
931         }
932         else {
933           TopoDS_Vertex VL = sae.LastVertex(edge);
934           gp_Pnt PL = BRep_Tool::Pnt(VL);
935           if( PL.Distance(PLoc) < tol ) {
936             double fp,lp;
937             Handle(Geom_Curve) C = BRep_Tool::Curve(edge,fp,lp);
938             gp_Pnt P1,P2;
939             C->D0(fp,P1);
940             if( P1.Distance(PLoc) < tol ) {
941               C->D0(fp+(lp-fp)/100,P2);
942             }
943             else {
944               C->D0(lp,P1);
945               C->D0(lp+(fp-lp)/100,P2);
946             }
947             PathNorm = gp_Vec(P2,P1);
948             break;
949           }
950         }
951       }
952       cout<<"VN("<<VN.X()<<","<<VN.Y()<<","<<VN.Z()<<")"<<endl;
953       cout<<"PathNorm("<<PathNorm.X()<<","<<PathNorm.Y()<<","<<PathNorm.Z()<<")"<<endl;
954       if(fabs(VN.Angle(PathNorm))>PI/2.) {
955         NeedReverse = true;
956         aShBase1.Reverse();
957       }
958     }
959     {
960       // second section
961       anExp.Init( aShBase2, TopAbs_FACE );
962       TopoDS_Shape aFace = anExp.Current();
963       TColgp_SequenceOfPnt aPnts;
964       double xc=0, yc=0, zc=0;
965       for ( anExp.Init( aFace, TopAbs_VERTEX ); anExp.More(); anExp.Next() ) {
966         TopoDS_Vertex V = TopoDS::Vertex(anExp.Current());
967         aPnts.Append(BRep_Tool::Pnt(V));
968         xc += aPnts.Last().X();
969         yc += aPnts.Last().Y();
970         zc += aPnts.Last().Z();
971       }
972       gp_Pnt PC( xc/aPnts.Length(), yc/aPnts.Length(), zc/aPnts.Length() );
973       gp_Vec V1(PC,aPnts.Value(1));
974       gp_Vec V2(PC,aPnts.Value(2));
975       gp_Vec VN = V1.Crossed(V2);
976       for(int ip=2; ip<aPnts.Length(); ip++) {
977         V1 = gp_Vec(PC,aPnts.Value(ip));
978         V2 = gp_Vec(PC,aPnts.Value(ip+1));
979         VN.Add(V1.Crossed(V2));
980       }
981       gp_Vec PathNorm;
982       gp_Pnt PLoc = BRep_Tool::Pnt(TopoDS::Vertex(VLocs(i+1)));
983       TopExp_Explorer WE;
984       for ( WE.Init( WPath, TopAbs_EDGE ); WE.More(); WE.Next() ) {
985         TopoDS_Edge edge = TopoDS::Edge(WE.Current());
986         double tol = BRep_Tool::Tolerance(edge);
987         TopoDS_Vertex VF = sae.FirstVertex(edge);
988         gp_Pnt PF = BRep_Tool::Pnt(VF);
989         if( PF.Distance(PLoc) < tol ) {
990           double fp,lp;
991           Handle(Geom_Curve) C = BRep_Tool::Curve(edge,fp,lp);
992           gp_Pnt P1,P2;
993           C->D0(fp,P1);
994           if( P1.Distance(PLoc) < tol ) {
995             C->D0(fp+(lp-fp)/100,P2);
996           }
997           else {
998             C->D0(lp,P1);
999             C->D0(lp+(fp-lp)/100,P2);
1000           }
1001           PathNorm = gp_Vec(P2,P1);
1002           break;
1003         }
1004         else {
1005           TopoDS_Vertex VL = sae.LastVertex(edge);
1006           gp_Pnt PL = BRep_Tool::Pnt(VL);
1007           if( PL.Distance(PLoc) < tol ) {
1008             double fp,lp;
1009             Handle(Geom_Curve) C = BRep_Tool::Curve(edge,fp,lp);
1010             gp_Pnt P1,P2;
1011             C->D0(fp,P1);
1012             if( P1.Distance(PLoc) < tol ) {
1013               C->D0(fp+(lp-fp)/100,P2);
1014             }
1015             else {
1016               C->D0(lp,P1);
1017               C->D0(lp+(fp-lp)/100,P2);
1018             }
1019             PathNorm = gp_Vec(P2,P1);
1020             break;
1021           }
1022         }
1023       }
1024       //cout<<"VN("<<VN.X()<<","<<VN.Y()<<","<<VN.Z()<<")"<<endl;
1025       //cout<<"PathNorm("<<PathNorm.X()<<","<<PathNorm.Y()<<","<<PathNorm.Z()<<")"<<endl;
1026       if(fabs(VN.Angle(PathNorm))>PI/2.)
1027         aShBase2.Reverse();
1028     }
1029     */
1030
1031     if(!CreateFewSolids) {
1032       // we can create only one solid
1033       TopoDS_Shape aWire1, aWire2;
1034       // prepare aWire1
1035       if(aType1==TopAbs_SHELL) {
1036         // create wire as boundary contour if shell is no closed
1037         // get free boundary shapes
1038         ShapeAnalysis_FreeBounds anAnalizer( aShBase1 );
1039         TopoDS_Compound aClosed = anAnalizer.GetClosedWires();
1040         //TopExp_Explorer anExp;
1041         Standard_Integer NbWires = 0;
1042         for ( anExp.Init( aClosed, TopAbs_WIRE ); anExp.More(); anExp.Next() ) {
1043           NbWires++;
1044           aWire1 = anExp.Current();
1045         }
1046         if(NbWires!=1) {
1047           // bad case
1048           if(aCI) delete aCI;
1049           Standard_ConstructionError::Raise("Bad shell is used as section ");
1050         }
1051       }
1052       else { // aType1==TopAbs_FACE
1053         TopExp_Explorer aExpW(aShBase1,TopAbs_WIRE);
1054         aWire1 = aExpW.Current();
1055       }
1056       // prepare aWire2
1057       if(aType2==TopAbs_SHELL) {
1058         // create wire as boundary contour if shell is no closed
1059         // get free boundary shapes
1060         ShapeAnalysis_FreeBounds anAnalizer( aShBase2 );
1061         TopoDS_Compound aClosed = anAnalizer.GetClosedWires();
1062         //TopExp_Explorer anExp;
1063         Standard_Integer NbWires = 0;
1064         for ( anExp.Init( aClosed, TopAbs_WIRE ); anExp.More(); anExp.Next() ) {
1065           NbWires++;
1066           aWire2 = anExp.Current();
1067         }
1068         if(NbWires!=1) {
1069           // bad case
1070           if(aCI) delete aCI;
1071           Standard_ConstructionError::Raise("Bad shell is used as section ");
1072         }
1073       }
1074       else { // aType2==TopAbs_FACE
1075         TopExp_Explorer aExpW(aShBase2,TopAbs_WIRE);
1076         aWire2 = aExpW.Current();
1077       }
1078       // make pipe using aWire1 and aWire2
1079       if( !aWire1.IsNull() && !aWire2.IsNull() ) {
1080         //BRepOffsetAPI_MakePipeShell aBuilder(aWirePath);
1081         BRepOffsetAPI_MakePipeShell aBuilder(WPath);
1082         aBuilder.Add(aWire1, TopoDS::Vertex(VLocs(i)),
1083                      aWithContact, aWithCorrect);
1084         aBuilder.Add(aWire2, TopoDS::Vertex(VLocs(i+1)),
1085                      aWithContact, aWithCorrect);
1086         if(!aBuilder.IsReady()) {
1087           if(aCI) delete aCI;
1088           Standard_ConstructionError::Raise("Invalid input data for building PIPE: bases are invalid");
1089         }
1090         aBuilder.Build();
1091         TopoDS_Shape aShape = aBuilder.Shape();
1092         TopoDS_Shell aShell;
1093         B.MakeShell(aShell);
1094         for ( anExp.Init( aShape, TopAbs_FACE ); anExp.More(); anExp.Next() ) {
1095           B.Add(aShell,anExp.Current());
1096         }
1097         for ( anExp.Init( aShBase1, TopAbs_FACE ); anExp.More(); anExp.Next() ) {
1098           B.Add(aShell,anExp.Current());
1099         }
1100         for ( anExp.Init( aShBase2, TopAbs_FACE ); anExp.More(); anExp.Next() ) {
1101           B.Add(aShell,anExp.Current());
1102         }
1103         // make sewing for this shell
1104         Handle(BRepBuilderAPI_Sewing) aSewing = new BRepBuilderAPI_Sewing;
1105         aSewing->SetTolerance(Precision::Confusion());
1106         aSewing->SetFaceMode(Standard_True);
1107         aSewing->SetFloatingEdgesMode(Standard_False);
1108         aSewing->SetNonManifoldMode(Standard_False);
1109         for ( anExp.Init( aShell, TopAbs_FACE ); anExp.More(); anExp.Next() ) {
1110           aSewing->Add(anExp.Current());
1111         }
1112         aSewing->Perform();
1113         const TopoDS_Shape aSewShape = aSewing->SewedShape();
1114         if( aSewShape.ShapeType() == TopAbs_SHELL ) {
1115           aShell = TopoDS::Shell(aSewShape);
1116           GProp_GProps aSystem;
1117           BRepGProp::VolumeProperties(aShell, aSystem);
1118           if(aSystem.Mass()<0) {
1119             aShell.Reverse();
1120           }
1121           if(BRep_Tool::IsClosed(aShell)) {
1122             TopoDS_Solid aSolid;
1123             B.MakeSolid(aSolid);
1124             B.Add(aSolid,aShell);
1125             B.Add(aComp,aSolid);
1126           }
1127           else {
1128             B.Add(aComp,aShell);
1129           }
1130         }
1131         else {
1132           B.Add(aComp,aShell);
1133         }
1134       }
1135     }
1136     else {
1137       // main block - creation few solids (for each pair of faces)
1138       TopTools_MapOfShape aFaces1,aFaces2;
1139       for ( anExp.Init( aShBase1, TopAbs_FACE ); anExp.More(); anExp.Next() ) {
1140         aFaces1.Add(anExp.Current());
1141       }
1142       for ( anExp.Init( aShBase2, TopAbs_FACE ); anExp.More(); anExp.Next() ) {
1143         aFaces2.Add(anExp.Current());
1144       }
1145       // creating map of edge faces
1146       TopTools_IndexedDataMapOfShapeListOfShape aMapEdgeFaces1;
1147       TopExp::MapShapesAndAncestors(aShBase1, TopAbs_EDGE, TopAbs_FACE, aMapEdgeFaces1);
1148       TopTools_IndexedDataMapOfShapeListOfShape aMapEdgeFaces2;
1149       TopExp::MapShapesAndAncestors(aShBase2, TopAbs_EDGE, TopAbs_FACE, aMapEdgeFaces2);
1150
1151       // constuct map face->face
1152       TopTools_IndexedDataMapOfShapeShape FF;
1153       TopoDS_Shape FS1,FS2;
1154       if(nbSubBases==0) {
1155         // find edge the most distant from location point
1156         // (this edge is not shared by two faces)
1157         double maxdist = 0.;
1158         TopoDS_Shape E1;
1159         TopoDS_Vertex V11,V21;
1160         for(j=1; j<=aMapEdgeFaces1.Extent(); j++) {
1161           TopoDS_Shape tmp = aMapEdgeFaces1.FindKey(j);
1162           const TopTools_ListOfShape& aList = aMapEdgeFaces1.FindFromKey(tmp);
1163           if(aList.Extent()>1)
1164             continue;
1165           TopExp_Explorer expv;
1166           expv.Init( tmp, TopAbs_VERTEX );
1167           TopoDS_Vertex V1 = TopoDS::Vertex(expv.Current());
1168           expv.Next();
1169           TopoDS_Vertex V2 = TopoDS::Vertex(expv.Current());
1170           gp_Pnt P1 = BRep_Tool::Pnt(V1);
1171           gp_Pnt P2 = BRep_Tool::Pnt(V2);
1172           double dist = PLocs.Value(i).Distance(P1) + PLocs.Value(i).Distance(P2);
1173           if(dist>maxdist) {
1174             E1 = tmp;
1175             V11 = V1;
1176             V21 = V2;
1177             TopTools_ListIteratorOfListOfShape anIter(aList);
1178             FS1 = anIter.Value();
1179             maxdist = dist;
1180           }
1181         }
1182         // main direction for comparing
1183         gp_Vec VM(PLocs.Value(i),PLocs.Value(i+1));
1184         // find corresponding edge from next section
1185         double minang = PI;
1186         gp_Pnt P11 = BRep_Tool::Pnt(V11);
1187         gp_Pnt P21 = BRep_Tool::Pnt(V21);
1188         TopoDS_Shape E2;
1189         TopoDS_Vertex V12,V22;
1190         for(j=1; j<=aMapEdgeFaces2.Extent(); j++) {
1191           TopoDS_Shape tmp = aMapEdgeFaces2.FindKey(j);
1192           const TopTools_ListOfShape& aList = aMapEdgeFaces2.FindFromKey(tmp);
1193           if(aList.Extent()>1)
1194             continue;
1195           TopExp_Explorer expv;
1196           expv.Init( tmp, TopAbs_VERTEX );
1197           TopoDS_Vertex V1tmp = TopoDS::Vertex(expv.Current());
1198           expv.Next();
1199           TopoDS_Vertex V2tmp = TopoDS::Vertex(expv.Current());
1200           gp_Pnt P1tmp = BRep_Tool::Pnt(V1tmp);
1201           gp_Pnt P2tmp = BRep_Tool::Pnt(V2tmp);
1202           double d1 = P1tmp.Distance(P11) + P2tmp.Distance(P21);
1203           double d2 = P1tmp.Distance(P21) + P2tmp.Distance(P11);
1204           TopoDS_Vertex V1,V2;
1205           gp_Pnt P1,P2;
1206           if(d1>d2) {
1207             V1 = V2tmp; P1 = P2tmp;
1208             V2 = V1tmp; P2 = P1tmp;
1209           }
1210           else {
1211             V1 = V1tmp; P1 = P1tmp;
1212             V2 = V2tmp; P2 = P2tmp;
1213           }
1214           gp_Vec Vec1(P11,P1);
1215           gp_Vec Vec2(P21,P2);
1216           double ang = fabs(Vec1.Angle(VM)) + fabs(Vec2.Angle(VM));
1217           if(ang<minang) {
1218             E2 = tmp;
1219             V12 = V1;
1220             V22 = V2;
1221             TopTools_ListIteratorOfListOfShape anIter(aList);
1222             FS2 = anIter.Value();
1223             minang = ang;
1224           }
1225         }
1226         // put all pairs to map FF
1227         FF.Add(FS1,FS2);
1228         FF.Add(E1,E2);
1229         FF.Add(V11,V12);
1230         FF.Add(V21,V22);
1231
1232         // add pairs of edges to FF
1233         bool stat =  FillForOtherEdges(FS1,E1,V11,FF);
1234         if( !stat ) {
1235           if(aCI) delete aCI;
1236           Standard_ConstructionError::Raise("FindForOtherEdges: Can not mapping other egdes");
1237         }
1238
1239       }
1240       else {
1241         { // 1 section
1242           Handle(Standard_Transient) anItem = aSubBasesObjs->Value(i);
1243           if(anItem.IsNull()) {
1244             if(aCI) delete aCI;
1245             Standard_ConstructionError::Raise("Invalid subbase shape");
1246           }
1247           Handle(GEOM_Function) aRefBase = Handle(GEOM_Function)::DownCast(anItem);
1248           if(aRefBase.IsNull()) {
1249             if(aCI) delete aCI;
1250             Standard_ConstructionError::Raise("Invalid subbase shape");
1251           }
1252           TopoDS_Shape aSh = aRefBase->GetValue();
1253           if(aSh.IsNull()) {
1254             if(aCI) delete aCI;
1255             Standard_ConstructionError::Raise("Invalid subbase shape");
1256           }
1257           if(aSh.ShapeType()!=TopAbs_FACE) {
1258             if(aCI) delete aCI;
1259             Standard_ConstructionError::Raise("Invalid subbase shape");
1260           }
1261           FS1 = aSh;
1262         }
1263         { // 2 section
1264           Handle(Standard_Transient) anItem = aSubBasesObjs->Value(i+1);
1265           if(anItem.IsNull()) {
1266             if(aCI) delete aCI;
1267             Standard_ConstructionError::Raise("Invalid subbase shape");
1268           }
1269           Handle(GEOM_Function) aRefBase = Handle(GEOM_Function)::DownCast(anItem);
1270           if(aRefBase.IsNull()) {
1271             if(aCI) delete aCI;
1272             Standard_ConstructionError::Raise("Invalid subbase shape");
1273           }
1274           TopoDS_Shape aSh = aRefBase->GetValue();
1275           if(aSh.IsNull()) {
1276             if(aCI) delete aCI;
1277             Standard_ConstructionError::Raise("Invalid subbase shape");
1278           }
1279           if(aSh.ShapeType()!=TopAbs_FACE) {
1280             if(aCI) delete aCI;
1281             Standard_ConstructionError::Raise("Invalid subbase shape");
1282           }
1283           FS2 = aSh;
1284         }
1285
1286         if( !aFaces1.Contains(FS1) || !aFaces2.Contains(FS2) ) {
1287           if(aCI) delete aCI;
1288           Standard_ConstructionError::Raise("Invalid subbase shape");
1289         }
1290
1291         FF.Add(FS1,FS2);
1292
1293         // add pairs of edges to FF
1294         bool stat =  FillCorrespondingEdges(FS1, FS2, TopoDS::Vertex(VLocs(i)),
1295                                             TopoDS::Vertex(VLocs(i+1)), WPath, FF);
1296         if( !stat ) {
1297           if(aCI) delete aCI;
1298           Standard_ConstructionError::Raise("Can not create correct pipe");
1299         }
1300       }
1301
1302       FindNextPairOfFaces(FS1, aMapEdgeFaces1, aMapEdgeFaces2, FF, aCI);
1303
1304       // make pipe for each pair of faces
1305       for(j=1; j<=FF.Extent(); j++) {
1306         TopoDS_Shape F1 = FF.FindKey(j);
1307         if( F1.ShapeType() != TopAbs_FACE )
1308           continue;
1309         TopoDS_Shape F2 = FF.FindFromIndex(j);
1310         TopExp_Explorer aExpW1(F1,TopAbs_WIRE);
1311         TopoDS_Wire aWire1 = TopoDS::Wire(aExpW1.Current());
1312         TopExp_Explorer aExpW2(F2,TopAbs_WIRE);
1313         TopoDS_Wire aWire2 = TopoDS::Wire(aExpW2.Current());
1314         // make pipe using aWire1 and aWire2
1315         if( !aWire1.IsNull() && !aWire2.IsNull() ) {
1316           BRepOffsetAPI_MakePipeShell aBuilder(WPath);
1317           aBuilder.Add(aWire1, TopoDS::Vertex(VLocs(i)),
1318                        aWithContact, aWithCorrect);
1319           aBuilder.Add(aWire2, TopoDS::Vertex(VLocs(i+1)),
1320                        aWithContact, aWithCorrect);
1321           if(!aBuilder.IsReady()) {
1322             if(aCI) delete aCI;
1323             Standard_ConstructionError::Raise("Invalid input data for building PIPE: bases are invalid");
1324           }
1325           aBuilder.Build();
1326           TopoDS_Shape aShape = aBuilder.Shape();
1327           TopoDS_Shell aShell;
1328           B.MakeShell(aShell);
1329           for ( anExp.Init( aShape, TopAbs_FACE ); anExp.More(); anExp.Next() ) {
1330             B.Add(aShell,anExp.Current());
1331           }
1332
1333           B.Add(aShell,F1);
1334           B.Add(aShell,F2);
1335           // make sewing for this shell
1336           Handle(BRepBuilderAPI_Sewing) aSewing = new BRepBuilderAPI_Sewing;
1337           aSewing->SetTolerance(Precision::Confusion());
1338           aSewing->SetFaceMode(Standard_True);
1339           aSewing->SetFloatingEdgesMode(Standard_False);
1340           aSewing->SetNonManifoldMode(Standard_False);
1341           for ( anExp.Init( aShell, TopAbs_FACE ); anExp.More(); anExp.Next() ) {
1342             aSewing->Add(anExp.Current());
1343           }
1344           aSewing->Perform();
1345           const TopoDS_Shape aSewShape = aSewing->SewedShape();
1346           if( aSewShape.ShapeType() == TopAbs_SHELL ) {
1347             aShell = TopoDS::Shell(aSewShape);
1348             GProp_GProps aSystem;
1349             BRepGProp::VolumeProperties(aShell, aSystem);
1350             if(aSystem.Mass()<0) {
1351               //cout<<"aSewShape is reversed"<<endl;
1352               aShell.Reverse();
1353             }
1354             if(BRep_Tool::IsClosed(aShell)) {
1355               TopoDS_Solid aSolid;
1356               B.MakeSolid(aSolid);
1357               B.Add(aSolid,aShell);
1358               B.Add(aComp,aSolid);
1359             }
1360             else {
1361               B.Add(aComp,aShell);
1362             }
1363           }
1364           else {
1365             B.Add(aComp,aShell);
1366           }
1367         }
1368       }
1369
1370     }
1371   }
1372
1373   //BRepTools::Write(aComp,"/dn02/users_Linux/skl/work/Bugs/14857/comp.brep");
1374   return aComp;
1375 }
1376
1377 //=======================================================================
1378 //function : CreatePipeShellsWithoutPath
1379 //purpose  : auxilary for Execute()
1380 //=======================================================================
1381 static TopoDS_Shape CreatePipeShellsWithoutPath(GEOMImpl_IPipe* aCI)
1382 {
1383   //cout<<"CreatePipeShellsWithoutPath"<<endl;
1384   int i,j;
1385   BRep_Builder B;
1386
1387   GEOMImpl_IPipeShellSect* aCIDS = (GEOMImpl_IPipeShellSect*)aCI;
1388   // shell sections
1389   Handle(TColStd_HSequenceOfTransient) aBasesObjs = aCIDS->GetBases();
1390   // vertex for recognition
1391   Handle(TColStd_HSequenceOfTransient) VObjs = aCIDS->GetLocations();
1392
1393   Standard_Integer nbBases = aBasesObjs->Length(),
1394     nbv = (VObjs.IsNull() ? 0 :VObjs->Length());
1395
1396   if( nbv != nbBases ) {
1397     if(aCI) delete aCI;
1398     Standard_ConstructionError::Raise("Number of shapes for recognition is invalid");
1399   }
1400
1401   TopTools_SequenceOfShape SecVs,Bases;
1402   for(i=1; i<=nbBases; i++) {
1403     // vertex
1404     Handle(Standard_Transient) anItem = VObjs->Value(i);
1405     if(anItem.IsNull())
1406       continue;
1407     Handle(GEOM_Function) aRef = Handle(GEOM_Function)::DownCast(anItem);
1408     TopoDS_Shape V = aRef->GetValue();
1409     if(V.IsNull() || V.ShapeType() != TopAbs_VERTEX)
1410       continue;
1411     SecVs.Append(V);
1412     // section
1413     anItem = aBasesObjs->Value(i);
1414     if(anItem.IsNull())
1415       continue;
1416     aRef = Handle(GEOM_Function)::DownCast(anItem);
1417     TopoDS_Shape aSh = aRef->GetValue();
1418     if(aSh.IsNull())
1419       continue;
1420     Bases.Append(aSh);
1421   }
1422   nbv = SecVs.Length();
1423   nbBases = Bases.Length();
1424   if( nbv != nbBases ) {
1425     if(aCI) delete aCI;
1426     Standard_ConstructionError::Raise("One of shapes for recognition is not a vertex");
1427   }
1428
1429   TopoDS_Compound aComp;
1430   B.MakeCompound(aComp);
1431
1432   for (i = 1; i < nbBases; i++) {
1433     MESSAGE ("Make pipe between sections "<<i<<" and "<<i+1);
1434     TopoDS_Shape aShBase1 = Bases.Value(i);
1435     TopoDS_Shape aShBase2 = Bases.Value(i+1);
1436     TopExp_Explorer anExp;
1437     Standard_Integer nbf1 = 0;
1438     for ( anExp.Init( aShBase1, TopAbs_FACE ); anExp.More(); anExp.Next() ) {
1439       nbf1++;
1440     }
1441     Standard_Integer nbf2 = 0;
1442     for ( anExp.Init( aShBase2, TopAbs_FACE ); anExp.More(); anExp.Next() ) {
1443       nbf2++;
1444     }
1445     //cout<<"nbf1="<<nbf1<<" nbf2="<<nbf2<<endl;
1446     if(nbf1!=nbf2) {
1447       if(aCI) delete aCI;
1448       Standard_ConstructionError::Raise("Different number of faces in the sections");
1449     }
1450
1451     TopTools_MapOfShape aFaces1,aFaces2;
1452     for ( anExp.Init( aShBase1, TopAbs_FACE ); anExp.More(); anExp.Next() ) {
1453       aFaces1.Add(anExp.Current());
1454     }
1455     for ( anExp.Init( aShBase2, TopAbs_FACE ); anExp.More(); anExp.Next() ) {
1456       aFaces2.Add(anExp.Current());
1457     }
1458
1459     // creating map of edge faces
1460     TopTools_IndexedDataMapOfShapeListOfShape aMapEdgeFaces1;
1461     TopExp::MapShapesAndAncestors(aShBase1, TopAbs_EDGE, TopAbs_FACE, aMapEdgeFaces1);
1462     TopTools_IndexedDataMapOfShapeListOfShape aMapEdgeFaces2;
1463     TopExp::MapShapesAndAncestors(aShBase2, TopAbs_EDGE, TopAbs_FACE, aMapEdgeFaces2);
1464
1465     // constuct map face->face (and subshapes)
1466     TopTools_IndexedDataMapOfShapeShape FF;
1467     //TopoDS_Shape FS1 = SecFs.Value(i), FS2 = SecFs.Value(i+1);
1468     TopoDS_Shape FS1, FS2;
1469     TopoDS_Vertex V1 = TopoDS::Vertex(SecVs(i));
1470     TopoDS_Vertex V2 = TopoDS::Vertex(SecVs(i+1));
1471     FindFirstPairFaces(aShBase1, aShBase2, V1, V2, FS1, FS2);
1472
1473     FF.Add(FS1,FS2);
1474     MESSAGE ("  first pair of corresponding faces is found");
1475
1476     // add pairs of edges and vertexes to FF
1477     bool stat =  FillCorrespondingEdges(FS1, FS2, V1, V2, FF);
1478     if( !stat ) {
1479       if(aCI) delete aCI;
1480       Standard_ConstructionError::Raise("Can not create correct pipe");
1481     }
1482     MESSAGE ("  correspondences for subshapes of first pair of faces is found");
1483
1484     FindNextPairOfFaces(FS1, aMapEdgeFaces1, aMapEdgeFaces2, FF, aCI);
1485     MESSAGE ("  other correspondences is found, make pipe for all pairs of faces");
1486
1487     // make pipe for each pair of faces
1488     // auxilary map vertex->edge for created pipe edges
1489     TopTools_IndexedDataMapOfShapeShape VPE;
1490     ShapeAnalysis_Edge sae;
1491     //cout<<"FF.Extent()="<<FF.Extent()<<endl;
1492     int nbff = 0;
1493     for(j=1; j<=FF.Extent(); j++) {
1494       TopoDS_Shape F1 = FF.FindKey(j);
1495       if( F1.ShapeType() != TopAbs_FACE )
1496         continue;
1497       TopoDS_Shape F2 = FF.FindFromIndex(j);
1498       nbff++;
1499
1500       //if(nbff!=3) continue;
1501
1502       MESSAGE ("    make pipe for "<<nbff<<" face");
1503
1504       Handle(Geom_Surface) S1 = BRep_Tool::Surface(TopoDS::Face(F1));
1505       if(S1->IsKind(STANDARD_TYPE(Geom_RectangularTrimmedSurface))) {
1506         Handle(Geom_RectangularTrimmedSurface) RTS =
1507           Handle(Geom_RectangularTrimmedSurface)::DownCast(S1);
1508         S1 = RTS->BasisSurface();
1509       }
1510       Handle(Geom_Plane) Pln1 = Handle(Geom_Plane)::DownCast(S1);
1511       if( Pln1.IsNull() ) {
1512         if(aCI) delete aCI;
1513         Standard_ConstructionError::Raise("Surface from face is not plane");
1514       }
1515       gp_Vec aDir1(Pln1->Axis().Direction());
1516
1517       Handle(Geom_Surface) S2 = BRep_Tool::Surface(TopoDS::Face(F2));
1518       if(S2->IsKind(STANDARD_TYPE(Geom_RectangularTrimmedSurface))) {
1519         Handle(Geom_RectangularTrimmedSurface) RTS =
1520           Handle(Geom_RectangularTrimmedSurface)::DownCast(S2);
1521         S2 = RTS->BasisSurface();
1522       }
1523       Handle(Geom_Plane) Pln2 =
1524           Handle(Geom_Plane)::DownCast(S2);
1525       if( Pln2.IsNull() ) {
1526         if(aCI) delete aCI;
1527         Standard_ConstructionError::Raise("Surface from face is not plane");
1528       }
1529       gp_Vec aDir2(Pln2->Axis().Direction());
1530
1531       gp_Pnt P1 = BRep_Tool::Pnt(TopoDS::Vertex(SecVs(i)));
1532       gp_Pnt P2 = BRep_Tool::Pnt(TopoDS::Vertex(SecVs(i+1)));
1533       gp_Vec aDir(P1,P2);
1534       if(fabs(aDir.Angle(aDir1))>PI/2.)
1535         aDir1.Reverse();
1536       if(fabs(aDir.Angle(aDir2))>PI/2.)
1537         aDir2.Reverse();
1538
1539       TopExp_Explorer anExpE(F1,TopAbs_EDGE);
1540       TopTools_SequenceOfShape aNewFs;
1541       //int nbee=0;
1542       for(; anExpE.More(); anExpE.Next()) {
1543         TopoDS_Edge E1 = TopoDS::Edge(anExpE.Current());
1544         //nbee++;
1545         if(!FF.Contains(E1))
1546           MESSAGE ("map FF not contains key E1");
1547
1548         if(VPE.Contains(E1)) {
1549           aNewFs.Append(VPE.FindFromKey(E1));
1550 #ifdef _DEBUG_
1551           MESSAGE ("    using existed face");
1552 #endif
1553           continue;
1554         }
1555
1556         TopoDS_Edge E3 = TopoDS::Edge(FF.FindFromKey(E1));
1557         TopoDS_Vertex V1 = sae.FirstVertex(E1);
1558         TopoDS_Vertex V2 = sae.LastVertex(E1);
1559         if(!FF.Contains(V1))
1560           MESSAGE ("map FF not contains key V1");
1561         if(!FF.Contains(V2))
1562           MESSAGE ("map FF not contains key V2");
1563         TopoDS_Vertex V3 = TopoDS::Vertex(FF.FindFromKey(V2));
1564         TopoDS_Vertex V4 = TopoDS::Vertex(FF.FindFromKey(V1));
1565         TopoDS_Vertex Vtmp = sae.FirstVertex(E3);
1566         if(Vtmp.IsSame(V4))
1567           E3.Reverse();
1568         gp_Pnt P1 = BRep_Tool::Pnt(V1);
1569         gp_Pnt P2 = BRep_Tool::Pnt(V2);
1570         gp_Pnt P3 = BRep_Tool::Pnt(V3);
1571         gp_Pnt P4 = BRep_Tool::Pnt(V4);
1572         // make E2
1573         TopoDS_Edge E2;
1574         Handle(Geom_BSplineCurve) C2;
1575         if(VPE.Contains(V2)) {
1576           E2 = TopoDS::Edge(VPE.FindFromKey(V2));
1577           double fp,lp;
1578           C2 = Handle(Geom_BSplineCurve)::DownCast(BRep_Tool::Curve(E2,fp,lp));
1579         }
1580         else {
1581           Handle(TColgp_HArray1OfPnt) HAP = new TColgp_HArray1OfPnt(1,2);
1582           HAP->SetValue(1,P2);
1583           HAP->SetValue(2,P3);
1584           GeomAPI_Interpolate anInt(HAP,Standard_False,1.e-7);
1585           anInt.Load(aDir1,aDir2);
1586           anInt.Perform();
1587           C2 = anInt.Curve();
1588           B.MakeEdge(E2,C2,1.e-7);
1589           B.Add(E2,TopoDS::Vertex(V2.Oriented(TopAbs_FORWARD)));
1590           B.Add(E2,TopoDS::Vertex(V3.Oriented(TopAbs_REVERSED)));
1591           VPE.Add(V2,E2);
1592         }
1593         // make E4
1594         TopoDS_Edge E4;
1595         Handle(Geom_BSplineCurve) C4;
1596         if(VPE.Contains(V1)) {
1597           E4 = TopoDS::Edge(VPE.FindFromKey(V1));
1598           double fp,lp;
1599           C4 = Handle(Geom_BSplineCurve)::DownCast(BRep_Tool::Curve(E4,fp,lp));
1600         }
1601         else {
1602           Handle(TColgp_HArray1OfPnt) HAP = new TColgp_HArray1OfPnt(1,2);
1603           HAP->SetValue(1,P1);
1604           HAP->SetValue(2,P4);
1605           GeomAPI_Interpolate anInt(HAP,Standard_False,1.e-7);
1606           anInt.Load(aDir1,aDir2);
1607           anInt.Perform();
1608           C4 = anInt.Curve();
1609           B.MakeEdge(E4,anInt.Curve(),1.e-7);
1610           B.Add(E4,TopoDS::Vertex(V1.Oriented(TopAbs_FORWARD)));
1611           B.Add(E4,TopoDS::Vertex(V4.Oriented(TopAbs_REVERSED)));
1612           VPE.Add(V1,E4);
1613         }
1614
1615         TopoDS_Wire W;
1616         B.MakeWire(W);
1617         B.Add(W,E1);
1618         B.Add(W,E2);
1619         B.Add(W,E3);
1620         B.Add(W,E4.Reversed());
1621         //cout<<"      wire for edge "<<nbee<<" is created"<<endl;
1622         //BRepTools::Write(W,"/dn02/users_Linux/skl/work/Bugs/14857/w.brep");
1623
1624         // make surface
1625
1626         double fp,lp;
1627         Handle(Geom_Curve) C1 = BRep_Tool::Curve(E1,fp,lp);
1628         //bool IsConicC1 = false;
1629         //if( C1->IsKind(STANDARD_TYPE(Geom_Conic)) ) {
1630         //  IsConicC1 = true;
1631         //  cout<<"C1 - Geom_Conic"<<endl;
1632         //}
1633         if( C1->IsKind(STANDARD_TYPE(Geom_Line)) || C1->IsKind(STANDARD_TYPE(Geom_Conic)) ) {
1634           C1 = new Geom_TrimmedCurve(C1,fp,lp);
1635         }
1636         //if(IsConicC1) {
1637         //  double tol = BRep_Tool::Tolerance(E1);
1638         //  GeomConvert_ApproxCurve ApxC1(C1,tol,GeomAbs_C1,10,5);
1639         //  C1 = ApxC1.Curve();
1640         //}
1641         Handle(Geom_Curve) C3 = BRep_Tool::Curve(E3,fp,lp);
1642         if( C3->IsKind(STANDARD_TYPE(Geom_Line)) || C3->IsKind(STANDARD_TYPE(Geom_Conic)) ) {
1643           C3 = new Geom_TrimmedCurve(C3,fp,lp);
1644         }
1645         //filebuf fic;
1646         //ostream os(&fic);
1647         //os.precision(15);
1648         Handle(Geom_BSplineCurve) CE1 =
1649           GeomConvert::CurveToBSplineCurve(C1,Convert_RationalC1);
1650         if(CE1->Degree()<3)
1651           CE1->IncreaseDegree(3);
1652         Handle(Geom_BSplineCurve) CE2 =
1653           GeomConvert::CurveToBSplineCurve(C2,Convert_RationalC1);
1654         if(CE2->Degree()<3)
1655           CE2->IncreaseDegree(3);
1656         Handle(Geom_BSplineCurve) CE3 =
1657           GeomConvert::CurveToBSplineCurve(C3,Convert_RationalC1);
1658         if(CE3->Degree()<3)
1659           CE3->IncreaseDegree(3);
1660         Handle(Geom_BSplineCurve) CE4 =
1661           GeomConvert::CurveToBSplineCurve(C4,Convert_RationalC1);
1662         if(CE4->Degree()<3)
1663           CE4->IncreaseDegree(3);
1664         //cout<<"CE1->Degree()="<<CE1->Degree()<<" CE2->Degree()="<<CE2->Degree()
1665         //    <<" CE3->Degree()="<<CE3->Degree()<<" CE4->Degree()="<<CE4->Degree()<<endl;
1666         //if(fic.open("/dn02/users_Linux/skl/work/Bugs/14857/ce1.brep",ios::out)) {
1667         //  os<<"DrawTrSurf_BSplineCurve"<<endl;
1668         //  GeomTools::Write(CE1,os);
1669         //  fic.close();
1670         //}
1671
1672         Handle(Geom_Surface) BS;
1673         try {
1674           GeomFill_BSplineCurves GF(CE1,CE2,CE3,CE4,GeomFill_CoonsStyle);
1675           //GeomFill_BSplineCurves GF(CE1,CE2,CE3,CE4,GeomFill_StretchStyle);
1676           BS = GF.Surface();
1677         }
1678         catch(...) {
1679           MESSAGE ("      can not create BSplineSurface - create Bezier");
1680           int NbP=26;
1681           TColgp_Array2OfPnt Points(1,NbP,1,NbP);
1682           double fp1,lp1,fp2,lp2;
1683           Handle(Geom_Curve) C1 = BRep_Tool::Curve(E1,fp1,lp1);
1684           Handle(Geom_Curve) C3 = BRep_Tool::Curve(E3,fp2,lp2);
1685           gp_Pnt P1C1,P2C1;
1686           C1->D0(fp1,P1C1);
1687           C1->D0(lp1,P2C1);
1688           gp_Pnt P1C3,P2C3;
1689           C3->D0(fp2,P1C3);
1690           C3->D0(lp2,P2C3);
1691           int n1,n2;
1692           double fp,lp;
1693           // get points from C1
1694           if(P1.Distance(P1C1)<1.e-6) {
1695             fp = fp1;
1696             lp = lp1;
1697           }
1698           else {
1699             fp = lp1;
1700             lp = fp1;
1701           }
1702           double step = (lp-fp)/(NbP-1);
1703           Points.SetValue(1,1,P1);
1704           double par = fp;
1705           for(n1=2; n1<NbP; n1++) {
1706             gp_Pnt P;
1707             par += step;
1708             C1->D0(par,P);
1709             Points.SetValue(1,n1,P);
1710           }
1711           Points.SetValue(1,NbP,P2);
1712           // get points from C3
1713           if(P4.Distance(P1C3)<1.e-6) {
1714             fp = fp2;
1715             lp = lp2;
1716           }
1717           else {
1718             fp = lp2;
1719             lp = fp2;
1720           }
1721           step = (lp-fp)/(NbP-1);
1722           Points.SetValue(NbP,1,P4);
1723           par = fp;
1724           for(n1=2; n1<NbP; n1++) {
1725             gp_Pnt P;
1726             par += step;
1727             C3->D0(par,P);
1728             Points.SetValue(NbP,n1,P);
1729           }
1730           Points.SetValue(NbP,NbP,P3);
1731           // create isolines and get points from them
1732           for(n1=1; n1<=NbP; n1++) {
1733             gp_Pnt PI1 = Points.Value(1,n1);
1734             gp_Pnt PI2 = Points.Value(NbP,n1);
1735             Handle(TColgp_HArray1OfPnt) HAP = new TColgp_HArray1OfPnt(1,2);
1736             HAP->SetValue(1,PI1);
1737             HAP->SetValue(2,PI2);
1738             GeomAPI_Interpolate anInt(HAP,Standard_False,1.e-7);
1739             anInt.Load(aDir1,aDir2);
1740             anInt.Perform();
1741             Handle(Geom_Curve) iso = anInt.Curve();
1742             fp = iso->FirstParameter();
1743             lp = iso->LastParameter();
1744             step = (lp-fp)/(NbP-1);
1745             par = fp;
1746             TopoDS_Compound VComp;
1747             B.MakeCompound(VComp);
1748             for(n2=2; n2<NbP; n2++) {
1749               gp_Pnt P;
1750               par += step;
1751               iso->D0(par,P);
1752               Points.SetValue(n2,n1,P);
1753             }
1754           }
1755           // create surface and face
1756           //Handle(Geom_BezierSurface) BS = new Geom_BezierSurface(Points);
1757           BS = new Geom_BezierSurface(Points);
1758         }
1759
1760         BRepBuilderAPI_MakeFace BB(BS,W);
1761         TopoDS_Face NewF = BB.Face();
1762         Handle(ShapeFix_Face) sff = new ShapeFix_Face(NewF);
1763         sff->Perform();
1764         sff->FixOrientation();
1765         TopoDS_Face FixedFace = sff->Face();
1766         aNewFs.Append(FixedFace);
1767         VPE.Add(E1,FixedFace);
1768         //cout<<"      face for edge "<<nbee<<" is created"<<endl;
1769         //BRepTools::Write(FixedFace,"/dn02/users_Linux/skl/work/Bugs/14857/f.brep");
1770       }
1771       // make shell
1772       TopoDS_Shell aShell;
1773       B.MakeShell(aShell);
1774       for(int nf=1; nf<=aNewFs.Length(); nf++) {
1775         B.Add(aShell,aNewFs(nf));
1776       }
1777       B.Add(aShell,F1);
1778       B.Add(aShell,F2);
1779
1780       // make sewing for this shell
1781       Handle(BRepBuilderAPI_Sewing) aSewing = new BRepBuilderAPI_Sewing;
1782       aSewing->SetTolerance(Precision::Confusion());
1783       aSewing->SetFaceMode(Standard_True);
1784       aSewing->SetFloatingEdgesMode(Standard_False);
1785       aSewing->SetNonManifoldMode(Standard_False);
1786       for ( anExp.Init( aShell, TopAbs_FACE ); anExp.More(); anExp.Next() ) {
1787         aSewing->Add(anExp.Current());
1788       }
1789       aSewing->Perform();
1790       MESSAGE ("    shell for face "<<nbff<<" is created");
1791       const TopoDS_Shape aSewShape = aSewing->SewedShape();
1792       //BRepTools::Write(aSewShape,"/dn02/users_Linux/skl/work/Bugs/14857/sew.brep");
1793       if( aSewShape.ShapeType() == TopAbs_SHELL ) {
1794         aShell = TopoDS::Shell(aSewShape);
1795         GProp_GProps aSystem;
1796         BRepGProp::VolumeProperties(aShell, aSystem);
1797         if(aSystem.Mass()<0) {
1798           //cout<<"aSewShape is reversed"<<endl;
1799           aShell.Reverse();
1800         }
1801         if(BRep_Tool::IsClosed(aShell)) {
1802           TopoDS_Solid aSolid;
1803           B.MakeSolid(aSolid);
1804           B.Add(aSolid,aShell);
1805           B.Add(aComp,aSolid);
1806           MESSAGE ("    solid for face "<<nbff<<" is created");
1807         }
1808         else {
1809           B.Add(aComp,aShell);
1810           MESSAGE ("    solid for face "<<nbff<<" is not created");
1811         }
1812       }
1813       else {
1814         B.Add(aComp,aShell);
1815         MESSAGE ("    solid for face "<<nbff<<" is not created");
1816       }
1817       //cout<<"    solid for face "<<nbff<<" is created"<<endl;
1818
1819       //Handle(ShapeFix_Shell) sfs = new ShapeFix_Shell(aShell);
1820       //sfs->Perform();
1821       //TopoDS_Shell FixedShell = sfs->Shell();
1822       /*
1823       GProp_GProps aSystem;
1824       BRepGProp::VolumeProperties(FixedShell, aSystem);
1825       if(aSystem.Mass()<0) {
1826         //cout<<"aSewShape is reversed"<<endl;
1827         FixedShell.Reverse();
1828       }
1829       if(BRep_Tool::IsClosed(FixedShell)) {
1830         TopoDS_Solid aSolid;
1831         B.MakeSolid(aSolid);
1832         B.Add(aSolid,aShell);
1833         B.Add(aComp,aSolid);
1834       }
1835       else {
1836         B.Add(aComp,FixedShell);
1837       }
1838       */
1839     }
1840   }
1841
1842   //BRepTools::Write(aComp,"/dn02/users_Linux/skl/work/Bugs/14857/comp.brep");
1843   return aComp;
1844 }
1845
1846 //=======================================================================
1847 //function : CreatePipeBiNormalAlongVector
1848 //purpose  : auxilary for Execute()
1849 //=======================================================================
1850 static TopoDS_Shape CreatePipeBiNormalAlongVector(const TopoDS_Wire& aWirePath,
1851                                                   GEOMImpl_IPipe* aCI)
1852 {
1853   GEOMImpl_IPipeBiNormal* aCIBN = (GEOMImpl_IPipeBiNormal*)aCI;
1854
1855   Handle(GEOM_Function) aRefBase = aCIBN->GetBase();
1856   Handle(GEOM_Function) aRefVec = aCIBN->GetVector();
1857   TopoDS_Shape aShapeBase = aRefBase->GetValue();
1858   TopoDS_Shape aShapeVec = aRefVec->GetValue();
1859
1860   if (aShapeBase.IsNull()) {
1861     if(aCIBN) delete aCIBN;
1862     Standard_NullObject::Raise("MakePipe aborted : null base argument");
1863   }
1864
1865   TopoDS_Shape aProf;
1866   if( aShapeBase.ShapeType() == TopAbs_VERTEX ) {
1867     aProf = aShapeBase;
1868   }
1869   else if( aShapeBase.ShapeType() == TopAbs_EDGE) {
1870     aProf = BRepBuilderAPI_MakeWire(TopoDS::Edge(aShapeBase)).Shape();
1871   }
1872   else if( aShapeBase.ShapeType() == TopAbs_WIRE) {
1873     aProf = aShapeBase;
1874   }
1875   else if( aShapeBase.ShapeType() == TopAbs_FACE) {
1876     TopExp_Explorer wexp(aShapeBase,TopAbs_WIRE);
1877     aProf = wexp.Current();
1878   }
1879   else {
1880     Standard_TypeMismatch::Raise
1881       ("MakePipe aborted : invalid type of base");
1882   }
1883   BRepOffsetAPI_MakePipeShell PipeBuilder(aWirePath);
1884   PipeBuilder.Add(aProf);
1885
1886   if (aShapeVec.IsNull()) {
1887     if(aCIBN) delete aCIBN;
1888     Standard_NullObject::Raise
1889       ("MakePipe aborted : null vector argument");
1890   }
1891   if (aShapeVec.ShapeType() != TopAbs_EDGE)
1892     Standard_TypeMismatch::Raise
1893       ("MakePipe aborted: invalid type of vector");
1894   TopoDS_Edge anEdge = TopoDS::Edge(aShapeVec);
1895   TopoDS_Vertex V1, V2;
1896   TopExp::Vertices(anEdge, V1, V2, Standard_True);
1897   if (V1.IsNull() || V2.IsNull())
1898     Standard_NullObject::Raise
1899       ("MakePipe aborted: vector is not defined");
1900   gp_Vec aVec(BRep_Tool::Pnt(V1), BRep_Tool::Pnt(V2));
1901   gp_Dir BiNormal(aVec);
1902   PipeBuilder.SetMode(BiNormal);
1903   PipeBuilder.Build();
1904   if( aShapeBase.ShapeType() == TopAbs_FACE) {
1905       PipeBuilder.MakeSolid();
1906   }
1907
1908   return PipeBuilder.Shape();
1909 }
1910
1911 //=======================================================================
1912 //function : Execute
1913 //purpose  :
1914 //=======================================================================
1915 Standard_Integer GEOMImpl_PipeDriver::Execute(TFunction_Logbook& log) const
1916 {
1917   //cout<<"PipeDriver::Execute"<<endl;
1918   if (Label().IsNull()) return 0;
1919   Handle(GEOM_Function) aFunction = GEOM_Function::GetFunction(Label());
1920   GEOMImpl_IPipe* aCI= 0;
1921   Standard_Integer aType = aFunction->GetType();
1922   if(aType == PIPE_BASE_PATH)
1923     aCI = new GEOMImpl_IPipe(aFunction);
1924   else if(aType == PIPE_DIFFERENT_SECTIONS)
1925     aCI = new GEOMImpl_IPipeDiffSect(aFunction);
1926   else if(aType == PIPE_SHELL_SECTIONS)
1927     aCI = new GEOMImpl_IPipeShellSect(aFunction);
1928   else if(aType == PIPE_SHELLS_WITHOUT_PATH)
1929     aCI = new GEOMImpl_IPipeShellSect(aFunction);
1930   else if(aType == PIPE_BI_NORMAL_ALONG_VECTOR)
1931     aCI = new GEOMImpl_IPipeBiNormal(aFunction);
1932   else
1933     return 0;
1934
1935   TopoDS_Wire aWirePath;
1936   if(aType != PIPE_SHELLS_WITHOUT_PATH) {
1937     // working with path
1938     Handle(GEOM_Function) aRefPath = aCI->GetPath();
1939     TopoDS_Shape aShapePath = aRefPath->GetValue();
1940
1941     if (aShapePath.IsNull()) {
1942       MESSAGE ("Driver : path is null");
1943       if(aCI) delete aCI;
1944       Standard_NullObject::Raise("MakePipe aborted : null path argument");
1945     }
1946
1947     // Get path contour
1948     if (aShapePath.ShapeType() == TopAbs_WIRE) {
1949       aWirePath = TopoDS::Wire(aShapePath);
1950     }
1951     else {
1952       if (aShapePath.ShapeType() == TopAbs_EDGE) {
1953         TopoDS_Edge anEdge = TopoDS::Edge(aShapePath);
1954         aWirePath = BRepBuilderAPI_MakeWire(anEdge);
1955       }
1956       else {
1957         if(aCI) delete aCI;
1958         Standard_TypeMismatch::Raise("MakePipe aborted : path shape is neither a wire nor an edge");
1959       }
1960     }
1961   }
1962
1963   TopoDS_Shape aShape;
1964
1965   if (aType == PIPE_BASE_PATH)
1966   {
1967     Handle(GEOM_Function) aRefBase = aCI->GetBase();
1968     TopoDS_Shape aShapeBase;
1969
1970     // Make copy to prevent modifying of base object 0020766 : EDF 1320
1971     BRepBuilderAPI_Copy Copy(aRefBase->GetValue());
1972     if( Copy.IsDone() )
1973       aShapeBase = Copy.Shape();
1974
1975     if (aShapeBase.IsNull()) {
1976       if(aCI) delete aCI;
1977       Standard_NullObject::Raise("MakePipe aborted : null base argument");
1978     }
1979
1980     // Make pipe
1981     aShape = BRepOffsetAPI_MakePipe(aWirePath, aShapeBase);
1982   }
1983
1984   //building pipe with different sections
1985   else if (aType == PIPE_DIFFERENT_SECTIONS) {
1986     GEOMImpl_IPipeDiffSect* aCIDS = (GEOMImpl_IPipeDiffSect*)aCI;
1987     //GEOMImpl_IPipeDiffSect* aCIDS = static_cast<GEOMImpl_IPipeDiffSect*>(aCI);
1988     //BRepOffsetAPI_MakePipeShell aBuilder(aWirePath);
1989     Handle(TColStd_HSequenceOfTransient) aBasesObjs = aCIDS->GetBases ();
1990     Handle(TColStd_HSequenceOfTransient) aLocObjs = aCIDS->GetLocations ();
1991     Standard_Boolean aWithContact = (aCIDS->GetWithContactMode());
1992     Standard_Boolean aWithCorrect = (aCIDS->GetWithCorrectionMode());
1993
1994     Standard_Integer i =1, nbBases = aBasesObjs->Length(),
1995       nbLocs = (aLocObjs.IsNull() ? 0 :aLocObjs->Length());
1996
1997     if(nbLocs && nbLocs != nbBases) {
1998       if(aCI) delete aCI;
1999       Standard_ConstructionError::Raise("Number of sections is not equal to number of locations ");
2000     }
2001     TopTools_SequenceOfShape aSeqBases;
2002     TopTools_SequenceOfShape aSeqLocs;
2003     TopTools_SequenceOfShape aSeqFaces;
2004     for (; i <= nbBases; i++) {
2005       Handle(Standard_Transient) anItem = aBasesObjs->Value(i);
2006       if(anItem.IsNull())
2007         continue;
2008       Handle(GEOM_Function) aRefBase = Handle(GEOM_Function)::DownCast(anItem);
2009       if(aRefBase.IsNull())
2010         continue;
2011
2012       if(aRefBase->GetValue().IsNull())
2013         continue;
2014
2015       // Make copy to prevent modifying of base object 0020766 : EDF 1320
2016       TopoDS_Shape aShapeBase;
2017       BRepBuilderAPI_Copy Copy(aRefBase->GetValue());
2018       if( Copy.IsDone() )
2019         aShapeBase = Copy.Shape();
2020
2021       TopAbs_ShapeEnum aTypeBase = aShapeBase.ShapeType();
2022
2023       //if for section was specified face with a few wires then a few
2024       //    pipes were build and make solid
2025       Standard_Boolean NeedCreateSolid = Standard_False;
2026       if(aTypeBase == TopAbs_SHELL) {
2027         // create wire as boundary contour if shell is no closed
2028         // get free boundary shapes
2029         ShapeAnalysis_FreeBounds anAnalizer( aShapeBase );
2030         TopoDS_Compound aClosed = anAnalizer.GetClosedWires();
2031         TopExp_Explorer anExp;
2032         TopoDS_Shape aWire;
2033         Standard_Integer NbWires = 0;
2034         for ( anExp.Init( aClosed, TopAbs_WIRE ); anExp.More(); anExp.Next() ) {
2035           NbWires++;
2036           aWire = anExp.Current();
2037         }
2038         if(NbWires!=1) {
2039           // bad case
2040           if(aCI) delete aCI;
2041           Standard_ConstructionError::Raise("Bad shell is used as section ");
2042         }
2043         NeedCreateSolid = Standard_True;
2044         aSeqFaces.Append(aShapeBase);
2045         aSeqBases.Append(aWire);
2046       }
2047       else if(aTypeBase == TopAbs_FACE) {
2048         NeedCreateSolid = Standard_True;
2049         //for case one path should be used other type function
2050         aSeqFaces.Append(aShapeBase);
2051         TopExp_Explorer aExpW(aShapeBase,TopAbs_WIRE);
2052         for (; aExpW.More(); aExpW.Next())
2053         {
2054           TopoDS_Shape aWireProf = aExpW.Current();
2055           aSeqBases.Append(aWireProf);
2056         }
2057       }
2058       else if(aTypeBase == TopAbs_WIRE || aTypeBase == TopAbs_VERTEX) {
2059         aSeqBases.Append(aShapeBase);
2060       }
2061       else if(aTypeBase == TopAbs_EDGE) {
2062         TopoDS_Edge anEdge = TopoDS::Edge(aShapeBase);
2063         TopoDS_Shape aWireProf = BRepBuilderAPI_MakeWire(anEdge);
2064         aSeqBases.Append(aWireProf);
2065       }
2066       if(nbLocs) {
2067         Handle(Standard_Transient) anItemLoc = aLocObjs->Value(i);
2068         if(anItemLoc.IsNull())
2069           continue;
2070         Handle(GEOM_Function) aRefLoc = Handle(GEOM_Function)::DownCast(anItemLoc);
2071         TopoDS_Shape aShapeLoc = aRefLoc->GetValue();
2072         if(aShapeLoc.IsNull() || aShapeLoc.ShapeType() != TopAbs_VERTEX)
2073           continue;
2074         aSeqLocs.Append(aShapeLoc);
2075       }
2076     }
2077
2078     nbLocs = aSeqLocs.Length();
2079
2080     // skl 02.05.2007
2081     TopTools_SequenceOfShape Edges;
2082     if(nbLocs>0) {
2083       // we have to check that each location shape is a vertex from
2084       // path and update aSeqLocs if it is needed (and possible)
2085       TColgp_SequenceOfPnt PLocs;
2086       for(i=1; i<=nbLocs; i++) {
2087         TopoDS_Vertex V = TopoDS::Vertex(aSeqLocs.Value(i));
2088         PLocs.Append(BRep_Tool::Pnt(V));
2089       }
2090       //TopTools_SequenceOfShape Edges;
2091       TopExp_Explorer anExp;
2092       for ( anExp.Init( aWirePath, TopAbs_EDGE ); anExp.More(); anExp.Next() ) {
2093         Edges.Append(anExp.Current());
2094       }
2095       int nbEdges = Edges.Length();
2096       ShapeAnalysis_Edge sae;
2097       TopoDS_Edge edge = TopoDS::Edge(Edges.First());
2098       double tol = BRep_Tool::Tolerance(edge);
2099       TopoDS_Vertex VF = sae.FirstVertex(edge);
2100       gp_Pnt PF = BRep_Tool::Pnt(VF);
2101       //cout<<"PF("<<PF.X()<<","<<PF.Y()<<","<<PF.Z()<<")"<<endl;
2102       if( PF.Distance(PLocs.First()) > tol ) {
2103         if(aCI) delete aCI;
2104         Standard_ConstructionError::Raise
2105           ("First location shapes is not coincided with first vertex of aWirePath");
2106       }
2107       aSeqLocs.ChangeValue(1) = VF;
2108       edge = TopoDS::Edge(Edges.Last());
2109       tol = BRep_Tool::Tolerance(edge);
2110       TopoDS_Vertex VL = sae.LastVertex(edge);
2111       gp_Pnt PL = BRep_Tool::Pnt(VL);
2112       if( PL.Distance(PLocs.Last()) > tol ) {
2113         if(aCI) delete aCI;
2114         Standard_ConstructionError::Raise
2115           ("Last location shapes is not coincided with last vertex of aWirePath");
2116       }
2117       aSeqLocs.ChangeValue(nbLocs) = VL;
2118       int jcurr = 2;
2119       for(i=1; i<=Edges.Length() && jcurr<nbLocs; i++) {
2120         TopoDS_Edge E = TopoDS::Edge(Edges.Value(i));
2121         tol = BRep_Tool::Tolerance(edge);
2122         TopoDS_Vertex V1 = sae.FirstVertex(E);
2123         TopoDS_Vertex V2 = sae.LastVertex(E);
2124         gp_Pnt P1 = BRep_Tool::Pnt(V1);
2125         gp_Pnt P2 = BRep_Tool::Pnt(V2);
2126         if( P2.Distance(PLocs.Value(jcurr)) < tol ) {
2127           aSeqLocs.ChangeValue(jcurr) = V2;
2128           jcurr++;
2129         }
2130         else {
2131           // find distance between E and aLocs(jcurr)
2132           double fp,lp;
2133           Handle(Geom_Curve) C = BRep_Tool::Curve(E,fp,lp);
2134           GeomAPI_ProjectPointOnCurve PPC (PLocs.Value(jcurr),C);
2135           if( PPC.NbPoints()>0 &&
2136               PLocs.Value(jcurr).Distance(PPC.Point(1)) < tol ) {
2137             double param = PPC.Parameter(1);
2138             gp_Pnt PC1;
2139             C->D0(param,PC1);
2140             // split current edge
2141             Handle(Geom_TrimmedCurve) tc1 = new Geom_TrimmedCurve(C,fp,param);
2142             Handle(Geom_TrimmedCurve) tc2 = new Geom_TrimmedCurve(C,param,lp);
2143             TopoDS_Edge E1,E2;
2144             BRep_Builder B;
2145             gp_Pnt Pfp;
2146             C->D0(fp,Pfp);
2147             if(Pfp.Distance(P1)<tol) {
2148               B.MakeEdge(E1,tc1,tol);
2149               B.Add(E1,V1);
2150               TopoDS_Shape tmpV = aSeqLocs.Value(jcurr).Oriented(TopAbs_REVERSED);
2151               B.Add(E1,TopoDS::Vertex(tmpV));
2152               B.MakeEdge(E2,tc2,tol);
2153               tmpV = aSeqLocs.Value(jcurr).Oriented(TopAbs_FORWARD);
2154               B.Add(E2,TopoDS::Vertex(tmpV));
2155               B.Add(E2,V2);
2156             }
2157             else {
2158               B.MakeEdge(E1,tc2,tol);
2159               TopoDS_Shape tmpV = aSeqLocs.Value(jcurr).Oriented(TopAbs_FORWARD);
2160               B.Add(E1,TopoDS::Vertex(tmpV));
2161               B.Add(E1,V1);
2162               E1.Reverse();
2163               B.MakeEdge(E2,tc1,tol);
2164               B.Add(E2,V2);
2165               tmpV = aSeqLocs.Value(jcurr).Oriented(TopAbs_REVERSED);
2166               B.Add(E2,TopoDS::Vertex(tmpV));
2167               E2.Reverse();
2168             }
2169             jcurr++;
2170             Edges.Remove(i);
2171             Edges.InsertAfter(i-1,E1);
2172             Edges.InsertAfter(i,E2);
2173           }
2174         }
2175       }
2176       if(nbEdges<Edges.Length()) {
2177         // one of edges was splitted => we have to update WirePath
2178         BRep_Builder B;
2179         TopoDS_Wire W;
2180         B.MakeWire(W);
2181         for(i=1; i<=Edges.Length(); i++) {
2182           B.Add(W,TopoDS::Edge(Edges.Value(i)));
2183         }
2184         aWirePath = W;
2185       }
2186     }
2187
2188     // check curvature of wire for condition that
2189     // max summary angle between directions along
2190     // wire path must be < 4*PI. If not - split wire
2191     // and seguences of shapes, perform pipe for each
2192     // and make sewing after that
2193     double fp,lp;
2194     Handle(Geom_Curve) C = BRep_Tool::Curve(TopoDS::Edge(Edges.Value(1)),fp,lp);
2195     gp_Pnt P1,P2;
2196     gp_Vec Vec1,Vec2;
2197     C->D1(fp,P1,Vec1);
2198     C->D1(lp,P2,Vec2);
2199     double SumAng = fabs(Vec1.Angle(Vec2));
2200     Vec1 = Vec2;
2201     P1 = P2;
2202     TColStd_SequenceOfInteger SplitEdgeNums,SplitLocNums;
2203     int LastLoc = 1;
2204     //cout<<"Edges.Length()="<<Edges.Length()<<endl;
2205     for(i=2; i<=Edges.Length(); i++) {
2206       TopoDS_Edge edge = TopoDS::Edge(Edges.Value(i));
2207       double tol = BRep_Tool::Tolerance(edge);
2208       Handle(Geom_Curve) C = BRep_Tool::Curve(edge,fp,lp);
2209       C->D1(lp,P2,Vec2);
2210       double ang = fabs(Vec1.Angle(Vec2));
2211       SumAng += ang;
2212       if(SumAng>4*PI) {
2213         SumAng = ang;
2214         SplitEdgeNums.Append(i-1);
2215         int j;
2216         for(j=LastLoc+1; j<=aSeqLocs.Length(); j++) {
2217           TopoDS_Vertex aVert = TopoDS::Vertex(aSeqLocs.Value(j));
2218           gp_Pnt P = BRep_Tool::Pnt(aVert);
2219           if( P1.Distance(P) < tol ) {
2220             SplitLocNums.Append(j);
2221             LastLoc = j;
2222             break;
2223           }
2224         }
2225       }
2226       Vec1 = Vec2;
2227       P1 = P2;
2228     }
2229
2230     //cout<<"SplitEdgeNums.Length()="<<SplitEdgeNums.Length()<<endl;
2231     //cout<<"SplitLocNums.Length()="<<SplitLocNums.Length()<<endl;
2232     if( SplitLocNums.Length()==SplitEdgeNums.Length() && SplitEdgeNums.Length()>0 ) {
2233       TopTools_SequenceOfShape aSeqRes;
2234       int nn, num1 = 1, num2 = 1;
2235       for(nn=1; nn<=SplitEdgeNums.Length(); nn++) {
2236         // create wirepath and sequences of shapes
2237         BRep_Builder B;
2238         TopoDS_Wire tmpW;
2239         B.MakeWire(tmpW);
2240         for(i=num1; i<=SplitEdgeNums.Value(nn); i++) {
2241           B.Add(tmpW,TopoDS::Edge(Edges.Value(i)));
2242         }
2243         num1 = SplitEdgeNums.Value(nn) + 1;
2244         TopTools_SequenceOfShape aTmpSeqBases;
2245         TopTools_SequenceOfShape aTmpSeqLocs;
2246         for(i=num2; i<=SplitLocNums.Value(nn); i++) {
2247           aTmpSeqBases.Append(aSeqBases.Value(i));
2248           aTmpSeqLocs.Append(aSeqLocs.Value(i));
2249         }
2250         num2 = SplitLocNums.Value(nn);
2251         // make pipe
2252         BRepOffsetAPI_MakePipeShell aBuilder(tmpW);
2253         Standard_Integer nbShapes = aTmpSeqBases.Length();
2254         for(i=1; i<=nbShapes; i++) {
2255           TopoDS_Shape aShapeLoc = aTmpSeqLocs.Value(i);
2256           TopoDS_Vertex aVert = TopoDS::Vertex(aShapeLoc);
2257           aBuilder.Add(aTmpSeqBases.Value(i), aVert, aWithContact, aWithCorrect);
2258         }
2259         if(!aBuilder.IsReady()) {
2260           if(aCI) delete aCI;
2261           Standard_ConstructionError::Raise("Invalid input data for building PIPE: bases are invalid");
2262         }
2263         aBuilder.Build();
2264         TopoDS_Shape resShape = aBuilder.Shape();
2265         aSeqRes.Append(resShape);
2266       }
2267       // create wirepath and sequences of shapes for last part
2268       BRep_Builder B;
2269       TopoDS_Wire tmpW;
2270       B.MakeWire(tmpW);
2271       for(i=num1; i<=Edges.Length(); i++) {
2272         B.Add(tmpW,TopoDS::Edge(Edges.Value(i)));
2273       }
2274       TopTools_SequenceOfShape aTmpSeqBases;
2275       TopTools_SequenceOfShape aTmpSeqLocs;
2276       for(i=num2; i<=aSeqLocs.Length(); i++) {
2277         aTmpSeqBases.Append(aSeqBases.Value(i));
2278         aTmpSeqLocs.Append(aSeqLocs.Value(i));
2279       }
2280       // make pipe for last part
2281       BRepOffsetAPI_MakePipeShell aBuilder(tmpW);
2282       Standard_Integer nbShapes = aTmpSeqBases.Length();
2283       for(i=1; i<=nbShapes; i++) {
2284         TopoDS_Shape aShapeLoc = aTmpSeqLocs.Value(i);
2285         TopoDS_Vertex aVert = TopoDS::Vertex(aShapeLoc);
2286         aBuilder.Add(aTmpSeqBases.Value(i), aVert, aWithContact, aWithCorrect);
2287       }
2288       if(!aBuilder.IsReady()) {
2289         if(aCI) delete aCI;
2290         Standard_ConstructionError::Raise("Invalid input data for building PIPE: bases are invalid");
2291       }
2292       aBuilder.Build();
2293       TopoDS_Shape resShape = aBuilder.Shape();
2294       aSeqRes.Append(resShape);
2295       // make sewing for result
2296       Handle(BRepBuilderAPI_Sewing) aSewing = new BRepBuilderAPI_Sewing;
2297       aSewing->SetTolerance(Precision::Confusion());
2298       aSewing->SetFaceMode(Standard_True);
2299       aSewing->SetFloatingEdgesMode(Standard_False);
2300       aSewing->SetNonManifoldMode(Standard_False);
2301       for(i=1; i<=aSeqRes.Length(); i++) {
2302         aSewing->Add(aSeqRes.Value(i));
2303       }
2304       aSewing->Perform();
2305       aShape = aSewing->SewedShape();
2306     }
2307     else {
2308       // old implementation without splitting
2309       BRepOffsetAPI_MakePipeShell aBuilder(aWirePath);
2310
2311       Standard_Integer nbShapes = aSeqBases.Length();
2312       Standard_Integer step = nbShapes/nbBases;
2313
2314       if (nbShapes < nbBases || fmod((double)nbShapes, (double)nbBases)) {
2315         if(aCI) delete aCI;
2316         Standard_ConstructionError::Raise("Invalid sections were specified for building pipe");
2317       }
2318       Standard_Integer ind =0;
2319       for (i = 1; i <= nbShapes && ind < nbShapes; i++) { //i+nbBases <= nbShapes
2320         TopTools_SequenceOfShape usedBases;
2321         Standard_Integer j = 1;
2322         for (; j <= nbBases; j++) {
2323           ind = i + (j-1)*step;
2324           TopoDS_Shape aWireProf = aSeqBases.Value(ind);
2325           usedBases.Append(aWireProf);
2326           if(nbLocs) {
2327             TopoDS_Shape aShapeLoc = aSeqLocs.Value(j);
2328             TopoDS_Vertex aVert = TopoDS::Vertex(aShapeLoc);
2329             aBuilder.Add(aWireProf,aVert,aWithContact,aWithCorrect);
2330           }
2331           else
2332             aBuilder.Add(aWireProf,aWithContact,aWithCorrect);
2333         }
2334         if(!aBuilder.IsReady()) {
2335           if(aCI) delete aCI;
2336           Standard_ConstructionError::Raise("Invalid input data for building PIPE: bases are invalid");
2337         }
2338         aBuilder.Build();
2339         aShape = aBuilder.Shape();
2340         aSeqFaces.Append(aShape);
2341         for( j = 1; j <=usedBases.Length(); j++)
2342           aBuilder.Delete(usedBases.Value(j));
2343       }
2344
2345       //for case if section is face
2346       if(aSeqFaces.Length() >1) {
2347         BRep_Builder aB;
2348         TopoDS_Compound aComp;
2349         aB.MakeCompound(aComp);
2350         for( i = 1; i <= aSeqFaces.Length(); i++)
2351           aB.Add(aComp,aSeqFaces.Value(i));
2352         aShape = aComp;
2353       }
2354     }
2355   }
2356
2357   //building pipe with shell sections
2358   else if (aType == PIPE_SHELL_SECTIONS) {
2359     aShape = CreatePipeForShellSections(aWirePath,aCI);
2360   }
2361
2362   //building pipe shell sections without path
2363   else if (aType == PIPE_SHELLS_WITHOUT_PATH) {
2364     aShape = CreatePipeShellsWithoutPath(aCI);
2365   }
2366
2367   //building a pipe with constant bi-normal along given vector
2368   else if (aType == PIPE_BI_NORMAL_ALONG_VECTOR) {
2369     aShape = CreatePipeBiNormalAlongVector(aWirePath, aCI);
2370   }
2371
2372   if (aCI) {
2373     delete aCI;
2374     aCI = 0;
2375   }
2376
2377   if (aShape.IsNull()) return 0;
2378
2379   BRepCheck_Analyzer ana (aShape, Standard_False);
2380   if (!ana.IsValid()) {
2381     ShapeFix_ShapeTolerance aSFT;
2382     aSFT.LimitTolerance(aShape,Precision::Confusion(),Precision::Confusion());
2383     Handle(ShapeFix_Shape) aSfs = new ShapeFix_Shape(aShape);
2384     aSfs->SetPrecision(Precision::Confusion());
2385     aSfs->Perform();
2386     aShape = aSfs->Shape();
2387
2388     ana.Init(aShape, Standard_False);
2389     if (!ana.IsValid())
2390       Standard_ConstructionError::Raise("Algorithm have produced an invalid shape result");
2391   }
2392
2393   // Glue (for bug 0020207)
2394   TopExp_Explorer anExpV (aShape, TopAbs_VERTEX);
2395   if (anExpV.More())
2396     aShape = GEOMImpl_GlueDriver::GlueFaces(aShape, Precision::Confusion(), Standard_True);
2397
2398   TopoDS_Shape aRes = GEOMImpl_IShapesOperations::CompsolidToCompound(aShape);
2399   aFunction->SetValue(aRes);
2400
2401   log.SetTouched(Label());
2402   return 1;
2403 }
2404
2405 //=======================================================================
2406 //function :  GEOMImpl_PipeDriver_Type_
2407 //purpose  :
2408 //=======================================================================
2409 Standard_EXPORT Handle_Standard_Type& GEOMImpl_PipeDriver_Type_()
2410 {
2411
2412   static Handle_Standard_Type aType1 = STANDARD_TYPE(TFunction_Driver);
2413   if ( aType1.IsNull()) aType1 = STANDARD_TYPE(TFunction_Driver);
2414   static Handle_Standard_Type aType2 = STANDARD_TYPE(MMgt_TShared);
2415   if ( aType2.IsNull()) aType2 = STANDARD_TYPE(MMgt_TShared);
2416   static Handle_Standard_Type aType3 = STANDARD_TYPE(Standard_Transient);
2417   if ( aType3.IsNull()) aType3 = STANDARD_TYPE(Standard_Transient);
2418
2419   static Handle_Standard_Transient _Ancestors[]= {aType1,aType2,aType3,NULL};
2420   static Handle_Standard_Type _aType = new Standard_Type("GEOMImpl_PipeDriver",
2421                                                          sizeof(GEOMImpl_PipeDriver),
2422                                                          1,
2423                                                          (Standard_Address)_Ancestors,
2424                                                          (Standard_Address)NULL);
2425
2426   return _aType;
2427 }
2428
2429 //=======================================================================
2430 //function : DownCast
2431 //purpose  :
2432 //=======================================================================
2433 const Handle(GEOMImpl_PipeDriver) Handle(GEOMImpl_PipeDriver)::DownCast(const Handle(Standard_Transient)& AnObject)
2434 {
2435   Handle(GEOMImpl_PipeDriver) _anOtherObject;
2436
2437   if (!AnObject.IsNull()) {
2438      if (AnObject->IsKind(STANDARD_TYPE(GEOMImpl_PipeDriver))) {
2439        _anOtherObject = Handle(GEOMImpl_PipeDriver)((Handle(GEOMImpl_PipeDriver)&)AnObject);
2440      }
2441   }
2442
2443   return _anOtherObject;
2444 }