File: itkResampleImageFilter.hxx

package info (click to toggle)
insighttoolkit5 5.4.3-5
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 704,384 kB
  • sloc: cpp: 783,592; ansic: 628,724; xml: 44,704; fortran: 34,250; python: 22,874; sh: 4,078; pascal: 2,636; lisp: 2,158; makefile: 464; yacc: 328; asm: 205; perl: 203; lex: 146; tcl: 132; javascript: 98; csh: 81
file content (678 lines) | stat: -rw-r--r-- 25,778 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
/*=========================================================================
 *
 *  Copyright NumFOCUS
 *
 *  Licensed under the Apache License, Version 2.0 (the "License");
 *  you may not use this file except in compliance with the License.
 *  You may obtain a copy of the License at
 *
 *         https://www.apache.org/licenses/LICENSE-2.0.txt
 *
 *  Unless required by applicable law or agreed to in writing, software
 *  distributed under the License is distributed on an "AS IS" BASIS,
 *  WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
 *  See the License for the specific language governing permissions and
 *  limitations under the License.
 *
 *=========================================================================*/
#ifndef itkResampleImageFilter_hxx
#define itkResampleImageFilter_hxx

#include "itkObjectFactory.h"
#include "itkIdentityTransform.h"
#include "itkTotalProgressReporter.h"
#include "itkImageRegionIteratorWithIndex.h"
#include "itkImageScanlineIterator.h"
#include "itkSpecialCoordinatesImage.h"
#include "itkDefaultConvertPixelTraits.h"
#include "itkImageAlgorithm.h"

#include <algorithm>   // For max.
#include <type_traits> // For is_same.

