ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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 * 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 itkDisplacementFieldJacobianDeterminantFilter_h
19#define itkDisplacementFieldJacobianDeterminantFilter_h
20
23#include "itkVector.h"
24#include "vnl/vnl_matrix.h"
25#include "vnl/vnl_det.h"
26
27namespace itk
28{
111template <typename TInputImage,
112 typename TRealType = float,
113 typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
115 : public ImageToImageFilter<TInputImage, TOutputImage>
116{
117public:
118 ITK_DISALLOW_COPY_AND_MOVE(DisplacementFieldJacobianDeterminantFilter);
119
125
127 itkNewMacro(Self);
128
130 itkOverrideGetNameOfClassMacro(DisplacementFieldJacobianDeterminantFilter);
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
161
170 void
172
179 void
181 itkGetConstMacro(UseImageSpacing, bool);
182 itkBooleanMacro(UseImageSpacing);
184
185#if !defined(ITK_FUTURE_LEGACY_REMOVE)
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
214 itkGetConstReferenceMacro(DerivativeWeights, WeightsType);
216
217protected:
220
224 void
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
258
261
265
266private:
267 bool m_UseImageSpacing{ true };
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
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
~DisplacementFieldJacobianDeterminantFilter() override=default
typename OutputImageType::RegionType OutputImageRegionType
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
virtual TRealType EvaluateAtNeighborhood(const ConstNeighborhoodIteratorType &it) const
Image< RealVectorType, TInputImage::ImageDimension > RealVectorImageType
void PrintSelf(std::ostream &os, Indent indent) const override
ConstNeighborhoodIterator< RealVectorImageType > ConstNeighborhoodIteratorType
void SetDerivativeWeights(const WeightsType &)
Simulate a standard C array with copy semantics.
typename OutputImageType::RegionType OutputImageRegionType
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
Implements transparent reference counting.
A templated class holding a n-Dimensional vector.
Definition itkVector.h:63
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned int ThreadIdType