ITK  4.6.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 >
51  public ImageToImageFilter< TInputImage, TOutputImage >
52 {
53 public:
59 
61  itkNewMacro(Self);
62 
65 
69  typedef typename InputImageType::Pointer InputImagePointer;
70  typedef typename OutputImageType::Pointer OutputImagePointer;
71  typedef typename InputImageType::RegionType InputImageRegionType;
72  typedef typename OutputImageType::RegionType OutputImageRegionType;
73 
75  itkStaticConstMacro(InputImageDimension, unsigned int,
76  TInputImage::ImageDimension);
77  itkStaticConstMacro(OutputImageDimension, unsigned int,
78  TOutputImage::ImageDimension);
80 
82  itkSetMacro(Spacing, double);
83  itkGetConstMacro(Spacing, double);
85 
87  itkSetMacro(Origin, double);
88  itkGetConstMacro(Origin, double);
90 
91 #ifdef ITK_USE_CONCEPT_CHECKING
92  // Begin concept checking
93  itkConceptMacro( InputConvertibleToOutputCheck,
94  ( Concept::Convertible< typename TInputImage::PixelType,
95  typename TOutputImage::PixelType > ) );
96  // End concept checking
97 #endif
98 
99 protected:
102  void PrintSelf(std::ostream & os, Indent indent) const;
103 
109  virtual void VerifyInputInformation();
110 
114  virtual void GenerateOutputInformation();
115 
120  virtual void GenerateInputRequestedRegion();
121 
125  virtual void ThreadedGenerateData(const OutputImageRegionType &
126  outputRegionForThread, ThreadIdType threadId);
127 
128 private:
129  JoinSeriesImageFilter(const Self &); //purposely not implemented
130  void operator=(const Self &); //purposely not implemented
131 
134  typedef unsigned int IndexValueType;
135 
136  double m_Spacing;
137  double m_Origin;
138 };
139 } // end namespace itk
140 
141 #ifndef ITK_MANUAL_INSTANTIATION
142 #include "itkJoinSeriesImageFilter.hxx"
143 #endif
144 
145 #endif
static const unsigned int OutputImageDimension
InputImageType::Pointer InputImagePointer
OutputImageType::Pointer OutputImagePointer
OutputImageType::RegionType OutputImageRegionType
Join N-D images into an (N+1)-D image.
Base class for all process objects that output image data.
Superclass::OutputImageType OutputImageType
InputImageType::RegionType InputImageRegionType
virtual void VerifyInputInformation()
virtual void GenerateInputRequestedRegion()
void operator=(const Self &)
void PrintSelf(std::ostream &os, Indent indent) const
ImageToImageFilter< TInputImage, TOutputImage > Superclass
Superclass::InputImageType InputImageType
static const unsigned int InputImageDimension
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
TOutputImage OutputImageType
virtual void GenerateOutputInformation()
#define itkConceptMacro(name, concept)
virtual void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
SmartPointer< const Self > ConstPointer