File: itkJoinSeriesImageFilter.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 (260 lines) | stat: -rw-r--r-- 8,979 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
/*=========================================================================
 *
 *  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 itkJoinSeriesImageFilter_hxx
#define itkJoinSeriesImageFilter_hxx

#include "itkProgressReporter.h"
#include "itkImageAlgorithm.h"
#include "itkTotalProgressReporter.h"

namespace itk
{
template <typename TInputImage, typename TOutputImage>
JoinSeriesImageFilter<TInputImage, TOutputImage>::JoinSeriesImageFilter()

{
  this->DynamicMultiThreadingOn();
  this->ThreaderUpdateProgressOff();
}

template <typename TInputImage, typename TOutputImage>
void
JoinSeriesImageFilter<TInputImage, TOutputImage>::VerifyInputInformation() ITKv5_CONST
{

  Superclass::VerifyInputInformation();

  // Verify that all input have the same number of components

  typename InputImageType::ConstPointer image = this->GetInput();

  if (image.IsNull())
  {
    itkExceptionMacro("Input not set as expected!");
  }

  const unsigned int numComponents = image->GetNumberOfComponentsPerPixel();

  typename Superclass::DataObjectPointerArraySizeType idx = 1;

  for (; idx < this->GetNumberOfIndexedInputs(); ++idx)
  {
    image = this->GetInput(idx);

    // If the input was not set it could still be null
    if (image.IsNull())
    {
      // An invalid requested region exception will be generated later.
      continue;
    }


    if (numComponents != image->GetNumberOfComponentsPerPixel())
    {
      itkExceptionMacro("Primary input has " << numComponents << " numberOfComponents "
                                             << "but input " << idx << " has " << image->GetNumberOfComponentsPerPixel()
                                             << '!');
    }
  }
}

template <typename TInputImage, typename TOutputImage>
void
JoinSeriesImageFilter<TInputImage, TOutputImage>::GenerateOutputInformation()
{
  // Do not call the superclass' implementation of this method since
  // this filter allows the input the output to be of different dimensions

  // Get pointers to the input and output
  typename Superclass::OutputImagePointer     outputPtr = this->GetOutput();
  typename Superclass::InputImageConstPointer inputPtr = this->GetInput();

  if (!outputPtr || !inputPtr)
  {
    return;
  }

  // Set the output image largest possible region.  Use a RegionCopier
  // so that the input and output images can be different dimensions.
  OutputImageRegionType outputLargestPossibleRegion;
  this->CallCopyInputRegionToOutputRegion(outputLargestPossibleRegion, inputPtr->GetLargestPossibleRegion());

  // for the new dimension, assuming the index has been set 0.
  outputLargestPossibleRegion.SetSize(static_cast<unsigned int>(InputImageDimension),
                                      static_cast<SizeValueType>(this->GetNumberOfIndexedInputs()));

  outputPtr->SetLargestPossibleRegion(outputLargestPossibleRegion);

  // Set the output spacing and origin
  const ImageBase<InputImageDimension> * phyData;

  phyData = dynamic_cast<const ImageBase<InputImageDimension> *>(this->GetInput());

  if (phyData)
  {
    // Copy what we can from the image from spacing and origin of the input
    // This logic needs to be augmented with logic that select which
    // dimensions to copy
    unsigned int                                 ii;
    const typename InputImageType::SpacingType & inputSpacing = inputPtr->GetSpacing();
    const typename InputImageType::PointType &   inputOrigin = inputPtr->GetOrigin();

    typename OutputImageType::SpacingType outputSpacing;
    typename OutputImageType::PointType   outputOrigin;

    // Copy the input to the output and fill the rest of the
    // output with zeros.
    for (ii = 0; ii < InputImageDimension; ++ii)
    {
      outputSpacing[ii] = inputSpacing[ii];
      outputOrigin[ii] = inputOrigin[ii];
    }
    for (; ii < OutputImageDimension; ++ii)
    {
      outputSpacing[ii] = 1.0;
      outputOrigin[ii] = 0.0;
    }

    // For the new dimension
    outputSpacing[InputImageDimension] = this->GetSpacing();
    outputOrigin[InputImageDimension] = this->GetOrigin();

    // Set the spacing and origin
    outputPtr->SetSpacing(outputSpacing);
    outputPtr->SetOrigin(outputOrigin);
    //
    // Copy the direction cosines from the input to the output.
    // On join, the output dim is always >= input dim
    using InputDirectionType = typename InputImageType::DirectionType;
    using OutputDirectionType = typename OutputImageType::DirectionType;
    InputDirectionType  inputDir = inputPtr->GetDirection();
    unsigned int        inputdim = InputImageType::GetImageDimension();
    unsigned int        outputdim = OutputImageType::GetImageDimension();
    OutputDirectionType outputDir = outputPtr->GetDirection();
    for (unsigned int i = 0; i < outputdim; ++i)
    {
      for (unsigned int j = 0; j < outputdim; ++j)
      {
        if (j < inputdim && i < inputdim)
        {
          outputDir[i][j] = inputDir[i][j];
        }
        else
        {
          outputDir[i][j] = i == j ? 1.0 : 0.0;
        }
      }
    }
    outputPtr->SetDirection(outputDir);
  }
  else
  {
    // Pointer could not be cast back down
    itkExceptionMacro("itk::JoinSeriesImageFilter::GenerateOutputInformation "
                      << "cannot cast input to " << typeid(ImageBase<InputImageDimension> *).name());
  }

  // Support VectorImages by setting number of components on output.
  const unsigned int numComponents = inputPtr->GetNumberOfComponentsPerPixel();
  if (numComponents != outputPtr->GetNumberOfComponentsPerPixel())
  {
    outputPtr->SetNumberOfComponentsPerPixel(numComponents);
  }
}

template <typename TInputImage, typename TOutputImage>
void
JoinSeriesImageFilter<TInputImage, TOutputImage>::GenerateInputRequestedRegion()
{
  Superclass::GenerateInputRequestedRegion();

  if (!this->GetOutput())
  {
    return;
  }
  OutputImageRegionType outputRegion = this->GetOutput()->GetRequestedRegion();
  IndexValueType        begin = outputRegion.GetIndex(InputImageDimension);
  IndexValueType        end = begin + outputRegion.GetSize(InputImageDimension);
  for (IndexValueType idx = 0; idx < this->GetNumberOfIndexedInputs(); ++idx)
  {
    InputImagePointer inputPtr = const_cast<InputImageType *>(this->GetInput(idx));
    if (!inputPtr)
    {
      // Because DataObject::PropagateRequestedRegion() allows only
      // InvalidRequestedRegionError, it's impossible to write simply:
      // itkExceptionMacro("Missing input " << idx);
      InvalidRequestedRegionError e(__FILE__, __LINE__);
      e.SetLocation(ITK_LOCATION);
      e.SetDescription("Missing input.");
      e.SetDataObject(this->GetOutput());
      throw e;
    }

    InputImageRegionType inputRegion;
    if (begin <= idx && idx < end)
    {
      this->CallCopyOutputRegionToInputRegion(inputRegion, outputRegion);
    }
    else
    {
      // Tell the pipeline that updating this input is unnecessary
      inputRegion = inputPtr->GetBufferedRegion();
    }
    inputPtr->SetRequestedRegion(inputRegion);
  }
}

template <typename TInputImage, typename TOutputImage>
void
JoinSeriesImageFilter<TInputImage, TOutputImage>::DynamicThreadedGenerateData(
  const OutputImageRegionType & outputRegionForThread)
{
  OutputImageRegionType outputRegion = outputRegionForThread;
  outputRegion.SetSize(InputImageDimension, 1);

  InputImageRegionType inputRegion;
  this->CallCopyOutputRegionToInputRegion(inputRegion, outputRegionForThread);

  const IndexValueType begin = outputRegionForThread.GetIndex(InputImageDimension);
  const IndexValueType end = begin + outputRegionForThread.GetSize(InputImageDimension);


  TOutputImage * output = this->GetOutput();

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

  for (IndexValueType idx = begin; idx < end; ++idx)
  {
    outputRegion.SetIndex(InputImageDimension, idx);
    ImageAlgorithm::Copy(this->GetInput(idx), output, inputRegion, outputRegion);
    progress.Completed(outputRegion.GetNumberOfPixels());
  }
}

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

  os << indent << "Spacing: " << m_Spacing << std::endl;
  os << indent << "Origin: " << m_Origin << std::endl;
}
} // end namespace itk

#endif