ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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
22#include "itkSpatialObject.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
32namespace itk
33{
52// Define how to print enumeration
53extern ITKDiffusionTensorImage_EXPORT std::ostream &
55
138
139template <typename TReferenceImagePixelType,
140 typename TGradientImagePixelType = TReferenceImagePixelType,
141 typename TTensorPixelType = double,
142 typename TMaskImageType = Image<unsigned char, 3>>
144 : public ImageToImageFilter<Image<TReferenceImagePixelType, 3>, Image<DiffusionTensor3D<TTensorPixelType>, 3>>
145{
146public:
152
154 itkNewMacro(Self);
155
157 itkOverrideGetNameOfClassMacro(DiffusionTensor3DReconstructionImageFilter);
158
159 using ReferencePixelType = TReferenceImagePixelType;
160
161 using GradientPixelType = TGradientImagePixelType;
162
164
168
170
172
174
177
183
186
188 using MaskImageType = TMaskImageType;
189
191 using TensorBasisMatrixType = vnl_matrix_fixed<double, 6, 6>;
192
193 using CoefficientMatrixType = vnl_matrix<double>;
194
196 using GradientDirectionType = vnl_vector_fixed<double, 3>;
197
200
202 void
204 const GradientImageType *
205 GetGradientImage(unsigned int index) const;
207
214 void
216
218 void
220 {
223 {
224 itkExceptionMacro(
225 "Cannot call both methods:AddGradientImage and SetGradientImage. Please call only one of them.");
226 }
227
228 this->ProcessObject::SetNthInput(0, referenceImage);
229
232 }
233
235 virtual ReferenceImageType *
237 {
238 return (static_cast<ReferenceImageType *>(this->ProcessObject::GetInput(0)));
239 }
240
242 virtual GradientDirectionType
243 GetGradientDirection(unsigned int idx) const
244 {
246 {
247 itkExceptionMacro("Gradient direction " << idx << " does not exist");
248 }
249 return m_GradientDirectionContainer->ElementAt(idx);
250 }
251
252
254 void
256
258 void
260
261
265 itkSetMacro(Threshold, ReferencePixelType);
266 itkGetConstMacro(Threshold, ReferencePixelType);
268
275 itkSetMacro(BValue, TTensorPixelType);
276#ifdef GetBValue
277# undef GetBValue
278#endif
279 itkGetConstReferenceMacro(BValue, TTensorPixelType);
281
282 itkConceptMacro(ReferenceEqualityComparableCheck, (Concept::EqualityComparable<ReferencePixelType>));
283 itkConceptMacro(TensorEqualityComparableCheck, (Concept::EqualityComparable<TensorPixelType>));
284 itkConceptMacro(GradientConvertibleToDoubleCheck, (Concept::Convertible<GradientPixelType, double>));
285 itkConceptMacro(DoubleConvertibleToTensorCheck, (Concept::Convertible<double, TensorPixelType>));
286 itkConceptMacro(GradientReferenceAdditiveOperatorsCheck,
288 itkConceptMacro(GradientReferenceAdditiveAndAssignOperatorsCheck,
290
291 itkConceptMacro(ReferenceOStreamWritableCheck, (Concept::OStreamWritable<ReferencePixelType>));
292 itkConceptMacro(TensorOStreamWritableCheck, (Concept::OStreamWritable<TensorPixelType>));
293
294protected:
297 void
298 PrintSelf(std::ostream & os, Indent indent) const override;
299
300 void
302
303 void
305
306 void
307 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
308
309 void
310 VerifyPreconditions() const override;
311
314#if !defined(ITK_LEGACY_REMOVE)
315 // We need to expose the enum values at the class level
316 // for backwards compatibility
317 static constexpr GradientImageTypeEnumeration GradientIsInASingleImage =
318 GradientImageTypeEnumeration::GradientIsInASingleImage;
319 static constexpr GradientImageTypeEnumeration GradientIsInManyImages =
320 GradientImageTypeEnumeration::GradientIsInManyImages;
321 static constexpr GradientImageTypeEnumeration Else = GradientImageTypeEnumeration::Else;
322#endif
323
324private:
325 /* Tensor basis coeffs */
327
329
332
335
338
341
343 TTensorPixelType m_BValue{};
344
347
350};
351} // namespace itk
352
353#ifndef ITK_MANUAL_INSTANTIATION
354# include "itkDiffusionTensor3DReconstructionImageFilter.hxx"
355#endif
356
357#endif
Contains all enum classes used by DiffusionTensor3DReconstructionImageFilter class.
void SetGradientImage(GradientDirectionContainerType *, const GradientImagesType *gradientImage)
typename OutputImageType::RegionType OutputImageRegionType
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
void SetMaskSpatialObject(MaskSpatialObjectType *maskSpatialObject)
virtual GradientDirectionType GetGradientDirection(unsigned int idx) const
ImageToImageFilter< Image< TReferenceImagePixelType, 3 >, Image< DiffusionTensor3D< TTensorPixelType >, 3 > > Superclass
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
~DiffusionTensor3DReconstructionImageFilter() override=default
void SetMaskImage(MaskImageType *maskImage)
DiffusionTensor3DReconstructionImageFilterEnums::GradientImageFormat GradientImageTypeEnumeration
const GradientImageType * GetGradientImage(unsigned int index) const
VectorContainer< unsigned int, GradientDirectionType > GradientDirectionContainerType
Represent a diffusion tensor as used in DTI images.
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
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.
#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