ITK  6.0.0
Insight Toolkit
itkSymmetricSecondRankTensor.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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  * https://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 VDimension = 3>
75 class ITK_TEMPLATE_EXPORT SymmetricSecondRankTensor : public FixedArray<TComponent, VDimension *(VDimension + 1) / 2>
76 {
77 public:
80  using Superclass = FixedArray<TComponent, VDimension *(VDimension + 1) / 2>;
81 
83  static constexpr unsigned int Dimension = VDimension;
84  static constexpr unsigned int InternalDimension = VDimension * (VDimension + 1) / 2;
85 
88 
91 
95 
97  using ComponentType = TComponent;
98  using typename Superclass::ValueType;
101 
104 
107 #ifdef ITK_FUTURE_LEGACY_REMOVE
108  SymmetricSecondRankTensor() = default;
109 #else
110  SymmetricSecondRankTensor() { this->Fill(0); }
111 #endif
112 
114  SymmetricSecondRankTensor(const ComponentType & r) { this->Fill(r); }
115 
117  template <typename TCoordRepB>
119  : BaseArray(pa)
120  {}
121 
122  using ComponentArrayType = ComponentType[Self::InternalDimension];
123 
126  : BaseArray(r)
127  {}
128 
130  template <typename TCoordRepB>
131  Self &
133  {
134  BaseArray::operator=(pa);
135  return *this;
136  }
140  Self &
141  operator=(const ComponentType & r);
142 
143  Self &
144  operator=(const ComponentArrayType r);
145 
148  Self
149  operator+(const Self & r) const;
150 
151  Self
152  operator-(const Self & r) const;
153 
154  const Self &
155  operator+=(const Self & r);
156 
157  const Self &
158  operator-=(const Self & r);
159 
161  Self operator*(const RealValueType & r) const;
162 
163  Self
164  operator/(const RealValueType & r) const;
165 
166  const Self &
167  operator*=(const RealValueType & r);
168 
169  const Self &
170  operator/=(const RealValueType & r);
171 
173  static unsigned int
175  {
176  return Self::InternalDimension;
177  }
178 
180  ComponentType
181  GetNthComponent(int c) const
182  {
183  return this->operator[](c);
184  }
185 
187  void
188  SetNthComponent(int c, const ComponentType & v)
189  {
190  this->operator[](c) = v;
191  }
192 
194  ValueType &
195  operator()(unsigned int row, unsigned int col);
196 
197  const ValueType &
198  operator()(unsigned int row, unsigned int col) const;
199 
202  void
203  SetIdentity();
204 
206  AccumulateValueType
207  GetTrace() const;
208 
210  void
211  ComputeEigenValues(EigenValuesArrayType & eigenValues) const;
212 
215  void
216  ComputeEigenAnalysis(EigenValuesArrayType & eigenValues, EigenVectorsMatrixType & eigenVectors) const;
217 
221  template <typename TMatrixValueType>
222  Self
223  Rotate(const Matrix<TMatrixValueType, VDimension, VDimension> & m) const;
224  template <typename TMatrixValueType>
225  Self
226  Rotate(const vnl_matrix_fixed<TMatrixValueType, VDimension, VDimension> & m) const
227  {
228  return this->Rotate(static_cast<Matrix<TMatrixValueType, VDimension, VDimension>>(m));
229  }
230  template <typename TMatrixValueType>
231  Self
232  Rotate(const vnl_matrix<TMatrixValueType> & m) const
233  {
234  return this->Rotate(static_cast<Matrix<TMatrixValueType>>(m));
235  }
239  MatrixType
240  PreMultiply(const MatrixType & m) const;
241 
243  MatrixType
244  PostMultiply(const MatrixType & m) const;
245 
246 private:
247 };
248 
251 using OutputStreamType = std::ostream;
252 using InputStreamType = std::istream;
253 
254 template <typename TComponent, unsigned int VDimension>
257 
258 template <typename TComponent, unsigned int VDimension>
261 
262 template <typename T>
263 inline void
265 {
266  a.swap(b);
267 }
268 } // end namespace itk
269 
271 
272 #ifndef ITK_MANUAL_INSTANTIATION
273 # include "itkSymmetricSecondRankTensor.hxx"
274 #endif
275 
276 #endif
itk::OutputStreamType
std::ostream OutputStreamType
Definition: itkSymmetricSecondRankTensor.h:251
itk::SymmetricSecondRankTensor::Rotate
Self Rotate(const vnl_matrix< TMatrixValueType > &m) const
Definition: itkSymmetricSecondRankTensor.h:232
itkMatrix.h
itk::operator*
CovariantVector< T, VVectorDimension > operator*(const T &scalar, const CovariantVector< T, VVectorDimension > &v)
Definition: itkCovariantVector.h:268
itk::SymmetricSecondRankTensor::GetNthComponent
ComponentType GetNthComponent(int c) const
Definition: itkSymmetricSecondRankTensor.h:181
itk::operator-
ConstNeighborhoodIterator< TImage > operator-(const ConstNeighborhoodIterator< TImage > &it, const typename ConstNeighborhoodIterator< TImage >::OffsetType &ind)
Definition: itkConstNeighborhoodIterator.h:672
itk::operator<<
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)
itkSymmetricEigenAnalysis.h
itk::SymmetricSecondRankTensor::SymmetricSecondRankTensor
SymmetricSecondRankTensor(const SymmetricSecondRankTensor< TCoordRepB, VDimension > &pa)
Definition: itkSymmetricSecondRankTensor.h:118
itkNumericTraitsTensorPixel.h
itk::SymmetricSecondRankTensor::Rotate
Self Rotate(const vnl_matrix_fixed< TMatrixValueType, VDimension, VDimension > &m) const
Definition: itkSymmetricSecondRankTensor.h:226
itk::SymmetricSecondRankTensor
Represent a symmetric tensor of second rank.
Definition: itkSymmetricSecondRankTensor.h:75
itk::SymmetricSecondRankTensor::SetNthComponent
void SetNthComponent(int c, const ComponentType &v)
Definition: itkSymmetricSecondRankTensor.h:188
itk::SymmetricSecondRankTensor< TComponent, 3 >::RealValueType
typename NumericTraits< ValueType >::RealType RealValueType
Definition: itkSymmetricSecondRankTensor.h:100
itk::SymmetricSecondRankTensor::operator=
Self & operator=(const SymmetricSecondRankTensor< TCoordRepB, VDimension > &pa)
Definition: itkSymmetricSecondRankTensor.h:132
itkIndent.h
itk::SymmetricSecondRankTensor::GetNumberOfComponents
static unsigned int GetNumberOfComponents()
Definition: itkSymmetricSecondRankTensor.h:174
itkFixedArray.h
itk::SymmetricSecondRankTensor< TComponent, 3 >::ComponentArrayType
ComponentType[Self::InternalDimension] ComponentArrayType
Definition: itkSymmetricSecondRankTensor.h:122
itk::SymmetricSecondRankTensor::SymmetricSecondRankTensor
SymmetricSecondRankTensor(const ComponentType &r)
Definition: itkSymmetricSecondRankTensor.h:114
itk::SymmetricSecondRankTensor< TComponent, 3 >::ComponentType
TComponent ComponentType
Definition: itkSymmetricSecondRankTensor.h:97
itk::FixedArray
Simulate a standard C array with copy semantics.
Definition: itkFixedArray.h:53
itk::Matrix
A templated class holding a M x N size Matrix.
Definition: itkMatrix.h:52
itk::SymmetricSecondRankTensor< TComponent, 3 >::AccumulateValueType
typename NumericTraits< ValueType >::RealType AccumulateValueType
Definition: itkSymmetricSecondRankTensor.h:99
itk::SymmetricEigenAnalysisFixedDimension
Definition: itkSymmetricEigenAnalysis.h:758
itk::swap
void swap(Array< T > &a, Array< T > &b) noexcept
Definition: itkArray.h:242
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::InputStreamType
std::istream InputStreamType
Definition: itkSymmetricSecondRankTensor.h:252
itk::FixedArray< TComponent, VDimension *(VDimension+1)/2 >::swap
void swap(FixedArray &other)
Definition: itkFixedArray.h:425
itk::operator+
ConstNeighborhoodIterator< TImage > operator+(const ConstNeighborhoodIterator< TImage > &it, const typename ConstNeighborhoodIterator< TImage >::OffsetType &ind)
Definition: itkConstNeighborhoodIterator.h:654
itk::NumericTraits::RealType
double RealType
Definition: itkNumericTraits.h:86
AddImageFilter
Definition: itkAddImageFilter.h:81
itk::SymmetricSecondRankTensor::SymmetricSecondRankTensor
SymmetricSecondRankTensor(const ComponentArrayType r)
Definition: itkSymmetricSecondRankTensor.h:125
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44
itk::operator>>
std::istream & operator>>(std::istream &is, Point< T, VPointDimension > &vct)
itk::FixedArray< TComponent, VDimension *(VDimension+1)/2 >::ValueType
TComponent ValueType
Definition: itkFixedArray.h:63