G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsCurveLoop.h
Go to the documentation of this file.
1 
15 #pragma once
16 
17 #include <iostream>
18 
19 # include <gsCore/gsCurve.h>
21 
22 namespace gismo
23 {
24 
35 template<class T>
37 {
38 
39 public:
40  typedef typename std::vector< gsCurve<T> * >::iterator iterator;
41  typedef typename std::vector< gsCurve<T> * >::const_iterator const_iterator;
42 
44  typedef memory::shared_ptr< gsCurveLoop > Ptr;
45 
47  typedef memory::unique_ptr< gsCurveLoop > uPtr;
48 
49 public:
52 
54  explicit gsCurveLoop( gsCurve<T> * cc ) { this->insertCurve(cc); }
55 
57  explicit gsCurveLoop(std::vector< gsCurve<T> *> const & curves ) : m_curves(curves) { }
58 
63  gsCurveLoop(const std::vector<gsVector3d<T> *> points3D, const std::vector<bool> isConvex, T margin, gsVector3d<T> *outNormal);
64 
68  gsCurveLoop(const std::vector<T> angles3D, const std::vector<T> lengths3D, const std::vector<bool> isConvex, T margin, const bool unitSquare=false);
69 
70  gsCurveLoop(const gsCurveLoop & other)
71  : m_curves( other.m_curves.size() )
72  {
73  cloneAll( other.m_curves.begin(), other.m_curves.end(),
74  this->m_curves.begin() );
75  }
76 
77  gsCurveLoop& operator=(const gsCurveLoop & other)
78  {
79  freeAll( m_curves );
80  m_curves.resize( other.m_curves.size() );
81  cloneAll( other.m_curves.begin(), other.m_curves.end(),
82  this->m_curves.begin() );
83  return *this;
84  }
85 
86  ~gsCurveLoop()
87  {
88  freeAll( m_curves );
89  }
90 
91  //GISMO_CLONE_FUNCTION(gsCurveLoop)
92  uPtr clone() const
93  {
94  return uPtr(new gsCurveLoop(*this));
95  }
96 
97 public:
98 
99  bool isInterior ( gsVector<T> const & p, const T& tol);
100 
105  bool isOn(gsMatrix<T> const &u,T & paramResult, T tol = 1e-3);
106 
109  bool is_ccw();
110 
113  void reverse();
114 
115  void translate(gsVector<T> const & v);
116 
117  void insertCurve( gsCurve<T> * c )
118  {
119  m_curves.push_back( c ) ;
120  }
121 
123  typename gsCurve<T>::uPtr singleCurve() const;
124 
125  int size() const { return m_curves.size(); }
126 
128  std::ostream &print(std::ostream &os) const
129  { os << "gsCurveLoop with "<< size() <<" curves :\n";
130  for (typename std::vector< gsCurve<T> *>::const_iterator it= m_curves.begin();
131  it!= m_curves.end(); it++ )
132  {
133  os <<"("<<*it<<") "<< **it ;
134  };
135  return os;
136  }
137 
140  gsMatrix<T> normal( int const & c, gsMatrix<T> const & u );
141 
145  uPtr split(int startIndex, int endIndex,
146  gsCurve<T> * newCurveThisFace, gsCurve<T> * newCurveNewFace);
147 
148 
151  // to do: replace return type with std::multimap<int,T>
152  std::vector<T> lineIntersections(int const & direction , T const & abscissa);
153 
154  std::vector< gsCurve<T> *> & curves() { return m_curves; }
155 
156  gsCurve<T> & curve(int i)
157  {
158  GISMO_ASSERT( i<int(m_curves.size()), "Curve does not exist" );
159  return *m_curves[i];
160  }
161 
162  const gsCurve<T> & curve(int i) const
163  {
164  GISMO_ASSERT( i<int(m_curves.size()), "Curve does not exist" );
165  return *m_curves[i];
166  }
167 
168  std::vector< gsCurve<T> *> const & curves() const
169  { return m_curves; }
170 
172  int numCurves() const { return m_curves.size(); }
173 
175  gsMatrix<T> sample(int npoints = 50, int numEndPoints=2) const;
176 
177  gsMatrix<T> getBoundingBox();
178 
183  static bool approximatingPolygon(const std::vector<T> &signedAngles, const std::vector<T> &lengths, T margin, gsMatrix<T> &result);
184 
186  std::vector<T> domainSizes() const;
187 
193  static void adjustPolygonToUnitSquare(gsMatrix<T> &corners, T const margin);
194 
198  gsMatrix<T> splitCurve(size_t curveId, T lengthRatio=.5);
199 
201  gsVector3d<T> initFrom3DPlaneFit(const std::vector<gsVector3d<T> *> points3D, T margin);
202 
205  bool initFrom3DByAngles(const std::vector<gsVector3d<T> *> points3D, const std::vector<bool> isConvex, T margin);
206 
208  void initFromIsConvex(const std::vector<bool> isConvex, T margin);
209 
211  void flip1(T minu = 0, T maxu = 1);
212 
213 private:
214  bool initFrom3DByAngles(const std::vector<T>& angles3D, const std::vector<T>& lengths3D, const std::vector<bool>& isConvex, T margin, bool unitsquare=false);
215 
216  bool parameterOf(gsMatrix<T> const &u, int i, T & result, T tol = 1e-5);
217 
218  void removeCurves(iterator begin, iterator end);
219 
220  void removeCurve(iterator it)
221  { removeCurves(it, it+1); }
222 
223  // Data members
224 
225  std::vector< gsCurve<T> *> m_curves; // CCW
226 
227 }; // class gsCurveLoop
228 
229 template<class T>
230 bool gsCurveLoop<T>::isInterior ( gsVector<T> const &, const T&)
231 {
233 }
234 
236 template<class T>
237 std::ostream &operator<<(std::ostream &os, const gsCurveLoop<T>& b)
238 {return b.print(os); }
239 
240 
241 } // namespace gismo
242 
243 #ifndef GISMO_BUILD_LIB
244 #include GISMO_HPP_HEADER(gsCurveLoop.hpp)
245 #endif
A fixed-size, statically allocated 3D vector.
Definition: gsVector.h:218
gsMatrix< T > normal(int const &c, gsMatrix< T > const &u)
Definition: gsCurveLoop.hpp:144
#define GISMO_NO_IMPLEMENTATION
Definition: gsDebug.h:129
Abstract base class representing a curve.
Definition: gsCurve.h:30
int numCurves() const
get Number of curves
Definition: gsCurveLoop.h:172
std::ostream & print(std::ostream &os) const
Prints the object as a string.
Definition: gsCurveLoop.h:128
#define GISMO_ASSERT(cond, message)
Definition: gsDebug.h:89
void cloneAll(It start, It end, ItOut out)
Clones all pointers in the range [start end) and stores new raw pointers in iterator out...
Definition: gsMemory.h:295
gsCurveLoop(std::vector< gsCurve< T > * > const &curves)
Constructor by a list of curves forming a loop.
Definition: gsCurveLoop.h:57
uPtr split(int startIndex, int endIndex, gsCurve< T > *newCurveThisFace, gsCurve< T > *newCurveNewFace)
Definition: gsCurveLoop.hpp:160
static bool approximatingPolygon(const std::vector< T > &signedAngles, const std::vector< T > &lengths, T margin, gsMatrix< T > &result)
Definition: gsCurveLoop.hpp:280
gsMatrix< T > sample(int npoints=50, int numEndPoints=2) const
Sample npoints uniformly distributed (in parameter domain) points on the loop.
Definition: gsCurveLoop.hpp:211
gsVector3d< T > initFrom3DPlaneFit(const std::vector< gsVector3d< T > * > points3D, T margin)
Initialize a curve loop from some 3d vertices by projecting them onto the plane of best fit and const...
Definition: gsCurveLoop.hpp:407
memory::unique_ptr< gsCurveLoop > uPtr
Unique pointer for gsCurveLoop.
Definition: gsCurveLoop.h:47
GISMO_DEPRECATED index_t direction(index_t s)
Returns the parametric direction that corresponds to side s.
Definition: gsBoundary.h:1048
void flip1(T minu=0, T maxu=1)
flip a curve loop in the u direction
Definition: gsCurveLoop.hpp:608
Provides forward declarations of types and structs.
void freeAll(It begin, It end)
Frees all pointers in the range [begin end)
Definition: gsMemory.h:312
static void adjustPolygonToUnitSquare(gsMatrix< T > &corners, T const margin)
Definition: gsCurveLoop.hpp:384
gsCurveLoop(gsCurve< T > *cc)
Constructor by one curve.
Definition: gsCurveLoop.h:54
gsMatrix< T > splitCurve(size_t curveId, T lengthRatio=.5)
Definition: gsCurveLoop.hpp:624
bool initFrom3DByAngles(const std::vector< gsVector3d< T > * > points3D, const std::vector< bool > isConvex, T margin)
Definition: gsCurveLoop.hpp:483
gsCurveLoop()
Default empty constructor.
Definition: gsCurveLoop.h:51
memory::shared_ptr< gsCurveLoop > Ptr
Shared pointer for gsCurveLoop.
Definition: gsCurveLoop.h:44
void initFromIsConvex(const std::vector< bool > isConvex, T margin)
Initialize a curve loop from a list of bools indicating whether the corners are convex or not...
Definition: gsCurveLoop.hpp:514
gsCurve< T >::uPtr singleCurve() const
Return the curve-loop as a single new B-Spline curve.
Definition: gsCurveLoop.hpp:127
Provides declaration of Curve abstract interface.
std::vector< T > domainSizes() const
return a vector containing whose nth entry is the length of the domin of the nth curve.
Definition: gsCurveLoop.hpp:369
A closed loop given by a collection of curves.
Definition: gsCurveLoop.h:36
std::vector< T > lineIntersections(int const &direction, T const &abscissa)
Definition: gsCurveLoop.hpp:192