SimpleITK  
Public Types | Public Member Functions | Private Types | Private Member Functions | Private Attributes | Friends | List of all members
itk::simple::JoinSeriesImageFilter Class Reference

Join N-D images into an (N+1)-D image. More...

#include <sitkJoinSeriesImageFilter.h>

+ Inheritance diagram for itk::simple::JoinSeriesImageFilter:
+ Collaboration diagram for itk::simple::JoinSeriesImageFilter:

Detailed Description

Join N-D images into an (N+1)-D image.

This filter is templated over the input image type and the output image type. The pixel type of them must be the same and the input dimension must be less than the output dimension. When the input images are N-dimensional, they are joined in order and the size of the N+1'th dimension of the output is same as the number of the inputs. The spacing and the origin (where the first input is placed) for the N+1'th dimension is specified in this filter. The output image information for the first N dimensions are taken from the first input. Note that all the inputs should have the same information.

Author
Hideaki Hiraki

Contributed in the users list https://public.kitware.com/pipermail/insight-users/2004-February/006542.html

See also
itk::simple::JoinSeries for the procedural interface

Definition at line 47 of file sitkJoinSeriesImageFilter.h.

Public Types

using PixelIDTypeList = NonLabelPixelIDTypeList
 
using Self = JoinSeriesImageFilter
 
- Public Types inherited from itk::simple::ImageFilter
using Self = ImageFilter
 
- Public Types inherited from itk::simple::ProcessObject
using Self = ProcessObject
 

Public Member Functions

Image Execute (const Image &image1)
 
Image Execute (const Image &image1, const Image &image2)
 
Image Execute (const Image &image1, const Image &image2, const Image &image3)
 
Image Execute (const Image &image1, const Image &image2, const Image &image3, const Image &image4)
 
Image Execute (const Image &image1, const Image &image2, const Image &image3, const Image &image4, const Image &image5)
 
Image Execute (const std::vector< Image > &images)
 
std::string GetName () const
 
double GetOrigin () const
 
double GetSpacing () const
 
 JoinSeriesImageFilter ()
 
SelfSetOrigin (double Origin)
 
SelfSetSpacing (double Spacing)
 
std::string ToString () const
 
virtual ~JoinSeriesImageFilter ()
 
- Public Member Functions inherited from itk::simple::ImageFilter
 ImageFilter ()
 
virtual ~ImageFilter ()=0
 
- Public Member Functions inherited from itk::simple::ProcessObject
virtual void Abort ()
 
virtual int AddCommand (itk::simple::EventEnum event, const std::function< void()> &func)
 Directly add a callback to observe an event. More...
 
virtual int AddCommand (itk::simple::EventEnum event, itk::simple::Command &cmd)
 Add a Command Object to observer the event. More...
 
virtual float GetProgress () const
 An Active Measurement of the progress of execution. More...
 
virtual bool HasCommand (itk::simple::EventEnum event) const
 Query of this object has any registered commands for event. More...
 
 ProcessObject ()
 
virtual void RemoveAllCommands ()
 Remove all registered commands. More...
 
virtual ~ProcessObject ()
 
virtual void DebugOn ()
 
virtual void DebugOff ()
 
virtual bool GetDebug () const
 
virtual void SetDebug (bool debugFlag)
 
virtual void SetNumberOfThreads (unsigned int n)
 
virtual unsigned int GetNumberOfThreads () const
 
virtual void SetNumberOfWorkUnits (unsigned int n)
 
virtual unsigned int GetNumberOfWorkUnits () const
 

Private Types

using MemberFunctionType = Image(Self::*)(const std::vector< Image > &)
 

Private Member Functions

template<class TImageType >
Image ExecuteInternal (const std::vector< Image > &images)
 

Private Attributes

std::unique_ptr< detail::MemberFunctionFactory< MemberFunctionType > > m_MemberFactory
 
double m_Origin {0.0}
 
double m_Spacing {1.0}
 

Friends

struct detail::MemberFunctionAddressor< MemberFunctionType >
 

