SimpleITK  
sitkProcessObject.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 * http://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 sitkProcessObject_h
19#define sitkProcessObject_h
20
21#include "sitkCommon.h"
22#include "sitkNonCopyable.h"
24#include "sitkEvent.h"
25#include "sitkImage.h"
26#include "sitkImageConvert.h"
27
28#include <iostream>
29#include <list>
30
31namespace itk
32{
33
34#ifndef SWIG
35
36template <typename T, unsigned int NVectorDimension>
37class Vector;
38
39class ProcessObject;
40class Command;
41class EventObject;
42#endif
43
44namespace simple
45{
46
47class Command;
48
49
55{
56public:
58
64
68 virtual ~ProcessObject();
69
70 // Print ourselves out
71 virtual std::string
72 ToString() const;
73
75 virtual std::string
76 GetName() const = 0;
77
84 virtual void
86 virtual void
89
93 virtual bool
94 GetDebug() const;
95 virtual void
96 SetDebug(bool debugFlag);
98
104 static void
106 static void
109
111 static bool
113 static void
114 SetGlobalDefaultDebug(bool debugFlag);
116
125 static void
127 static void
129 static void
131 static bool
134
147 static double
149 static void
151
152 static double
154 static void
157
158
183 static bool
184 SetGlobalDefaultThreader(const std::string & threader);
185 static std::string
188
202 static void
204 static unsigned int
207
222 virtual void
223 SetNumberOfThreads(unsigned int n);
224 virtual unsigned int
227
238 virtual void
239 SetNumberOfWorkUnits(unsigned int n);
240 virtual unsigned int
243
244
271 virtual int
273
274#ifndef SWIG
280 virtual int
281 AddCommand(itk::simple::EventEnum event, const std::function<void()> & func);
282#endif
283
289 virtual void
291
293 virtual bool
295
296
306 virtual float
307 GetProgress() const;
308
324 virtual void
326
327protected:
328#ifndef SWIG
329
331 {
335
336 // set to max if currently not registered
337 unsigned long m_ITKTag;
338
339 inline bool
340 operator==(const EventCommand & o) const
341 {
342 return m_Command == o.m_Command;
343 }
344 inline bool
345 operator<(const EventCommand & o) const
346 {
347 return m_Command < o.m_Command;
348 }
349 };
350
351 // method called before filter update to set parameters and
352 // connect commands.
353 virtual void
355
356 // overridable method to add a command, the return value is
357 // placed in the m_ITKTag of the EventCommand object.
358 virtual unsigned long
360
361 // overridable method to remove a command
362 virtual void
364
365 // Create an ITK EventObject from the SimpleITK enumerated type.
366 static const itk::EventObject &
368
369 // returns the current active process, if no active process then
370 // an exception is throw.
371 virtual itk::ProcessObject *
373
374 // overridable callback when the active process has completed
375 virtual void
377
379 // method call by command when it's deleted, maintains internal
380 // references between command and process objects.
381 virtual void
383#endif
384
385
386 template <class TImageType>
387 static typename TImageType::ConstPointer
389 {
390 typename TImageType::ConstPointer itkImage = dynamic_cast<const TImageType *>(img.GetITKBase());
391
392 if (itkImage.IsNull())
393 {
394 sitkExceptionMacro("Failure to convert SimpleITK image of dimension: "
395 << img.GetDimension() << " and pixel type: \"" << img.GetPixelIDTypeAsString()
396 << "\" to ITK image of dimension: " << TImageType::GetImageDimension() << " and pixel type: \""
398 }
399 return itkImage;
400 }
401
402 template <class TImageType>
403 static Image
404 CastITKToImage(TImageType * img)
405 {
406 return Image(img);
407 }
408
409#ifndef SWIG
410 template <class TPixelType,
411 unsigned int VImageDimension,
412 unsigned int VLength,
413 template <typename, unsigned int> class TVector>
414 static Image
415 CastITKToImage(itk::Image<TVector<TPixelType, VLength>, VImageDimension> * img)
416 {
417 // The implementation defined int sitkImageConvert.hxx needs
418 // to be manually included when this method is used.
419 auto out = GetVectorImageFromImage(img, true);
420
421 return Image(out.GetPointer());
422 }
423
424 template <unsigned int VImageDimension, unsigned int VLength, template <unsigned int> class TVector>
425 static Image
426 CastITKToImage(itk::Image<TVector<VLength>, VImageDimension> * img)
427 {
428 // The implementation defined int sitkImageConvert.hxx needs
429 // to be manually included when this method is used.
430 auto out = GetVectorImageFromImage(img, true);
431
432 return Image(out.GetPointer());
433 }
434
435#endif
436
444 template <typename T>
445 static std::ostream &
446 ToStringHelper(std::ostream & os, const T & v)
447 {
448 os << v;
449 return os;
450 }
451 static std::ostream &
452 ToStringHelper(std::ostream & os, const char & v);
453 static std::ostream &
454 ToStringHelper(std::ostream & os, const signed char & v);
455 static std::ostream &
456 ToStringHelper(std::ostream & os, const unsigned char & v);
458
459private:
460 // Add command to active process object, the EventCommand's
461 // ITKTag must be unset as max or else an exception is
462 // thrown. The EventCommand's ITKTag is updated to the command
463 // registered to ITK's ProcessObject. It's assumed that there is
464 // an current active process
465 unsigned long
467
468 // Remove the command from the active processes. Its is assumed
469 // that an active process exists. The tag is set to max after it
470 // is removed.
471 void
473
475
476 unsigned int m_NumberOfThreads;
478
479 std::list<EventCommand> m_Commands;
480
482
483 //
485};
486
487
488} // namespace simple
489} // namespace itk
490#endif
An implementation of the Command design pattern for callback.
Definition sitkCommand.h:45
The Image class for SimpleITK.
Definition sitkImage.h:77
unsigned int GetDimension() const
itk::DataObject * GetITKBase()
std::string GetPixelIDTypeAsString() const
unsigned long AddObserverToActiveProcessObject(EventCommand &e)
static void GlobalWarningDisplayOn()
virtual unsigned int GetNumberOfWorkUnits() const
static void SetGlobalDefaultDebug(bool debugFlag)
virtual int AddCommand(itk::simple::EventEnum event, itk::simple::Command &cmd)
Add a Command Object to observer the event.
static unsigned int GetGlobalDefaultNumberOfThreads()
Set/Get the default threader used for process objects.
virtual void PreUpdate(itk::ProcessObject *p)
static Image CastITKToImage(itk::Image< TVector< VLength >, VImageDimension > *img)
static void SetGlobalDefaultNumberOfThreads(unsigned int n)
virtual void RemoveITKObserver(EventCommand &e)
static std::ostream & ToStringHelper(std::ostream &os, const unsigned char &v)
virtual float GetProgress() const
An Active Measurement of the progress of execution.
static std::ostream & ToStringHelper(std::ostream &os, const char &v)
static void GlobalDefaultDebugOn()
static std::ostream & ToStringHelper(std::ostream &os, const signed char &v)
static double GetGlobalDefaultDirectionTolerance()
Access the global tolerance to determine congruent spaces.
virtual int AddCommand(itk::simple::EventEnum event, const std::function< void()> &func)
Directly add a callback to observe an event.
void RemoveObserverFromActiveProcessObject(EventCommand &e)
std::list< EventCommand > m_Commands
static void SetGlobalDefaultCoordinateTolerance(double)
Access the global tolerance to determine congruent spaces.
static TImageType::ConstPointer CastImageToITK(const Image &img)
static std::ostream & ToStringHelper(std::ostream &os, const T &v)
static bool SetGlobalDefaultThreader(const std::string &threader)
Set/Get the default threader used for process objects.
itk::ProcessObject * m_ActiveProcess
static bool GetGlobalDefaultDebug()
virtual bool HasCommand(itk::simple::EventEnum event) const
Query of this object has any registered commands for event.
virtual void OnActiveProcessDelete()
virtual void SetDebug(bool debugFlag)
static bool GetGlobalWarningDisplay()
virtual unsigned int GetNumberOfThreads() const
static void SetGlobalWarningDisplay(bool flag)
virtual void onCommandDelete(const itk::simple::Command *cmd) noexcept
virtual void SetNumberOfThreads(unsigned int n)
static void GlobalWarningDisplayOff()
virtual std::string ToString() const
virtual void SetNumberOfWorkUnits(unsigned int n)
friend class itk::simple::Command
virtual void RemoveAllCommands()
Remove all registered commands.
static const itk::EventObject & GetITKEventObject(EventEnum e)
virtual unsigned long AddITKObserver(const itk::EventObject &, itk::Command *)
virtual std::string GetName() const =0
static void SetGlobalDefaultDirectionTolerance(double)
Access the global tolerance to determine congruent spaces.
virtual bool GetDebug() const
static void GlobalDefaultDebugOff()
static Image CastITKToImage(TImageType *img)
virtual itk::ProcessObject * GetActiveProcess()
static double GetGlobalDefaultCoordinateTolerance()
Access the global tolerance to determine congruent spaces.
static std::string GetGlobalDefaultThreader()
Set/Get the default threader used for process objects.
static Image CastITKToImage(itk::Image< TVector< TPixelType, VLength >, VImageDimension > *img)
const std::string SITKCommon_EXPORT GetPixelIDValueAsString(PixelIDValueType type)
SITKCommon_HIDDEN itk::VectorImage< TPixelType, NImageDimension >::Pointer GetVectorImageFromImage(itk::Image< itk::Vector< TPixelType, NLength >, NImageDimension > *img, bool transferOwnership=false)
Utility methods to convert between itk image types efficiently by sharing the buffer between the inpu...
EventEnum
Events which can be observed from ProcessObject.
Definition sitkEvent.h:32
#define SITKCommon_EXPORT
Definition sitkCommon.h:41
#define sitkExceptionMacro(x)
Definition sitkMacro.h:65
bool operator==(const EventCommand &o) const
bool operator<(const EventCommand &o) const
EventCommand(EventEnum e, Command *c)