ITK  5.4.0
Insight Toolkit
itkTestingExtractSliceImageFilter.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 * https://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"
22#include "itkImageSource.h"
24
25namespace itk
26{
27namespace Testing
28{
34{
35public:
41 {
42 DIRECTIONCOLLAPSETOUNKOWN = 0,
43 DIRECTIONCOLLAPSETOIDENTITY = 1,
44 DIRECTIONCOLLAPSETOSUBMATRIX = 2,
45 DIRECTIONCOLLAPSETOGUESS = 3
46 };
47};
48// Define how to print enumeration
49extern std::ostream &
51
107template <typename TInputImage, typename TOutputImage>
108class ITK_TEMPLATE_EXPORT ExtractSliceImageFilter : public ImageSource<TOutputImage>
109{
110public:
111 ITK_DISALLOW_COPY_AND_MOVE(ExtractSliceImageFilter);
112
118
120 itkNewMacro(Self);
121
123 itkOverrideGetNameOfClassMacro(ExtractSliceImageFilter);
124
126 using InputImageType = TInputImage;
127 using OutputImageType = TOutputImage;
128
132
134 using OutputImagePixelType = typename TOutputImage::PixelType;
135 using InputImagePixelType = typename TInputImage::PixelType;
136
142
148#if !defined(ITK_LEGACY_REMOVE)
149 // We need to expose the enum values at the class level
150 // for backwards compatibility
151 static constexpr DIRECTIONCOLLAPSESTRATEGY DIRECTIONCOLLAPSETOUNKOWN =
152 DIRECTIONCOLLAPSESTRATEGY::DIRECTIONCOLLAPSETOUNKOWN;
153 static constexpr DIRECTIONCOLLAPSESTRATEGY DIRECTIONCOLLAPSETOIDENTITY =
154 DIRECTIONCOLLAPSESTRATEGY::DIRECTIONCOLLAPSETOIDENTITY;
155 static constexpr DIRECTIONCOLLAPSESTRATEGY DIRECTIONCOLLAPSETOSUBMATRIX =
156 DIRECTIONCOLLAPSESTRATEGY::DIRECTIONCOLLAPSETOSUBMATRIX;
157 static constexpr DIRECTIONCOLLAPSESTRATEGY DIRECTIONCOLLAPSETOGUESS =
158 DIRECTIONCOLLAPSESTRATEGY::DIRECTIONCOLLAPSETOGUESS;
159#endif
160
185 void
187 {
188 switch (choosenStrategy)
189 {
190 case TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOGUESS:
191 case TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOIDENTITY:
192 case TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOSUBMATRIX:
193 break;
194 case TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOUNKOWN:
195 default:
196 itkExceptionMacro("Invalid Strategy Chosen for itk::ExtractSliceImageFilter");
197 }
200 this->m_DirectionCollaspeStrategy = choosenStrategy;
201 this->Modified();
202 }
203
212 DIRECTIONCOLLAPSESTRATEGY
213 GetDirectionCollapseToStrategy() const { return this->m_DirectionCollaspeStrategy; }
214
216 void
218 {
219 this->SetDirectionCollapseToStrategy(TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOGUESS);
220 }
221
223 void
225 {
226 this->SetDirectionCollapseToStrategy(TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOIDENTITY);
227 }
228
230 void
232 {
233 this->SetDirectionCollapseToStrategy(TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOSUBMATRIX);
234 }
235
236
238 static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
239 static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
240
243
249 void
251 itkGetConstMacro(ExtractionRegion, InputImageRegionType);
255 using Superclass::SetInput;
256 virtual void
257 SetInput(const TInputImage * input);
258 const TInputImage *
259 GetInput() const;
262#ifdef ITK_USE_CONCEPT_CHECKING
263 // Begin concept checking
265 // End concept checking
266#endif
267
268protected:
270 ~ExtractSliceImageFilter() override = default;
271 void
272 PrintSelf(std::ostream & os, Indent indent) const override;
273
282 void
284
294 virtual void
296
305 void
306 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
307
308 InputImageRegionType m_ExtractionRegion{};
309
310 OutputImageRegionType m_OutputImageRegion{};
311
312private:
313 DIRECTIONCOLLAPSESTRATEGY m_DirectionCollaspeStrategy{
314 TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOUNKOWN
315 };
316};
317} // end namespace Testing
318} // end namespace itk
319
320#ifndef ITK_MANUAL_INSTANTIATION
321# include "itkTestingExtractSliceImageFilter.hxx"
322#endif
323
324#endif
Base class for all process objects that output image data.
typename OutputImageType::PixelType OutputImagePixelType
typename OutputImageType::RegionType OutputImageRegionType
TOutputImage OutputImageType
A special variation of ImageRegionCopier for when the output image has fewer dimensions than the inpu...
Control indentation during Print() invocation.
Definition: itkIndent.h:50
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Contains all enum classes used by the ExtractSliceImageFilterEnums class.
Decrease the image size by cropping the image to the selected region bounds.
~ExtractSliceImageFilter() override=default
void SetExtractionRegion(InputImageRegionType extractRegion)
void PrintSelf(std::ostream &os, Indent indent) const override
void SetDirectionCollapseToStrategy(const DIRECTIONCOLLAPSESTRATEGY choosenStrategy)
virtual void CallCopyOutputRegionToInputRegion(InputImageRegionType &destRegion, const OutputImageRegionType &srcRegion)
virtual void SetInput(const TInputImage *input)
DIRECTIONCOLLAPSESTRATEGY GetDirectionCollapseToStrategy() const
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
const TInputImage * GetInput() const
#define itkConceptMacro(name, concept)
std::ostream & operator<<(std::ostream &out, const ExtractSliceImageFilterEnums::TestExtractSliceImageFilterCollapseStrategy value)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....