Additional Inherited Members

- Static Public Member Functions inherited from itk::simple::ProcessObject
static bool GetGlobalDefaultDebug ()
 
static void GlobalDefaultDebugOff ()
 
static void GlobalDefaultDebugOn ()
 
static void SetGlobalDefaultDebug (bool debugFlag)
 
static void GlobalWarningDisplayOn ()
 
static void GlobalWarningDisplayOff ()
 
static void SetGlobalWarningDisplay (bool flag)
 
static bool GetGlobalWarningDisplay ()
 
static double GetGlobalDefaultCoordinateTolerance ()
 Access the global tolerance to determine congruent spaces. More...
 
static void SetGlobalDefaultCoordinateTolerance (double)
 Access the global tolerance to determine congruent spaces. More...
 
static double GetGlobalDefaultDirectionTolerance ()
 Access the global tolerance to determine congruent spaces. More...
 
static void SetGlobalDefaultDirectionTolerance (double)
 Access the global tolerance to determine congruent spaces. More...
 
static bool SetGlobalDefaultThreader (const std::string &threader)
 Set/Get the default threader used for process objects. More...
 
static std::string GetGlobalDefaultThreader ()
 Set/Get the default threader used for process objects. More...
 
static void SetGlobalDefaultNumberOfThreads (unsigned int n)
 
static unsigned int GetGlobalDefaultNumberOfThreads ()
 Set/Get the default threader used for process objects. More...
 
- Protected Member Functions inherited from itk::simple::ImageFilter
void CheckImageMatchingDimension (const Image &image1, const Image &image2, const std::string &image2Name)
 
void CheckImageMatchingPixelType (const Image &image1, const Image &image2, const std::string &image2Name)
 
void CheckImageMatchingSize (const Image &image1, const Image &image2, const std::string &image2Name)
 
- Protected Member Functions inherited from itk::simple::ProcessObject
virtual unsigned long AddITKObserver (const itk::EventObject &, itk::Command *)
 
virtual itk::ProcessObjectGetActiveProcess ()
 
virtual void OnActiveProcessDelete ()
 
virtual void onCommandDelete (const itk::simple::Command *cmd) noexcept
 
virtual void PreUpdate (itk::ProcessObject *p)
 
virtual void RemoveITKObserver (EventCommand &e)
 
- Protected Member Functions inherited from itk::simple::NonCopyable
 NonCopyable ()=default
 
 NonCopyable (const NonCopyable &)=delete
 
NonCopyableoperator= (const NonCopyable &)=delete
 
- Static Protected Member Functions inherited from itk::simple::ImageFilter
template<class TImageType >
static void FixNonZeroIndex (TImageType *img)
 
- Static Protected Member Functions inherited from itk::simple::ProcessObject
template<class TImageType >
static TImageType::ConstPointer CastImageToITK (const Image &img)
 
template<class TPixelType , unsigned int VImageDimension, unsigned int VLength, template< typename, unsigned int > class TVector>
static Image CastITKToImage (itk::Image< TVector< TPixelType, VLength >, VImageDimension > *img)
 
template<unsigned int VImageDimension, unsigned int VLength, template< unsigned int > class TVector>
static Image CastITKToImage (itk::Image< TVector< VLength >, VImageDimension > *img)
 
template<class TImageType >
static Image CastITKToImage (TImageType *img)
 
static const itk::EventObjectGetITKEventObject (EventEnum e)
 
template<typename T >
static std::ostream & ToStringHelper (std::ostream &os, const T &v)
 
static std::ostream & ToStringHelper (std::ostream &os, const char &v)
 
static std::ostream & ToStringHelper (std::ostream &os, const signed char &v)
 
static std::ostream & ToStringHelper (std::ostream &os, const unsigned char &v)
 

Member Typedef Documentation

◆ MemberFunctionType

using itk::simple::JoinSeriesImageFilter::MemberFunctionType = Image (Self::*)( const std::vector<Image> & )
private

Setup for member function dispatching

Definition at line 108 of file sitkJoinSeriesImageFilter.h.

