ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkBSplineInterpolateImageFunction.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/*=========================================================================
19 *
20 * Portions of this file are subject to the VTK Toolkit Version 3 copyright.
21 *
22 * Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
23 *
24 * For complete copyright, license and disclaimer of warranty information
25 * please refer to the NOTICE file at the top of the ITK source tree.
26 *
27 *=========================================================================*/
28#ifndef itkBSplineInterpolateImageFunction_h
29#define itkBSplineInterpolateImageFunction_h
30
32#include "vnl/vnl_matrix.h"
33
35#include "itkConceptChecking.h"
36#include "itkCovariantVector.h"
37
38#include <memory> // For unique_ptr.
39#include <vector>
40
41namespace itk
42{
68template <typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
69class ITK_TEMPLATE_EXPORT BSplineInterpolateImageFunction : public InterpolateImageFunction<TImageType, TCoordinate>
70{
71public:
72 ITK_DISALLOW_COPY_AND_MOVE(BSplineInterpolateImageFunction);
73
79
81 itkOverrideGetNameOfClassMacro(BSplineInterpolateImageFunction);
82
84 itkNewMacro(Self);
85
87 using typename Superclass::OutputType;
88
90 using typename Superclass::InputImageType;
91
93 static constexpr unsigned int ImageDimension = Superclass::ImageDimension;
94
96 using typename Superclass::IndexType;
97
99 using typename Superclass::SizeType;
100
102 using typename Superclass::ContinuousIndexType;
103
105 using typename Superclass::PointType;
106
109
111 using CoefficientDataType = TCoefficientType;
113
117
120
130 Evaluate(const PointType & point) const override
131 {
132 const ContinuousIndexType index =
133 this->GetInputImage()->template TransformPhysicalPointToContinuousIndex<TCoordinate>(point);
134 // No thread info passed in, so call method that doesn't need thread ID.
135 return (this->EvaluateAtContinuousIndex(index));
136 }
137
138
139 virtual OutputType
140 Evaluate(const PointType & point, ThreadIdType threadId) const
141 {
142 const ContinuousIndexType index =
143 this->GetInputImage()->template TransformPhysicalPointToContinuousIndex<TCoordinate>(point);
144 return (this->EvaluateAtContinuousIndex(index, threadId));
145 }
146
147 OutputType
148 EvaluateAtContinuousIndex(const ContinuousIndexType & index) const override
149 {
150 // Don't know thread information, make evaluateIndex, weights on the stack.
151 // Slower, but safer.
152 vnl_matrix<long> evaluateIndex(ImageDimension, (m_SplineOrder + 1));
153 vnl_matrix<double> weights(ImageDimension, (m_SplineOrder + 1));
154
155 // Pass evaluateIndex, weights by reference. They're only good as long
156 // as this method is in scope.
157 return this->EvaluateAtContinuousIndexInternal(index, evaluateIndex, weights);
158 }
159
160 virtual OutputType
162 {
163 // Pass evaluateIndex, weights by reference. Different threadIDs get different instances.
165 }
166
167 CovariantVectorType
169 {
170 const ContinuousIndexType index =
171 this->GetInputImage()->template TransformPhysicalPointToContinuousIndex<TCoordinate>(point);
172
173 // No thread info passed in, so call method that doesn't need thread ID.
174 return (this->EvaluateDerivativeAtContinuousIndex(index));
175 }
176
177 CovariantVectorType
179 {
180 const ContinuousIndexType index =
181 this->GetInputImage()->template TransformPhysicalPointToContinuousIndex<TCoordinate>(point);
182 return (this->EvaluateDerivativeAtContinuousIndex(index, threadId));
183 }
184
185 CovariantVectorType
187 {
188 // Don't know thread information, make evaluateIndex, weights,
189 // weightsDerivative
190 // on the stack.
191 // Slower, but safer.
192 vnl_matrix<long> evaluateIndex(ImageDimension, (m_SplineOrder + 1));
193 vnl_matrix<double> weights(ImageDimension, (m_SplineOrder + 1));
194 vnl_matrix<double> weightsDerivative(ImageDimension, (m_SplineOrder + 1));
195
196 // Pass evaluateIndex, weights, weightsDerivative by reference. They're only
197 // good
198 // as long as this method is in scope.
199 return this->EvaluateDerivativeAtContinuousIndexInternal(x, evaluateIndex, weights, weightsDerivative);
200 }
201
202 CovariantVectorType
208
209 void
211 {
212 const ContinuousIndexType index =
213 this->GetInputImage()->template TransformPhysicalPointToContinuousIndex<TCoordinate>(point);
214
215 // No thread info passed in, so call method that doesn't need thread ID.
216 this->EvaluateValueAndDerivativeAtContinuousIndex(index, value, deriv);
217 }
218
219 void
221 OutputType & value,
222 CovariantVectorType & deriv,
223 ThreadIdType threadId) const
224 {
225 const ContinuousIndexType index =
226 this->GetInputImage()->template TransformPhysicalPointToContinuousIndex<TCoordinate>(point);
227 this->EvaluateValueAndDerivativeAtContinuousIndex(index, value, deriv, threadId);
228 }
229
230 void
232 OutputType & value,
233 CovariantVectorType & deriv) const
234 {
235 // Don't know thread information, make evaluateIndex, weights,
236 // weightsDerivative
237 // on the stack.
238 // Slower, but safer.
239 vnl_matrix<long> evaluateIndex(ImageDimension, (m_SplineOrder + 1));
240 vnl_matrix<double> weights(ImageDimension, (m_SplineOrder + 1));
241 vnl_matrix<double> weightsDerivative(ImageDimension, (m_SplineOrder + 1));
242
243 // Pass evaluateIndex, weights, weightsDerivative by reference. They're only
244 // good
245 // as long as this method is in scope.
247 x, value, deriv, evaluateIndex, weights, weightsDerivative);
248 }
249
250 void
252 OutputType & value,
253 CovariantVectorType & derivativeValue,
254 ThreadIdType threadId) const
255 {
257 value,
258 derivativeValue,
259 m_ThreadedEvaluateIndex[threadId],
260 m_ThreadedWeights[threadId],
262 }
263
266 void
267 SetSplineOrder(unsigned int SplineOrder);
268
269 itkGetConstMacro(SplineOrder, unsigned int);
270
271 void
273
274 itkGetConstMacro(NumberOfWorkUnits, ThreadIdType);
275
277 void
278 SetInputImage(const TImageType * inputData) override;
279
290 itkSetMacro(UseImageDirection, bool);
291 itkGetConstMacro(UseImageDirection, bool);
292 itkBooleanMacro(UseImageDirection);
294
296 GetRadius() const override
297 {
298 return SizeType::Filled(m_SplineOrder + 1);
299 }
300
301protected:
320 virtual OutputType
322 vnl_matrix<long> & evaluateIndex,
323 vnl_matrix<double> & weights) const;
324
325 virtual void
327 OutputType & value,
328 CovariantVectorType & derivativeValue,
329 vnl_matrix<long> & evaluateIndex,
330 vnl_matrix<double> & weights,
331 vnl_matrix<double> & weightsDerivative) const;
332
333 virtual CovariantVectorType
335 vnl_matrix<long> & evaluateIndex,
336 vnl_matrix<double> & weights,
337 vnl_matrix<double> & weightsDerivative) const;
338
341 void
342 PrintSelf(std::ostream & os, Indent indent) const override;
343
344 // These are needed by the smoothing spline routine.
345 // temp storage for processing of Coefficients
346 std::vector<CoefficientDataType> m_Scratch{};
347 // Image size
348 typename TImageType::SizeType m_DataLength{};
349 // User specified spline order (3rd or cubic is the default)
350 unsigned int m_SplineOrder{};
351
352 // Spline coefficients
354
355private:
357 void
359 const vnl_matrix<long> & EvaluateIndex,
360 vnl_matrix<double> & weights,
361 unsigned int splineOrder) const;
362
364 void
366 const vnl_matrix<long> & EvaluateIndex,
367 vnl_matrix<double> & weights,
368 unsigned int splineOrder) const;
369
372 void
374
376 void
377 DetermineRegionOfSupport(vnl_matrix<long> & evaluateIndex,
378 const ContinuousIndexType & x,
379 unsigned int splineOrder) const;
380
383 void
384 ApplyMirrorBoundaryConditions(vnl_matrix<long> & evaluateIndex, unsigned int splineOrder) const;
385
386 Iterator m_CIterator{}; // Iterator for
387 // traversing spline
388 // coefficients.
389 unsigned long m_MaxNumberInterpolationPoints{}; // number of
390 // neighborhood
391 // points used for
392 // interpolation
393 std::vector<IndexType> m_PointsToIndex{}; // Preallocation of
394 // interpolation
395 // neighborhood
396 // indices
397
399
400 // flag to take or not the image direction into account when computing the
401 // derivatives.
403
405 std::unique_ptr<vnl_matrix<long>[]> m_ThreadedEvaluateIndex;
406 std::unique_ptr<vnl_matrix<double>[]> m_ThreadedWeights;
407 std::unique_ptr<vnl_matrix<double>[]> m_ThreadedWeightsDerivative;
408};
409} // namespace itk
410
411#ifndef ITK_MANUAL_INSTANTIATION
412# include "itkBSplineInterpolateImageFunction.hxx"
413#endif
414
415#endif
Calculates the B-Spline coefficients of an image. Spline order may be from 0 to 5.
typename CoefficientFilter::Pointer CoefficientFilterPointer
virtual CovariantVectorType EvaluateDerivativeAtContinuousIndexInternal(const ContinuousIndexType &x, vnl_matrix< long > &evaluateIndex, vnl_matrix< double > &weights, vnl_matrix< double > &weightsDerivative) const
OutputType Evaluate(const PointType &point) const override
virtual OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &x, ThreadIdType threadId) const
ImageLinearIteratorWithIndex< TImageType > Iterator
CovariantVector< OutputType, Self::ImageDimension > CovariantVectorType
void EvaluateValueAndDerivativeAtContinuousIndex(const ContinuousIndexType &x, OutputType &value, CovariantVectorType &derivativeValue, ThreadIdType threadId) const
virtual void EvaluateValueAndDerivativeAtContinuousIndexInternal(const ContinuousIndexType &x, OutputType &value, CovariantVectorType &derivativeValue, vnl_matrix< long > &evaluateIndex, vnl_matrix< double > &weights, vnl_matrix< double > &weightsDerivative) const
void PrintSelf(std::ostream &os, Indent indent) const override
void DetermineRegionOfSupport(vnl_matrix< long > &evaluateIndex, const ContinuousIndexType &x, unsigned int splineOrder) const
virtual OutputType EvaluateAtContinuousIndexInternal(const ContinuousIndexType &x, vnl_matrix< long > &evaluateIndex, vnl_matrix< double > &weights) const
void SetInputImage(const TImageType *inputData) override
void EvaluateValueAndDerivative(const PointType &point, OutputType &value, CovariantVectorType &deriv) const
void EvaluateValueAndDerivative(const PointType &point, OutputType &value, CovariantVectorType &deriv, ThreadIdType threadId) const
void SetSplineOrder(unsigned int SplineOrder)
~BSplineInterpolateImageFunction() override=default
CovariantVectorType EvaluateDerivative(const PointType &point) const
void SetDerivativeWeights(const ContinuousIndexType &x, const vnl_matrix< long > &EvaluateIndex, vnl_matrix< double > &weights, unsigned int splineOrder) const
BSplineDecompositionImageFilter< TImageType, CoefficientImageType > CoefficientFilter
void SetNumberOfWorkUnits(ThreadIdType numWorkUnits)
CovariantVectorType EvaluateDerivativeAtContinuousIndex(const ContinuousIndexType &x, ThreadIdType threadId) const
void ApplyMirrorBoundaryConditions(vnl_matrix< long > &evaluateIndex, unsigned int splineOrder) const
void EvaluateValueAndDerivativeAtContinuousIndex(const ContinuousIndexType &x, OutputType &value, CovariantVectorType &deriv) const
InterpolateImageFunction< TImageType, TCoordinate > Superclass
virtual OutputType Evaluate(const PointType &point, ThreadIdType threadId) const
CovariantVectorType EvaluateDerivative(const PointType &point, ThreadIdType threadId) const
Image< CoefficientDataType, Self::ImageDimension > CoefficientImageType
void SetInterpolationWeights(const ContinuousIndexType &x, const vnl_matrix< long > &EvaluateIndex, vnl_matrix< double > &weights, unsigned int splineOrder) const
A templated class holding a n-Dimensional covariant vector.
const InputImageType * GetInputImage() const
A multi-dimensional image iterator that visits image pixels within a region in a "scan-line" order.
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
Point< TCoordinate, Self::ImageDimension > PointType
ContinuousIndex< TCoordinate, Self::ImageDimension > ContinuousIndexType
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned int ThreadIdType
*par Constraints *The filter image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating point(float or double) or *a user-defined "real" numerical type with arithmetic operations defined *sufficient to compute derivatives. **\par Performance *This filter will automatically multithread if run with *SetUsePrincipleComponents