Loading [MathJax]/extensions/tex2jax.js
SimpleITK  2.5.0.dev
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
182 static bool
183 SetGlobalDefaultThreader(const std::string & threader);
184 static std::string
187
201 static void
203 static unsigned int
206
221 virtual void
222 SetNumberOfThreads(unsigned int n);
223 virtual unsigned int
226
237 virtual void
238 SetNumberOfWorkUnits(unsigned int n);
239 virtual unsigned int
242
243
270 virtual int
272
273#ifndef SWIG
279 virtual int
280 AddCommand(itk::simple::EventEnum event, const std::function<void()> & func);
281#endif
282
288 virtual void
290
292 virtual bool
294
295
305 virtual float
306 GetProgress() const;
307
323 virtual void
325
326protected:
327#ifndef SWIG
328
330 {
334
335 // set to max if currently not registered
336 unsigned long m_ITKTag;
337
338 inline bool
339 operator==(const EventCommand & o) const
340 {
341 return m_Command == o.m_Command;
342 }
343 inline bool
344 operator<(const EventCommand & o) const
345 {
346 return m_Command < o.m_Command;
347 }
348 };
349
350 // method called before filter update to set parameters and
351 // connect commands.
352 virtual void
354
355 // overridable method to add a command, the return value is
356 // placed in the m_ITKTag of the EventCommand object.
357 virtual unsigned long
359
360 // overridable method to remove a command
361 virtual void
363
364 // Create an ITK EventObject from the SimpleITK enumerated type.
365 static const itk::EventObject &
367
368 // returns the current active process, if no active process then
369 // an exception is throw.
370 virtual itk::ProcessObject *
372
373 // overridable callback when the active process has completed
374 virtual void
376
378 // method call by command when it's deleted, maintains internal
379 // references between command and process objects.
380 virtual void
382#endif
383
384
385 template <class TImageType>
386 static typename TImageType::ConstPointer
388 {
389 typename TImageType::ConstPointer itkImage = dynamic_cast<const TImageType *>(img.GetITKBase());
390
391 if (itkImage.IsNull())
392 {
393 sitkExceptionMacro("Failure to convert SimpleITK image of dimension: "
394 << img.GetDimension() << " and pixel type: \"" << img.GetPixelIDTypeAsString()
395 << "\" to ITK image of dimension: " << TImageType::GetImageDimension() << " and pixel type: \""
397 }
398 return itkImage;
399 }
400
401 template <class TImageType>
402 static Image
403 CastITKToImage(TImageType * img)
404 {
405 return Image(img);
406 }
407
408#ifndef SWIG
409 template <class TPixelType,
410 unsigned int VImageDimension,
411 unsigned int VLength,
412 template <typename, unsigned int> class TVector>
413 static Image
414 CastITKToImage(itk::Image<TVector<TPixelType, VLength>, VImageDimension> * img)
415 {
416 // The implementation defined int sitkImageConvert.hxx needs
417 // to be manually included when this method is used.
418 auto out = GetVectorImageFromImage(img, true);
419
420 return Image(out.GetPointer());
421 }
422
423 template <unsigned int VImageDimension, unsigned int VLength, template <unsigned int> class TVector>
424 static Image
425 CastITKToImage(itk::Image<TVector<VLength>, VImageDimension> * img)
426 {
427 // The implementation defined int sitkImageConvert.hxx needs
428 // to be manually included when this method is used.
429 auto out = GetVectorImageFromImage(img, true);
430
431 return Image(out.GetPointer());
432 }
433
434#endif
435
443 template <typename T>
444 static std::ostream &
445 ToStringHelper(std::ostream & os, const T & v)
446 {
447 os << v;
448 return os;
449 }
450 static std::ostream &
451 ToStringHelper(std::ostream & os, const char & v);
452 static std::ostream &
453 ToStringHelper(std::ostream & os, const signed char & v);
454 static std::ostream &
455 ToStringHelper(std::ostream & os, const unsigned char & v);
457
458private:
459 // Add command to active process object, the EventCommand's
460 // ITKTag must be unset as max or else an exception is
461 // thrown. The EventCommand's ITKTag is updated to the command
462 // registered to ITK's ProcessObject. It's assumed that there is
463 // an current active process
464 unsigned long
466
467 // Remove the command from the active processes. Its is assumed
468 // that an active process exists. The tag is set to max after it
469 // is removed.
470 void
472
474
475 unsigned int m_NumberOfThreads;
477
478 std::list<EventCommand> m_Commands;
479
481
482 //
484};
485
486
487} // namespace simple
488} // namespace itk
489#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:70
bool operator==(const EventCommand &o) const
bool operator<(const EventCommand &o) const
EventCommand(EventEnum e, Command *c)