]> Creatis software - clitk.git/blobdiff - itk/clitkExtractSliceFilter.txx
Allow to display in all directions for images with size 2 and 3 in 3rd direction
[clitk.git] / itk / clitkExtractSliceFilter.txx
index 40f7982b71e62445477f9b102f21f76a848f3f44..af008ec3e9379e6d8f60bfd8cc51123285e26a64 100644 (file)
@@ -3,7 +3,7 @@
 
   Authors belong to: 
   - University of LYON              http://www.universite-lyon.fr/
-  - Léon Bérard cancer center       http://oncora1.lyon.fnclcc.fr
+  - Léon Bérard cancer center       http://www.centreleonberard.fr
   - CREATIS CNRS laboratory         http://www.creatis.insa-lyon.fr
 
   This software is distributed WITHOUT ANY WARRANTY; without even
 
   - BSD        See included LICENSE.txt file
   - CeCILL-B   http://www.cecill.info/licences/Licence_CeCILL-B_V1-en.html
-  ======================================================================-====*/
+  ===========================================================================**/
 
 // clitk
 #include "clitkCommon.h"
 
+// itk
+#include <itkExtractImageFilter.h>
 
 //--------------------------------------------------------------------
 template <class ImageType>
 clitk::ExtractSliceFilter<ImageType>::
 ExtractSliceFilter():
   clitk::FilterBase(),
-  itk::ImageToImageFilter<ImageType, ImageType>()
+  Superclass()
 {
   this->SetNumberOfRequiredInputs(1);
   SetDirection(2);
@@ -37,35 +39,41 @@ ExtractSliceFilter():
 template <class ImageType>
 void 
 clitk::ExtractSliceFilter<ImageType>::
-SetInput(const ImageType * image) {
+SetInput(const ImageType * image) 
+{
   // Process object is not const-correct so the const casting is required.
   this->SetNthInput(0, const_cast<ImageType *>(image));
 }
 //--------------------------------------------------------------------
   
 
+//--------------------------------------------------------------------
+template <class ImageType>
+void
+clitk::ExtractSliceFilter<ImageType>::
+GetOutputSlices(std::vector<typename SliceType::Pointer> & o)
+{
+  o.clear();
+  for(unsigned int i=0; i<this->GetNumberOfOutputs(); i++) {
+    o.push_back(this->GetOutput(i));
+  }
+}
+//--------------------------------------------------------------------
+  
+
 //--------------------------------------------------------------------
 template <class ImageType>
 void 
 clitk::ExtractSliceFilter<ImageType>::
-GenerateOutputInformation() 
-  DD("GenerateOutputInformation");
+GenerateOutputInformation() 
+{ 
   ImagePointer input = dynamic_cast<ImageType*>(itk::ProcessObject::GetInput(0));
-  //   ImagePointer outputImage = this->GetOutput(0);
-  //   outputImage->SetRegions(input->GetLargestPossibleRegion());
-  
-  output = this->GetOutput(0);
   
-  // create vector
-  typename SliceType::RegionType SliceRegionType;
-  typename SliceType::SizeType SliceSizeType;
-  typename SliceType::IndexType SliceIndexType;
-  // SliceRegionType region;
-
-  // create region
-  // loop ExtractImageFilter with region updated, push_back
-
-
+  // Create region to extract
+  m_region = input->GetLargestPossibleRegion();
+  m_size   = m_region.GetSize();
+  m_index  = m_region.GetIndex();
+  m_NumberOfSlices = m_region.GetSize()[GetDirection()];
 }
 //--------------------------------------------------------------------
 
@@ -75,9 +83,8 @@ template <class ImageType>
 void 
 clitk::ExtractSliceFilter<ImageType>::
 GenerateInputRequestedRegion() {
-  DD("GenerateInputRequestedRegion");
   // Call default
-  itk::ImageToImageFilter<ImageType, ImageType>::GenerateInputRequestedRegion();
+  Superclass::GenerateInputRequestedRegion();
   // Get input pointers and set requested region to common region
   ImagePointer input = dynamic_cast<ImageType*>(itk::ProcessObject::GetInput(0));
   input->SetRequestedRegion(input->GetLargestPossibleRegion());
@@ -89,17 +96,32 @@ template <class ImageType>
 void 
 clitk::ExtractSliceFilter<ImageType>::
 GenerateData() {
-  DD("GenerateData");
-
+  //--------------------------------------------------------------------
   // Get input pointer
   input = dynamic_cast<ImageType*>(itk::ProcessObject::GetInput(0));
 
-
+  //--------------------------------------------------------------------
+  // Create region to extract in 3D (and 2D = not used)
+  m_size[GetDirection()] = 0;
+  m_region.SetSize(m_size);
+  int start = m_index[GetDirection()];
+  this->SetNumberOfIndexedInputs(m_NumberOfSlices);
 
   //--------------------------------------------------------------------
-  //--------------------------------------------------------------------  
-  // Final Step -> set output
-  //this->SetNthOutput(0, working_image);
+  // loop ExtractImageFilter with region updated, push_back
+  typedef itk::ExtractImageFilter<ImageType, SliceType> ExtractImageFilterType;
+  typename ExtractImageFilterType::Pointer extract = ExtractImageFilterType::New();
+  extract->SetInput(input);
+  for(int i=0; i<m_NumberOfSlices; i++) {
+    extract = ExtractImageFilterType::New();
+    extract->SetInput(input);
+    m_index[GetDirection()] = start + i;
+    m_region.SetIndex(m_index);
+    extract->SetExtractionRegion(m_region);
+    extract->SetDirectionCollapseToSubmatrix();
+    extract->Update();
+    this->SetNthOutput(i, extract->GetOutput());
+  }
   return;
 }
 //--------------------------------------------------------------------