ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkJoinSeriesImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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 {
49 template< typename TInputImage, typename TOutputImage >
50 class ITK_TEMPLATE_EXPORT JoinSeriesImageFilter:
51  public ImageToImageFilter< TInputImage, TOutputImage >
52 {
53 public:
54  ITK_DISALLOW_COPY_AND_ASSIGN(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,
93  ( Concept::Convertible< typename TInputImage::PixelType,
94  typename TOutputImage::PixelType > ) );
95  // End concept checking
96 #endif
97 
98 protected:
100  ~JoinSeriesImageFilter() override = default;
101  void PrintSelf(std::ostream & os, Indent indent) const override;
102 
108  void VerifyInputInformation() ITKv5_CONST override;
109 
113  void GenerateOutputInformation() override;
114 
119  void GenerateInputRequestedRegion() override;
120 
124  void DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
125 
126 
127 private:
130  using IndexValueType = unsigned int;
131 
132  double m_Spacing{ 1.0 };
133  double m_Origin{ 0.0 };
134 };
135 } // end namespace itk
136 
137 #ifndef ITK_MANUAL_INSTANTIATION
138 #include "itkJoinSeriesImageFilter.hxx"
139 #endif
140 
141 #endif
typename OutputImageType::Pointer OutputImagePointer
Join N-D images into an (N+1)-D image.
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Base class for all process objects that output image data.
typename InputImageType::Pointer InputImagePointer
typename OutputImageType::RegionType OutputImageRegionType
TOutputImage OutputImageType
typename InputImageType::RegionType InputImageRegionType
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
#define itkConceptMacro(name, concept)