namespace itk
{

template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::
  ResampleImageFilter()
  : m_Extrapolator(nullptr)
  , m_OutputSpacing(MakeFilled<SpacingType>(1.0))
  , m_OutputOrigin()

{

  m_Size.Fill(0);
  m_OutputStartIndex.Fill(0);

  m_OutputDirection.SetIdentity();

  // Pipeline input configuration

  // implicit:
  // #0 "Primary" required

  // #1 "ReferenceImage" optional
  Self::AddRequiredInputName("ReferenceImage", 1);
  Self::RemoveRequiredInputName("ReferenceImage");

  // "Transform" required ( not numbered )
  Self::AddRequiredInputName("Transform");
  this->InitializeTransform();

  m_Interpolator = dynamic_cast<InterpolatorType *>(LinearInterpolatorType::New().GetPointer());

  m_DefaultPixelValue = NumericTraits<PixelType>::ZeroValue(m_DefaultPixelValue);
  this->DynamicMultiThreadingOn();
}

template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
void
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::
  InitializeTransform()
{
  using IdentityTransformType =
    Transform<TTransformPrecisionType, Self::OutputImageDimension, Self::OutputImageDimension>;
  typename IdentityTransformType::Pointer defaultTransform =
    IdentityTransform<TTransformPrecisionType, OutputImageDimension>::New();
  if (InputImageDimension == OutputImageDimension)
  {
    using DecoratorType = DataObjectDecorator<IdentityTransformType>;
    auto decoratedInput = DecoratorType::New();
    decoratedInput->Set(defaultTransform);
    this->ProcessObject::SetInput(
      "Transform", const_cast<DataObjectDecorator<IdentityTransformType> *>(decoratedInput.GetPointer()));
  }
  else
  {
    // Initialize with rectangular identity?
  }
  this->Modified();
}

template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
void
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::
  VerifyPreconditions() ITKv5_CONST
{
  this->Superclass::VerifyPreconditions();
  const ReferenceImageBaseType * const referenceImage = this->GetReferenceImage();
  if (this->m_Size[0] == 0 && referenceImage && !m_UseReferenceImage)
  {
    itkExceptionMacro("Output image size is zero in all dimensions.  Consider using UseReferenceImageOn()."
                      "or SetUseReferenceImage(true) to define the resample output from the ReferenceImage.");
  }
}

template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
void
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::SetOutputSpacing(
  const double * spacing)
{
  SpacingType s;
  for (unsigned int i = 0; i < TOutputImage::ImageDimension; ++i)
  {
    s[i] = static_cast<typename SpacingType::ValueType>(spacing[i]);
  }
  this->SetOutputSpacing(s);
}

template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
void
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::SetOutputOrigin(
  const double * origin)
{
  this->SetOutputOrigin(OriginPointType(origin));
}

template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
void
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::
  SetOutputParametersFromImage(const ImageBaseType * image)
{
  this->SetOutputOrigin(image->GetOrigin());
  this->SetOutputSpacing(image->GetSpacing());
  this->SetOutputDirection(image->GetDirection());
  this->SetOutputStartIndex(image->GetLargestPossibleRegion().GetIndex());
  this->SetSize(image->GetLargestPossibleRegion().GetSize());
}

template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
void
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::
  BeforeThreadedGenerateData()
{
  m_Interpolator->SetInputImage(this->GetInput());

  // Connect input image to extrapolator
  if (!m_Extrapolator.IsNull())
  {
    m_Extrapolator->SetInputImage(this->GetInput());
  }

  unsigned int nComponents = DefaultConvertPixelTraits<PixelType>::GetNumberOfComponents(m_DefaultPixelValue);

  if (nComponents == 0)
  {
    PixelComponentType tempZeroComponents{ 0 };
    PixelComponentType zeroComponent = NumericTraits<PixelComponentType>::ZeroValue(tempZeroComponents);
    nComponents = this->GetInput()->GetNumberOfComponentsPerPixel();
    NumericTraits<PixelType>::SetLength(m_DefaultPixelValue, nComponents);
    for (unsigned int n = 0; n < nComponents; ++n)
    {
      PixelConvertType::SetNthComponent(n, m_DefaultPixelValue, zeroComponent);
    }
  }
}

template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
void
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::
  AfterThreadedGenerateData()
{
  // Disconnect input image from the interpolator
  m_Interpolator->SetInputImage(nullptr);
  if (!m_Extrapolator.IsNull())
  {
    // Disconnect input image from the extrapolator
    m_Extrapolator->SetInputImage(nullptr);
  }
}

template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
void
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::
  DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread)
{
  // Check whether the input or the output is a
  // SpecialCoordinatesImage. If either are, then we cannot use the
  // fast path since index mapping will definitely not be linear.
  using OutputSpecialCoordinatesImageType = SpecialCoordinatesImage<PixelType, OutputImageDimension>;
  using InputSpecialCoordinatesImageType = SpecialCoordinatesImage<InputPixelType, InputImageDimension>;

  if (outputRegionForThread.GetNumberOfPixels() == 0)
  {
    return;
  }

  const bool isSpecialCoordinatesImage =
    ((dynamic_cast<const InputSpecialCoordinatesImageType *>(this->GetInput()) != nullptr) ||
     (dynamic_cast<const OutputSpecialCoordinatesImageType *>(this->GetOutput()) != nullptr));


  // Check whether we can use a fast path for resampling. Fast path
  // can be used if the transformation is linear. Transform respond
  // to the IsLinear() call.
  if (!isSpecialCoordinatesImage &&
      this->GetTransform()->GetTransformCategory() == TransformType::TransformCategoryEnum::Linear)
  {
    this->LinearThreadedGenerateData(outputRegionForThread);
    return;
  }

  // Otherwise, we use the normal method where the transform is called
  // for computing the transformation of every point.
  this->NonlinearThreadedGenerateData(outputRegionForThread);
}


#ifndef ITK_LEGACY_REMOVE
template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
typename ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::PixelType
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::
  CastPixelWithBoundsChecking(const InterpolatorOutputType value,
                              const ComponentType          minComponent,
                              const ComponentType          maxComponent) const
{
  const unsigned int nComponents = InterpolatorConvertType::GetNumberOfComponents(value);
  PixelType          outputValue;

  NumericTraits<PixelType>::SetLength(outputValue, nComponents);

  for (unsigned int n = 0; n < nComponents; ++n)
  {
    ComponentType component = InterpolatorConvertType::GetNthComponent(n, value);

    if (component < minComponent)
    {
      PixelConvertType::SetNthComponent(n, outputValue, static_cast<PixelComponentType>(minComponent));
    }
    else if (component > maxComponent)
    {
      PixelConvertType::SetNthComponent(n, outputValue, static_cast<PixelComponentType>(maxComponent));
    }
    else
    {
      PixelConvertType::SetNthComponent(n, outputValue, static_cast<PixelComponentType>(component));
    }
  }

  return outputValue;
}
#endif // ITK_LEGACY_REMOVE


template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
auto
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::
  CastComponentWithBoundsChecking(const PixelComponentType value) -> PixelComponentType
{
  // Just return the argument. In this case, there is no need to cast or clamp its value.
  return value;
}


