/Users/enrico/Ricerca/develop/PINS/pins-mechatronix/LibSources/submodules/Clothoids/src/Clothoids/Circle.hxx Source File

Clothoids: /Users/enrico/Ricerca/develop/PINS/pins-mechatronix/LibSources/submodules/Clothoids/src/Clothoids/Circle.hxx Source File
Clothoids
Circle.hxx
1/*--------------------------------------------------------------------------*\
2 | |
3 | Copyright (C) 2017 |
4 | |
5 | , __ , __ |
6 | /|/ \ /|/ \ |
7 | | __/ _ ,_ | __/ _ ,_ |
8 | | \|/ / | | | | \|/ / | | | |
9 | |(__/|__/ |_/ \_/|/|(__/|__/ |_/ \_/|/ |
10 | /| /| |
11 | \| \| |
12 | |
13 | Enrico Bertolazzi |
14 | Dipartimento di Ingegneria Industriale |
15 | Universita` degli Studi di Trento |
16 | email: enrico.bertolazzi@unitn.it |
17 | |
18\*--------------------------------------------------------------------------*/
19
23
24namespace G2lib {
25
26 /*\
27 | ____ _ _ _
28 | / ___(_)_ __ ___| | ___ / \ _ __ ___
29 | | | | | '__/ __| |/ _ \ / _ \ | '__/ __|
30 | | |___| | | | (__| | __// ___ \| | | (__
31 | \____|_|_| \___|_|\___/_/ \_\_| \___|
32 \*/
33
37 class CircleArc : public BaseCurve {
38
39 friend class Biarc;
40
41 real_type m_x0{0};
42 real_type m_y0{0};
43 real_type m_theta0{0};
44 real_type m_c0{1};
45 real_type m_s0{0};
46 real_type m_k{0};
47
48 real_type m_L{0};
49
50 public:
51
52 #include "BaseCurve_using.hxx"
53
57 CircleArc() = delete;
58 CircleArc( string const & name ) : BaseCurve( name ) {};
59
60 void setup( GenericContainer const & gc ) override;
61
65 CircleArc( CircleArc const & s ) : BaseCurve( s.name() )
66 { this->copy(s); }
67
78 explicit
80 real_type x0,
81 real_type y0,
82 real_type theta0,
83 real_type k,
84 real_type L,
85 string const & name
86 )
87 : BaseCurve( name )
88 , m_x0(x0)
89 , m_y0(y0)
90 , m_theta0(theta0)
91 , m_c0(cos(theta0))
92 , m_s0(sin(theta0))
93 , m_k(k)
94 , m_L(L)
95 {}
96
101 explicit
102 CircleArc( LineSegment const & LS )
103 : BaseCurve( LS.name() )
104 , m_x0(LS.x_begin())
105 , m_y0(LS.y_begin())
106 , m_theta0(LS.m_theta0)
107 , m_c0(LS.m_c0)
108 , m_s0(LS.m_s0)
109 , m_k(0)
110 , m_L(LS.length())
111 {}
112
116 void
117 copy( CircleArc const & c ) {
118 m_x0 = c.m_x0;
119 m_y0 = c.m_y0;
120 m_theta0 = c.m_theta0;
121 m_c0 = c.m_c0;
122 m_s0 = c.m_s0;
123 m_k = c.m_k;
124 m_L = c.m_L;
125 }
126
130 explicit
131 CircleArc( BaseCurve const * pC );
132
133 CurveType type() const override { return CurveType::CIRCLE; }
134
138 CircleArc const &
140 { this->copy(s); return *this; }
141
151 void
153 real_type x0,
154 real_type y0,
155 real_type theta0,
156 real_type k,
157 real_type L
158 ) {
159 m_x0 = x0;
160 m_y0 = y0;
161 m_theta0 = theta0;
162 m_k = k;
163 m_L = L;
164 }
165
176 bool
177 build_G1(
178 real_type x0,
179 real_type y0,
180 real_type theta0,
181 real_type x1,
182 real_type y1
183 );
184
196 bool
197 build_3P(
198 real_type x0,
199 real_type y0,
200 real_type x1,
201 real_type y1,
202 real_type x2,
203 real_type y2
204 );
205
210 void build( LineSegment const & );
211 void build( CircleArc const & );
212 void build( Biarc const & );
213 void build( ClothoidCurve const & );
214 void build( PolyLine const & );
215 void build( BiarcList const & );
216 void build( ClothoidList const & );
217 void build( Dubins const & );
218 void build( Dubins3p const & );
219
220 // . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
221
233 bool
235 real_type & x0, real_type & y0,
236 real_type & x1, real_type & y1,
237 real_type & x2, real_type & y2
238 ) const;
239
253 bool
255 real_type offs,
256 real_type & x0, real_type & y0,
257 real_type & x1, real_type & y1,
258 real_type & x2, real_type & y2
259 ) const;
260
274 bool
276 real_type offs,
277 real_type & x0, real_type & y0,
278 real_type & x1, real_type & y1,
279 real_type & x2, real_type & y2
280 ) const {
281 return this->bbTriangle_ISO( -offs, x0, y0, x1, y1, x2, y2 );
282 }
283
293 bool
295 real_type p0[],
296 real_type p1[],
297 real_type p2[]
298 ) const {
299 return bbTriangle( p0[0], p0[1], p1[0], p1[1], p2[0], p2[1] );
300 }
301
312 bool
314 real_type offs,
315 real_type p0[],
316 real_type p1[],
317 real_type p2[]
318 ) const {
319 return bbTriangle_ISO( offs, p0[0], p0[1], p1[0], p1[1], p2[0], p2[1] );
320 }
321
332 bool
334 real_type offs,
335 real_type p0[],
336 real_type p1[],
337 real_type p2[]
338 ) const {
339 return bbTriangle_SAE( offs, p0[0], p0[1], p1[0], p1[1], p2[0], p2[1] );
340 }
341
351 bool
353 Triangle2D & t,
354 real_type ss0 = 0,
355 real_type ss1 = 0,
356 integer icurve = 0
357 ) const {
358 real_type p0[2], p1[2], p2[2];
359 bool ok = bbTriangle( p0, p1, p2 );
360 if ( ok ) t.build( p0, p1, p2, ss0, ss1, icurve );
361 return ok;
362 }
363
375 bool
377 real_type offs,
378 Triangle2D & t,
379 real_type ss0 = 0,
380 real_type ss1 = 0,
381 integer icurve = 0
382 ) const {
383 real_type p0[2], p1[2], p2[2];
384 bool ok = bbTriangle_ISO( offs, p0, p1, p2 );
385 if ( ok ) t.build( p0, p1, p2, ss0, ss1, icurve );
386 return ok;
387 }
388
400 bool
402 real_type offs,
403 Triangle2D & t,
404 real_type ss0 = 0,
405 real_type ss1 = 0,
406 integer icurve = 0
407 ) const {
408 return this->bbTriangle_ISO( -offs, t, ss0, ss1, icurve );
409 }
410
419 void
421 vector<Triangle2D> & tvec,
422 real_type max_angle = Utils::m_pi/18,
423 real_type max_size = 1e100,
424 integer icurve = 0
425 ) const override; // 10 degree
426
436 void
438 real_type offs,
439 vector<Triangle2D> & tvec,
440 real_type max_angle = Utils::m_pi/18,
441 real_type max_size = 1e100,
442 integer icurve = 0
443 ) const override; // 10 degree
444
454 void
456 real_type offs,
457 vector<Triangle2D> & tvec,
458 real_type max_angle = Utils::m_pi/18,
459 real_type max_size = 1e100,
460 integer icurve = 0
461 ) const override {
462 this->bb_triangles_ISO( -offs, tvec, max_angle, max_size, icurve );
463 }
464
465 // . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
466
467 void
468 bbox(
469 real_type & xmin,
470 real_type & ymin,
471 real_type & xmax,
472 real_type & ymax
473 ) const override;
474
475 void
476 bbox_ISO(
477 real_type offs,
478 real_type & xmin,
479 real_type & ymin,
480 real_type & xmax,
481 real_type & ymax
482 ) const override;
483
484 // . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
485
487 length() const override
488 { return m_L; }
489
491 length_ISO( real_type offs ) const override
492 { return m_L*(1+m_k*offs); }
493
494 real_type theta_begin() const override { return m_theta0; }
495 real_type kappa_begin() const override { return m_k; }
496 real_type kappa_end() const override { return m_k; }
497 real_type x_begin() const override { return m_x0; }
498 real_type y_begin() const override { return m_y0; }
499 real_type tx_begin() const override { return m_c0; }
500 real_type ty_begin() const override { return m_s0; }
501 real_type nx_begin_ISO() const override { return m_s0; }
502 real_type ny_begin_ISO() const override { return -m_c0; }
503
504 // . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
505
506 real_type theta ( real_type s ) const override { return m_theta0 + s*m_k; }
507 real_type theta_D ( real_type ) const override { return m_k; }
508 real_type theta_DD ( real_type ) const override { return 0; }
509 real_type theta_DDD( real_type ) const override { return 0; }
510
511 // . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
512
513 void
515 real_type s,
516 real_type & th,
518 real_type & x,
519 real_type & y
520 ) const override {
521 eval( s, x, y );
522 th = m_theta0 + s*m_k;
523 kappa = m_k;
524 }
525
526 // . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
527
528 real_type X( real_type s ) const override;
529 real_type Y( real_type s ) const override;
530
531 real_type X_D( real_type ) const override;
532 real_type Y_D( real_type ) const override;
533
534 real_type X_DD( real_type ) const override;
535 real_type Y_DD( real_type ) const override;
536
537 real_type X_DDD( real_type ) const override;
538 real_type Y_DDD( real_type ) const override;
539
540 // . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
541
542 void
543 eval(
544 real_type s,
545 real_type & x,
546 real_type & y
547 ) const override;
548
549 void
550 eval_D(
551 real_type,
552 real_type & x_D,
553 real_type & y_D
554 ) const override;
555
556 void
557 eval_DD(
558 real_type,
559 real_type & x_DD,
560 real_type & y_DD
561 ) const override;
562
563 void
564 eval_DDD(
565 real_type,
566 real_type & x_DDD,
567 real_type & y_DDD
568 ) const override;
569
570 /*\
571 | _____ _ _ _
572 | |_ _| __ _ _ __ __| | | \ | |
573 | | | / _` | '_ \ / _` | | \| |
574 | | | | (_| | | | | (_| | | |\ |
575 | |_| \__,_|_| |_|\__,_| |_| \_|
576 \*/
577
578 real_type tx( real_type s ) const override { return cos(theta(s)); }
579 real_type ty( real_type s ) const override { return sin(theta(s)); }
580
581 real_type tx_D( real_type s ) const override { return -sin(theta(s))*m_k; }
582 real_type ty_D( real_type s ) const override { return cos(theta(s))*m_k; }
583
584 real_type tx_DD( real_type s ) const override { return -cos(theta(s))*m_k*m_k; }
585 real_type ty_DD( real_type s ) const override { return -sin(theta(s))*m_k*m_k; }
586
587 real_type tx_DDD( real_type s ) const override { return sin(theta(s))*m_k*m_k*m_k; }
588 real_type ty_DDD( real_type s ) const override { return -cos(theta(s))*m_k*m_k*m_k; }
589
590 // . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
591
592 void tg( real_type s, real_type & tx, real_type & ty ) const override;
593 void tg_D( real_type s, real_type & tx_D, real_type & ty_D ) const override;
594 void tg_DD( real_type s, real_type & tx_DD, real_type & ty_DD ) const override;
595 void tg_DDD( real_type s, real_type & tx_DDD, real_type & ty_DDD ) const override;
596
597 /*\
598 | _ __
599 | | |_ _ __ __ _ _ __ ___ / _| ___ _ __ _ __ ___
600 | | __| '__/ _` | '_ \/ __| |_ / _ \| '__| '_ ` _ \
601 | | |_| | | (_| | | | \__ \ _| (_) | | | | | | | |
602 | \__|_| \__,_|_| |_|___/_| \___/|_| |_| |_| |_|
603 \*/
604
605 void
607 { m_x0 += tx; m_y0 += ty; }
608
609 void rotate( real_type angle, real_type cx, real_type cy ) override;
610 void reverse() override;
611
612 void
613 change_origin( real_type newx0, real_type newy0 ) override
614 { m_x0 = newx0; m_y0 = newy0; }
615
616 void scale( real_type s ) override;
617 void trim( real_type s_begin, real_type s_end ) override;
618
619 /*\
620 | _ _ ____ _ _
621 | ___| | ___ ___ ___ ___| |_| _ \ ___ (_)_ __ | |_
622 | / __| |/ _ \/ __|/ _ \/ __| __| |_) / _ \| | '_ \| __|
623 | | (__| | (_) \__ \ __/\__ \ |_| __/ (_) | | | | | |_
624 | \___|_|\___/|___/\___||___/\__|_| \___/|_|_| |_|\__|
625 |
626 \*/
627
628 integer
630 real_type qx,
631 real_type qy,
632 real_type & x,
633 real_type & y,
634 real_type & s,
635 real_type & t,
636 real_type & dst
637 ) const override;
638
639 integer
641 real_type qx,
642 real_type qy,
643 real_type offs,
644 real_type & x,
645 real_type & y,
646 real_type & s,
647 real_type & t,
648 real_type & dst
649 ) const override;
650
651 string info() const;
652
653 void
654 info( ostream_type & stream ) const override
655 { stream << this->info(); }
656
657 friend
659 operator << ( ostream_type & stream, CircleArc const & bi );
660
661 // . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
662 // . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
663 // . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
664
665 /*\
666 | _ _ _ _
667 | ___ ___ | | (_)___(_) ___ _ __
668 | / __/ _ \| | | / __| |/ _ \| '_ \
669 | | (_| (_) | | | \__ \ | (_) | | | |
670 | \___\___/|_|_|_|___/_|\___/|_| |_|
671 \*/
672
676 bool
677 collision( CircleArc const & ) const;
678
686 bool
688 real_type offs,
689 CircleArc const & C,
690 real_type offs_obj
691 ) const;
692
693 bool
694 collision( BaseCurve const * pC ) const override;
695
696 bool
698 real_type offs,
699 BaseCurve const * pC,
700 real_type offs_C
701 ) const override;
702
703 /*\
704 | _ _ _
705 | (_)_ __ | |_ ___ _ __ ___ ___ ___| |_
706 | | | '_ \| __/ _ \ '__/ __|/ _ \/ __| __|
707 | | | | | | || __/ | \__ \ __/ (__| |_
708 | |_|_| |_|\__\___|_| |___/\___|\___|\__|
709 \*/
710
717 void
718 intersect(
719 CircleArc const & obj,
720 IntersectList & ilist
721 ) const;
722
731 void
733 real_type offs,
734 CircleArc const & C,
735 real_type offs_obj,
736 IntersectList & ilist
737 ) const;
738
739 void
740 intersect(
741 BaseCurve const * pC,
742 IntersectList & ilist
743 ) const override;
744
745 void
747 real_type offs,
748 BaseCurve const * pC,
749 real_type offs_LS,
750 IntersectList & ilist
751 ) const override;
752
757 real_type sin_theta0() const { return sin(m_theta0); }
758
763 real_type cos_theta0() const { return cos(m_theta0); }
764
768 real_type curvature() const { return m_k; }
769
774 real_type len_tolerance( real_type tol ) const;
775
779 real_type delta_theta() const { return std::abs(m_L*m_k); }
780
786
795 theta_min_max( real_type & thMin, real_type & thMax ) const;
796
804 void
806
813 void
814 center( real_type & cx, real_type & cy ) const;
815
819 real_type ray() const { return 1/std::abs(m_k); }
820
821 /*\
822 | _ _ _ _ ____ ____ ____
823 | | \ | | | | | _ \| __ ) ___|
824 | | \| | | | | |_) | _ \___ \
825 | | |\ | |_| | _ <| |_) |__) |
826 | |_| \_|\___/|_| \_\____/____/
827 \*/
828
835 void
836 paramNURBS( integer & n_knots, integer & n_pnts ) const;
837
844 void
845 toNURBS( real_type knots[], real_type Poly[][3] ) const;
846
847 friend class ClothoidCurve;
848
849#ifdef CLOTHOIDS_BACK_COMPATIBILITY
850#include "Circle_compatibility.hxx"
851#endif
852
853 };
854
855}
856
Definition BaseCurve.hxx:192
real_type kappa(real_type s) const
Definition BaseCurve.hxx:560
Definition Biarc.hxx:39
Definition BiarcList.hxx:42
Definition Circle.hxx:37
real_type tx_DDD(real_type s) const override
Definition Circle.hxx:587
void eval(real_type s, real_type &x, real_type &y) const override
Definition Circle.cc:294
void bbox(real_type &xmin, real_type &ymin, real_type &xmax, real_type &ymax) const override
Definition Circle.cc:557
void bb_triangles_SAE(real_type offs, vector< Triangle2D > &tvec, real_type max_angle=Utils::m_pi/18, real_type max_size=1e100, integer icurve=0) const override
Definition Circle.hxx:455
real_type ny_begin_ISO() const override
Definition Circle.hxx:502
friend ostream_type & operator<<(ostream_type &stream, CircleArc const &bi)
Definition Circle.cc:995
bool bbTriangle(real_type p0[], real_type p1[], real_type p2[]) const
Definition Circle.hxx:294
bool bbTriangle(Triangle2D &t, real_type ss0=0, real_type ss1=0, integer icurve=0) const
Definition Circle.hxx:352
bool bbTriangle_SAE(real_type offs, Triangle2D &t, real_type ss0=0, real_type ss1=0, integer icurve=0) const
Definition Circle.hxx:401
integer closest_point_ISO(real_type qx, real_type qy, real_type &x, real_type &y, real_type &s, real_type &t, real_type &dst) const override
Definition Circle.cc:808
real_type kappa_end() const override
Definition Circle.hxx:496
void tg_DD(real_type s, real_type &tx_DD, real_type &ty_DD) const override
Definition Circle.cc:268
bool bbTriangle_SAE(real_type offs, real_type &x0, real_type &y0, real_type &x1, real_type &y1, real_type &x2, real_type &y2) const
Definition Circle.hxx:275
void tg_DDD(real_type s, real_type &tx_DDD, real_type &ty_DDD) const override
Definition Circle.cc:280
real_type delta_theta() const
Definition Circle.hxx:779
void bbox_ISO(real_type offs, real_type &xmin, real_type &ymin, real_type &xmax, real_type &ymax) const override
Definition Circle.cc:579
bool bbTriangle_ISO(real_type offs, real_type &x0, real_type &y0, real_type &x1, real_type &y1, real_type &x2, real_type &y2) const
Definition Circle.cc:456
void trim(real_type s_begin, real_type s_end) override
Definition Circle.cc:349
real_type theta_min_max(real_type &thMin, real_type &thMax) const
Definition Circle.cc:182
bool bbTriangle_ISO(real_type offs, Triangle2D &t, real_type ss0=0, real_type ss1=0, integer icurve=0) const
Definition Circle.hxx:376
real_type ty(real_type s) const override
Definition Circle.hxx:579
real_type Y_DDD(real_type) const override
Definition Circle.cc:231
real_type X_D(real_type) const override
Definition Circle.cc:198
void tg(real_type s, real_type &tx, real_type &ty) const override
Definition Circle.cc:246
CircleArc(real_type x0, real_type y0, real_type theta0, real_type k, real_type L, string const &name)
Definition Circle.hxx:79
real_type ty_DD(real_type s) const override
Definition Circle.hxx:585
CircleArc const & operator=(CircleArc const &s)
Definition Circle.hxx:139
CircleArc()=delete
bool collision_ISO(real_type offs, CircleArc const &C, real_type offs_obj) const
Definition Circle.cc:647
void info(ostream_type &stream) const override
Definition Circle.hxx:654
void bb_triangles_ISO(real_type offs, vector< Triangle2D > &tvec, real_type max_angle=Utils::m_pi/18, real_type max_size=1e100, integer icurve=0) const override
Definition Circle.cc:518
real_type tx(real_type s) const override
Definition Circle.hxx:578
void intersect_ISO(real_type offs, CircleArc const &C, real_type offs_obj, IntersectList &ilist) const
Definition Circle.cc:727
real_type Y_DD(real_type) const override
Definition Circle.cc:226
bool build_G1(real_type x0, real_type y0, real_type theta0, real_type x1, real_type y1)
Definition Circle.cc:119
real_type X(real_type s) const override
Definition Circle.cc:192
void rotate(real_type angle, real_type cx, real_type cy) override
Definition Circle.cc:368
real_type theta(real_type s) const override
Definition Circle.hxx:506
real_type Y(real_type s) const override
Definition Circle.cc:215
real_type X_DDD(real_type) const override
Definition Circle.cc:208
real_type curvature() const
Definition Circle.hxx:768
real_type y_begin() const override
Definition Circle.hxx:498
void translate(real_type tx, real_type ty) override
translate curve by
Definition Circle.hxx:606
real_type theta_total_variation() const
Definition Circle.hxx:785
void copy(CircleArc const &c)
Definition Circle.hxx:117
void paramNURBS(integer &n_knots, integer &n_pnts) const
Definition Circle.cc:898
void evaluate(real_type s, real_type &th, real_type &kappa, real_type &x, real_type &y) const override
Definition Circle.hxx:514
CircleArc(CircleArc const &s)
Definition Circle.hxx:65
CircleArc(LineSegment const &LS)
Definition Circle.hxx:102
bool bbTriangle(real_type &x0, real_type &y0, real_type &x1, real_type &y1, real_type &x2, real_type &y2) const
Definition Circle.cc:431
void center(real_type &cx, real_type &cy) const
Definition Circle.cc:409
real_type theta_D(real_type) const override
Definition Circle.hxx:507
real_type nx_begin_ISO() const override
Definition Circle.hxx:501
bool bbTriangle_ISO(real_type offs, real_type p0[], real_type p1[], real_type p2[]) const
Definition Circle.hxx:313
void eval_D(real_type, real_type &x_D, real_type &y_D) const override
Definition Circle.cc:309
real_type len_tolerance(real_type tol) const
Definition Circle.cc:969
real_type length() const override
Definition Circle.hxx:487
void bb_triangles(vector< Triangle2D > &tvec, real_type max_angle=Utils::m_pi/18, real_type max_size=1e100, integer icurve=0) const override
Definition Circle.cc:481
real_type ty_begin() const override
Definition Circle.hxx:500
void build(real_type x0, real_type y0, real_type theta0, real_type k, real_type L)
Definition Circle.hxx:152
real_type ray() const
Definition Circle.hxx:819
void reverse() override
Definition Circle.cc:393
void eval_DDD(real_type, real_type &x_DDD, real_type &y_DDD) const override
Definition Circle.cc:335
real_type tx_begin() const override
Definition Circle.hxx:499
real_type ty_D(real_type s) const override
Definition Circle.hxx:582
bool bbTriangle_SAE(real_type offs, real_type p0[], real_type p1[], real_type p2[]) const
Definition Circle.hxx:333
real_type Y_D(real_type) const override
Definition Circle.cc:221
bool collision(CircleArc const &) const
Definition Circle.cc:610
real_type sin_theta0() const
Definition Circle.hxx:757
real_type theta_DD(real_type) const override
Definition Circle.hxx:508
CurveType type() const override
Definition Circle.hxx:133
real_type x_begin() const override
Definition Circle.hxx:497
void change_origin(real_type newx0, real_type newy0) override
Definition Circle.hxx:613
void tg_D(real_type s, real_type &tx_D, real_type &ty_D) const override
Definition Circle.cc:257
real_type theta_DDD(real_type) const override
Definition Circle.hxx:509
void intersect(CircleArc const &obj, IntersectList &ilist) const
Definition Circle.cc:703
real_type kappa_begin() const override
Definition Circle.hxx:495
void eval_DD(real_type, real_type &x_DD, real_type &y_DD) const override
Definition Circle.cc:322
real_type ty_DDD(real_type s) const override
Definition Circle.hxx:588
bool build_3P(real_type x0, real_type y0, real_type x1, real_type y1, real_type x2, real_type y2)
Definition Circle.cc:146
void toNURBS(real_type knots[], real_type Poly[][3]) const
Definition Circle.cc:912
real_type tx_DD(real_type s) const override
Definition Circle.hxx:584
real_type cos_theta0() const
Definition Circle.hxx:763
real_type X_DD(real_type) const override
Definition Circle.cc:203
real_type tx_D(real_type s) const override
Definition Circle.hxx:581
real_type length_ISO(real_type offs) const override
Definition Circle.hxx:491
real_type theta_begin() const override
Definition Circle.hxx:494
void change_curvilinear_origin(real_type s0, real_type newL)
Definition Circle.cc:419
void scale(real_type s) override
Definition Circle.cc:385
Definition Clothoid.hxx:48
Definition ClothoidList.hxx:861
Definition Dubins3p.hxx:78
Definition Dubins.hxx:74
Definition Line.hxx:37
Class to manage a collection of straight segment.
Definition PolyLine.hxx:42
Definition Triangle2D.hxx:37
Definition BBox.cc:42
GC_namespace::GenericContainer GenericContainer
Generic container object.
Definition Clothoids.hh:84
std::basic_ostream< char > ostream_type
output streaming
Definition Clothoids.hh:78
std::vector< Ipair > IntersectList
Vector of pair of two real number.
Definition BaseCurve.hxx:36
enum class CurveType :integer { LINE, POLYLINE, CIRCLE, BIARC, BIARC_LIST, CLOTHOID, CLOTHOID_LIST, DUBINS, DUBINS3P } CurveType
Definition Clothoids.hh:89
double real_type
real type number
Definition Clothoids.hh:79
int integer
integer type number
Definition Clothoids.hh:80