ITK  6.0.0
Insight Toolkit
itkDiffusionTensor3DReconstructionImageFilter.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 #ifndef itkDiffusionTensor3DReconstructionImageFilter_h
19 #define itkDiffusionTensor3DReconstructionImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include "itkSpatialObject.h"
23 #include "itkDiffusionTensor3D.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"
28 #include "itkVectorContainer.h"
29 #include "itkVectorImage.h"
30 #include "ITKDiffusionTensorImageExport.h"
31 
32 namespace itk
33 {
39 {
40 public:
45  enum class GradientImageFormat : uint8_t
46  {
49  Else
50  };
51 };
52 // Define how to print enumeration
53 extern ITKDiffusionTensorImage_EXPORT std::ostream &
55 
145 template <typename TReferenceImagePixelType,
146  typename TGradientImagePixelType = TReferenceImagePixelType,
147  typename TTensorPixelType = double,
148  typename TMaskImageType = Image<unsigned char, 3>>
150  : public ImageToImageFilter<Image<TReferenceImagePixelType, 3>, Image<DiffusionTensor3D<TTensorPixelType>, 3>>
151 {
152 public:
156  using Superclass =
158 
160  itkNewMacro(Self);
161 
163  itkOverrideGetNameOfClassMacro(DiffusionTensor3DReconstructionImageFilter);
164 
165  using ReferencePixelType = TReferenceImagePixelType;
166 
167  using GradientPixelType = TGradientImagePixelType;
168 
170 
173  using ReferenceImageType = typename Superclass::InputImageType;
174 
176 
178 
179  using typename Superclass::OutputImageRegionType;
180 
183 
189 
192 
194  using MaskImageType = TMaskImageType;
195 
197  using TensorBasisMatrixType = vnl_matrix_fixed<double, 6, 6>;
198 
199  using CoefficientMatrixType = vnl_matrix<double>;
200 
202  using GradientDirectionType = vnl_vector_fixed<double, 3>;
203 
206 
208  void
209  AddGradientImage(const GradientDirectionType &, const GradientImageType * gradientImage);
210  const GradientImageType *
211  GetGradientImage(unsigned int index) const;
220  void
221  SetGradientImage(GradientDirectionContainerType *, const GradientImagesType * gradientImage);
222 
224  void
226  {
227  if (m_GradientImageTypeEnumeration ==
229  {
230  itkExceptionMacro(
231  "Cannot call both methods:AddGradientImage and SetGradientImage. Please call only one of them.");
232  }
233 
234  this->ProcessObject::SetNthInput(0, referenceImage);
235 
236  m_GradientImageTypeEnumeration =
238  }
239 
241  virtual ReferenceImageType *
243  {
244  return (static_cast<ReferenceImageType *>(this->ProcessObject::GetInput(0)));
245  }
246 
248  virtual GradientDirectionType
249  GetGradientDirection(unsigned int idx) const
250  {
251  if (idx >= m_NumberOfGradientDirections)
252  {
253  itkExceptionMacro("Gradient direction " << idx << " does not exist");
254  }
255  return m_GradientDirectionContainer->ElementAt(idx);
256  }
260  void
261  SetMaskImage(MaskImageType * maskImage);
262 
264  void
265  SetMaskSpatialObject(MaskSpatialObjectType * maskSpatialObject);
266 
267 
271  itkSetMacro(Threshold, ReferencePixelType);
272  itkGetConstMacro(Threshold, ReferencePixelType);
281  itkSetMacro(BValue, TTensorPixelType);
282 #ifdef GetBValue
283 # undef GetBValue
284 #endif
285  itkGetConstReferenceMacro(BValue, TTensorPixelType);
288 #ifdef ITK_USE_CONCEPT_CHECKING
289  // Begin concept checking
290  itkConceptMacro(ReferenceEqualityComparableCheck, (Concept::EqualityComparable<ReferencePixelType>));
291  itkConceptMacro(TensorEqualityComparableCheck, (Concept::EqualityComparable<TensorPixelType>));
292  itkConceptMacro(GradientConvertibleToDoubleCheck, (Concept::Convertible<GradientPixelType, double>));
293  itkConceptMacro(DoubleConvertibleToTensorCheck, (Concept::Convertible<double, TensorPixelType>));
294  itkConceptMacro(GradientReferenceAdditiveOperatorsCheck,
296  itkConceptMacro(GradientReferenceAdditiveAndAssignOperatorsCheck,
298 
299  itkConceptMacro(ReferenceOStreamWritableCheck, (Concept::OStreamWritable<ReferencePixelType>));
300  itkConceptMacro(TensorOStreamWritableCheck, (Concept::OStreamWritable<TensorPixelType>));
301  // End concept checking
302 #endif
303 
304 protected:
306  ~DiffusionTensor3DReconstructionImageFilter() override = default;
307  void
308  PrintSelf(std::ostream & os, Indent indent) const override;
309 
310  void
311  ComputeTensorBasis();
312 
313  void
314  BeforeThreadedGenerateData() override;
315 
316  void
317  DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
318 
319  void
320  VerifyPreconditions() const override;
321 
324 #if !defined(ITK_LEGACY_REMOVE)
325  // We need to expose the enum values at the class level
326  // for backwards compatibility
327  static constexpr GradientImageTypeEnumeration GradientIsInASingleImage =
328  GradientImageTypeEnumeration::GradientIsInASingleImage;
329  static constexpr GradientImageTypeEnumeration GradientIsInManyImages =
330  GradientImageTypeEnumeration::GradientIsInManyImages;
331  static constexpr GradientImageTypeEnumeration Else = GradientImageTypeEnumeration::Else;
332 #endif
333 
334 private:
335  /* Tensor basis coeffs */
336  TensorBasisMatrixType m_TensorBasis{};
337 
339 
341  GradientDirectionContainerType::Pointer m_GradientDirectionContainer{};
342 
344  unsigned int m_NumberOfGradientDirections{};
345 
347  unsigned int m_NumberOfBaselineImages{};
348 
350  ReferencePixelType m_Threshold{};
351 
353  TTensorPixelType m_BValue{};
354 
356  GradientImageTypeEnumeration m_GradientImageTypeEnumeration{};
357 
359  bool m_MaskImagePresent{};
360 };
361 } // namespace itk
362 
363 #ifndef ITK_MANUAL_INSTANTIATION
364 # include "itkDiffusionTensor3DReconstructionImageFilter.hxx"
365 #endif
366 
367 #endif
itk::DiffusionTensor3DReconstructionImageFilter::ReferenceImageType
typename Superclass::InputImageType ReferenceImageType
Definition: itkDiffusionTensor3DReconstructionImageFilter.h:173
itk::DiffusionTensor3DReconstructionImageFilter::TensorBasisMatrixType
vnl_matrix_fixed< double, 6, 6 > TensorBasisMatrixType
Definition: itkDiffusionTensor3DReconstructionImageFilter.h:197
itk::Concept::OStreamWritable
Definition: itkConceptChecking.h:639
itk::detail::VectorContainer
Define a front-end to the STL "vector" container that conforms to the IndexedContainerInterface.
Definition: itkVectorContainer.h:51
itk::DiffusionTensor3DReconstructionImageFilterEnums
Contains all enum classes used by DiffusionTensor3DReconstructionImageFilter class.
Definition: itkDiffusionTensor3DReconstructionImageFilter.h:38
itk::DiffusionTensor3DReconstructionImageFilter::CoefficientMatrixType
vnl_matrix< double > CoefficientMatrixType
Definition: itkDiffusionTensor3DReconstructionImageFilter.h:199
itk::DiffusionTensor3DReconstructionImageFilter::GradientPixelType
TGradientImagePixelType GradientPixelType
Definition: itkDiffusionTensor3DReconstructionImageFilter.h:167
itkSpatialObject.h
itk::VectorImage
Templated n-dimensional vector image class.
Definition: itkImageAlgorithm.h:29
itkDiffusionTensor3D.h
itk::DiffusionTensor3D
Represent a diffusion tensor as used in DTI images.
Definition: itkDiffusionTensor3D.h:79
itk::operator<<
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::DiffusionTensor3DReconstructionImageFilterEnums::GradientImageFormat
GradientImageFormat
Definition: itkDiffusionTensor3DReconstructionImageFilter.h:45
itk::DiffusionTensor3DReconstructionImageFilter::ReferencePixelType
TReferenceImagePixelType ReferencePixelType
Definition: itkDiffusionTensor3DReconstructionImageFilter.h:165
itk::Concept::AdditiveAndAssignOperators
Definition: itkConceptChecking.h:393
itkVectorImage.h
itk::DiffusionTensor3DReconstructionImageFilter::GetGradientDirection
virtual GradientDirectionType GetGradientDirection(unsigned int idx) const
Definition: itkDiffusionTensor3DReconstructionImageFilter.h:249
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::DiffusionTensor3DReconstructionImageFilterEnums::GradientImageFormat::GradientIsInManyImages
itk::DiffusionTensor3DReconstructionImageFilterEnums::GradientImageFormat::Else
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itk::DiffusionTensor3DReconstructionImageFilter::SetReferenceImage
void SetReferenceImage(ReferenceImageType *referenceImage)
Definition: itkDiffusionTensor3DReconstructionImageFilter.h:225
itk::DiffusionTensor3DReconstructionImageFilter::MaskImageType
TMaskImageType MaskImageType
Definition: itkDiffusionTensor3DReconstructionImageFilter.h:194
itk::SpatialObject
Implementation of the composite pattern.
Definition: itkSpatialObject.h:58
itk::DiffusionTensor3DReconstructionImageFilter::GradientDirectionType
vnl_vector_fixed< double, 3 > GradientDirectionType
Definition: itkDiffusionTensor3DReconstructionImageFilter.h:202
itkImageToImageFilter.h
itkVectorContainer.h
itk::Concept::AdditiveOperators
Definition: itkConceptChecking.h:361
itk::DiffusionTensor3DReconstructionImageFilterEnums::GradientImageFormat::GradientIsInASingleImage
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ProcessObject::SetNthInput
virtual void SetNthInput(DataObjectPointerArraySizeType idx, DataObject *input)
itk::Concept::Convertible
Definition: itkConceptChecking.h:217
GradientImageFormat
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::ProcessObject::GetInput
DataObject * GetInput(const DataObjectIdentifierType &key)
Return an input.
itk::Concept::EqualityComparable
Definition: itkConceptChecking.h:307
itk::DiffusionTensor3DReconstructionImageFilter::GetReferenceImage
virtual ReferenceImageType * GetReferenceImage()
Definition: itkDiffusionTensor3DReconstructionImageFilter.h:242
itk::DiffusionTensor3DReconstructionImageFilter
This class takes as input one or more reference image (acquired in the absence of diffusion sensitizi...
Definition: itkDiffusionTensor3DReconstructionImageFilter.h:149