ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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{
47// Define how to print enumeration
48extern std::ostream &
50
104
105template <typename TInputImage, typename TOutputImage>
106class ITK_TEMPLATE_EXPORT ExtractSliceImageFilter : public ImageSource<TOutputImage>
107{
108public:
109 ITK_DISALLOW_COPY_AND_MOVE(ExtractSliceImageFilter);
110
116
118 itkNewMacro(Self);
119
121 itkOverrideGetNameOfClassMacro(ExtractSliceImageFilter);
122
124 using InputImageType = TInputImage;
125 using OutputImageType = TOutputImage;
126
128 using OutputImageRegionType = typename TOutputImage::RegionType;
129 using InputImageRegionType = typename TInputImage::RegionType;
130
132 using OutputImagePixelType = typename TOutputImage::PixelType;
133 using InputImagePixelType = typename TInputImage::PixelType;
134
136 using OutputImageIndexType = typename TOutputImage::IndexType;
137 using InputImageIndexType = typename TInputImage::IndexType;
138 using OutputImageSizeType = typename TOutputImage::SizeType;
139 using InputImageSizeType = typename TInputImage::SizeType;
140
146#if !defined(ITK_LEGACY_REMOVE)
147 // We need to expose the enum values at the class level
148 // for backwards compatibility
149 static constexpr DIRECTIONCOLLAPSESTRATEGY DIRECTIONCOLLAPSETOUNKOWN =
150 DIRECTIONCOLLAPSESTRATEGY::DIRECTIONCOLLAPSETOUNKOWN;
151 static constexpr DIRECTIONCOLLAPSESTRATEGY DIRECTIONCOLLAPSETOIDENTITY =
152 DIRECTIONCOLLAPSESTRATEGY::DIRECTIONCOLLAPSETOIDENTITY;
153 static constexpr DIRECTIONCOLLAPSESTRATEGY DIRECTIONCOLLAPSETOSUBMATRIX =
154 DIRECTIONCOLLAPSESTRATEGY::DIRECTIONCOLLAPSETOSUBMATRIX;
155 static constexpr DIRECTIONCOLLAPSESTRATEGY DIRECTIONCOLLAPSETOGUESS =
156 DIRECTIONCOLLAPSESTRATEGY::DIRECTIONCOLLAPSETOGUESS;
157#endif
158
183 void
185 {
186 switch (choosenStrategy)
187 {
188 case TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOGUESS:
189 case TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOIDENTITY:
190 case TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOSUBMATRIX:
191 break;
192 case TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOUNKOWN:
193 default:
194 itkExceptionMacro("Invalid Strategy Chosen for itk::ExtractSliceImageFilter");
195 }
196
197 this->m_DirectionCollaspeStrategy = choosenStrategy;
198 this->Modified();
199 }
200
205
209 DIRECTIONCOLLAPSESTRATEGY
211
213 void
215 {
216 this->SetDirectionCollapseToStrategy(TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOGUESS);
217 }
218
220 void
222 {
223 this->SetDirectionCollapseToStrategy(TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOIDENTITY);
224 }
225
227 void
229 {
230 this->SetDirectionCollapseToStrategy(TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOSUBMATRIX);
231 }
232
233
235 static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
236 static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
237
240
246 void
248 itkGetConstMacro(ExtractionRegion, InputImageRegionType);
249
252 virtual void
253 SetInput(const TInputImage * input);
254 const TInputImage *
255 GetInput() const;
256
258
259protected:
261 ~ExtractSliceImageFilter() override = default;
262 void
263 PrintSelf(std::ostream & os, Indent indent) const override;
264
273 void
275
285 virtual void
287
296 void
297 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
298
300
302
303private:
305 TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOUNKOWN
306 };
307};
308} // end namespace Testing
309} // end namespace itk
310
311#ifndef ITK_MANUAL_INSTANTIATION
312# include "itkTestingExtractSliceImageFilter.hxx"
313#endif
314
315#endif
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
virtual void Modified() const
virtual void SetInput(const DataObjectIdentifierType &key, DataObject *input)
Protected method for setting indexed and named inputs.
Implements transparent reference counting.
Contains all enum classes used by the ExtractSliceImageFilterEnums class.
~ExtractSliceImageFilter() override=default
ImageToImageFilterDetail::ExtractImageFilterRegionCopier< Self::InputImageDimension, Self::OutputImageDimension > ExtractSliceImageFilterRegionCopierType
void SetExtractionRegion(InputImageRegionType extractRegion)
TestExtractSliceImageFilterCollapseStrategyEnum DIRECTIONCOLLAPSESTRATEGY
void PrintSelf(std::ostream &os, Indent indent) const override
void SetDirectionCollapseToStrategy(const DIRECTIONCOLLAPSESTRATEGY choosenStrategy)
virtual void CallCopyOutputRegionToInputRegion(InputImageRegionType &destRegion, const OutputImageRegionType &srcRegion)
ExtractSliceImageFilterEnums::TestExtractSliceImageFilterCollapseStrategy TestExtractSliceImageFilterCollapseStrategyEnum
virtual void SetInput(const TInputImage *input)
DIRECTIONCOLLAPSESTRATEGY GetDirectionCollapseToStrategy() const
typename TOutputImage::RegionType OutputImageRegionType
TestExtractSliceImageFilterCollapseStrategyEnum DirectionCollaspeStrategyEnum
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....