template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
template <typename TComponent>
auto
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::
  CastComponentWithBoundsChecking(const TComponent value) -> PixelComponentType
{
  static_assert(std::is_same_v<TComponent, ComponentType>, "TComponent should just be the same as the ComponentType!");
  static_assert(!std::is_same_v<TComponent, PixelComponentType>,
                "For PixelComponentType there is a more appropriate overload, that should be called instead!");

  // Retrieve minimum and maximum values at compile-time:
  constexpr auto minPixelComponent = NumericTraits<PixelComponentType>::NonpositiveMin();
  constexpr auto maxPixelComponent = NumericTraits<PixelComponentType>::max();
  constexpr auto minComponent = static_cast<ComponentType>(minPixelComponent);
  constexpr auto maxComponent = static_cast<ComponentType>(maxPixelComponent);

  // Clamp the value between minPixelComponent and maxPixelComponent:
  return (value <= minComponent) ? minPixelComponent
                                 : (value >= maxComponent) ? maxPixelComponent : static_cast<PixelComponentType>(value);
}


template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
auto
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::
  CastPixelWithBoundsChecking(const ComponentType value) -> PixelType
{
  return CastComponentWithBoundsChecking(value);
}


template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
template <typename TPixel>
auto
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::
  CastPixelWithBoundsChecking(const TPixel value) -> PixelType
{
  static_assert(std::is_same_v<TPixel, InterpolatorOutputType>,
                "TPixel should just be the same as the InterpolatorOutputType!");
  static_assert(!std::is_same_v<TPixel, ComponentType>,
                "For ComponentType there is a more efficient overload, that should be called instead!");

  const unsigned int nComponents = InterpolatorConvertType::GetNumberOfComponents(value);
  PixelType          outputValue;

  NumericTraits<PixelType>::SetLength(outputValue, nComponents);

  for (unsigned int n = 0; n < nComponents; ++n)
  {
    const ComponentType component = InterpolatorConvertType::GetNthComponent(n, value);
    PixelConvertType::SetNthComponent(n, outputValue, Self::CastComponentWithBoundsChecking(component));
  }

  return outputValue;
}


template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
void
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::
  NonlinearThreadedGenerateData(const OutputImageRegionType & outputRegionForThread)
{
  OutputImageType *      outputPtr = this->GetOutput();
  const InputImageType * inputPtr = this->GetInput();
  const TransformType *  transformPtr = this->GetTransform();

  TotalProgressReporter progress(this, outputPtr->GetRequestedRegion().GetNumberOfPixels());

  // Honor the SpecialCoordinatesImage isInside value returned
  // by TransformPhysicalPointToContinuousIndex
  using InputSpecialCoordinatesImageType = SpecialCoordinatesImage<InputPixelType, InputImageDimension>;
  const bool isSpecialCoordinatesImage = (dynamic_cast<const InputSpecialCoordinatesImageType *>(inputPtr) != nullptr);


  // Create an iterator that will walk the output region for this thread.
  using OutputIterator = ImageRegionIteratorWithIndex<TOutputImage>;

  using OutputType = typename InterpolatorType::OutputType;

  // Walk the output region
  for (OutputIterator outIt(outputPtr, outputRegionForThread); !outIt.IsAtEnd(); ++outIt)
  {
    // Determine the index of the current output pixel

    OutputPointType outputPoint; // Coordinates of current output pixel
    outputPtr->TransformIndexToPhysicalPoint(outIt.GetIndex(), outputPoint);

    // Compute corresponding input pixel position
    const InputPointType inputPoint = transformPtr->TransformPoint(outputPoint);

    ContinuousInputIndexType inputIndex;
    const bool               isInsideInput = inputPtr->TransformPhysicalPointToContinuousIndex(inputPoint, inputIndex);

    OutputType value;
    // Evaluate input at right position and copy to the output
    if (m_Interpolator->IsInsideBuffer(inputIndex) && (!isSpecialCoordinatesImage || isInsideInput))
    {
      value = m_Interpolator->EvaluateAtContinuousIndex(inputIndex);
      outIt.Set(Self::CastPixelWithBoundsChecking(value));
    }
    else
    {
      if (m_Extrapolator.IsNull())
      {
        outIt.Set(m_DefaultPixelValue); // default background value
      }
      else
      {
        value = m_Extrapolator->EvaluateAtContinuousIndex(inputIndex);
        outIt.Set(Self::CastPixelWithBoundsChecking(value));
      }
    }
    progress.CompletedPixel();
  }
}

