ITK  5.4.0
Insight Toolkit
itkProcessObject.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 /*=========================================================================
19  *
20  * Portions of this file are subject to the VTK Toolkit Version 3 copyright.
21  *
22  * Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
23  *
24  * For complete copyright, license and disclaimer of warranty information
25  * please refer to the NOTICE file at the top of the ITK source tree.
26  *
27  *=========================================================================*/
28 #ifndef itkProcessObject_h
29 #define itkProcessObject_h
30 
31 #include "itkDataObject.h"
32 #include "itkObjectFactory.h"
33 #include "itkNumericTraits.h"
34 #include "itkThreadSupport.h"
35 #include "itkIntTypes.h"
36 #include <vector>
37 #include <map>
38 #include <set>
39 #include <algorithm>
40 #include <thread>
41 
42 namespace itk
43 {
44 
45 class MultiThreaderBase;
46 
139 class ITKCommon_EXPORT ProcessObject : public Object
140 {
141 public:
142  ITK_DISALLOW_COPY_AND_MOVE(ProcessObject);
143 
149 
151  itkOverrideGetNameOfClassMacro(ProcessObject);
152 
155 
157  // using ConstDataObjectPointerArray = std::vector< const DataObject * >;
158 
160  using DataObjectPointerArray = std::vector<DataObjectPointer>;
161 
163 
165  using NameArray = std::vector<DataObjectIdentifierType>;
166 
169 
175  NameArray
176  GetInputNames() const;
177 
179  NameArray
180  GetRequiredInputNames() const;
181 
188  GetInputs();
189 
191  bool
192  HasInput(const DataObjectIdentifierType & key) const;
193 
194  using DataObjectPointerArraySizeType = DataObjectPointerArray::size_type;
195 
207  GetNumberOfInputs() const;
208 
211  GetNumberOfOutputs() const;
212 
218  NameArray
219  GetOutputNames() const;
220 
227  GetOutputs();
228 
230  bool
231  HasOutput(const DataObjectIdentifierType & key) const;
232 
235  GetIndexedInputs();
236 
251  GetNumberOfIndexedInputs() const;
252 
263  GetNumberOfValidRequiredInputs() const;
264 
267  GetIndexedOutputs();
268 
271  GetNumberOfIndexedOutputs() const;
272 
288  virtual DataObjectPointer
289  MakeOutput(DataObjectPointerArraySizeType idx);
290 
296  itkSetMacro(AbortGenerateData, bool);
297 
299  itkGetConstReferenceMacro(AbortGenerateData, bool);
300 
302  itkBooleanMacro(AbortGenerateData);
303 
312 #if !defined(ITK_LEGACY_REMOVE)
313  void
314  SetProgress(float progress)
315  {
316  m_Progress = progressFloatToFixed(progress);
317  }
318 #endif
319 
325  virtual float
326  GetProgress() const
327  {
328  return progressFixedToFloat(m_Progress);
329  }
330 
340  void
341  UpdateProgress(float progress);
342 
350  void
351  IncrementProgress(float increment);
352 
374  virtual void
375  Update();
376 
386  virtual void
387  UpdateLargestPossibleRegion();
388 
404  virtual void
405  UpdateOutputInformation();
406 
409  virtual void
410  PropagateRequestedRegion(DataObject * output);
411 
413  virtual void
414  UpdateOutputData(DataObject * output);
415 
425  virtual void
427  {}
428 
435  virtual void
436  ResetPipeline();
437 
454  virtual DataObjectPointer
455  MakeOutput(const DataObjectIdentifierType &);
456 
462  virtual void
463  SetReleaseDataFlag(bool val);
464  virtual bool
465  GetReleaseDataFlag() const;
466  void
468  {
469  this->SetReleaseDataFlag(true);
470  }
471  void
473  {
474  this->SetReleaseDataFlag(false);
475  }
490  itkSetMacro(ReleaseDataBeforeUpdateFlag, bool);
491  itkGetConstReferenceMacro(ReleaseDataBeforeUpdateFlag, bool);
492  itkBooleanMacro(ReleaseDataBeforeUpdateFlag);
496  itkSetClampMacro(NumberOfWorkUnits, ThreadIdType, 1, ITK_MAX_THREADS);
497  itkGetConstReferenceMacro(NumberOfWorkUnits, ThreadIdType);
500 #if !defined(ITK_LEGACY_REMOVE) || defined(ITKV4_COMPATIBILITY)
501  itkLegacyMacro(void SetNumberOfThreads(ThreadIdType count)) { this->SetNumberOfWorkUnits(count); }
502 
503  itkLegacyMacro(ThreadIdType GetNumberOfThreads() const) { return this->GetNumberOfWorkUnits(); }
504 #endif // !ITK_LEGACY_REMOVE
505 
507  MultiThreaderType *
509  {
510  return m_MultiThreader;
511  }
512 
514  void
515  SetMultiThreader(MultiThreaderType * threader);
516 
523  virtual void
524  PrepareOutputs();
525 
526 protected:
527  ProcessObject();
528  ~ProcessObject() override;
529 
530  void
531  PrintSelf(std::ostream & os, Indent indent) const override;
532 
533  //
534  // Input Methods
535  //
536 
541  DataObject *
542  GetInput(const DataObjectIdentifierType & key);
543  const DataObject *
544  GetInput(const DataObjectIdentifierType & key) const;
545 
548  DataObject *
550  {
551  return idx < m_IndexedInputs.size() ? m_IndexedInputs[idx]->second.GetPointer() : nullptr;
552  }
553  const DataObject *
555  {
556  return idx < m_IndexedInputs.size() ? m_IndexedInputs[idx]->second.GetPointer() : nullptr;
557  }
566  virtual void
567  SetInput(const DataObjectIdentifierType & key, DataObject * input);
568  virtual void
569  SetNthInput(DataObjectPointerArraySizeType idx, DataObject * input);
570 
572  virtual void
573  AddInput(DataObject * input);
574 
583  virtual void
584  PushBackInput(const DataObject * input);
585  virtual void
586  PopBackInput();
587  virtual void
588  PushFrontInput(const DataObject * input);
589  virtual void
590  PopFrontInput();
591 
598  virtual void
599  RemoveInput(const DataObjectIdentifierType & key);
600  virtual void RemoveInput(DataObjectPointerArraySizeType);
601 
603  DataObject *
605  {
606  return m_IndexedInputs[0]->second;
607  }
608  const DataObject *
610  {
611  return m_IndexedInputs[0]->second;
612  }
616  virtual void
617  SetPrimaryInputName(const DataObjectIdentifierType & key);
618  virtual const char *
620  {
621  return this->m_IndexedInputs[0]->first.c_str();
622  }
626  virtual void
627  SetPrimaryInput(DataObject * object);
628 
635  void
636  SetNumberOfIndexedInputs(DataObjectPointerArraySizeType num);
637 
647  virtual void SetNumberOfRequiredInputs(DataObjectPointerArraySizeType);
648  itkGetConstReferenceMacro(NumberOfRequiredInputs, DataObjectPointerArraySizeType);
649 
650 
655  bool
656  RemoveRequiredInputName(const DataObjectIdentifierType &);
657 
659  bool
660  IsRequiredInputName(const DataObjectIdentifierType &) const;
661 
666  void
667  SetRequiredInputNames(const NameArray &);
668 
679  bool
680  AddRequiredInputName(const DataObjectIdentifierType &);
681  bool
682  AddRequiredInputName(const DataObjectIdentifierType &, DataObjectPointerArraySizeType idx);
683 
694  void
695  AddOptionalInputName(const DataObjectIdentifierType &);
696  void
697  AddOptionalInputName(const DataObjectIdentifierType &, DataObjectPointerArraySizeType idx);
698 
699 
700  //
701  // Output Methods
702  //
703 
705  DataObject *
706  GetOutput(const DataObjectIdentifierType & key);
707  const DataObject *
708  GetOutput(const DataObjectIdentifierType & key) const;
712  virtual void
713  SetPrimaryOutputName(const DataObjectIdentifierType & key);
714  virtual const char *
716  {
717  return this->m_IndexedOutputs[0]->first.c_str();
718  }
722  DataObject *
723  GetOutput(DataObjectPointerArraySizeType i);
724  const DataObject *
725  GetOutput(DataObjectPointerArraySizeType i) const;
729  virtual void
730  SetOutput(const DataObjectIdentifierType & name, DataObject * output);
731 
733  virtual void
734  RemoveOutput(const DataObjectIdentifierType & key);
735 
737  DataObject *
739  {
740  return m_IndexedOutputs[0]->second;
741  }
742  const DataObject *
744  {
745  return m_IndexedOutputs[0]->second;
746  }
750  virtual void
751  SetPrimaryOutput(DataObject * object);
752 
755  virtual void
756  SetNthOutput(DataObjectPointerArraySizeType idx, DataObject * output);
757 
758  virtual void
759  AddOutput(DataObject * output);
760 
761  virtual void
762  RemoveOutput(DataObjectPointerArraySizeType idx);
763 
764  itkSetMacro(NumberOfRequiredOutputs, DataObjectPointerArraySizeType);
765  itkGetConstReferenceMacro(NumberOfRequiredOutputs, DataObjectPointerArraySizeType);
766 
768  void
769  SetNumberOfIndexedOutputs(DataObjectPointerArraySizeType num);
770 
771 
772  DataObjectIdentifierType
773  MakeNameFromInputIndex(DataObjectPointerArraySizeType idx) const;
774  DataObjectIdentifierType
775  MakeNameFromOutputIndex(DataObjectPointerArraySizeType idx) const;
776  DataObjectPointerArraySizeType
777  MakeIndexFromInputName(const DataObjectIdentifierType & name) const;
778  DataObjectPointerArraySizeType
779  MakeIndexFromOutputName(const DataObjectIdentifierType & name) const;
780  bool
781  IsIndexedInputName(const DataObjectIdentifierType &) const;
782  bool
783  IsIndexedOutputName(const DataObjectIdentifierType &) const;
784 
785  //
786  // Pipeline Methods
787  //
788 
800  virtual void
801  VerifyPreconditions() ITKv5_CONST;
802 
813  virtual void
814  VerifyInputInformation() ITKv5_CONST;
815 
829  virtual void
830  GenerateInputRequestedRegion();
831 
843  virtual void
844  GenerateOutputRequestedRegion(DataObject * output);
845 
856  virtual void
857  GenerateOutputInformation();
858 
860  virtual void
861  GenerateData()
862  {}
863 
868  virtual void
869  PropagateResetPipeline();
870 
882  virtual void
883  ReleaseInputs();
884 
893  virtual void
894  CacheInputReleaseDataFlags();
895 
899  virtual void
900  RestoreInputReleaseDataFlags();
901 
906  itkGetConstMacro(ThreaderUpdateProgress, bool);
907  itkBooleanMacro(ThreaderUpdateProgress);
908  virtual void
909  SetThreaderUpdateProgress(bool arg);
916  static inline constexpr float
917  progressFixedToFloat(uint32_t fixed)
918  {
919  return static_cast<double>(fixed) / static_cast<double>(std::numeric_limits<uint32_t>::max());
920  };
927  static inline uint32_t
929  {
930  if (f <= 0.0f)
931  {
932  return 0;
933  }
934  if (f >= 1.0f)
935  {
936  return std::numeric_limits<uint32_t>::max();
937  }
938  double temp = static_cast<double>(f) * std::numeric_limits<uint32_t>::max();
939  return static_cast<uint32_t>(temp);
940  };
948  bool m_Updating{};
949 
951  TimeStamp m_OutputInformationMTime{};
952 
953 private:
954  DataObjectIdentifierType MakeNameFromIndex(DataObjectPointerArraySizeType) const;
955  DataObjectPointerArraySizeType
956  MakeIndexFromName(const DataObjectIdentifierType &) const;
957 
959  using DataObjectPointerMap = std::map<DataObjectIdentifierType, DataObjectPointer>;
960 
961 
964  DataObjectPointerMap m_Outputs{};
965 
966  std::vector<DataObjectPointerMap::iterator> m_IndexedInputs{};
967  std::vector<DataObjectPointerMap::iterator> m_IndexedOutputs{};
968 
970  std::map<DataObjectIdentifierType, bool> m_CachedInputReleaseDataFlags{};
971 
972  DataObjectPointerArraySizeType m_NumberOfRequiredInputs{};
973  DataObjectPointerArraySizeType m_NumberOfRequiredOutputs{};
974 
976  using NameSet = std::set<DataObjectIdentifierType>;
977 
979  NameSet m_RequiredInputNames{};
980 
982  bool m_AbortGenerateData{};
983  std::atomic<uint32_t> m_Progress{};
984 
985 
986  std::thread::id m_UpdateThreadID{};
987 
991  ThreadIdType m_NumberOfWorkUnits{};
992 
993  bool m_ThreaderUpdateProgress{ true };
994 
996  bool m_ReleaseDataBeforeUpdateFlag{};
997 
999  friend class DataObject;
1000 
1001  friend class ProgressReporter;
1003 
1007 
1008  friend class DataObjectIterator;
1011 
1012  friend class TestProcessObject;
1013 };
1014 } // end namespace itk
1015 
1016 #endif
itk::MultiThreaderBase
A class for performing multithreaded execution.
Definition: itkMultiThreaderBase.h:106
itkObjectFactory.h
itk::ProcessObject::DataObjectPointerArray
std::vector< DataObjectPointer > DataObjectPointerArray
Definition: itkProcessObject.h:160
itk::InputDataObjectConstIterator
A forward iterator over inputs of a ProcessObject.
Definition: itkInputDataObjectConstIterator.h:30
itk::ProcessObject::NameArray
std::vector< DataObjectIdentifierType > NameArray
Definition: itkProcessObject.h:165
itk::ProcessObject::GetPrimaryInput
const DataObject * GetPrimaryInput() const
Definition: itkProcessObject.h:609
itk::ProcessObject::DataObjectIdentifierType
DataObject::DataObjectIdentifierType DataObjectIdentifierType
Definition: itkProcessObject.h:162
itk::ProcessObject::GetInput
DataObject * GetInput(DataObjectPointerArraySizeType idx)
Definition: itkProcessObject.h:549
itk::ProcessObject::m_MultiThreader
itk::SmartPointer< MultiThreaderType > m_MultiThreader
Definition: itkProcessObject.h:990
itk::DataObjectIterator
A forward iterator over the DataObject of a ProcessObject.
Definition: itkDataObjectIterator.h:30
itk::ProcessObject::progressFloatToFixed
static uint32_t progressFloatToFixed(float f)
Definition: itkProcessObject.h:928
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::ProcessObject::DataObjectPointerArraySizeType
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
Definition: itkProcessObject.h:194
itk::ITK_MAX_THREADS
constexpr vcl_size_t ITK_MAX_THREADS
Definition: itkThreadSupport.h:79
itkDataObject.h
itk::ProcessObject::GetProgress
virtual float GetProgress() const
Get the execution progress of a process object.
Definition: itkProcessObject.h:326
itk::ProcessObject::DataObjectPointerMap
std::map< DataObjectIdentifierType, DataObjectPointer > DataObjectPointerMap
Definition: itkProcessObject.h:959
itk::ThreadIdType
unsigned int ThreadIdType
Definition: itkIntTypes.h:99
itk::ProcessObject::NameSet
std::set< DataObjectIdentifierType > NameSet
Definition: itkProcessObject.h:976
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itkThreadSupport.h
itk::InputDataObjectIterator
A forward iterator over inputs of a ProcessObject.
Definition: itkInputDataObjectIterator.h:30
itk::ProcessObject::ReleaseDataFlagOn
void ReleaseDataFlagOn()
Definition: itkProcessObject.h:467
itk::TimeStamp
Generate a unique, increasing time value.
Definition: itkTimeStamp.h:60
itk::ProcessObject::GetPrimaryOutput
DataObject * GetPrimaryOutput()
Definition: itkProcessObject.h:738
itkIntTypes.h
itk::DataObject::DataObjectIdentifierType
std::string DataObjectIdentifierType
Definition: itkDataObject.h:304
itk::DataObjectConstIterator
A forward iterator over the DataObject of a ProcessObject.
Definition: itkDataObjectConstIterator.h:30
itk::ProcessObject::GetInput
const DataObject * GetInput(DataObjectPointerArraySizeType idx) const
Definition: itkProcessObject.h:554
itk::ProcessObject::ReleaseDataFlagOff
void ReleaseDataFlagOff()
Definition: itkProcessObject.h:472
itk::OutputDataObjectConstIterator
A forward iterator over outputs of a ProcessObject.
Definition: itkOutputDataObjectConstIterator.h:30
itk::ProcessObject::GetPrimaryOutput
const DataObject * GetPrimaryOutput() const
Definition: itkProcessObject.h:743
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ProcessObject::EnlargeOutputRequestedRegion
virtual void EnlargeOutputRequestedRegion(DataObject *)
Definition: itkProcessObject.h:426
itk::TotalProgressReporter
A progress reporter for concurrent threads.
Definition: itkTotalProgressReporter.h:40
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::OutputDataObjectIterator
A forward iterator over outputs of a ProcessObject.
Definition: itkOutputDataObjectIterator.h:30
itk::ProcessObject
class ITK_FORWARD_EXPORT ProcessObject
Definition: itkDataObject.h:41
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:61
itkNumericTraits.h
itk::ProcessObject::GetMultiThreader
MultiThreaderType * GetMultiThreader() const
Definition: itkProcessObject.h:508
itk::ProgressReporter
Implements progress tracking for a filter.
Definition: itkProgressReporter.h:60
itk::ProcessObject::GetPrimaryInputName
virtual const char * GetPrimaryInputName() const
Definition: itkProcessObject.h:619
itk::ProcessObject::progressFixedToFloat
static constexpr float progressFixedToFloat(uint32_t fixed)
Definition: itkProcessObject.h:917
itk::DataObject::Pointer
SmartPointer< Self > Pointer
Definition: itkDataObject.h:301
itk::ProcessObject::GetPrimaryInput
DataObject * GetPrimaryInput()
Definition: itkProcessObject.h:604
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293
itk::ProcessObject::GetPrimaryOutputName
virtual const char * GetPrimaryOutputName() const
Definition: itkProcessObject.h:715