ITK  4.6.0
Insight Segmentation and Registration Toolkit
itkSymmetricSecondRankTensor.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 __itkSymmetricSecondRankTensor_h
19 #define __itkSymmetricSecondRankTensor_h
20 
21 // Undefine an eventual SymmetricSecondRankTensor macro
22 #ifdef SymmetricSecondRankTensor
23 #undef SymmetricSecondRankTensor
24 #endif
25 
26 #include "itkIndent.h"
27 #include "itkFixedArray.h"
28 #include "itkMatrix.h"
30 
31 namespace itk
32 {
74 template< typename TComponent, unsigned int NDimension = 3 >
76  FixedArray< TComponent, NDimension *( NDimension + 1 ) / 2 >
77 {
78 public:
81  typedef FixedArray< TComponent, NDimension *( NDimension + 1 ) / 2 > Superclass;
82 
84  itkStaticConstMacro(Dimension, unsigned int, NDimension);
85  itkStaticConstMacro( InternalDimension, unsigned int, ( NDimension * ( NDimension + 1 ) / 2 ) );
87 
89  typedef FixedArray< TComponent,
90  itkGetStaticConstMacro(InternalDimension) > BaseArray;
91 
94 
98 
100  typedef TComponent ComponentType;
104 
107 
110 
111  SymmetricSecondRankTensor (const ComponentType & r) { this->Fill(r); }
112 
114  template< typename TCoordRepB >
116  BaseArray(pa) {}
117 
118  typedef ComponentType ComponentArrayType[itkGetStaticConstMacro(InternalDimension)];
119 
122 
124  template< typename TCoordRepB >
126  {
128  return *this;
129  }
131 
133  Self & operator=(const ComponentType & r);
134 
135  Self & operator=(const ComponentArrayType r);
136 
139  Self operator+(const Self & vec) const;
140 
141  Self operator-(const Self & vec) const;
142 
143  const Self & operator+=(const Self & vec);
144 
145  const Self & operator-=(const Self & vec);
146 
148  Self operator *(const RealValueType & scalar) const;
149 
150  Self operator/(const RealValueType & scalar) const;
151 
152  const Self & operator*=(const RealValueType & scalar);
153 
154  const Self & operator/=(const RealValueType & scalar);
155 
157  static unsigned int GetNumberOfComponents()
158  {
159  return itkGetStaticConstMacro(InternalDimension);
160  }
161 
163  ComponentType GetNthComponent(int c) const { return this->operator[](c); }
164 
166  void SetNthComponent(int c, const ComponentType & v) { this->operator[](c) = v; }
167 
169  ValueType & operator()(unsigned int row, unsigned int col);
170 
171  const ValueType & operator()(unsigned int row, unsigned int col) const;
172 
175  void SetIdentity();
176 
179 
181  void ComputeEigenValues(EigenValuesArrayType & eigenValues) const;
182 
185  void ComputeEigenAnalysis(EigenValuesArrayType & eigenValues,
186  EigenVectorsMatrixType & eigenVectors) const;
187 
191  template<typename TMatrixValueType>
193  template<typename TMatrixValueType>
194  Self Rotate( const vnl_matrix_fixed<TMatrixValueType, NDimension, NDimension> & m) const
195  {
196  return this->Rotate( static_cast<Matrix<TMatrixValueType, NDimension, NDimension> >(m) );
197  }
198  template<typename TMatrixValueType>
199  Self Rotate( const vnl_matrix<TMatrixValueType> & m) const
200  {
201  return this->Rotate( static_cast<Matrix<TMatrixValueType> >(m) );
202  }
204 
206  MatrixType PreMultiply(const MatrixType & m) const;
207 
209  MatrixType PostMultiply(const MatrixType & m) const;
210 
211 private:
212 };
213 
216 typedef std::ostream OutputStreamType;
217 typedef std::istream InputStreamType;
218 
219 template< typename TComponent, unsigned int NDimension >
222 
223 template< typename TComponent, unsigned int NDimension >
226 } // end namespace itk
227 
229 
230 #ifndef ITK_MANUAL_INSTANTIATION
231 #include "itkSymmetricSecondRankTensor.hxx"
232 #endif
233 
234 #endif
Self operator-(const Self &vec) const
Self operator+(const Self &vec) const
A templated class holding a M x N size Matrix.
Definition: itkMatrix.h:46
Self & operator=(const SymmetricSecondRankTensor< TCoordRepB, NDimension > &pa)
ValueType & operator()(unsigned int row, unsigned int col)
Find Eigen values of a real 2D symmetric matrix. It serves as a thread-safe alternative to the class:...
FixedArray< TComponent, NDimension > EigenValuesArrayType
Self operator/(const RealValueType &scalar) const
std::istream InputStreamType
Self Rotate(const Matrix< TMatrixValueType, NDimension, NDimension > &m) const
Self operator*(const RealValueType &scalar) const
Represent a symmetric tensor of second rank.
std::istream & operator>>(std::istream &is, Point< T, NPointDimension > &v)
SymmetricSecondRankTensor(const ComponentArrayType r)
NumericTraits< ValueType >::RealType AccumulateValueType
void SetNthComponent(int c, const ComponentType &v)
const Self & operator/=(const RealValueType &scalar)
std::ostream & operator<<(std::ostream &os, const Array< TValue > &arr)
Definition: itkArray.h:184
const Self & operator+=(const Self &vec)
std::ostream OutputStreamType
AccumulateValueType GetTrace() const
Simulate a standard C array with copy semnatics.
Definition: itkFixedArray.h:50
SymmetricSecondRankTensor(const ComponentType &r)
static const unsigned int InternalDimension
ComponentType GetNthComponent(int c) const
FixedArray< TComponent, itkGetStaticConstMacro(InternalDimension) > BaseArray
void ComputeEigenValues(EigenValuesArrayType &eigenValues) const
MatrixType PreMultiply(const MatrixType &m) const
Self Rotate(const vnl_matrix_fixed< TMatrixValueType, NDimension, NDimension > &m) const
Matrix< TComponent, NDimension, NDimension > EigenVectorsMatrixType
FixedArray & operator=(const FixedArray< TFixedArrayValueType, VLength > &r)
ComponentType ComponentArrayType[itkGetStaticConstMacro(InternalDimension)]
NumericTraits< ValueType >::RealType RealValueType
const Self & operator-=(const Self &vec)
const Self & operator*=(const RealValueType &scalar)
Matrix< TComponent, NDimension, NDimension > MatrixType
void ComputeEigenAnalysis(EigenValuesArrayType &eigenValues, EigenVectorsMatrixType &eigenVectors) const
SymmetricEigenAnalysis< MatrixType, EigenValuesArrayType, EigenVectorsMatrixType > SymmetricEigenAnalysisType
MatrixType PostMultiply(const MatrixType &m) const
SymmetricSecondRankTensor(const SymmetricSecondRankTensor< TCoordRepB, NDimension > &pa)
Self Rotate(const vnl_matrix< TMatrixValueType > &m) const