1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166
|
/*=========================================================================
*
* Copyright NumFOCUS
*
* Licensed under the Apache License, Version 2.0 (the "License");
* you may not use this file except in compliance with the License.
* You may obtain a copy of the License at
*
* https://www.apache.org/licenses/LICENSE-2.0.txt
*
* Unless required by applicable law or agreed to in writing, software
* distributed under the License is distributed on an "AS IS" BASIS,
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
* See the License for the specific language governing permissions and
* limitations under the License.
*
*=========================================================================*/
#ifndef itkTileImageFilter_h
#define itkTileImageFilter_h
#include "itkImageToImageFilter.h"
#include "itkFixedArray.h"
namespace itk
{
/**
* \class TileImageFilter
* \brief Tile multiple input images into a single output image.
*
* This filter will tile multiple images using a user-specified
* layout. The tile sizes will be large enough to accommodate the
* largest image for each tile. The layout is specified with the
* SetLayout method. The layout has the same dimension as the output
* image. If all entries of the layout are positive, the tiled output
* will contain the exact number of tiles. If the layout contains a 0
* in the last dimension, the filter will compute a size that will
* accommodate all of the images. Empty tiles are filled with the
* value specified with the SetDefault value method. The input images
* must have a dimension less than or equal to the output image. The
* output image have a larger dimension than the input images. This
* filter can be used to create a volume from a series of inputs by
* specifying a layout of 1,1,0.
* \ingroup ITKImageGrid
*
* \sphinx
* \sphinxexample{Filtering/ImageGrid/Stack2DImagesInto3DImage,Stack 2D Images Into 3D Image}
* \sphinxexample{Filtering/ImageGrid/TileImagesSideBySide,Tile Images Side By Side}
* \endsphinx
*/
template <typename TInputImage, typename TOutputImage>
class ITK_TEMPLATE_EXPORT TileImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
{
public:
ITK_DISALLOW_COPY_AND_MOVE(TileImageFilter);
/** Standard Self type alias */
using Self = TileImageFilter;
using Superclass = ImageToImageFilter<TInputImage, TOutputImage>;
using Pointer = SmartPointer<Self>;
using ConstPointer = SmartPointer<const Self>;
/** Method for creation through the object factory. */
itkNewMacro(Self);
/** \see LightObject::GetNameOfClass() */
itkOverrideGetNameOfClassMacro(TileImageFilter);
/** Image pixel value type alias. */
using InputPixelType = typename TInputImage::PixelType;
using OutputPixelType = typename TOutputImage::PixelType;
using OutputPixelComponentType = typename NumericTraits<OutputPixelType>::ValueType;
/** Image related type alias. */
using InputImagePointer = typename TInputImage::Pointer;
using OutputImagePointer = typename TOutputImage::Pointer;
using InputSizeType = typename TInputImage::SizeType;
using InputIndexType = typename TInputImage::IndexType;
using InputImageRegionType = typename TInputImage::RegionType;
using OutputSizeType = typename TOutputImage::SizeType;
using OutputIndexType = typename TOutputImage::IndexType;
using OutputImageRegionType = typename TOutputImage::RegionType;
/** Image related type alias. */
static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
/**
* \class TileInfo
* Define a tile structure
* \ingroup ITKImageGrid
*/
class TileInfo
{
public:
int m_ImageNumber{ -1 };
OutputImageRegionType m_Region;
TileInfo() = default;
};
using TileImageType = Image<TileInfo, Self::OutputImageDimension>;
/** LayoutArray type. */
using LayoutArrayType = FixedArray<unsigned int, Self::OutputImageDimension>;
/** Set/Get the layout of the tiles. If the last Layout value is 0,
* the filter will compute a value that will accommodate all of the
* images. */
itkSetMacro(Layout, LayoutArrayType);
itkGetConstMacro(Layout, LayoutArrayType);
/** Set the pixel value for locations that are not covered by an
* input image. The default default pixel value is Zero. */
itkSetMacro(DefaultPixelValue, OutputPixelType);
/** Get the pixel value for locations that are not covered by an
* input image. */
itkGetConstMacro(DefaultPixelValue, OutputPixelType);
#ifdef ITK_USE_CONCEPT_CHECKING
// Begin concept checking
itkConceptMacro(OutputEqualityComparableCheck, (Concept::EqualityComparable<OutputPixelType>));
itkConceptMacro(SameTypeCheck, (Concept::SameType<InputPixelType, OutputPixelType>));
itkConceptMacro(OutputOStreamWritableCheck, (Concept::OStreamWritable<OutputPixelType>));
// End concept checking
#endif
protected:
TileImageFilter();
~TileImageFilter() override = default;
void
PrintSelf(std::ostream & os, Indent indent) const override;
void
GenerateInputRequestedRegion() override;
void
GenerateOutputInformation() override;
void
GenerateData() override;
/** Override VerifyInputInformation() since this filter's inputs do
* not need to occupy the same physical space.
*
* \sa ProcessObject::VerifyInputInformation
*/
void
VerifyInputInformation() ITKv5_CONST override;
private:
typename TileImageType::Pointer m_TileImage{};
OutputPixelType m_DefaultPixelValue{};
LayoutArrayType m_Layout{};
}; // end of class
} // end namespace itk
#ifndef ITK_MANUAL_INSTANTIATION
# include "itkTileImageFilter.hxx"
#endif
#endif
|