ANIMA  4.0
animaImageDataSplitter.hxx
Go to the documentation of this file.
1 #pragma once
3 
4 #include <itkImageRegionConstIteratorWithIndex.h>
5 #include <itkImageRegionConstIterator.h>
6 #include <itkImageRegionIterator.h>
7 
8 #include <itkMacro.h>
9 
10 namespace anima
11 {
12 template <typename TInputImage> ImageDataSplitter<TInputImage>::ImageDataSplitter()
13 {
14  m_NbImages = 0;
15  m_NeedsUpdate = true;
16 
17  m_NbBlocks.Fill (0);
18  m_Block.Fill (0);
19  m_Margin.Fill (0);
20 
21  for (unsigned int i = 0;i < m_GlobalRegionOfInterest.GetImageDimension();++i)
22  {
23  m_GlobalRegionOfInterest.SetIndex(i,0);
24  m_GlobalRegionOfInterest.SetSize(i,0);
25 
26  m_BlockRegion.SetIndex(i,0);
27  m_BlockRegion.SetSize(i,0);
28 
29  m_BlockRegionWithMargin.SetIndex(i,0);
30  m_BlockRegionWithMargin.SetSize(i,0);
31  }
32 
33  m_Images.clear();
34  m_FileNames.clear();
35 }
36 
37 template <typename TInputImage> ImageDataSplitter<TInputImage>::~ImageDataSplitter()
38 {
39  m_Images.clear();
40  m_FileNames.clear();
41 }
42 
43 template <typename TInputImage> void ImageDataSplitter<TInputImage>::SetUniqueFileName(std::string &inputFileName)
44 {
45  m_FileNames.clear();
46  m_FileNames.push_back(inputFileName);
47  m_NbImages = 1;
48 }
49 
50 template <typename TInputImage> void ImageDataSplitter<TInputImage>::SetFileNames(std::string &inputFileList)
51 {
52  std::ifstream fileIn(inputFileList.c_str());
53 
54  if (!fileIn.is_open())
55  {
56  std::string error("Invalid input file list. ");
57  error += inputFileList;
58  error += " could not be opened...";
59 
60  throw itk::ExceptionObject(__FILE__, __LINE__,error,ITK_LOCATION);
61  }
62 
63  m_FileNames.clear();
64 
65  char tmpStr[2048];
66  while (!fileIn.eof())
67  {
68  fileIn.getline(tmpStr,2048);
69 
70  if (strcmp(tmpStr,"") != 0)
71  m_FileNames.push_back(std::string(tmpStr));
72  }
73 
74  fileIn.close();
75  m_NeedsUpdate = true;
76  m_NbImages = m_FileNames.size();
77 }
78 
79 template <typename TInputImage> void ImageDataSplitter<TInputImage>::SetComputationMask(MaskImageType::Pointer &maskIm)
80 {
81  m_MaskImage = maskIm;
82  typedef itk::ImageRegionConstIteratorWithIndex< MaskImageType > MaskIteratorType;
83 
84  MaskIteratorType maskIterator(m_MaskImage,m_MaskImage->GetLargestPossibleRegion());
85 
86  bool foundOneNonNull = false;
87  while (!maskIterator.IsAtEnd())
88  {
89  if (maskIterator.Get() == 0)
90  {
91  ++maskIterator;
92  continue;
93  }
94 
95  if (foundOneNonNull)
96  {
97  for (unsigned int i = 0;i < MaskImageType::GetImageDimension();++i)
98  {
99  if (m_GlobalRegionOfInterest.GetIndex()[i] > maskIterator.GetIndex()[i])
100  {
101  unsigned int newSize = m_GlobalRegionOfInterest.GetSize()[i] - maskIterator.GetIndex()[i] + m_GlobalRegionOfInterest.GetIndex()[i];
102  m_GlobalRegionOfInterest.SetSize(i,newSize);
103  m_GlobalRegionOfInterest.SetIndex(i,maskIterator.GetIndex()[i]);
104  }
105  }
106 
107  for (unsigned int i = 0;i < MaskImageType::GetImageDimension();++i)
108  {
109  unsigned int tmpNum = m_GlobalRegionOfInterest.GetIndex()[i] + m_GlobalRegionOfInterest.GetSize()[i];
110  if (tmpNum <= maskIterator.GetIndex()[i])
111  m_GlobalRegionOfInterest.SetSize(i,maskIterator.GetIndex()[i] - m_GlobalRegionOfInterest.GetIndex()[i] + 1);
112  }
113  }
114  else
115  {
116  m_GlobalRegionOfInterest.SetIndex(maskIterator.GetIndex());
117  for (unsigned int i = 0;i < m_GlobalRegionOfInterest.GetImageDimension();++i)
118  m_GlobalRegionOfInterest.SetSize(i,1);
119 
120  foundOneNonNull = true;
121  }
122 
123  ++maskIterator;
124  }
125 
126  m_NeedsUpdate = true;
127 }
128 
129 template <typename TInputImage> void ImageDataSplitter<TInputImage>::SetNumberOfBlocks(TInputIndexType &bNumBlocks)
130 {
131  if (m_NbBlocks != bNumBlocks)
132  {
133  m_NeedsUpdate = true;
134  m_NbBlocks = bNumBlocks;
135  }
136 }
137 
138 template <typename TInputImage> void ImageDataSplitter<TInputImage>::SetBlockIndex(TInputIndexType &bIndex)
139 {
140  if (m_Block != bIndex)
141  {
142  m_NeedsUpdate = true;
143  m_Block = bIndex;
144  }
145 }
146 
147 template <typename TInputImage> void ImageDataSplitter<TInputImage>::SetBlockMargin(TInputIndexType &bMargin)
148 {
149  if (m_Margin != bMargin)
150  {
151  m_NeedsUpdate = true;
152  m_Margin = bMargin;
153  }
154 }
155 
156 template <typename TInputImage> typename TInputImage::RegionType ImageDataSplitter<TInputImage>::GetSpecificBlockRegion(TInputIndexType &block)
157 {
158  TInputRegionType tmpBlock;
159 
160  for (unsigned int i = 0;i < m_GlobalRegionOfInterest.GetImageDimension();++i)
161  {
162  tmpBlock.SetIndex(i,m_GlobalRegionOfInterest.GetIndex(i) + m_GlobalRegionOfInterest.GetSize(i)*block[i]/m_NbBlocks[i]);
163  unsigned int tmpMax = m_GlobalRegionOfInterest.GetIndex(i) + m_GlobalRegionOfInterest.GetSize(i)*(1+block[i])/m_NbBlocks[i] - tmpBlock.GetIndex(i);
164  tmpBlock.SetSize(i,tmpMax);
165  }
166 
167  return tmpBlock;
168 }
169 
170 template <typename TInputImage> bool ImageDataSplitter<TInputImage>::EmptyMask(TInputIndexType &bIndex)
171 {
172  if (!m_MaskImage)
173  throw itk::ExceptionObject(__FILE__, __LINE__,"No mask input. This is required. Exiting...",ITK_LOCATION);
174 
175  TInputRegionType tmpBlock = this->GetSpecificBlockRegion(bIndex);
176  itk::ImageRegionConstIteratorWithIndex < MaskImageType > maskIt(m_MaskImage,tmpBlock);
177 
178  bool hasNonNullValue = false;
179  while (!maskIt.IsAtEnd())
180  {
181  if (maskIt.Get() != 0)
182  {
183  hasNonNullValue = true;
184  break;
185  }
186 
187  ++maskIt;
188  }
189 
190  return (!hasNonNullValue);
191 }
192 
193 template <typename TInputImage> void ImageDataSplitter<TInputImage>::Update()
194 {
195  if (!m_NeedsUpdate)
196  return;
197 
198  m_Images.clear();
199 
200  if (!m_MaskImage)
201  throw itk::ExceptionObject(__FILE__, __LINE__,"No mask input. This is required. Exiting...",ITK_LOCATION);
202 
203  for (unsigned int i = 0;i < m_GlobalRegionOfInterest.GetImageDimension();++i)
204  {
205  m_BlockRegion.SetIndex(i,m_GlobalRegionOfInterest.GetIndex(i) + m_GlobalRegionOfInterest.GetSize(i)*m_Block[i]/m_NbBlocks[i]);
206  unsigned int tmpMax = m_GlobalRegionOfInterest.GetIndex(i) + m_GlobalRegionOfInterest.GetSize(i)*(1+m_Block[i])/m_NbBlocks[i] - m_BlockRegion.GetIndex(i);
207  m_BlockRegion.SetSize(i,tmpMax);
208  }
209 
210  m_BlockRegionWithMargin = m_BlockRegion;
211  for (unsigned int i = 0;i < m_GlobalRegionOfInterest.GetImageDimension();++i)
212  {
213  if (m_Margin[i] != 0)
214  {
215  unsigned int tmpMin = m_BlockRegion.GetIndex()[i] - m_Margin[i];
216  if (tmpMin < m_MaskImage->GetLargestPossibleRegion().GetIndex()[i])
217  tmpMin = m_MaskImage->GetLargestPossibleRegion().GetIndex()[i];
218  unsigned int tmpMax = m_BlockRegion.GetIndex()[i] + m_BlockRegion.GetSize()[i] + m_Margin[i] - 1;
219  if (tmpMax >= m_MaskImage->GetLargestPossibleRegion().GetIndex()[i] + m_MaskImage->GetLargestPossibleRegion().GetSize()[i])
220  tmpMax = m_MaskImage->GetLargestPossibleRegion().GetIndex()[i] + m_MaskImage->GetLargestPossibleRegion().GetSize()[i] - 1;
221 
222  m_BlockRegionWithMargin.SetIndex(i,tmpMin);
223  m_BlockRegionWithMargin.SetSize(i,tmpMax - tmpMin + 1);
224  }
225  }
226 
227  MaskImageType::RegionType tmpRegion = m_BlockRegion;
228  for (unsigned int i = 0;i < MaskImageType::GetImageDimension();++i)
229  tmpRegion.SetIndex(i,0);
230 
231  m_SmallMask = MaskImageType::New();
232  m_SmallMask->Initialize();
233  m_SmallMask->SetRegions(tmpRegion);
234  m_SmallMask->SetOrigin(m_MaskImage->GetOrigin());
235  m_SmallMask->SetDirection(m_MaskImage->GetDirection());
236  m_SmallMask->SetSpacing(m_MaskImage->GetSpacing());
237  m_SmallMask->Allocate();
238 
239  itk::ImageRegionIterator <MaskImageType> smallIt(m_SmallMask,tmpRegion);
240  itk::ImageRegionConstIterator <MaskImageType> maskIt(m_MaskImage,m_BlockRegion);
241 
242  while (!maskIt.IsAtEnd())
243  {
244  smallIt.Set(maskIt.Get());
245  ++smallIt;
246  ++maskIt;
247  }
248 
249  for (unsigned int i = 0;i < MaskImageType::GetImageDimension();++i)
250  tmpRegion.SetSize(i,m_BlockRegionWithMargin.GetSize()[i]);
251 
252  m_SmallMaskWithMargin = MaskImageType::New();
253  m_SmallMaskWithMargin->Initialize();
254  m_SmallMaskWithMargin->SetRegions(tmpRegion);
255  m_SmallMaskWithMargin->SetOrigin(m_MaskImage->GetOrigin());
256  m_SmallMaskWithMargin->SetDirection(m_MaskImage->GetDirection());
257  m_SmallMaskWithMargin->SetSpacing(m_MaskImage->GetSpacing());
258  m_SmallMaskWithMargin->Allocate();
259 
260  itk::ImageRegionIterator <MaskImageType> smallItWM(m_SmallMaskWithMargin,tmpRegion);
261  itk::ImageRegionConstIterator <MaskImageType> maskItWM(m_MaskImage,m_BlockRegionWithMargin);
262 
263  while (!maskItWM.IsAtEnd())
264  {
265  smallItWM.Set(maskItWM.Get());
266  ++smallItWM;
267  ++maskItWM;
268  }
269 
270  for (unsigned int i = 0;i < m_FileNames.size();++i)
271  {
272  std::cout << "Processing image file " << m_FileNames[i] << "..." << std::endl;
273  InputReaderPointer tmpImReader = InputReaderType::New();
274  tmpImReader->SetFileName(m_FileNames[i]);
275  tmpImReader->Update();
276 
277  m_Images.push_back(TInputImage::New());
278  m_Images[i]->Initialize();
279  m_Images[i]->SetRegions(tmpRegion);
280  m_Images[i]->SetOrigin(m_MaskImage->GetOrigin());
281  m_Images[i]->SetDirection(m_MaskImage->GetDirection());
282  m_Images[i]->SetSpacing(m_MaskImage->GetSpacing());
283 
284  m_Images[i]->SetNumberOfComponentsPerPixel(tmpImReader->GetOutput()->GetNumberOfComponentsPerPixel());
285 
286  m_Images[i]->Allocate();
287 
288  itk::ImageRegionIterator <TInputImage> cropImIt(m_Images[i],tmpRegion);
289  itk::ImageRegionConstIterator <TInputImage> reImIt(tmpImReader->GetOutput(),m_BlockRegionWithMargin);
290 
291  while (!cropImIt.IsAtEnd())
292  {
293  cropImIt.Set(reImIt.Get());
294 
295  ++cropImIt;
296  ++reImIt;
297  }
298  }
299 
300  m_NeedsUpdate = false;
301 }
302 
303 template <typename TInputImage> typename TInputImage::RegionType ImageDataSplitter<TInputImage>::GetBlockRegionInsideMargin()
304 {
305  if (m_NeedsUpdate)
306  this->Update();
307 
308  TInputRegionType bRegInsideMargin = m_BlockRegion;
309 
310  for (unsigned int i = 0;i < TInputImage::GetImageDimension();++i)
311  bRegInsideMargin.SetIndex(i,m_BlockRegion.GetIndex()[i] - m_BlockRegionWithMargin.GetIndex()[i]);
312 
313  return bRegInsideMargin;
314 }
315 
316 } // end namespace anima
void SetUniqueFileName(std::string &inputFileName)
TInputRegionType GetBlockRegionInsideMargin()
TInputImage::RegionType TInputRegionType
void SetFileNames(std::string &inputFileList)
void SetComputationMask(MaskImageType::Pointer &maskIm)
void SetBlockMargin(TInputIndexType &bMargin)
void SetBlockIndex(TInputIndexType &bIndex)
void SetNumberOfBlocks(TInputIndexType &bNumBlocks)
InputReaderType::Pointer InputReaderPointer
TInputRegionType GetSpecificBlockRegion(TInputIndexType &block)
bool EmptyMask(TInputIndexType &bIndex)
TInputImage::IndexType TInputIndexType