ITK  5.2.0
Insight Toolkit
itkJoinSeriesImageFilter.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  * 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 itkJoinSeriesImageFilter_h
19 #define itkJoinSeriesImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 
23 namespace itk
24 {
50 template <typename TInputImage, typename TOutputImage>
51 class ITK_TEMPLATE_EXPORT JoinSeriesImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
52 {
53 public:
54  ITK_DISALLOW_COPY_AND_MOVE(JoinSeriesImageFilter);
55 
61 
63  itkNewMacro(Self);
64 
67 
69  using InputImageType = typename Superclass::InputImageType;
70  using OutputImageType = typename Superclass::OutputImageType;
71  using InputImagePointer = typename InputImageType::Pointer;
72  using OutputImagePointer = typename OutputImageType::Pointer;
75 
77  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
78  static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
79 
81  itkSetMacro(Spacing, double);
82  itkGetConstMacro(Spacing, double);
84 
86  itkSetMacro(Origin, double);
87  itkGetConstMacro(Origin, double);
89 
90 #ifdef ITK_USE_CONCEPT_CHECKING
91  // Begin concept checking
92  itkConceptMacro(InputConvertibleToOutputCheck,
94  // End concept checking
95 #endif
96 
97 protected:
99  ~JoinSeriesImageFilter() override = default;
100  void
101  PrintSelf(std::ostream & os, Indent indent) const override;
102 
108  void
109  VerifyInputInformation() ITKv5_CONST override;
110 
114  void
115  GenerateOutputInformation() override;
116 
121  void
122  GenerateInputRequestedRegion() override;
123 
127  void
128  DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
129 
130 
131 private:
134  using IndexValueType = unsigned int;
135 
136  double m_Spacing{ 1.0 };
137  double m_Origin{ 0.0 };
138 };
139 } // end namespace itk
140 
141 #ifndef ITK_MANUAL_INSTANTIATION
142 # include "itkJoinSeriesImageFilter.hxx"
143 #endif
144 
145 #endif
itk::ImageSource::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkImageSource.h:91
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::JoinSeriesImageFilter::IndexValueType
unsigned int IndexValueType
Definition: itkJoinSeriesImageFilter.h:134
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::ImageToImageFilter::InputImagePointer
typename InputImageType::Pointer InputImagePointer
Definition: itkImageToImageFilter.h:130
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itkImageToImageFilter.h
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk::JoinSeriesImageFilter
Join N-D images into an (N+1)-D image.
Definition: itkJoinSeriesImageFilter.h:51
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:138
itk::Concept::Convertible
Definition: itkConceptChecking.h:216
itk::ImageToImageFilter::InputImageRegionType
typename InputImageType::RegionType InputImageRegionType
Definition: itkImageToImageFilter.h:132
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90