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
105
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
129 using OutputImageRegionType = typename TOutputImage::RegionType;
130 using InputImageRegionType = typename TInputImage::RegionType;
131
133 using OutputImagePixelType = typename TOutputImage::PixelType;
134 using InputImagePixelType = typename TInputImage::PixelType;
135
137 using OutputImageIndexType = typename TOutputImage::IndexType;
138 using InputImageIndexType = typename TInputImage::IndexType;
139 using OutputImageSizeType = typename TOutputImage::SizeType;
140 using InputImageSizeType = typename TInputImage::SizeType;
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 }
198
199 this->m_DirectionCollaspeStrategy = choosenStrategy;
200 this->Modified();
201 }
202
207
211 DIRECTIONCOLLAPSESTRATEGY
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);
252
255 virtual void
256 SetInput(const TInputImage * input);
257 const TInputImage *
258 GetInput() const;
260
262
263protected:
265 ~ExtractSliceImageFilter() override = default;
266 void
267 PrintSelf(std::ostream & os, Indent indent) const override;
268
277 void
279
289 virtual void
291
300 void
301 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
302
304
306
307private:
309 TestExtractSliceImageFilterCollapseStrategyEnum::DIRECTIONCOLLAPSETOUNKOWN
310 };
311};
312} // end namespace Testing
313} // end namespace itk
314
315#ifndef ITK_MANUAL_INSTANTIATION
316# include "itkTestingExtractSliceImageFilter.hxx"
317#endif
318
319#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....