ANIMA  4.0
animaLowMemLocalPatchMeanDistanceBridge.cxx
Go to the documentation of this file.
2 #include <itkImageFileReader.h>
3 #include <itkImageFileWriter.h>
4 
5 namespace anima
6 {
7 
9 {
10  m_NbSplits = 2;
11  m_NumThreads = itk::MultiThreaderBase::GetGlobalDefaultNumberOfThreads();
12 
13  m_DatabaseImages = new ImageSplitterType;
14  m_ComputationMask = NULL;
15 
16  m_PatchHalfSize = 1;
17 }
18 
20 {
21  if (m_DatabaseImages)
22  delete m_DatabaseImages;
23 }
24 
26 {
27  itk::ImageFileReader <MaskImageType>::Pointer tmpRead = itk::ImageFileReader <MaskImageType>::New();
28  tmpRead->SetFileName(cMask);
29  tmpRead->Update();
30 
31  m_ComputationMask = tmpRead->GetOutput();
32 
33  m_DatabaseImages->SetComputationMask(m_ComputationMask);
34 }
35 
36 void LowMemoryLocalPatchMeanDistanceBridge::Update(int specificSplitToDo, bool genOutputDescriptionData)
37 {
38  if (!m_ComputationMask)
39  itkExceptionMacro("No computation mask... Exiting...");
40 
42  for (unsigned int i = 0;i < MaskImageType::GetImageDimension();++i)
43  tmpInd[i] = m_NbSplits;
44 
46 
47  for (unsigned int i = 0;i < MaskImageType::GetImageDimension();++i)
48  marginBlock[i] = m_PatchHalfSize;
49 
50  m_DatabaseImages->SetBlockMargin(marginBlock);
51  m_DatabaseImages->SetNumberOfBlocks(tmpInd);
52 
53  std::vector < ImageSplitterType::TInputIndexType > splitIndexesToProcess;
54 
55  if (specificSplitToDo != -1)
56  {
57  tmpInd[0] = (unsigned int)floor((double)(specificSplitToDo/(m_NbSplits*m_NbSplits)));
58  unsigned int tmpVal = specificSplitToDo - tmpInd[0]*m_NbSplits*m_NbSplits;
59  tmpInd[1] = (unsigned int)floor((double)(tmpVal/m_NbSplits));
60  tmpInd[2] = tmpVal - tmpInd[1]*m_NbSplits;
61 
62  if (!m_DatabaseImages->EmptyMask(tmpInd))
63  splitIndexesToProcess.push_back(tmpInd);
64  }
65  else
66  {
67  for (unsigned int i = 0;i < m_NbSplits;++i)
68  {
69  tmpInd[0] = i;
70  for (unsigned int j = 0;j < m_NbSplits;++j)
71  {
72  tmpInd[1] = j;
73  for (unsigned int k = 0;k < m_NbSplits;++k)
74  {
75  tmpInd[2] = k;
76 
77  if (!m_DatabaseImages->EmptyMask(tmpInd))
78  splitIndexesToProcess.push_back(tmpInd);
79  }
80  }
81  }
82  }
83 
84  for (unsigned int i = 0;i < splitIndexesToProcess.size();++i)
85  {
86  std::cout << "Processing block : " << splitIndexesToProcess[i][0] << " "
87  << splitIndexesToProcess[i][1] << " " << splitIndexesToProcess[i][2] << std::endl;
88 
89  m_DatabaseImages->SetBlockIndex(splitIndexesToProcess[i]);
90  m_DatabaseImages->Update();
91 
93 
94  for (unsigned int j = 0;j < m_DatabaseImages->GetNbImages();++j)
95  mainFilter->SetInput(j,m_DatabaseImages->GetOutput(j));
96 
97  mainFilter->SetComputationMask(m_DatabaseImages->GetSmallMaskWithMargin());
98  mainFilter->SetNumberOfWorkUnits(m_NumThreads);
99 
100  mainFilter->SetPatchHalfSize(m_PatchHalfSize);
101 
102  mainFilter->Update();
103 
104  std::cout << "Results computed... Writing output parcel..." << std::endl;
105 
106  char numSplit[2048];
107  sprintf(numSplit,"_%ld_%ld_%ld.nrrd",splitIndexesToProcess[i][0],splitIndexesToProcess[i][1],splitIndexesToProcess[i][2]);
108 
109  //Write outputs
110  std::string outputMeanName = m_OutputMeanName + numSplit;
111  std::string outputStdName = m_OutputStdName + numSplit;
112  this->BuildAndWrite(mainFilter->GetOutput(0),outputMeanName,m_DatabaseImages->GetBlockRegionInsideMargin());
113 
114  if (m_OutputStdName != "")
115  this->BuildAndWrite(mainFilter->GetOutput(1),outputStdName,m_DatabaseImages->GetBlockRegionInsideMargin());
116  }
117 
118  if (genOutputDescriptionData)
119  {
120  std::vector <std::ofstream> mainOutFiles;
121 
122  std::string tmpOutStdName = m_OutputStdName + ".txt";
123  std::ofstream tmpFileStdOut;
124 
125  if (m_OutputStdName != "")
126  tmpFileStdOut.open(tmpOutStdName.c_str());
127 
128  std::string tmpOutMeanName = m_OutputMeanName + ".txt";
129  std::ofstream tmpFileMeanOut(tmpOutMeanName.c_str());
130 
131  for (unsigned int i = 0;i < m_NbSplits;++i)
132  {
133  tmpInd[0] = i;
134  for (unsigned int j = 0;j < m_NbSplits;++j)
135  {
136  tmpInd[1] = j;
137  for (unsigned int k = 0;k < m_NbSplits;++k)
138  {
139  tmpInd[2] = k;
140 
141  if (!m_DatabaseImages->EmptyMask(tmpInd))
142  {
143  OutputImageType::RegionType tmpBlRegion = m_DatabaseImages->GetSpecificBlockRegion(tmpInd);
144  char numSplit[2048];
145  sprintf(numSplit,"_%ld_%ld_%ld.nrrd",tmpInd[0],tmpInd[1],tmpInd[2]);
146 
147  if (tmpFileStdOut.is_open())
148  {
149  tmpFileStdOut << "<BLOCK>" << std::endl;
150  tmpFileStdOut << "BLOCK_FILE=" << m_OutputStdName + numSplit << std::endl;
151  tmpFileStdOut << "STARTING_INDEX=" << tmpBlRegion.GetIndex()[0] << " "
152  << tmpBlRegion.GetIndex()[1] << " " << tmpBlRegion.GetIndex()[2] << std::endl;
153  tmpFileStdOut << "</BLOCK>" << std::endl;
154  }
155 
156  tmpFileMeanOut << "<BLOCK>" << std::endl;
157  tmpFileMeanOut << "BLOCK_FILE=" << m_OutputMeanName + numSplit << std::endl;
158  tmpFileMeanOut << "STARTING_INDEX=" << tmpBlRegion.GetIndex()[0] << " "
159  << tmpBlRegion.GetIndex()[1] << " " << tmpBlRegion.GetIndex()[2] << std::endl;
160  tmpFileMeanOut << "</BLOCK>" << std::endl;
161  }
162  }
163  }
164  }
165 
166  tmpFileMeanOut.close();
167  tmpFileStdOut.close();
168  }
169 }
170 
172  OutputImageType::RegionType finalROI)
173 {
174  OutputImageType::RegionType tmpRegion = finalROI;
175  for (unsigned int i = 0;i < OutputImageType::GetImageDimension();++i)
176  tmpRegion.SetIndex(i,0);
177 
178  OutputImageType::Pointer tmpRes = OutputImageType::New();
179  tmpRes->Initialize();
180 
181  tmpRes->SetOrigin(m_ComputationMask->GetOrigin());
182  tmpRes->SetRegions(tmpRegion);
183  tmpRes->SetDirection(m_ComputationMask->GetDirection());
184  tmpRes->SetSpacing(m_ComputationMask->GetSpacing());
185 
186  tmpRes->Allocate();
187 
188  itk::ImageRegionIterator <OutputImageType> tmpImIt (tmpIm,finalROI);
189  itk::ImageRegionIterator <OutputImageType> tmpResIt (tmpRes,tmpRegion);
190 
191  while (!tmpImIt.IsAtEnd())
192  {
193  tmpResIt.Set(tmpImIt.Get());
194 
195  ++tmpImIt;
196  ++tmpResIt;
197  }
198 
199  itk::ImageFileWriter <OutputImageType>::Pointer outWriter = itk::ImageFileWriter <OutputImageType>::New();
200  outWriter->SetInput(tmpRes);
201  outWriter->SetFileName(resName);
202  outWriter->SetUseCompression(true);
203 
204  outWriter->Update();
205 }
206 
207 } // end namesapce itk
TInputRegionType GetBlockRegionInsideMargin()
MaskImageType * GetSmallMaskWithMargin()
void SetComputationMask(MaskImageType::Pointer &maskIm)
void BuildAndWrite(OutputImageType *tmpIm, std::string resName, OutputImageType::RegionType finalROI)
void SetBlockMargin(TInputIndexType &bMargin)
void SetBlockIndex(TInputIndexType &bIndex)
void SetNumberOfBlocks(TInputIndexType &bNumBlocks)
LocalPatchMeanDistanceImageFilter< double >::OutputImageType OutputImageType
TInputImage * GetOutput(unsigned int i)
void Update(int specificSplitToDo=-1, bool genOutputDescriptionData=false)
TInputRegionType GetSpecificBlockRegion(TInputIndexType &block)
bool EmptyMask(TInputIndexType &bIndex)
anima::ImageDataSplitter< InputImageType > ImageSplitterType
TInputImage::IndexType TInputIndexType