ITK  6.0.0
Insight Toolkit
itkPipelineMonitorImageFilter.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 itkPipelineMonitorImageFilter_h
19 #define itkPipelineMonitorImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 
23 namespace itk
24 {
25 
67 template <typename TImageType>
68 class ITK_TEMPLATE_EXPORT PipelineMonitorImageFilter : public ImageToImageFilter<TImageType, TImageType>
69 {
70 public:
71  ITK_DISALLOW_COPY_AND_MOVE(PipelineMonitorImageFilter);
72 
77 
78  using PointType = typename TImageType::PointType;
80  using SpacingType = typename TImageType::SpacingType;
83  using ImageRegionType = typename Superclass::InputImageRegionType;
84 
85  using RegionVectorType = std::vector<typename TImageType::RegionType>;
86 
88  itkNewMacro(Self);
89 
91  itkOverrideGetNameOfClassMacro(PipelineMonitorImageFilter);
92 
102  itkSetMacro(ClearPipelineOnGenerateOutputInformation, bool);
103  itkGetMacro(ClearPipelineOnGenerateOutputInformation, bool);
104  itkBooleanMacro(ClearPipelineOnGenerateOutputInformation);
113  bool
114  VerifyAllInputCanStream(int expectedNumber);
115 
116 
120  bool
121  VerifyAllInputCanNotStream();
122 
126  bool
127  VerifyAllNoUpdate();
128 
129  bool
130  VerifyDownStreamFilterExecutedPropagation();
131 
141  bool
142  VerifyInputFilterExecutedStreaming(int expectedNumber);
143 
148  bool
149  VerifyInputFilterMatchedUpdateOutputInformation();
150 
152  bool
153  VerifyInputFilterBufferedRequestedRegions();
154 
155  bool
156  VerifyInputFilterMatchedRequestedRegions();
157 
158  bool
159  VerifyInputFilterRequestedLargestRegion();
160 
161 
162  unsigned int
164  {
165  return m_NumberOfUpdates;
166  }
167  RegionVectorType
169  {
170  return m_OutputRequestedRegions;
171  }
172  RegionVectorType
174  {
175  return m_InputRequestedRegions;
176  }
177  RegionVectorType
179  {
180  return m_UpdatedBufferedRegions;
181  }
182  RegionVectorType
184  {
185  return m_UpdatedRequestedRegions;
186  }
187 
188  itkGetConstMacro(UpdatedOutputOrigin, PointType);
189  itkGetConstMacro(UpdatedOutputDirection, DirectionType);
190  itkGetConstMacro(UpdatedOutputSpacing, SpacingType);
191  itkGetConstMacro(UpdatedOutputLargestPossibleRegion, ImageRegionType);
192 
195  void
196  ClearPipelineSavedInformation();
197 
198 
202  void
203  GenerateOutputInformation() override;
204  void
205  PropagateRequestedRegion(DataObject * output) override;
206  void
207  EnlargeOutputRequestedRegion(DataObject * output) override;
208  void
209  GenerateInputRequestedRegion() override;
210  void
211  GenerateData() override;
214 protected:
216 
217  // ~PipelineMonitorImageFilter() { } default implementation OK
218 
219  void
220  PrintSelf(std::ostream & os, Indent indent) const override;
221 
222 private:
223  bool m_ClearPipelineOnGenerateOutputInformation{};
224 
225  unsigned int m_NumberOfUpdates{};
226 
227  unsigned int m_NumberOfClearPipeline{};
228 
229  RegionVectorType m_OutputRequestedRegions{};
230  RegionVectorType m_InputRequestedRegions{};
231  RegionVectorType m_UpdatedBufferedRegions{};
232  RegionVectorType m_UpdatedRequestedRegions{};
233 
234  PointType m_UpdatedOutputOrigin{};
235  DirectionType m_UpdatedOutputDirection{};
236  SpacingType m_UpdatedOutputSpacing{};
237  ImageRegionType m_UpdatedOutputLargestPossibleRegion{};
238 };
239 } // end namespace itk
240 
241 #ifndef ITK_MANUAL_INSTANTIATION
242 # include "itkPipelineMonitorImageFilter.hxx"
243 #endif
244 
245 #endif // itkPipelineMonitorImageFilter_hxx
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::PipelineMonitorImageFilter::InputImageConstPointer
typename TImageType::ConstPointer InputImageConstPointer
Definition: itkPipelineMonitorImageFilter.h:82
itk::PipelineMonitorImageFilter::GetNumberOfUpdates
unsigned int GetNumberOfUpdates() const
Definition: itkPipelineMonitorImageFilter.h:163
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::GTest::TypedefsAndConstructors::Dimension2::DirectionType
ImageBaseType::DirectionType DirectionType
Definition: itkGTestTypedefsAndConstructors.h:52
itk::PipelineMonitorImageFilter::GetUpdatedRequestedRegions
RegionVectorType GetUpdatedRequestedRegions() const
Definition: itkPipelineMonitorImageFilter.h:183
itk::PipelineMonitorImageFilter::GetUpdatedBufferedRegions
RegionVectorType GetUpdatedBufferedRegions() const
Definition: itkPipelineMonitorImageFilter.h:178
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::PipelineMonitorImageFilter
Enables monitoring, recording and debugging of the pipeline execution and information exchange.
Definition: itkPipelineMonitorImageFilter.h:68
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itk::PipelineMonitorImageFilter::RegionVectorType
std::vector< typename TImageType::RegionType > RegionVectorType
Definition: itkPipelineMonitorImageFilter.h:85
itk::PipelineMonitorImageFilter::PointType
typename TImageType::PointType PointType
Definition: itkPipelineMonitorImageFilter.h:78
itkImageToImageFilter.h
itk::PipelineMonitorImageFilter::GetOutputRequestedRegions
RegionVectorType GetOutputRequestedRegions() const
Definition: itkPipelineMonitorImageFilter.h:168
itk::PipelineMonitorImageFilter::GetInputRequestedRegions
RegionVectorType GetInputRequestedRegions() const
Definition: itkPipelineMonitorImageFilter.h:173
itk::PipelineMonitorImageFilter::DirectionType
typename TImageType::DirectionType DirectionType
Definition: itkPipelineMonitorImageFilter.h:79
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::PipelineMonitorImageFilter::SpacingType
typename TImageType::SpacingType SpacingType
Definition: itkPipelineMonitorImageFilter.h:80
itk::PipelineMonitorImageFilter::InputImagePointer
typename TImageType::Pointer InputImagePointer
Definition: itkPipelineMonitorImageFilter.h:81
itk::PipelineMonitorImageFilter::ImageRegionType
typename Superclass::InputImageRegionType ImageRegionType
Definition: itkPipelineMonitorImageFilter.h:83
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293