G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsVector.h
Go to the documentation of this file.
1
14# pragma once
15
16// Assumes that Eigen library has been already included
17
18namespace gismo
19{
20
34template<class T, int _Rows, int _Options>
35class gsVector : public gsMatrix<T, _Rows, 1, _Options>
36//class gsVector : public gsMatrix<T, _Rows, (_Rows!=-1 ? 1 : -1), _Options>
37{
38public:
40 //typedef gsMatrix<T,_Rows,(_Rows!=-1 ? 1 : -1), _Options> gsBase;
41
42 // Base is the single-column dense matrix class of Eigen
43 typedef typename gsBase::Base Base;
44
45 // Self type
47
48 typedef typename gsEigen::aligned_allocator<Self> aalloc;
49
50 // The type of the coefficients of the matrix
51 typedef T Scalar_t;
52
53 // Type pointing to a block of the vector
54 typedef typename gsEigen::Block<Base> Block;
55
56 // Type pointing to a block view of the vector
58
60 typedef memory::shared_ptr< gsVector > Ptr;
61
63 typedef memory::unique_ptr< gsVector > uPtr;
64
65 // Type for copying a vector as a permutation matrix
66 typedef gsEigen::PermutationMatrix<_Rows,Base::SizeAtCompileTime,index_t> Permutation;
67
68 // Type for treating a vector as a permutation matrix
69 typedef gsEigen::PermutationWrapper<Base> PermutationWrap;
70
71 typedef gsEigen::Ref<Base> Ref;
72
73 typedef const gsEigen::Ref<const Base> ConstRef;
74
75 // Type for a vector of dimension one less
76 typedef gsMatrix< T, ChangeDim<_Rows, -1>::D, ColMajor> Projection_t;
77
78public:
79
80 typedef T * iterator;
81
82 typedef const T * const_iterator;
83
84 T * begin()
85 { return this->data(); }
86
87 const T * begin() const
88 { return this->data(); }
89
90 T * end()
91 { return this->data() + this->size(); }
92
93 const T * end() const
94 { return this->data() + this->size(); }
95
96public:
97
98 gsVector() ;
99
100 gsVector(const Base& a) ;
101
102 // implicitly deleted in C++11
103 //gsVector(const gsVector& a) : gsBase(a) { }
104
105 explicit gsVector(index_t dimension) ;
106
107 // To enable pybind11 in gsPointLoads
108 explicit gsVector(index_t _rows, index_t _cols)
109 : gsVector(_rows)
110 {
111 GISMO_ASSERT(1==_cols,"Columns should be 1");
112 }
113
114 inline operator Ref () { return Ref(*this); }
115
116 inline operator const ConstRef () { return ConstRef(*this); }
117
118 void clear() { this->resize(0); }
119
120 // This constructor allows constructing a gsVector from Eigen expressions
121 template<typename OtherDerived>
122 gsVector(const gsEigen::EigenBase<OtherDerived>& other) : gsBase(other) { }
123
124 // This constructor allows constructing a gsVector from Eigen expressions
125 template<typename OtherDerived>
126 gsVector(const gsEigen::MatrixBase<OtherDerived>& other) : gsBase(other) { }
127
128 // This constructor allows constructing a gsVector from Eigen expressions
129 template<typename OtherDerived>
130 gsVector(const gsEigen::ReturnByValue<OtherDerived>& other) : gsBase(other) { }
131
132 static gsVector<T,2> vec( T x, T y)
133 {
134 return typename gsVector<T,2>::Base(x, y);
135 }
136
137 static gsVector<T,3> vec( T x, T y, T z)
138 {
139 return typename gsVector<T,3>::Base(x, y, z);
140 }
141
142 template<class iterator> void
143 assign(iterator from, iterator to)
144 {
145 this->resize(std::distance(from,to));
146 T * a = this->data();
147 for(iterator it = from; it!=to; ++it) *(a++) = *it;
148 }
149
150/*
151 // Using the assignment operators of Eigen
152 // Note: using Base::operator=; is ambiguous in MSVC
153#ifdef _MSC_VER
154 template <class EigenExpr>
155 gsVector& operator= (const EigenExpr & other)
156 {
157 this->Base::operator=(other);
158 return *this;
159 }
160#else
161 using Base::operator=;
162#endif
163*/
164#if !EIGEN_HAS_RVALUE_REFERENCES
165 gsVector & operator=(typename gsEigen::internal::conditional<
166 -1==_Rows,gsVector, const gsVector &>::type other)
167 {
168 if (-1==_Rows)
169 this->swap(other);
170 else
171 this->Base::operator=(other);
172 return *this;
173 }
174#endif
175
177 inline T at(index_t i) const { return *(this->data()+i);}
178
180 inline T & at(index_t i) { return *(this->data()+i);}
181
183 inline T last() const { return *(this->data()+this->size()-1);}
184
186 inline T & last() { return *(this->data()+this->size()-1);}
187
190 {
191 return BlockView(*this, rowSizes);
192 }
193
194 PermutationWrap asPermutation() const { return PermutationWrap(*this);}
195
198 void removeElement(const index_t i )
199 {
200 GISMO_ASSERT( i < this->size(), "Invalid vector element." );
201 const T * ce = this->data() + this->size();
202 for ( T * c = this->data()+i+1; c!= ce; ++c ) *(c-1) = *c;
203 this->conservativeResize(this->size()-1,gsEigen::NoChange);
204 }
205
206}; // class gsVector
207
208
209
210
217template<class T>
218class gsVector3d : public gsEigen::Matrix<T,3,1>
219{
220public:
221 typedef T scalar_t;
222 typedef gsEigen::Matrix<T,3,1> Base ;
223
225 typedef memory::shared_ptr< gsVector3d > Ptr;
226
227public:
228
229 gsVector3d();
230
231 gsVector3d(scalar_t x, scalar_t y, scalar_t z = 0 );
232
233 // implcitly declared deleted in C++11
234 //gsVector3d(const gsVector3d& a) : Base(a) { }
235
236 gsVector3d(const Base& a) ;
237
239 template<typename OtherDerived>
240 gsVector3d(const gsEigen::MatrixBase<OtherDerived>& other) : Base(other) { }
241
242 T angle(const gsVector3d<T> & other)
243 {
244 return math::acos(this->normalized().dot(other.normalized()));
245 }
246
247 /*
248 // Using the assignment operators of Eigen
249 // Note: using Base::operator=; is ambiguous in MSVC
250#ifdef _MSC_VER
251 template <class EigenExpr>
252 gsVector3d& operator= (const EigenExpr & other)
253 {
254 this->Base::operator=(other);
255 return *this;
256 }
257#else
258 using Base::operator=;
259#endif
260 */
261
262 // Copy constructor
263 gsVector3d(const gsVector3d & other) : Base(other) { }
264
265 // implicitly deleted in C++11
266 gsVector3d & operator=(const gsVector3d & other)
267 {
268 this->Base::operator=(other);
269 return *this;
270 }
271
272 inline T at (index_t i) const { return (*this)(i,0); }
273 inline T & at (index_t i) { return (*this)(i,0); }
274 inline T x () const { return (*this)(0); }
275 inline T & x () { return (*this)(0); }
276 inline T y () const { return (*this)(1); }
277 inline T & y () { return (*this)(1); }
278 inline T z () const { return (*this)(2); }
279 inline T & z () { return (*this)(2); }
280
281}; // class gsVector3d
282
283
284template<class T, int _Rows, int _Options> inline
285gsVector<T,_Rows,_Options>::gsVector() : gsBase() { }
286
287template<class T, int _Rows, int _Options> inline
288gsVector<T,_Rows,_Options>::gsVector(const Base& a): gsBase(a) { }
289
290template<class T, int _Rows, int _Options> inline
291gsVector<T,_Rows,_Options>::gsVector(index_t dimension): gsBase(dimension,1) { }
292
293template<class T> inline
294gsVector3d<T>::gsVector3d() : Base() { }
295
296template<class T> inline
297gsVector3d<T>::gsVector3d(scalar_t x, scalar_t y,scalar_t z )
298 {
299 (*this)(0,0)=x;
300 (*this)(1,0)=y;
301 (*this)(2,0)=z;
302 }
303
304template<class T> inline
305gsVector3d<T>::gsVector3d(const Base& a): Base(a) { }
306
307// template<class T>
308// template<typename OtherDerived> inline
309// gsVector3d<T>::gsVector3d(const gsEigen::MatrixBase<OtherDerived>& other) : Base(other) { }
310
311
312// template<class T> inline
313// inline T gsVector3d<T>::x () const { return (*this)(0); }
314// template<class T> inline
315// inline T & gsVector3d<T>::x () { return (*this)(0); }
316// template<class T> inline
317// inline T gsVector3d<T>::y () const { return (*this)(1); }
318// template<class T> inline
319// inline T & gsVector3d<T>::y () { return (*this)(1); }
320// template<class T> inline
321// inline T gsVector3d<T>::z () const { return (*this)(2); }
322// template<class T> inline
323// inline T & gsVector3d<T>::z () { return (*this)(2); }
324
325#ifdef GISMO_WITH_PYBIND11
326
330 template<typename T>
331 void pybind11_init_gsVector(pybind11::module &m, const std::string & typestr)
332 {
333 using Class = gsVector<T>;
334 std::string pyclass_name = std::string("gsVector") + typestr;
335 pybind11::class_<Class>(m, pyclass_name.c_str(), pybind11::buffer_protocol(), pybind11::dynamic_attr())
336 // Constructors
337 .def(pybind11::init<>())
338 .def(pybind11::init<index_t, index_t>())
339 // Member functions
340 .def("size", &Class::size)
341 .def("rows", &Class::rows)
342 // .def("transpose", &Class::transpose)
343 ;
344 }
345
346#endif // GISMO_WITH_PYBIND11
347
348
349
350} // namespace gismo
Represents a block-view of the given matrix.
Definition gsMatrixBlockView.h:32
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
A fixed-size, statically allocated 3D vector.
Definition gsVector.h:219
gsVector3d(const gsEigen::MatrixBase< OtherDerived > &other)
This constructor allows constructing a gsVector3d from Eigen expressions.
Definition gsVector.h:240
memory::shared_ptr< gsVector3d > Ptr
Shared pointer for gsVector3d.
Definition gsVector.h:225
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
T at(index_t i) const
Returns the i-th element of the vector.
Definition gsVector.h:177
memory::shared_ptr< gsVector > Ptr
Shared pointer for gsVector.
Definition gsVector.h:60
T & at(index_t i)
Returns the i-th element of the vector.
Definition gsVector.h:180
T last() const
Returns the last (bottom) element of the vector.
Definition gsVector.h:183
void removeElement(const index_t i)
Definition gsVector.h:198
memory::unique_ptr< gsVector > uPtr
Unique pointer for gsVector.
Definition gsVector.h:63
BlockView blockView(const gsVector< index_t > &rowSizes)
Return a row-block view of the vector with rowSizes.
Definition gsVector.h:189
T & last()
Returns the last (bottom) element of the vector.
Definition gsVector.h:186
#define index_t
Definition gsConfig.h:32
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
The G+Smo namespace, containing all definitions for the library.