ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkImageBase.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 itkImageBase_h
29#define itkImageBase_h
30
31#include "itkDataObject.h"
32
33#include "itkImageRegion.h"
34#include "itkMatrix.h"
35#include "itkObjectFactory.h"
36#include "itkOffset.h"
37#include "itkFixedArray.h"
38#include "itkImageHelper.h"
39#include "itkFloatTypes.h"
40
41#include <vxl_version.h>
42#include "vnl/vnl_matrix_fixed.hxx" // Get the templates
43
44namespace itk
45{
46
52inline constexpr double DefaultImageCoordinateTolerance = 1e-6;
53inline constexpr double DefaultImageDirectionTolerance = 1e-6;
54
55
113template <unsigned int VImageDimension = 2>
114class ITK_TEMPLATE_EXPORT ImageBase : public DataObject
115{
116public:
117 ITK_DISALLOW_COPY_AND_MOVE(ImageBase);
118
124
126 itkNewMacro(Self);
127
129 itkOverrideGetNameOfClassMacro(ImageBase);
130
132 using ImageDimensionType = unsigned int;
133
138 static constexpr ImageDimensionType ImageDimension = VImageDimension;
139
143
162
172
174 void
175 Initialize() override;
176
178 static constexpr unsigned int
180 {
181 return VImageDimension;
182 }
183
188 itkSetMacro(Origin, PointType);
189 virtual void
190 SetOrigin(const double origin[VImageDimension]);
191 virtual void
192 SetOrigin(const float origin[VImageDimension]);
194
220 virtual void
221 SetDirection(const DirectionType & direction);
222
226 itkGetConstReferenceMacro(Direction, DirectionType);
227
231 itkGetConstReferenceMacro(InverseDirection, DirectionType);
232
237 itkGetConstReferenceMacro(Spacing, SpacingType);
238
243 itkGetConstReferenceMacro(Origin, PointType);
244
252 virtual void
253 Allocate(bool itkNotUsed(initialize) = false)
254 {}
255
258 void
260 {
261 return this->Allocate(true);
262 }
263
270 virtual void
272
279 virtual const RegionType &
281 {
283 }
284
288 virtual void
290
294 virtual const RegionType &
296 {
297 return m_BufferedRegion;
298 }
299
307 virtual void
309
317 void
318 SetRequestedRegion(const DataObject * data) override;
319
324 virtual const RegionType &
326 {
327 return m_RequestedRegion;
328 }
329
333 virtual void
334 SetRegions(const RegionType & region)
335 {
336 this->SetLargestPossibleRegion(region);
337 this->SetBufferedRegion(region);
338 this->SetRequestedRegion(region);
339 }
340
341
342 virtual void
343 SetRegions(const SizeType & size)
344 {
345 RegionType region;
346 region.SetSize(size);
347
348 this->Self::SetRegions(region);
349 }
350
361 const OffsetValueType *
363 {
364 return m_OffsetTable;
365 }
366
367
374 inline OffsetValueType
375 ComputeOffset(const IndexType & ind) const
376 {
377 OffsetValueType offset = 0;
378
380 this->GetBufferedRegion().GetIndex(), ind, m_OffsetTable, offset);
381 return offset;
382 /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
383 * Leaving here for documentation purposes
384 * OffsetValueType ComputeOffset(const IndexType & ind) const
385 * {
386 * // need to add bounds checking for the region/buffer?
387 * OffsetValueType offset = 0;
388 * const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
389 * // data is arranged as [][][][slice][row][col]
390 * // with Index[0] = col, Index[1] = row, Index[2] = slice
391 * for ( int i = VImageDimension - 1; i > 0; i-- )
392 * {
393 * offset += ( ind[i] - bufferedRegionIndex[i] ) * m_OffsetTable[i];
394 * }
395 * offset += ( ind[0] - bufferedRegionIndex[0] );
396 * return offset;
397 * }
398 */
399 }
400
408 inline IndexType
410 {
411 IndexType index;
412 const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
414
416 return index;
417 /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
418 * Leaving here for documentation purposes
419 * IndexType ComputeIndex(OffsetValueType offset) const
420 * {
421 * IndexType index;
422 * const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
423 * for ( int i = VImageDimension - 1; i > 0; i-- )
424 * {
425 * index[i] = static_cast< IndexValueType >( offset / m_OffsetTable[i] );
426 * offset -= ( index[i] * m_OffsetTable[i] );
427 * index[i] += bufferedRegionIndex[i];
428 * }
429 * index[0] = bufferedRegionIndex[0] + static_cast< IndexValueType >( offset );
430 * return index;
431 * }
432 */
433 }
434
441 virtual void
442 SetSpacing(const SpacingType & spacing);
443 virtual void
444 SetSpacing(const double spacing[VImageDimension]);
445 virtual void
446 SetSpacing(const float spacing[VImageDimension]);
448
464 template <typename TCoordinate>
465 [[nodiscard]] IndexType
467 {
468 IndexType index;
470
471 for (unsigned int i = 0; i < VImageDimension; ++i)
472 {
473 TCoordinate sum{};
474 for (unsigned int j = 0; j < VImageDimension; ++j)
475 {
476 sum += this->m_PhysicalPointToIndex[i][j] * (point[j] - this->m_Origin[j]);
477 }
479 }
480 return index;
481 }
482
491 template <typename TCoordinate>
492 ITK_NODISCARD("Call the overload which has the point as the only parameter and returns the index")
493 bool TransformPhysicalPointToIndex(const Point<TCoordinate, VImageDimension> & point, IndexType & index) const
494 {
496
497 // Now, check to see if the index is within allowed bounds
498 const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
499 return isInside;
500 }
501
502
517 template <typename TIndexRep, typename TCoordinate>
520 {
523
524 for (unsigned int k = 0; k < VImageDimension; ++k)
525 {
526 cvector[k] = point[k] - this->m_Origin[k];
527 }
528 cvector = m_PhysicalPointToIndex * cvector;
529 for (unsigned int i = 0; i < VImageDimension; ++i)
530 {
531 index[i] = static_cast<TIndexRep>(cvector[i]);
532 }
533 return index;
534 }
535
544 template <typename TCoordinate, typename TIndexRep>
545 ITK_NODISCARD("Call the overload which has the point as the only parameter and returns the index")
546 bool TransformPhysicalPointToContinuousIndex(const Point<TCoordinate, VImageDimension> & point,
547 ContinuousIndex<TIndexRep, VImageDimension> & index) const
548 {
550
551 // Now, check to see if the index is within allowed bounds
552 const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
553 return isInside;
554 }
555
560 template <typename TCoordinate, typename TIndexRep>
561 void
564 {
565 for (unsigned int r = 0; r < VImageDimension; ++r)
566 {
567 TCoordinate sum{};
568 for (unsigned int c = 0; c < VImageDimension; ++c)
569 {
570 sum += this->m_IndexToPhysicalPoint(r, c) * index[c];
571 }
572 point[r] = sum + this->m_Origin[r];
573 }
574 }
575
576
581 template <typename TCoordinate, typename TIndexRep>
589
590
596 template <typename TCoordinate>
597 void
599 {
600 for (unsigned int i = 0; i < VImageDimension; ++i)
601 {
602 point[i] = this->m_Origin[i];
603 for (unsigned int j = 0; j < VImageDimension; ++j)
604 {
605 point[i] += m_IndexToPhysicalPoint[i][j] * index[j];
606 }
607 }
608 }
609
610
616 template <typename TCoordinate>
624
625
640 template <typename TCoordinate>
641 void
643 FixedArray<TCoordinate, VImageDimension> & outputGradient) const
644 {
645 const DirectionType & direction = this->GetDirection();
646
647 itkAssertInDebugAndIgnoreInReleaseMacro(inputGradient.GetDataPointer() != outputGradient.GetDataPointer());
648
649 for (unsigned int i = 0; i < VImageDimension; ++i)
650 {
651 using CoordSumType = typename NumericTraits<TCoordinate>::AccumulateType;
652 CoordSumType sum{};
653 for (unsigned int j = 0; j < VImageDimension; ++j)
654 {
655 sum += direction[i][j] * inputGradient[j];
656 }
657 outputGradient[i] = static_cast<TCoordinate>(sum);
658 }
659 }
660
667 template <typename TVector>
668 [[nodiscard]] TVector
669 TransformLocalVectorToPhysicalVector(const TVector & inputGradient) const
670 {
671 TVector outputGradient;
672 TransformLocalVectorToPhysicalVector(inputGradient, outputGradient);
673 return outputGradient;
674 }
675
676
677
690 template <typename TCoordinate>
691 void
693 FixedArray<TCoordinate, VImageDimension> & outputGradient) const
694 {
695 const DirectionType & inverseDirection = this->GetInverseDirection();
696
697 itkAssertInDebugAndIgnoreInReleaseMacro(inputGradient.GetDataPointer() != outputGradient.GetDataPointer());
698
699 for (unsigned int i = 0; i < VImageDimension; ++i)
700 {
701 using CoordSumType = typename NumericTraits<TCoordinate>::AccumulateType;
702 CoordSumType sum{};
703 for (unsigned int j = 0; j < VImageDimension; ++j)
704 {
705 sum += inverseDirection[i][j] * inputGradient[j];
706 }
707 outputGradient[i] = static_cast<TCoordinate>(sum);
708 }
709 }
710
717 template <typename TVector>
718 [[nodiscard]] TVector
719 TransformPhysicalVectorToLocalVector(const TVector & inputGradient) const
720 {
721 TVector outputGradient;
722 TransformPhysicalVectorToLocalVector(inputGradient, outputGradient);
723 return outputGradient;
724 }
725
726
736 void
737 CopyInformation(const DataObject * data) override;
738
749 virtual void
750 Graft(const Self * image);
751
759 void
761
769 void
771
775 void
777
787 bool
789
798 bool
800
804 bool
805 IsCongruentImageGeometry(const ImageBase * otherImage, double coordinateTolerance, double directionTolerance) const;
806
814 bool
816 double coordinateTolerance = DefaultImageCoordinateTolerance,
817 double directionTolerance = DefaultImageDirectionTolerance) const;
818
837 virtual unsigned int
839 virtual void
842
843protected:
844 ImageBase() = default;
845 ~ImageBase() override = default;
846 void
847 PrintSelf(std::ostream & os, Indent indent) const override;
848
853 void
855
861 virtual void
863
864protected:
873
879
884 virtual void
886
898 FastComputeOffset(const IndexType & ind) const
899 {
900 OffsetValueType offset = 0;
902 Self::GetBufferedRegion().GetIndex(), ind, m_OffsetTable, offset);
903 return offset;
904 }
905
906
918 IndexType
920 {
921 IndexType index;
922 const IndexType & bufferedRegionIndex = Self::GetBufferedRegion().GetIndex();
924 return index;
925 }
926
927
928 void
929 Graft(const DataObject * data) override;
930
931private:
932 OffsetValueType m_OffsetTable[VImageDimension + 1]{};
933
937};
938} // end namespace itk
939
940#ifndef ITK_MANUAL_INSTANTIATION
941# include "itkImageBase.hxx"
942#endif
943
944#endif
A templated class holding a point in n-Dimensional image space.
Simulate a standard C array with copy semantics.
ValueType * GetDataPointer()
void TransformPhysicalVectorToLocalVector(const FixedArray< TCoordinate, VImageDimension > &inputGradient, FixedArray< TCoordinate, VImageDimension > &outputGradient) const
IndexType FastComputeIndex(OffsetValueType offset) const
Vector< SpacingValueType, VImageDimension > SpacingType
virtual void SetRequestedRegion(const RegionType &region)
IndexType ComputeIndex(OffsetValueType offset) const
virtual void InitializeBufferedRegion()
void SetRequestedRegionToLargestPossibleRegion() override
ImageRegion< VImageDimension > RegionType
ContinuousIndex< TIndexRep, VImageDimension > TransformPhysicalPointToContinuousIndex(const Point< TCoordinate, VImageDimension > &point) const
Returns the continuous index from a physical point.
bool VerifyRequestedRegion() override
SpacePrecisionType SpacingValueType
TVector TransformPhysicalVectorToLocalVector(const TVector &inputGradient) const
IndexType TransformPhysicalPointToIndex(const Point< TCoordinate, VImageDimension > &point) const
virtual void SetDirection(const DirectionType &direction)
bool IsSameImageGeometryAs(const ImageBase *otherImage, double coordinateTolerance=DefaultImageCoordinateTolerance, double directionTolerance=DefaultImageDirectionTolerance) const
virtual const DirectionType & GetDirection() const
void AllocateInitialized()
void PrintSelf(std::ostream &os, Indent indent) const override
virtual const RegionType & GetBufferedRegion() const
bool RequestedRegionIsOutsideOfTheBufferedRegion() override
virtual void SetRegions(const SizeType &size)
virtual void SetBufferedRegion(const RegionType &region)
Index< VImageDimension > IndexType
Offset< VImageDimension > OffsetType
virtual void SetSpacing(const SpacingType &spacing)
void UpdateOutputData() override
virtual void SetNumberOfComponentsPerPixel(unsigned int)
OffsetValueType m_OffsetTable[VImageDimension+1]
OffsetValueType FastComputeOffset(const IndexType &ind) const
virtual void Allocate(bool initialize=false)
TVector TransformLocalVectorToPhysicalVector(const TVector &inputGradient) const
virtual void SetOrigin(const float origin[VImageDimension])
typename OffsetType::OffsetValueType OffsetValueType
typename IndexType::IndexValueType IndexValueType
void UpdateOutputInformation() override
virtual void Graft(const Self *image)
void TransformContinuousIndexToPhysicalPoint(const ContinuousIndex< TIndexRep, VImageDimension > &index, Point< TCoordinate, VImageDimension > &point) const
SmartPointer< const Self > ConstPointer
virtual unsigned int GetNumberOfComponentsPerPixel() const
virtual void SetSpacing(const float spacing[VImageDimension])
const OffsetValueType * GetOffsetTable() const
virtual const RegionType & GetRequestedRegion() const
virtual void ComputeIndexToPhysicalPointMatrices()
void SetRequestedRegion(const DataObject *data) override
static constexpr ImageDimensionType ImageDimension
OffsetValueType ComputeOffset(const IndexType &ind) const
virtual void SetRegions(const RegionType &region)
Matrix< SpacePrecisionType, VImageDimension, VImageDimension > DirectionType
void CopyInformation(const DataObject *data) override
void TransformIndexToPhysicalPoint(const IndexType &index, Point< TCoordinate, VImageDimension > &point) const
virtual void SetOrigin(const double origin[VImageDimension])
Size< VImageDimension > SizeType
virtual void SetSpacing(const double spacing[VImageDimension])
ImageBase()=default
Point< PointValueType, VImageDimension > PointType
Point< TCoordinate, VImageDimension > TransformContinuousIndexToPhysicalPoint(const ContinuousIndex< TIndexRep, VImageDimension > &index) const
virtual void SetLargestPossibleRegion(const RegionType &region)
virtual const RegionType & GetLargestPossibleRegion() const
Point< TCoordinate, VImageDimension > TransformIndexToPhysicalPoint(const IndexType &index) const
void ComputeOffsetTable()
virtual const DirectionType & GetInverseDirection() const
~ImageBase() override=default
void Initialize() override
bool IsCongruentImageGeometry(const ImageBase *otherImage, double coordinateTolerance, double directionTolerance) const
typename SizeType::SizeValueType SizeValueType
void TransformLocalVectorToPhysicalVector(const FixedArray< TCoordinate, VImageDimension > &inputGradient, FixedArray< TCoordinate, VImageDimension > &outputGradient) const
void Graft(const DataObject *data) override
static constexpr unsigned int GetImageDimension()
static void ComputeIndex(const IndexType &bufferedRegionIndex, OffsetValueType offset, const OffsetValueType offsetTable[], IndexType &index)
static void ComputeOffset(const IndexType &bufferedRegionIndex, const IndexType &index, const OffsetValueType offsetTable[], OffsetValueType &offset)
An image region represents a structured region of data.
void SetSize(const SizeType &size)
Vector< SpacingValueType, VImageDimension > SpacingType
Control indentation during Print() invocation.
Definition itkIndent.h:50
A templated class holding a M x N size Matrix.
Definition itkMatrix.h:53
A templated class holding a geometric point in n-Dimensional space.
Definition itkPoint.h:54
Implements transparent reference counting.
Point< PointValueType, VImageDimension > PointType
typename OffsetType::OffsetValueType OffsetValueType
A templated class holding a n-Dimensional vector.
Definition itkVector.h:63
TInput RoundHalfIntegerUp(TInput x) template< typename TReturn
Round towards nearest integer (This is a synonym for RoundHalfIntegerUp)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
constexpr TContainer MakeFilled(typename TContainer::const_reference value)
constexpr double DefaultImageDirectionTolerance
constexpr double DefaultImageCoordinateTolerance
*par Constraints *The filter image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating point(float or double) or *a user-defined "real" numerical type with arithmetic operations defined *sufficient to compute derivatives. **\par Performance *This filter will automatically multithread if run with *SetUsePrincipleComponents
long OffsetValueType
Definition itkIntTypes.h:97
double SpacePrecisionType
Represent a n-dimensional index in a n-dimensional image.
Definition itkIndex.h:69
itk::IndexValueType IndexValueType
Definition itkIndex.h:79
Represent a n-dimensional offset between two n-dimensional indexes of n-dimensional image.
Definition itkOffset.h:67
itk::OffsetValueType OffsetValueType
Definition itkOffset.h:77
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition itkSize.h:70
itk::SizeValueType SizeValueType
Definition itkSize.h:80