18#ifndef itkDiffusionTensor3DReconstructionImageFilter_h
19#define itkDiffusionTensor3DReconstructionImageFilter_h
24#include "vnl/vnl_matrix.h"
25#include "vnl/vnl_vector_fixed.h"
26#include "vnl/vnl_matrix_fixed.h"
27#include "vnl/algo/vnl_svd.h"
30#include "ITKDiffusionTensorImageExport.h"
53extern ITKDiffusionTensorImage_EXPORT std::ostream &
139template <
typename TReferenceImagePixelType,
140 typename TGradientImagePixelType = TReferenceImagePixelType,
141 typename TTensorPixelType = double,
144 :
public ImageToImageFilter<Image<TReferenceImagePixelType, 3>, Image<DiffusionTensor3D<TTensorPixelType>, 3>>
225 "Cannot call both methods:AddGradientImage and SetGradientImage. Please call only one of them.");
235 virtual ReferenceImageType *
242 virtual GradientDirectionType
247 itkExceptionMacro(
"Gradient direction " << idx <<
" does not exist");
275 itkSetMacro(BValue, TTensorPixelType);
279 itkGetConstReferenceMacro(BValue, TTensorPixelType);
314#if !defined(ITK_LEGACY_REMOVE)
318 GradientImageTypeEnumeration::GradientIsInASingleImage;
320 GradientImageTypeEnumeration::GradientIsInManyImages;
353#ifndef ITK_MANUAL_INSTANTIATION
354# include "itkDiffusionTensor3DReconstructionImageFilter.hxx"
Contains all enum classes used by DiffusionTensor3DReconstructionImageFilter class.
@ GradientIsInASingleImage
VectorImage< GradientPixelType, 3 > GradientImagesType
TTensorPixelType m_BValue
void SetReferenceImage(ReferenceImageType *referenceImage)
vnl_matrix< double > CoefficientMatrixType
Image< GradientPixelType, 3 > GradientImageType
void SetGradientImage(GradientDirectionContainerType *, const GradientImagesType *gradientImage)
DiffusionTensor3DReconstructionImageFilter Self
typename OutputImageType::RegionType OutputImageRegionType
SmartPointer< const Self > ConstPointer
void AddGradientImage(const GradientDirectionType &, const GradientImageType *gradientImage)
void VerifyPreconditions() const override
Verifies that the process object has been configured correctly, that all required inputs are set,...
void PrintSelf(std::ostream &os, Indent indent) const override
SmartPointer< Self > Pointer
vnl_matrix_fixed< double, 6, 6 > TensorBasisMatrixType
TensorImageType OutputImageType
void SetMaskSpatialObject(MaskSpatialObjectType *maskSpatialObject)
virtual ReferenceImageType * GetReferenceImage()
typename Superclass::InputImageType ReferenceImageType
GradientImageTypeEnumeration m_GradientImageTypeEnumeration
DiffusionTensor3DReconstructionImageFilter()
virtual GradientDirectionType GetGradientDirection(unsigned int idx) const
ReferencePixelType m_Threshold
vnl_vector_fixed< double, 3 > GradientDirectionType
unsigned int m_NumberOfBaselineImages
ImageToImageFilter< Image< TReferenceImagePixelType, 3 >, Image< DiffusionTensor3D< TTensorPixelType >, 3 > > Superclass
GradientDirectionContainerType::Pointer m_GradientDirectionContainer
void BeforeThreadedGenerateData() override
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
~DiffusionTensor3DReconstructionImageFilter() override=default
DiffusionTensor3D< TTensorPixelType > TensorPixelType
void SetMaskImage(MaskImageType *maskImage)
CoefficientMatrixType m_BMatrix
TReferenceImagePixelType ReferencePixelType
SpatialObject< 3 > MaskSpatialObjectType
TensorBasisMatrixType m_TensorBasis
TGradientImagePixelType GradientPixelType
void ComputeTensorBasis()
TMaskImageType MaskImageType
DiffusionTensor3DReconstructionImageFilterEnums::GradientImageFormat GradientImageTypeEnumeration
const GradientImageType * GetGradientImage(unsigned int index) const
VectorContainer< unsigned int, GradientDirectionType > GradientDirectionContainerType
unsigned int m_NumberOfGradientDirections
Image< TensorPixelType, 3 > TensorImageType
Represent a diffusion tensor as used in DTI images.
Image< TReferenceImagePixelType, 3 > InputImageType
typename OutputImageType::RegionType OutputImageRegionType
Templated n-dimensional image class.
Control indentation during Print() invocation.
virtual void SetNthInput(DataObjectPointerArraySizeType idx, DataObject *input)
DataObject * GetInput(const DataObjectIdentifierType &key)
Return an input.
Implements transparent reference counting.
Implementation of the composite pattern.
Templated n-dimensional vector image class.
SmartPointer< Self > Pointer
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)
detail::VectorContainer< std::conditional_t< std::is_void_v< T2 >, SizeValueType, T1 >, std::conditional_t< std::is_void_v< T2 >, T1, T2 > > VectorContainer