ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkDiffusionTensor3DReconstructionImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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  * http://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 
31 namespace itk
32 {
122 template< typename TReferenceImagePixelType,
123  typename TGradientImagePixelType = TReferenceImagePixelType,
124  typename TTensorPixelType = double,
125  typename TMaskImageType = Image<unsigned char, 3> >
127  public ImageToImageFilter< Image< TReferenceImagePixelType, 3 >,
128  Image< DiffusionTensor3D< TTensorPixelType >, 3 > >
129 {
130 public:
131 
138 
140  itkNewMacro(Self);
141 
145 
146  typedef TReferenceImagePixelType ReferencePixelType;
147 
148  typedef TGradientImagePixelType GradientPixelType;
149 
151 
154  typedef typename Superclass::InputImageType ReferenceImageType;
155 
157 
159 
160  typedef typename Superclass::OutputImageRegionType
162 
165 
171 
174 
176  typedef TMaskImageType MaskImageType;
177 
179  typedef vnl_matrix_fixed< double, 6, 6 > TensorBasisMatrixType;
180 
181  typedef vnl_matrix< double > CoefficientMatrixType;
182 
184  typedef vnl_vector_fixed< double, 3 > GradientDirectionType;
185 
187  typedef VectorContainer< unsigned int,
189 
191  void AddGradientImage(const GradientDirectionType &, const GradientImageType *image);
192  const GradientImageType *GetGradientImage(unsigned index) const;
194 
201  void SetGradientImage(GradientDirectionContainerType *,
202  const GradientImagesType *image);
203 
205  void SetReferenceImage(ReferenceImageType *referenceImage)
206  {
207  if ( m_GradientImageTypeEnumeration == GradientIsInASingleImage )
208  {
209  itkExceptionMacro(<< "Cannot call both methods:"
210  << "AddGradientImage and SetGradientImage. Please call only one of them.");
211  }
212 
213  this->ProcessObject::SetNthInput(0, referenceImage);
214 
215  m_GradientImageTypeEnumeration = GradientIsInManyImages;
216  }
217 
220  { return ( static_cast< ReferenceImageType * >( this->ProcessObject::GetInput(0) ) ); }
221 
223  virtual GradientDirectionType GetGradientDirection(unsigned int idx) const
224  {
225  if ( idx >= m_NumberOfGradientDirections )
226  {
227  itkExceptionMacro(<< "Gradient direction " << idx << "does not exist");
228  }
229  return m_GradientDirectionContainer->ElementAt(idx + 1);
230  }
232 
234  void SetMaskImage(MaskImageType *maskImage);
235 
237  void SetMaskSpatialObject(MaskSpatialObjectType *maskSpatialObject);
238 
239 
243  itkSetMacro(Threshold, ReferencePixelType);
244  itkGetConstMacro(Threshold, ReferencePixelType);
246 
253  itkSetMacro(BValue, TTensorPixelType);
254 #ifdef GetBValue
255 #undef GetBValue
256 #endif
257  itkGetConstReferenceMacro(BValue, TTensorPixelType);
259 
260 #ifdef ITK_USE_CONCEPT_CHECKING
261  // Begin concept checking
262  itkConceptMacro( ReferenceEqualityComparableCheck,
264  itkConceptMacro( TensorEqualityComparableCheck,
266  itkConceptMacro( GradientConvertibleToDoubleCheck,
268  itkConceptMacro( DoubleConvertibleToTensorCheck,
270  itkConceptMacro( GradientReferenceAdditiveOperatorsCheck,
271  ( Concept::AdditiveOperators< GradientPixelType, GradientPixelType,
272  ReferencePixelType > ) );
273  itkConceptMacro( GradientReferenceAdditiveAndAssignOperatorsCheck,
274  ( Concept::AdditiveAndAssignOperators< GradientPixelType,
275  ReferencePixelType > ) );
276 
277  itkConceptMacro( ReferenceOStreamWritableCheck,
279  itkConceptMacro( TensorOStreamWritableCheck,
281  // End concept checking
282 #endif
283 
284 protected:
287  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
288 
289  void ComputeTensorBasis();
290 
291  void BeforeThreadedGenerateData() ITK_OVERRIDE;
292 
293  void ThreadedGenerateData(const
294  OutputImageRegionType & outputRegionForThread, ThreadIdType) ITK_OVERRIDE;
295 
296  void VerifyPreconditions() ITK_OVERRIDE;
297 
300  typedef enum {
301  GradientIsInASingleImage = 1,
303  Else
304  } GradientImageTypeEnumeration;
305 
306 private:
307 
308  /* Tensor basis coeffs */
310 
312 
315 
318 
321 
324 
326  TTensorPixelType m_BValue;
327 
330 
333 };
334 }
335 
336 #ifndef ITK_MANUAL_INSTANTIATION
337 #include "itkDiffusionTensor3DReconstructionImageFilter.hxx"
338 #endif
339 
340 #endif
Light weight base class for most itk classes.
virtual GradientDirectionType GetGradientDirection(unsigned int idx) const
Templated n-dimensional vector image class.
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
DataObject * GetInput(const DataObjectIdentifierType &key)
Return an input.
This class takes as input one or more reference image (acquired in the absence of diffusion sensitizi...
Base class for filters that take an image as input and produce an image as output.
Define a front-end to the STL &quot;vector&quot; container that conforms to the IndexedContainerInterface.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
virtual void SetNthInput(DataObjectPointerArraySizeType num, DataObject *input)
#define itkConceptMacro(name, concept)
Represent a diffusion tensor as used in DTI images.
ImageToImageFilter< Image< TReferenceImagePixelType, 3 >, Image< DiffusionTensor3D< TTensorPixelType >, 3 > > Superclass
Templated n-dimensional image class.
Definition: itkImage.h:75
VectorContainer< unsigned int, GradientDirectionType > GradientDirectionContainerType