ANIMA  4.0
animaFibersApplyTransformSerie.cxx
Go to the documentation of this file.
1 #include <tclap/CmdLine.h>
3 #include <animaShapesReader.h>
4 #include <animaShapesWriter.h>
5 
6 #include <vtkPolyData.h>
7 
9 {
10  typedef itk::Image <unsigned short, 3>::PointType PointType;
11  PointType pointPositionIn, pointPositionOut;
12  double pointPositionVTK[3];
13 
14  for (unsigned int i = 0;i < dataPoints->GetNumberOfPoints();++i)
15  {
16  for (unsigned int k = 0; k < 3; ++k)
17  pointPositionIn[k] = dataPoints->GetPoint(i)[k];
18 
19  pointPositionOut = transform->TransformPoint(pointPositionIn);
20 
21  for (unsigned int k = 0; k < 3; ++k)
22  pointPositionVTK[k] = pointPositionOut[k];
23 
24  dataPoints->SetPoint(i,pointPositionVTK);
25  }
26 }
27 
28 int main(int ac, const char** av)
29 {
30  TCLAP::CmdLine cmd("INRIA / IRISA - VisAGeS/Empenn Team", ' ',ANIMA_VERSION);
31 
32  TCLAP::ValueArg<std::string> inArg("i","input","input tracks file",true,"","input tracks",cmd);
33  TCLAP::ValueArg<std::string> outArg("o","output","output tracks name",true,"","output tracks",cmd);
34  TCLAP::ValueArg<std::string> trArg("t","trsf","Transformations XML list",true,"","transformations list",cmd);
35 
36  TCLAP::ValueArg<unsigned int> expOrderArg("e","exp-order","Order of field exponentiation approximation (in between 0 and 1, default: 0)",false,0,"exponentiation order",cmd);
37  TCLAP::SwitchArg invertArg("I","invert","Invert the transformation series",cmd,false);
38  TCLAP::ValueArg<unsigned int> nbpArg("p","numberofthreads","Number of threads to run on (default: all cores)",false,itk::MultiThreaderBase::GetGlobalDefaultNumberOfThreads(),"number of threads",cmd);
39 
40  try
41  {
42  cmd.parse(ac,av);
43  }
44  catch (TCLAP::ArgException& e)
45  {
46  std::cerr << "Error: " << e.error() << "for argument " << e.argId() << std::endl;
47  return EXIT_FAILURE;
48  }
49 
50  typedef anima::TransformSeriesReader <double, 3> TransformSeriesReaderType;
51  TransformSeriesReaderType trsfReader;
52  trsfReader.SetInput(trArg.getValue());
53  trsfReader.SetInvertTransform(!invertArg.isSet());
54  trsfReader.SetExponentiationOrder(expOrderArg.getValue());
55  trsfReader.SetNumberOfWorkUnits(nbpArg.getValue());
56  trsfReader.Update();
57 
58  anima::ShapesReader trackReader;
59  trackReader.SetFileName(inArg.getValue());
60  trackReader.Update();
61 
62  vtkSmartPointer <vtkPolyData> tracks = trackReader.GetOutput();
63  ApplyTransformToTracks(tracks->GetPoints(), trsfReader.GetOutputTransform());
64 
65  anima::ShapesWriter writer;
66  writer.SetInputData(tracks);
67  writer.SetFileName(outArg.getValue());
68  std::cout << "Writing tracks: " << outArg.getValue() << std::endl;
69  writer.Update();
70 
71  return EXIT_SUCCESS;
72 }
void SetInputData(vtkPolyData *data)
itk::CompositeTransform< TScalarType, NDimensions > OutputTransformType
void SetFileName(std::string &name)
void ApplyTransformToTracks(vtkPoints *dataPoints, anima::TransformSeriesReader< double, 3 >::OutputTransformType *transform)
void SetInput(std::string const &trListName)
vtkPolyData * GetOutput()
int main(int ac, const char **av)