template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
void
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::
  LinearThreadedGenerateData(const OutputImageRegionType & outputRegionForThread)
{
  OutputImageType *      outputPtr = this->GetOutput();
  const InputImageType * inputPtr = this->GetInput();
  const TransformType *  transformPtr = this->GetTransform();

  TotalProgressReporter progress(this, outputPtr->GetRequestedRegion().GetNumberOfPixels());

  const OutputImageRegionType & largestPossibleRegion = outputPtr->GetLargestPossibleRegion();

  const auto firstIndexValueOfLargestPossibleRegion = largestPossibleRegion.GetIndex(0);
  const auto firstSizeValueOfLargestPossibleRegion = static_cast<double>(largestPossibleRegion.GetSize(0));

  // Cache information from the superclass
  PixelType defaultValue = this->GetDefaultPixelValue();

  // As we walk across a scan line in the output image, we trace
  // an oriented/scaled/translated line in the input image. Each scan
  // line has a starting and ending point. Since all transforms
  // are linear, the path between the points is linear and can be
  // defined by interpolating between the two points. By using
  // interpolation we avoid accumulation errors, and by using the
  // whole scan line from the largest possible region we make the
  // computation independent for each point and independent of the
  // region we are processing which makes the method independent of
  // how the whole image is split for processing ( threading,
  // streaming, etc ).
  //

  const auto transformIndex = [outputPtr, transformPtr, inputPtr](const IndexType & index) {
    return inputPtr->template TransformPhysicalPointToContinuousIndex<TInterpolatorPrecisionType>(
      transformPtr->TransformPoint(outputPtr->template TransformIndexToPhysicalPoint<double>(index)));
  };

  // Create an iterator that will walk the output region for this thread.
  for (ImageScanlineIterator outIt(outputPtr, outputRegionForThread); !outIt.IsAtEnd(); outIt.NextLine())
  {
    // Determine the continuous index of the first and end pixel of output
    // scan line when mapped to the input coordinate frame.

    IndexType index = outIt.GetIndex();
    index[0] = firstIndexValueOfLargestPossibleRegion;

    const ContinuousInputIndexType startIndex = transformIndex(index);
    index[0] += firstSizeValueOfLargestPossibleRegion;
    const auto vectorFromStartIndex = transformIndex(index) - startIndex;

    IndexValueType scanlineIndex = outIt.GetIndex()[0];


    while (!outIt.IsAtEndOfLine())
    {

      // Perform linear interpolation from startIndex, along vectorFromStartIndex
      const double alpha =
        (scanlineIndex - firstIndexValueOfLargestPossibleRegion) / firstSizeValueOfLargestPossibleRegion;

      ContinuousInputIndexType inputIndex(startIndex);
      for (unsigned int i = 0; i < InputImageDimension; ++i)
      {
        inputIndex[i] += alpha * vectorFromStartIndex[i];
      }

      // Evaluate input at right position and copy to the output
      if (m_Interpolator->IsInsideBuffer(inputIndex))
      {
        outIt.Set(Self::CastPixelWithBoundsChecking(m_Interpolator->EvaluateAtContinuousIndex(inputIndex)));
      }
      else
      {
        if (m_Extrapolator.IsNull())
        {
          outIt.Set(defaultValue); // default background value
        }
        else
        {
          outIt.Set(Self::CastPixelWithBoundsChecking(m_Extrapolator->EvaluateAtContinuousIndex(inputIndex)));
        }
      }

      ++outIt;
      ++scanlineIndex;
    }
    progress.Completed(outputRegionForThread.GetSize()[0]);
  }
}

