ITK  6.0.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:
40 {
41 DIRECTIONCOLLAPSETOUNKOWN = 0,
42 DIRECTIONCOLLAPSETOIDENTITY = 1,
43 DIRECTIONCOLLAPSETOSUBMATRIX = 2,
44 DIRECTIONCOLLAPSETOGUESS = 3
45 };
46};
47// Define how to print enumeration
48extern std::ostream &
50
106template <typename TInputImage, typename TOutputImage>
107class ITK_TEMPLATE_EXPORT ExtractSliceImageFilter : public ImageSource<TOutputImage>
108{
109public:
110 ITK_DISALLOW_COPY_AND_MOVE(ExtractSliceImageFilter);
111
117
119 itkNewMacro(Self);
120
122 itkOverrideGetNameOfClassMacro(ExtractSliceImageFilter);
123
125 using InputImageType = TInputImage;
126 using OutputImageType = TOutputImage;
127
131
133 using OutputImagePixelType = typename TOutputImage::PixelType;
134 using InputImagePixelType = typename TInputImage::PixelType;
135
141
147#if !defined(ITK_LEGACY_REMOVE)
148 // We need to expose the enum values at the class level
149 // for backwards compatibility
150 static constexpr DIRECTIONCOLLAPSESTRATEGY DIRECTIONCOLLAPSETOUNKOWN =
151 DIRECTIONCOLLAPSESTRATEGY::DIRECTIONCOLLAPSETOUNKOWN;
152 static constexpr DIRECTIONCOLLAPSESTRATEGY DIRECTIONCOLLAPSETOIDENTITY =
153 DIRECTIONCOLLAPSESTRATEGY::DIRECTIONCOLLAPSETOIDENTITY;
154 static constexpr DIRECTIONCOLLAPSESTRATEGY DIRECTIONCOLLAPSETOSUBMATRIX =
155 DIRECTIONCOLLAPSESTRATEGY::DIRECTIONCOLLAPSETOSUBMATRIX;
156 static constexpr DIRECTIONCOLLAPSESTRATEGY DIRECTIONCOLLAPSETOGUESS =
157 DIRECTIONCOLLAPSESTRATEGY::DIRECTIONCOLLAPSETOGUESS;
158#endif
159
184 void
186 {
187 switch (choosenStrategy)
188 {
189 case TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOGUESS:
190 case TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOIDENTITY:
191 case TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOSUBMATRIX:
192 break;
193 case TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOUNKOWN:
194 default:
195 itkExceptionMacro("Invalid Strategy Chosen for itk::ExtractSliceImageFilter");
196 }
199 this->m_DirectionCollaspeStrategy = choosenStrategy;
200 this->Modified();
201 }
202
211 DIRECTIONCOLLAPSESTRATEGY
212 GetDirectionCollapseToStrategy() const { return this->m_DirectionCollaspeStrategy; }
213
215 void
217 {
218 this->SetDirectionCollapseToStrategy(TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOGUESS);
219 }
220
222 void
224 {
225 this->SetDirectionCollapseToStrategy(TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOIDENTITY);
226 }
227
229 void
231 {
232 this->SetDirectionCollapseToStrategy(TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOSUBMATRIX);
233 }
234
235
237 static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
238 static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
239
242
248 void
250 itkGetConstMacro(ExtractionRegion, InputImageRegionType);
254 using Superclass::SetInput;
255 virtual void
256 SetInput(const TInputImage * input);
257 const TInputImage *
258 GetInput() const;
261#ifdef ITK_USE_CONCEPT_CHECKING
262 // Begin concept checking
264 // End concept checking
265#endif
266
267protected:
269 ~ExtractSliceImageFilter() override = default;
270 void
271 PrintSelf(std::ostream & os, Indent indent) const override;
272
281 void
283
293 virtual void
295
304 void
305 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
306
307 InputImageRegionType m_ExtractionRegion{};
308
309 OutputImageRegionType m_OutputImageRegion{};
310
311private:
312 DIRECTIONCOLLAPSESTRATEGY m_DirectionCollaspeStrategy{
313 TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOUNKOWN
314 };
315};
316} // end namespace Testing
317} // end namespace itk
318
319#ifndef ITK_MANUAL_INSTANTIATION
320# include "itkTestingExtractSliceImageFilter.hxx"
321#endif
322
323#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....