ITK  4.3.0
Insight Segmentation and Registration Toolkit
itkTestingExtractSliceImageFilter.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 __itkTestingExtractSliceImageFilter_h
19 #define __itkTestingExtractSliceImageFilter_h
20 
21 #include "itkSmartPointer.h"
23 
24 namespace itk
25 {
26 namespace Testing
27 {
81 template< class TInputImage, class TOutputImage >
82 class ITK_EXPORT ExtractSliceImageFilter:
83  public ImageSource< TOutputImage >
84 {
85 public:
91 
93  itkNewMacro(Self);
94 
97 
99  typedef TInputImage InputImageType;
100  typedef TOutputImage OutputImageType;
101 
103  typedef typename TOutputImage::RegionType OutputImageRegionType;
104  typedef typename TInputImage::RegionType InputImageRegionType;
105 
107  typedef typename TOutputImage::PixelType OutputImagePixelType;
108  typedef typename TInputImage::PixelType InputImagePixelType;
109 
111  typedef typename TOutputImage::IndexType OutputImageIndexType;
112  typedef typename TInputImage::IndexType InputImageIndexType;
113  typedef typename TOutputImage::SizeType OutputImageSizeType;
114  typedef typename TInputImage::SizeType InputImageSizeType;
115 
117  DIRECTIONCOLLAPSETOUNKOWN=0,
118  DIRECTIONCOLLAPSETOIDENTITY=1,
119  DIRECTIONCOLLAPSETOSUBMATRIX=2,
120  DIRECTIONCOLLAPSETOGUESS=3
121  } DIRECTIONCOLLAPSESTRATEGY;
122 
123 
148  void SetDirectionCollapseToStrategy(const DIRECTIONCOLLAPSESTRATEGY choosenStrategy)
149  {
150  switch(choosenStrategy)
151  {
152  case DIRECTIONCOLLAPSETOGUESS:
153  case DIRECTIONCOLLAPSETOIDENTITY:
154  case DIRECTIONCOLLAPSETOSUBMATRIX:
155  break;
156  case DIRECTIONCOLLAPSETOUNKOWN:
157  default:
158  itkExceptionMacro( << "Invalid Strategy Chosen for itk::ExtractSliceImageFilter" );
159  }
161 
162  this->m_DirectionCollaspeStrategy=choosenStrategy;
163  this->Modified();
164  }
165 
174  DIRECTIONCOLLAPSESTRATEGY GetDirectionCollapseToStrategy() const
175  {
176  return this->m_DirectionCollaspeStrategy;
177  }
178 
180  void SetDirectionCollapseToGuess()
181  {
182  this->SetDirectionCollapseToStrategy(DIRECTIONCOLLAPSETOGUESS);
183  }
184 
186  void SetDirectionCollapseToIdentity()
187  {
188  this->SetDirectionCollapseToStrategy(DIRECTIONCOLLAPSETOIDENTITY);
189  }
190 
192  void SetDirectionCollapseToSubmatrix()
193  {
194  this->SetDirectionCollapseToStrategy(DIRECTIONCOLLAPSETOSUBMATRIX);
195  }
196 
197 
199  itkStaticConstMacro(InputImageDimension, unsigned int,
200  TInputImage::ImageDimension);
201  itkStaticConstMacro(OutputImageDimension, unsigned int,
202  TOutputImage::ImageDimension);
204 
206  itkGetStaticConstMacro(InputImageDimension),
207  itkGetStaticConstMacro(OutputImageDimension) > ExtractSliceImageFilterRegionCopierType;
208 
214  void SetExtractionRegion(InputImageRegionType extractRegion);
215  itkGetConstMacro(ExtractionRegion, InputImageRegionType);
217 
219  using Superclass::SetInput;
220  virtual void SetInput(const TInputImage *image);
221  const TInputImage * GetInput(void) const;
223 
224 #ifdef ITK_USE_CONCEPT_CHECKING
225 
226  itkConceptMacro( InputCovertibleToOutputCheck,
228 
230 #endif
231 
232 protected:
235  void PrintSelf(std::ostream & os, Indent indent) const;
236 
245  virtual void GenerateOutputInformation();
246 
256  virtual void CallCopyOutputRegionToInputRegion(InputImageRegionType & destRegion,
257  const OutputImageRegionType & srcRegion);
258 
267  void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
268  ThreadIdType threadId);
269 
271 
273 
274 private:
275  ExtractSliceImageFilter(const Self &); //purposely not implemented
276  void operator=(const Self &); //purposely not implemented
277 
279 };
280 } // end namespace Testing
281 } // end namespace itk
282 
283 #ifndef ITK_MANUAL_INSTANTIATION
284 #include "itkTestingExtractSliceImageFilter.hxx"
285 #endif
286 
287 #endif
288