◆ PixelIDTypeList

Define the pixels types supported by this filter

Definition at line 62 of file sitkJoinSeriesImageFilter.h.

◆ Self

Definition at line 51 of file sitkJoinSeriesImageFilter.h.

Constructor & Destructor Documentation

◆ ~JoinSeriesImageFilter()

virtual itk::simple::JoinSeriesImageFilter::~JoinSeriesImageFilter ( )
virtual

Destructor

◆ JoinSeriesImageFilter()

itk::simple::JoinSeriesImageFilter::JoinSeriesImageFilter ( )

Default Constructor that takes no arguments and initializes default parameters

Member Function Documentation

◆ Execute() [1/6]

Image itk::simple::JoinSeriesImageFilter::Execute ( const Image image1)

◆ Execute() [2/6]

Image itk::simple::JoinSeriesImageFilter::Execute ( const Image image1,
const Image image2 
)

◆ Execute() [3/6]

Image itk::simple::JoinSeriesImageFilter::Execute ( const Image image1,
const Image image2,
const Image image3 
)

◆ Execute() [4/6]

Image itk::simple::JoinSeriesImageFilter::Execute ( const Image image1,
const Image image2,
const Image image3,
const Image image4 
)

◆ Execute() [5/6]

Image itk::simple::JoinSeriesImageFilter::Execute ( const Image image1,
const Image image2,
const Image image3,
const Image image4,
const Image image5 
)

◆ Execute() [6/6]

Image itk::simple::JoinSeriesImageFilter::Execute ( const std::vector< Image > &  images)

Execute the filter on the input images

◆ ExecuteInternal()

template<class TImageType >
Image itk::simple::JoinSeriesImageFilter::ExecuteInternal ( const std::vector< Image > &  images)
private

◆ GetName()

std::string itk::simple::JoinSeriesImageFilter::GetName ( ) const
inlinevirtual

Name of this class

Implements itk::simple::ProcessObject.

Definition at line 88 of file sitkJoinSeriesImageFilter.h.

◆ GetOrigin()

double itk::simple::JoinSeriesImageFilter::GetOrigin ( ) const
inline

Set/Get origin of the new dimension

Definition at line 75 of file sitkJoinSeriesImageFilter.h.

◆ GetSpacing()

double itk::simple::JoinSeriesImageFilter::GetSpacing ( ) const
inline

Set/Get spacing of the new dimension

Definition at line 85 of file sitkJoinSeriesImageFilter.h.

◆ SetOrigin()

Self& itk::simple::JoinSeriesImageFilter::SetOrigin ( double  Origin)
inline

Set/Get origin of the new dimension

Definition at line 70 of file sitkJoinSeriesImageFilter.h.

◆ SetSpacing()

Self& itk::simple::JoinSeriesImageFilter::SetSpacing ( double  Spacing)
inline

Set/Get spacing of the new dimension

Definition at line 80 of file sitkJoinSeriesImageFilter.h.

◆ ToString()

std::string itk::simple::JoinSeriesImageFilter::ToString ( ) const
virtual

Print ourselves out

Reimplemented from itk::simple::ProcessObject.

Friends And Related Function Documentation

◆ detail::MemberFunctionAddressor< MemberFunctionType >

Definition at line 113 of file sitkJoinSeriesImageFilter.h.

Member Data Documentation

◆ m_MemberFactory

std::unique_ptr<detail::MemberFunctionFactory<MemberFunctionType> > itk::simple::JoinSeriesImageFilter::m_MemberFactory
private

Definition at line 115 of file sitkJoinSeriesImageFilter.h.

◆ m_Origin

double itk::simple::JoinSeriesImageFilter::m_Origin {0.0}
private

Definition at line 118 of file sitkJoinSeriesImageFilter.h.

◆ m_Spacing

double itk::simple::JoinSeriesImageFilter::m_Spacing {1.0}
private

Definition at line 120 of file sitkJoinSeriesImageFilter.h.


The documentation for this class was generated from the following file: