ITK  5.2.0
Insight Toolkit
itkDisplacementFieldJacobianDeterminantFilter.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  * 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 itkDisplacementFieldJacobianDeterminantFilter_h
19 #define itkDisplacementFieldJacobianDeterminantFilter_h
20 
22 #include "itkImageToImageFilter.h"
23 #include "itkVector.h"
24 #include "vnl/vnl_matrix.h"
25 #include "vnl/vnl_det.h"
26 
27 namespace itk
28 {
111 template <typename TInputImage,
112  typename TRealType = float,
113  typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
115  : public ImageToImageFilter<TInputImage, TOutputImage>
116 {
117 public:
118  ITK_DISALLOW_COPY_AND_MOVE(DisplacementFieldJacobianDeterminantFilter);
119 
125 
127  itkNewMacro(Self);
128 
131 
134  using OutputPixelType = typename TOutputImage::PixelType;
135  using InputPixelType = typename TInputImage::PixelType;
136 
138  using InputImageType = TInputImage;
139  using OutputImageType = TOutputImage;
140  using InputImagePointer = typename InputImageType::Pointer;
141  using OutputImagePointer = typename OutputImageType::Pointer;
142 
144  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
145 
147  static constexpr unsigned int VectorDimension = InputPixelType::Dimension;
148 
150  using RealType = TRealType;
153 
158 
160  using OutputImageRegionType = typename Superclass::OutputImageRegionType;
161 
170  void
171  GenerateInputRequestedRegion() override;
172 
179  void
180  SetUseImageSpacing(bool);
181  itkGetConstMacro(UseImageSpacing, bool);
182  itkBooleanMacro(UseImageSpacing);
184 
185 #if !defined(ITK_FUTURE_LEGACY_REMOVE)
186 
191  void
192  SetUseImageSpacingOn()
193  {
194  this->SetUseImageSpacing(true);
195  }
196 
201  void
202  SetUseImageSpacingOff()
203  {
204  this->SetUseImageSpacing(false);
205  }
206 #endif
207 
209 
212  void
213  SetDerivativeWeights(const WeightsType &);
214  itkGetConstReferenceMacro(DerivativeWeights, WeightsType);
216 
217 protected:
219  ~DisplacementFieldJacobianDeterminantFilter() override = default;
220 
224  void
225  BeforeThreadedGenerateData() override;
226 
239  void
240  DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
241 
242 
243  void
244  PrintSelf(std::ostream & os, Indent indent) const override;
245 
246  using ImageBaseType = typename InputImageType::Superclass;
247 
249  itkGetConstObjectMacro(RealValuedInputImage, ImageBaseType);
250 
252  itkGetConstReferenceMacro(NeighborhoodRadius, RadiusType);
253  itkSetMacro(NeighborhoodRadius, RadiusType);
255 
256  virtual TRealType
257  EvaluateAtNeighborhood(const ConstNeighborhoodIteratorType & it) const;
258 
261 
265 
266 private:
268 
270 
271  typename ImageBaseType::ConstPointer m_RealValuedInputImage;
272 
274 };
275 } // end namespace itk
276 
277 #ifndef ITK_MANUAL_INSTANTIATION
278 # include "itkDisplacementFieldJacobianDeterminantFilter.hxx"
279 #endif
280 
281 #endif
itk::DisplacementFieldJacobianDeterminantFilter::InputPixelType
typename TInputImage::PixelType InputPixelType
Definition: itkDisplacementFieldJacobianDeterminantFilter.h:135
itk::ConstNeighborhoodIterator::RadiusType
typename Superclass::RadiusType RadiusType
Definition: itkConstNeighborhoodIterator.h:71
itk::ImageSource::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkImageSource.h:91
itkNeighborhoodIterator.h
itk::Vector
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::DisplacementFieldJacobianDeterminantFilter::m_DerivativeWeights
WeightsType m_DerivativeWeights
Definition: itkDisplacementFieldJacobianDeterminantFilter.h:260
itk::DisplacementFieldJacobianDeterminantFilter::RealType
TRealType RealType
Definition: itkDisplacementFieldJacobianDeterminantFilter.h:150
itk::ThreadIdType
unsigned int ThreadIdType
Definition: itkIntTypes.h:99
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::DisplacementFieldJacobianDeterminantFilter::m_NeighborhoodRadius
RadiusType m_NeighborhoodRadius
Definition: itkDisplacementFieldJacobianDeterminantFilter.h:273
itk::DisplacementFieldJacobianDeterminantFilter::OutputPixelType
typename TOutputImage::PixelType OutputPixelType
Definition: itkDisplacementFieldJacobianDeterminantFilter.h:134
itk::ImageToImageFilter::InputImagePointer
typename InputImageType::Pointer InputImagePointer
Definition: itkImageToImageFilter.h:130
itk::DisplacementFieldJacobianDeterminantFilter::m_RealValuedInputImage
ImageBaseType::ConstPointer m_RealValuedInputImage
Definition: itkDisplacementFieldJacobianDeterminantFilter.h:271
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itkImageToImageFilter.h
itk::FixedArray< TRealType, ImageDimension >
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ConstNeighborhoodIterator
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
Definition: itkConstNeighborhoodIterator.h:51
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:138
itkVector.h
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::DisplacementFieldJacobianDeterminantFilter
Computes a scalar image from a vector image (e.g., deformation field) input, where each output scalar...
Definition: itkDisplacementFieldJacobianDeterminantFilter.h:114
itk::DisplacementFieldJacobianDeterminantFilter::RadiusType
typename ConstNeighborhoodIteratorType::RadiusType RadiusType
Definition: itkDisplacementFieldJacobianDeterminantFilter.h:157
itk::DisplacementFieldJacobianDeterminantFilter::ImageBaseType
typename InputImageType::Superclass ImageBaseType
Definition: itkDisplacementFieldJacobianDeterminantFilter.h:246
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44
itk::DisplacementFieldJacobianDeterminantFilter::m_RequestedNumberOfThreads
ThreadIdType m_RequestedNumberOfThreads
Definition: itkDisplacementFieldJacobianDeterminantFilter.h:269
itk::DisplacementFieldJacobianDeterminantFilter::m_UseImageSpacing
bool m_UseImageSpacing
Definition: itkDisplacementFieldJacobianDeterminantFilter.h:267
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90
itk::DisplacementFieldJacobianDeterminantFilter::m_HalfDerivativeWeights
WeightsType m_HalfDerivativeWeights
Definition: itkDisplacementFieldJacobianDeterminantFilter.h:264