ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkStrainImageFilter.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 itkStrainImageFilter_h
19 #define itkStrainImageFilter_h
20 
21 #include "itkCovariantVector.h"
22 #include "itkImageToImageFilter.h"
25 
26 namespace itk
27 {
28 
62 template< typename TInputImage, typename TOperatorValueType=float, typename TOutputValueType=float >
63 class StrainImageFilter : public
64  ImageToImageFilter< TInputImage,
65  Image< SymmetricSecondRankTensor< TOutputValueType, TInputImage::ImageDimension >,
66  TInputImage::ImageDimension > >
67 {
68 public:
69  ITK_DISALLOW_COPY_AND_ASSIGN(StrainImageFilter);
70 
72  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
73 
74  using InputImageType = TInputImage;
78 
82 
85 
86 
91 
94 
96  itkNewMacro( Self );
97 
99  itkTypeMacro( StrainImageFilter, ImageToImageFilter );
100 
103  itkSetObjectMacro( GradientFilter, GradientFilterType );
104  itkGetConstObjectMacro( GradientFilter, GradientFilterType );
106 
111  itkSetObjectMacro( VectorGradientFilter, VectorGradientFilterType );
112  itkGetConstObjectMacro( VectorGradientFilter, VectorGradientFilterType );
114 
121 
122  itkSetMacro( StrainForm, StrainFormType );
123  itkGetConstMacro( StrainForm, StrainFormType );
124 
125 protected:
127 
129 
130  void BeforeThreadedGenerateData() override;
131 
132  void DynamicThreadedGenerateData( const OutputRegionType& outputRegion ) override;
133 
135 
136  void PrintSelf ( std::ostream& os, Indent indent ) const override;
137 
138 private:
140 
142 
144 
146 };
147 
148 } // end namespace itk
149 
150 #ifndef ITK_MANUAL_INSTANTIATION
151 #include "itkStrainImageFilter.hxx"
152 #endif
153 
154 #endif
Extract components of an Image with multi-component pixels.
Light weight base class for most itk classes.
Represent a symmetric tensor of second rank.
Generate a strain field image from a displacement field image.
typename OutputImageType::RegionType OutputRegionType
InputComponentsImageFilterType::Pointer m_InputComponentsFilter
static constexpr unsigned int ImageDimension
VectorGradientFilterType::Pointer m_VectorGradientFilter
typename Superclass::RegionType RegionType
Definition: itkImage.h:139
void BeforeThreadedGenerateData() override
GradientFilterType::Pointer m_GradientFilter
void DynamicThreadedGenerateData(const OutputRegionType &outputRegion) override
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
void PrintSelf(std::ostream &os, Indent indent) const override
A templated class holding a n-Dimensional covariant vector.
Templated n-dimensional image class.
Definition: itkImage.h:75