ITK  4.3.0
Insight Segmentation and Registration Toolkit
itkVector.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef __itkVector_h
19 #define __itkVector_h
20 
21 #include "itkFixedArray.h"
22 
23 #include "vnl/vnl_vector_ref.h" // GetVnlVector method return
24 
25 namespace itk
26 {
61 template< class T, unsigned int NVectorDimension = 3 >
62 class Vector:public FixedArray< T, NVectorDimension >
63 {
64 public:
66  typedef Vector Self;
68 
71  typedef T ValueType;
73 
75  itkStaticConstMacro(Dimension, unsigned int, NVectorDimension);
76 
78  typedef Self VectorType;
79 
81  typedef T ComponentType;
82 
85 
87  static unsigned int GetVectorDimension() { return NVectorDimension; }
88 
90  void SetVnlVector(const vnl_vector< T > &);
91 
93  vnl_vector_ref< T > GetVnlVector(void);
94 
96  vnl_vector< T > GetVnlVector(void) const;
97 
100  itkLegacyMacro(void Set_vnl_vector(const vnl_vector< T > &));
101 
104  itkLegacyMacro(vnl_vector_ref< T > Get_vnl_vector(void));
105 
108  itkLegacyMacro(vnl_vector< T > Get_vnl_vector(void) const);
109 
112  Vector(const ValueType & r);
113 
115  template< class TVectorValueType >
119 
121  template< class TVectorValueType >
123  {
125  return *this;
126  }
128 
129  Vector & operator=(const ValueType r[NVectorDimension]);
130 
132  template< class Tt >
133  inline const Self & operator*=(const Tt & value)
134  {
135  for ( unsigned int i = 0; i < NVectorDimension; i++ )
136  {
137  ( *this )[i] = static_cast< ValueType >( ( *this )[i] * value );
138  }
139  return *this;
140  }
142 
144  template< class Tt >
145  inline const Self & operator/=(const Tt & value)
146  {
147  for ( unsigned int i = 0; i < NVectorDimension; i++ )
148  {
149  ( *this )[i] = static_cast< ValueType >( ( *this )[i] / value );
150  }
151  return *this;
152  }
154 
156  const Self & operator+=(const Self & vec);
157 
159  const Self & operator-=(const Self & vec);
160 
163  Self operator-() const;
164 
166  Self operator+(const Self & vec) const;
167 
169  Self operator-(const Self & vec) const;
170 
173  ValueType operator *(const Self & vec) const;
174 
177  inline Self operator*(const ValueType & value) const
178  {
179  Self result;
180 
181  for ( unsigned int i = 0; i < NVectorDimension; i++ )
182  {
183  result[i] = static_cast< ValueType >( ( *this )[i] * value );
184  }
185  return result;
186  }
187 
190  template< class Tt >
191  inline Self operator/(const Tt & value) const
192  {
193  Self result;
194 
195  for ( unsigned int i = 0; i < NVectorDimension; i++ )
196  {
197  result[i] = static_cast< ValueType >( ( *this )[i] / value );
198  }
199  return result;
200  }
201 
206  bool operator==(const Self & v) const
207  { return Superclass::operator==(v); }
208  bool operator!=(const Self & v) const
209  { return !operator==(v); }
211 
213  RealValueType GetNorm(void) const;
214 
216  RealValueType GetSquaredNorm(void) const;
217 
219  static unsigned int GetNumberOfComponents() { return NVectorDimension; }
220 
222  void Normalize(void);
223 
224  void SetNthComponent(int c, const ComponentType & v)
225  { this->operator[](c) = v; }
226 
229  template< typename TCoordRepB >
231  {
232  for ( unsigned int i = 0; i < NVectorDimension; i++ )
233  {
234  ( *this )[i] = static_cast< T >( pa[i] );
235  }
236  }
238 
239  template<typename TCoordRepB>
241  {
243  for (unsigned int i = 0; i < NVectorDimension; i++)
244  {
245  r[i] = static_cast<TCoordRepB> ((*this)[i]);
246  }
247  return r;
248  }
249 
250 };
251 
254 template< class T, unsigned int NVectorDimension >
255 inline
256 Vector< T, NVectorDimension >
257 operator*(const T & scalar, const Vector< T, NVectorDimension > & v)
258 {
259  return v * scalar;
260 }
261 
262 template< class T, unsigned int NVectorDimension >
263 std::ostream & operator<<(std::ostream & os,
264  const Vector< T, NVectorDimension > & v);
265 
266 template< class T, unsigned int NVectorDimension >
267 std::istream & operator>>(std::istream & is,
268  Vector< T, NVectorDimension > & v);
269 
270 ITKCommon_EXPORT Vector< double, 3 > CrossProduct(const Vector< double, 3 > &,
271  const Vector< double, 3 > &);
272 
273 ITKCommon_EXPORT Vector< float, 3 > CrossProduct(const Vector< float, 3 > &,
274  const Vector< float, 3 > &);
275 
276 ITKCommon_EXPORT Vector< int, 3 > CrossProduct(const Vector< int, 3 > &,
277  const Vector< int, 3 > &);
278 } // end namespace itk
279 
280 #ifndef ITK_MANUAL_INSTANTIATION
281 #include "itkVector.hxx"
282 #endif
283 
284 #endif
285