ANIMA  4.0
animaGMMT2RelaxometryEstimation.cxx
Go to the documentation of this file.
1 #include <iostream>
2 #include <tclap/CmdLine.h>
3 
6 #include <itkTimeProbe.h>
7 
8 #include <itkCommand.h>
9 
10 //Update progression of the process
11 void eventCallback (itk::Object* caller, const itk::EventObject& event, void* clientData)
12 {
13  itk::ProcessObject * processObject = (itk::ProcessObject*) caller;
14  std::cout<<"\033[K\rProgression: "<<(int)(processObject->GetProgress() * 100)<<"%"<<std::flush;
15 }
16 
17 int main(int argc, char **argv)
18 {
19  TCLAP::CmdLine cmd("INRIA / IRISA - VisAGeS/Empenn Team", ' ',ANIMA_VERSION);
20 
21  TCLAP::ValueArg<std::string> t2Arg("i","t2","T2 relaxometry images (list or 4D volume)",true,"","T2 relaxometry images",cmd);
22  TCLAP::ValueArg<std::string> maskArg("m","maskname","Computation mask",false,"","computation mask",cmd);
23 
24  TCLAP::ValueArg<std::string> t1MapArg("","t1","T1 map",false,"","T1 map",cmd);
25  TCLAP::ValueArg<std::string> gaussMeansArg("g","gauss-mean","Text file with Gaussian means",false,"","Gaussian means",cmd);
26  TCLAP::ValueArg<std::string> gaussVarsArg("G","gauss-var","Text file with Gaussian variances",false,"","Gaussian variances",cmd);
27 
28  TCLAP::ValueArg<std::string> resWeightsArg("O","out-weights","Result weights image",false,"","result weights image",cmd);
29  TCLAP::ValueArg<std::string> resM0Arg("","out-m0","Result M0 image",false,"","result M0 image",cmd);
30  TCLAP::ValueArg<std::string> resMWFArg("o","out-mwf","Result MWF image",true,"","result MWF image",cmd);
31  TCLAP::ValueArg<std::string> resB1Arg("","out-b1","Result B1 image",false,"","result B1 image",cmd);
32  TCLAP::ValueArg<std::string> resSigmaSqArg("","out-sig","Result sigma square image",false,"","result sigma square image",cmd);
33 
34  TCLAP::ValueArg<double> echoSpacingArg("e","echo-spacing","Spacing between two successive echoes (default: 10)",false,10,"Spacing between echoes",cmd);
35  TCLAP::ValueArg<double> excitationT2FlipAngleArg("","t2-ex-flip","Excitation flip angle for T2 (in degrees, default: 90)",false,90,"T2 excitation flip angle",cmd);
36  TCLAP::ValueArg<double> t2FlipAngleArg("","t2-flip","All flip angles for T2 (in degrees, default: 180)",false,180,"T2 flip angle",cmd);
37  TCLAP::ValueArg<double> backgroundSignalThresholdArg("t","signal-thr","Background signal threshold (default: 10)",false,10,"Background signal threshold",cmd);
38 
39  TCLAP::ValueArg<double> gaussianToleranceApproxArg("","g-tol","Gaussian approximation tolerance (default: 1.0e-8)",false,1.0e-8,"Gaussian approximation tolerance",cmd);
40 
41  TCLAP::ValueArg<unsigned int> nbpArg("T","numberofthreads","Number of threads to run on (default : all cores)",false,itk::MultiThreaderBase::GetGlobalDefaultNumberOfThreads(),"number of threads",cmd);
42 
43  try
44  {
45  cmd.parse(argc,argv);
46  }
47  catch (TCLAP::ArgException& e)
48  {
49  std::cerr << "Error: " << e.error() << "for argument " << e.argId() << std::endl;
50  return EXIT_FAILURE;
51  }
52 
53  typedef itk::Image <double,3> InputImageType;
54  typedef InputImageType OutputImageType;
56  typedef FilterType::OutputImageType OutputImageType;
57  typedef FilterType::VectorOutputImageType VectorOutputImageType;
58  typedef itk::ImageFileReader <InputImageType> InputImageReaderType;
59 
60  FilterType::Pointer mainFilter = FilterType::New();
61 
62  unsigned int numInputs = anima::setMultipleImageFilterInputsFromFileName<InputImageType,FilterType>(t2Arg.getValue(),mainFilter);
63 
64  mainFilter->SetEchoSpacing(echoSpacingArg.getValue());
65  mainFilter->SetT2FlipAngles(t2FlipAngleArg.getValue() * M_PI / 180.0,numInputs);
66  mainFilter->SetT2ExcitationFlipAngle(excitationT2FlipAngleArg.getValue() * M_PI / 180.0);
67 
68  mainFilter->SetGaussianIntegralTolerance(gaussianToleranceApproxArg.getValue());
69 
70  if (gaussMeansArg.getValue() != "")
71  mainFilter->SetGaussianMeans(gaussMeansArg.getValue());
72 
73  if (gaussVarsArg.getValue() != "")
74  mainFilter->SetGaussianVariances(gaussVarsArg.getValue());
75 
76  if (t1MapArg.getValue() != "")
77  {
78  InputImageReaderType::Pointer t1MapRead = InputImageReaderType::New();
79  t1MapRead->SetFileName(t1MapArg.getValue().c_str());
80  t1MapRead->Update();
81 
82  mainFilter->SetT1Map(t1MapRead->GetOutput());
83  }
84 
85  if (maskArg.getValue() != "")
86  {
87  typedef itk::ImageFileReader < itk::Image <unsigned char, 3> > itkMaskReader;
88  itkMaskReader::Pointer maskRead = itkMaskReader::New();
89  maskRead->SetFileName(maskArg.getValue().c_str());
90  maskRead->Update();
91 
92  mainFilter->SetComputationMask(maskRead->GetOutput());
93  }
94 
95  mainFilter->SetAverageSignalThreshold(backgroundSignalThresholdArg.getValue());
96  mainFilter->SetNumberOfWorkUnits(nbpArg.getValue());
97 
98  itk::CStyleCommand::Pointer callback = itk::CStyleCommand::New();
99  callback->SetCallback(eventCallback);
100  mainFilter->AddObserver(itk::ProgressEvent(), callback );
101 
102  itk::TimeProbe tmpTime;
103  tmpTime.Start();
104 
105  try
106  {
107  mainFilter->Update();
108  }
109  catch (itk::ExceptionObject &e)
110  {
111  std::cerr << e << std::endl;
112  return EXIT_FAILURE;
113  }
114 
115  tmpTime.Stop();
116 
117  std::cout << "\nRegular estimation computation time: " << tmpTime.GetTotal() << std::endl;
118 
119  anima::writeImage<OutputImageType> (resMWFArg.getValue(),mainFilter->GetMWFOutputImage());
120 
121  if (resWeightsArg.getValue() != "")
122  anima::writeImage<VectorOutputImageType> (resWeightsArg.getValue(),mainFilter->GetWeightsImage());
123 
124  if (resM0Arg.getValue() != "")
125  anima::writeImage<OutputImageType> (resM0Arg.getValue(),mainFilter->GetM0OutputImage());
126 
127  if (resB1Arg.getValue() != "")
128  anima::writeImage<OutputImageType> (resB1Arg.getValue(),mainFilter->GetB1OutputImage());
129 
130  if (resSigmaSqArg.getValue() != "")
131  anima::writeImage<OutputImageType> (resSigmaSqArg.getValue(),mainFilter->GetSigmaSquareOutputImage());
132 
133  return EXIT_SUCCESS;
134 }
void eventCallback(itk::Object *caller, const itk::EventObject &event, void *clientData)
int main(int argc, char **argv)