template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
void
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::
  GenerateInputRequestedRegion()
{
  if (!m_Interpolator)
  {
    itkExceptionMacro("Interpolator not set");
  }

  // Get pointers to the input and output
  auto * input = const_cast<InputImageType *>(this->GetInput());

  // Some interpolators need to look at their images in GetRadius()
  m_Interpolator->SetInputImage(input);

  // Check whether the input or the output is a
  // SpecialCoordinatesImage. If either are, then we cannot use the
  // fast path since index mapping will definitely not be linear.
  using OutputSpecialCoordinatesImageType = SpecialCoordinatesImage<PixelType, OutputImageDimension>;
  using InputSpecialCoordinatesImageType = SpecialCoordinatesImage<InputPixelType, InputImageDimension>;

  const bool isSpecialCoordinatesImage =
    ((dynamic_cast<const InputSpecialCoordinatesImageType *>(this->GetInput()) != nullptr) ||
     (dynamic_cast<const OutputSpecialCoordinatesImageType *>(this->GetOutput()) != nullptr));

  const OutputImageType * output = this->GetOutput();
  // Get the input transform
  const TransformType * transform = this->GetTransform();

  // Check whether we can use upstream streaming for resampling. Upstream streaming
  // can be used if the transformation is linear. Transform respond
  // to the IsLinear() call.
  if (!isSpecialCoordinatesImage && transform->GetTransformCategory() == TransformType::TransformCategoryEnum::Linear)
  {
    typename TInputImage::RegionType inputRequestedRegion;
    inputRequestedRegion = ImageAlgorithm::EnlargeRegionOverBox(output->GetRequestedRegion(), output, input, transform);

    const typename TInputImage::RegionType inputLargestRegion(input->GetLargestPossibleRegion());
    if (inputLargestRegion.IsInside(inputRequestedRegion.GetIndex()) ||
        inputLargestRegion.IsInside(inputRequestedRegion.GetUpperIndex()))
    {
      // Input requested region is partially outside the largest possible region.
      //   or
      // Input requested region is completely inside the largest possible region.
      inputRequestedRegion.PadByRadius(m_Interpolator->GetRadius());
      inputRequestedRegion.Crop(inputLargestRegion);
      input->SetRequestedRegion(inputRequestedRegion);
    }
    else if (inputRequestedRegion.IsInside(inputLargestRegion))
    {
      // Input requested region completely surrounds the largest possible region.
      input->SetRequestedRegion(inputLargestRegion);
    }
    else
    {
      // Input requested region is completely outside the largest possible region. Do not set the requested region in
      // this case.
    }
    return;
  }

  // Otherwise, determining the actual input region is non-trivial, especially
  // when we cannot assume anything about the transform being used.
  // So we do the easy thing and request the entire input image.
  //
  input->SetRequestedRegionToLargestPossibleRegion();
}

template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
void
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::
  GenerateOutputInformation()
{
  // Call the superclass' implementation of this method
  if (InputImageDimension == OutputImageDimension)
  {
    Superclass::GenerateOutputInformation();
  }
  // Get pointers to the input and output
  OutputImageType * outputPtr = this->GetOutput();

  const ReferenceImageBaseType * const referenceImage = this->GetReferenceImage();

  // Set the size of the output region
  if (m_UseReferenceImage && referenceImage)
  {
    outputPtr->SetLargestPossibleRegion(referenceImage->GetLargestPossibleRegion());
  }
  else
  {
    const typename TOutputImage::RegionType outputLargestPossibleRegion(m_OutputStartIndex, m_Size);
    outputPtr->SetLargestPossibleRegion(outputLargestPossibleRegion);
  }

  // Set spacing and origin
  if (m_UseReferenceImage && referenceImage)
  {
    outputPtr->SetSpacing(referenceImage->GetSpacing());
    outputPtr->SetOrigin(referenceImage->GetOrigin());
    outputPtr->SetDirection(referenceImage->GetDirection());
  }
  else
  {
    outputPtr->SetSpacing(m_OutputSpacing);
    outputPtr->SetOrigin(m_OutputOrigin);
    outputPtr->SetDirection(m_OutputDirection);
  }
}

template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
ModifiedTimeType
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::GetMTime() const
{
  ModifiedTimeType latestTime = Object::GetMTime();

  if (m_Interpolator)
  {
    latestTime = std::max(latestTime, m_Interpolator->GetMTime());
  }

  return latestTime;
}

template <typename TInputImage,
          typename TOutputImage,
          typename TInterpolatorPrecisionType,
          typename TTransformPrecisionType>
void
ResampleImageFilter<TInputImage, TOutputImage, TInterpolatorPrecisionType, TTransformPrecisionType>::PrintSelf(
  std::ostream & os,
  Indent         indent) const
{
  Superclass::PrintSelf(os, indent);

  os << indent
     << "DefaultPixelValue: " << static_cast<typename NumericTraits<PixelType>::PrintType>(m_DefaultPixelValue)
     << std::endl;
  os << indent << "Size: " << m_Size << std::endl;
  os << indent << "OutputStartIndex: " << m_OutputStartIndex << std::endl;
  os << indent << "OutputSpacing: " << m_OutputSpacing << std::endl;
  os << indent << "OutputOrigin: " << m_OutputOrigin << std::endl;
  os << indent << "OutputDirection: " << m_OutputDirection << std::endl;
  os << indent << "Transform: " << this->GetTransform() << std::endl;
  os << indent << "Interpolator: " << m_Interpolator.GetPointer() << std::endl;
  os << indent << "Extrapolator: " << m_Extrapolator.GetPointer() << std::endl;
  os << indent << "UseReferenceImage: " << (m_UseReferenceImage ? "On" : "Off") << std::endl;
}
} // end namespace itk

#endif