ITK  4.8.0
Insight Segmentation and Registration Toolkit
itkDisplacementFieldJacobianDeterminantFilter.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 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,
114  TInputImage::ImageDimension >
115  >
117  public ImageToImageFilter< TInputImage, TOutputImage >
118 {
119 public:
125 
127  itkNewMacro(Self);
128 
131 
134  typedef typename TOutputImage::PixelType OutputPixelType;
135  typedef typename TInputImage::PixelType InputPixelType;
136 
138  typedef TInputImage InputImageType;
139  typedef TOutputImage OutputImageType;
140  typedef typename InputImageType::Pointer InputImagePointer;
141  typedef typename OutputImageType::Pointer OutputImagePointer;
142 
144  itkStaticConstMacro(ImageDimension, unsigned int,
145  TOutputImage::ImageDimension);
146 
148  itkStaticConstMacro(VectorDimension, unsigned int,
149  InputPixelType::Dimension);
150 
152  typedef TRealType RealType;
157 
162 
165 
174  virtual void GenerateInputRequestedRegion() ITK_OVERRIDE;
175 
181  { this->SetUseImageSpacing(true); }
182 
187  { this->SetUseImageSpacing(false); }
188 
191  void SetUseImageSpacing(bool);
192 
193  itkGetConstMacro(UseImageSpacing, bool);
194 
196 
199  void SetDerivativeWeights(const WeightsType &);
200  itkGetConstReferenceMacro(DerivativeWeights, WeightsType);
202 
203 protected:
206 
210  void BeforeThreadedGenerateData() ITK_OVERRIDE;
211 
224  void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
225  ThreadIdType threadId) ITK_OVERRIDE;
226 
227  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
228 
230 
232  itkGetConstObjectMacro(RealValuedInputImage, ImageBaseType);
233 
235  itkGetConstReferenceMacro(NeighborhoodRadius, RadiusType);
236  itkSetMacro(NeighborhoodRadius, RadiusType);
238 
239  virtual TRealType EvaluateAtNeighborhood(const ConstNeighborhoodIteratorType & it) const;
240 
243 
247 
248 private:
250 
252 
253  typename ImageBaseType::ConstPointer m_RealValuedInputImage;
254 
255  DisplacementFieldJacobianDeterminantFilter(const Self &); //purposely not
256  // implemented
257  void operator=(const Self &); //purposely not
258 
259  // implemented
260 
262 };
263 } // end namespace itk
264 
265 #ifndef ITK_MANUAL_INSTANTIATION
266 #include "itkDisplacementFieldJacobianDeterminantFilter.hxx"
267 #endif
268 
269 #endif
Image< RealVectorType, TInputImage::ImageDimension > RealVectorImageType
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
virtual TRealType EvaluateAtNeighborhood(const ConstNeighborhoodIteratorType &it) const
void PrintSelf(std::ostream &os, Indent indent) const override
Base class for all process objects that output image data.
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
virtual void GenerateInputRequestedRegion() override
void SetDerivativeWeights(const WeightsType &)
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId) override
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
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
Superclass::OutputImageRegionType OutputImageRegionType
Templated n-dimensional image class.
Definition: itkImage.h:75
Computes a scalar image from a vector image (e.g., deformation field) input, where each output scalar...
ConstNeighborhoodIterator< RealVectorImageType > ConstNeighborhoodIteratorType