22 #include "itkPSMProcrustesRegistration.h" 24 #include "itkImageFileReader.h" 25 #include "itkPSMEntropyModelFilter.h" 26 #include "itkPSMProjectReader.h" 27 #include "itkPSMParticleSystem.h" 28 #include "itkPSMRegionDomain.h" 29 #include "vnl/vnl_matrix_fixed.h" 30 #include "vnl/vnl_vector_fixed.h" 31 #include "itkCommand.h" 40 typedef Command Superclass;
41 typedef SmartPointer< Self > Pointer;
42 typedef SmartPointer< const Self > ConstPointer;
44 typedef Image<float, 2> ImageType;
54 virtual void Execute(Object *caller,
const EventObject &)
59 if (this->procrustesRegistration->GetProcrustesInterval() != 0)
61 this->m_ProcrustesCounter++;
63 if (this->m_ProcrustesCounter >= (
int)this->procrustesRegistration->GetProcrustesInterval())
66 this->m_ProcrustesCounter = 0;
68 std::cout <<
"Run Procrustes Registration" << std::endl;
76 std::cout <<
" Eigenmode variances: ";
81 std::cout << std::endl;
82 std::cout <<
" Regularization = " << o->GetRegularizationConstant() << std::endl;
84 virtual void Execute(
const Object *,
const EventObject &)
86 std::cout <<
"SHOULDN'T BE HERE" << std::endl;
89 { procrustesRegistration = p; }
94 m_ProcrustesCounter = 0;
98 int m_ProcrustesCounter;
100 void operator=(
const Self &);
116 typedef T ObjectType;
123 std::vector<ObjectType> &GetOutput()
128 void SetFileName(
const char *fn)
130 void SetFileName(
const std::string &fn)
132 const std::string& GetFileName()
const 133 {
return m_FileName; }
142 std::ifstream in( m_FileName.c_str(), std::ios::binary );
146 std::cerr <<
"Could not open filename " << m_FileName << std::endl;
151 in.read(reinterpret_cast<char *>(&N),
sizeof(
int));
153 int sz =
sizeof(ObjectType);
155 for (
unsigned int i = 0; i < (
unsigned int)N; i++)
158 in.read(reinterpret_cast<char *>(&q), sz);
159 m_Output.push_back(q);
169 void operator=(
const Self&);
171 std::vector<ObjectType> m_Output;
172 std::string m_FileName;
177 int itkPSMProcrustesRegistration2DTest(
int argc,
char* argv[] )
180 std::string errstring =
"";
181 std::string output_path =
"";
182 std::string input_path_prefix =
"";
187 std::cout <<
"Wrong number of arguments. \nUse: " 188 <<
"itkPSMProcrustesRegistration2DTest parameter_file transforms_file [output_path] [input_path]\n" 189 <<
"See itk::PSMParameterFileReader for documentation on the parameter file format.\n" 190 <<
" Note that input_path will be prefixed to any file names and paths in the xml parameter file.\n" 197 output_path = std::string(argv[3]);
202 input_path_prefix = std::string(argv[4]);
205 typedef itk::Image<float, 2> ImageType;
210 itk::PSMProjectReader::Pointer xmlreader =
211 itk::PSMProjectReader::New();
212 xmlreader->SetFileName(argv[1]);
215 itk::PSMProject::Pointer project = xmlreader->GetOutput();
218 itk::PSMEntropyModelFilter<ImageType>::Pointer P
223 itk::MyPSMProcrustesIterationCommand::Pointer mycommand
224 = itk::MyPSMProcrustesIterationCommand::New();
225 P->AddObserver(itk::IterationEvent(), mycommand);
228 itk::PSMProcrustesRegistration<2>::Pointer procrustesRegistration
231 mycommand->SetPSMProcrustesRegistration( procrustesRegistration );
233 const std::vector<std::string> &dt_files = project->GetDistanceTransforms();
234 itk::ImageFileReader<ImageType>::Pointer reader =
235 itk::ImageFileReader<ImageType>::New();
237 std::cout <<
"Reading distance transforms ..." << std::endl;
238 for (
unsigned int i = 0; i < dt_files.size(); i++)
240 reader->SetFileName(input_path_prefix + dt_files[i]);
243 std::cout <<
" " << dt_files[i] << std::endl;
246 for(
unsigned int i = 0; i < 103; i++)
250 std::cout <<
"Done!" << std::endl;
251 std::cout <<
"Number of inputs: " << P->GetNumberOfInputs() << std::endl;
254 const std::vector<std::string> &pt_files = project->GetModel(std::string(
"initialization"));
255 std::vector<itk::PSMEntropyModelFilter<ImageType>::PointType> c;
256 std::cout <<
"Reading the initial model correspondences ..." << std::endl;
257 unsigned int numOfPoints;
258 for (
unsigned int i = 0; i < pt_files.size(); i++)
263 std::ifstream in( (input_path_prefix + pt_files[0]).c_str() );
266 errstring +=
"Could not open point file for input.";
276 for (
unsigned int d = 0; d < 2; d++)
285 std::cout <<
"Read " << numOfPoints-1 <<
" points. " << std::endl;
289 for(
unsigned int i = 0; i < 100; i++)
294 std::cout <<
"Done!" << std::endl;
298 transform_reader.SetFileName(argv[2]);
299 transform_reader.Update();
301 std::cout <<
"Reading transforms." << std::endl;
303 for (
unsigned int i = 0; i < P->GetParticleSystem()->GetNumberOfDomains(); i++)
305 for(
unsigned int j = 0; j < numOfPoints; j++)
309 itk::PSMParticleSystem<2>::Pointer PS = P->GetParticleSystem();
316 PS->SetPosition( trPoint, j, i);
321 double regularization_initial = 100.0f;
322 double regularization_final = 5.0f;
323 double regularization_decayspan = 2000.0f;
324 double tolerance = 1.0e-8;
325 unsigned int maximum_iterations = 200000;
326 unsigned int procrustes_interval = 1;
327 if ( project->HasOptimizationAttribute(
"regularization_initial") )
329 regularization_initial = project->GetOptimizationAttribute(
"regularization_initial");
331 if ( project->HasOptimizationAttribute(
"regularization_final") )
333 regularization_final = project->GetOptimizationAttribute(
"regularization_final");
335 if ( project->HasOptimizationAttribute(
"regularization_decayspan") )
337 regularization_decayspan = project->GetOptimizationAttribute(
"regularization_decayspan");
339 if ( project->HasOptimizationAttribute(
"tolerance") )
341 tolerance = project->GetOptimizationAttribute(
"tolerance");
343 if ( project->HasOptimizationAttribute(
"maximum_iterations") )
346 =
static_cast<unsigned int>(project->GetOptimizationAttribute(
"maximum_iterations"));
348 if ( project->HasOptimizationAttribute(
"procrustes_interval") )
351 =
static_cast<unsigned int>(project->GetOptimizationAttribute(
"procrustes_interval"));
358 std::cout <<
"Optimization parameters: " << std::endl;
359 std::cout <<
" regularization_initial = " << regularization_initial << std::endl;
360 std::cout <<
" regularization_final = " << regularization_final << std::endl;
361 std::cout <<
" regularization_decayspan = " << regularization_decayspan << std::endl;
362 std::cout <<
" tolerance = " << tolerance << std::endl;
363 std::cout <<
" maximum_iterations = " << maximum_iterations << std::endl;
364 std::cout <<
" procrustes_interval = " << procrustes_interval << std::endl;
369 P->SetRegularizationFinal(regularization_final);
370 P->SetRegularizationDecaySpan(regularization_decayspan);
376 errstring +=
"Optimization did not converge based on tolerance criteria.\n";
381 std::string output_transform_file =
"output_transforms_PSMProcrustesRegistration2DTest.txt";
382 std::string out_file = output_path + output_transform_file;
383 std::ofstream out(out_file.c_str());
384 for (
unsigned int d = 0; d < P->GetParticleSystem()->GetNumberOfDomains(); d++)
388 errstring +=
"Could not open file for output: ";
392 out << P->GetParticleSystem()->GetTransform(d);
399 const std::vector<std::string> &out_files = project->GetModel(std::string(
"optimized"));
401 for (
unsigned int d = 0; d < P->GetParticleSystem()->GetNumberOfDomains(); d++)
404 std::ostringstream ss;
406 std::string fname = output_path + out_files[0] +
"_" + ss.str() +
".lpts";
407 std::ofstream out( fname.c_str() );
410 errstring +=
"Could not open point file for output: ";
414 for (
unsigned int j = 0; j < P->GetParticleSystem()->GetNumberOfParticles(d); j++)
419 out << P->GetParticleSystem()->GetPosition(j,d)[0] <<
" " << P->GetParticleSystem()->GetPosition(j,d)[1] <<
" " << 0.0;
426 catch(itk::ExceptionObject &e)
428 errstring =
"ITK exception with description: " + std::string(e.GetDescription())
429 + std::string(
"\n at location:") + std::string(e.GetLocation())
430 + std::string(
"\n in file:") + std::string(e.GetFile());
435 errstring =
"Unknown exception thrown";
441 std::cout <<
"All tests passed" << std::endl;
446 std::cout <<
"Test failed with the following error:" << std::endl;
447 std::cout << errstring << std::endl;
void SetMaximumNumberOfIterations(const std::vector< unsigned int > &n)
unsigned int GetNumberOfElapsedIterations() const
void SetProcrustesInterval(int i)
PointType & GetPosition(unsigned long int k, unsigned int d=0)
void SetInputCorrespondencePoints(unsigned int index, const std::vector< PointType > &corr)
ImageType::PointType PointType
const std::vector< double > & GetShapePCAVariances() const
void SetRegularizationInitial(const std::vector< double > &v)
const std::vector< ObjectType > & GetOutput() const
void SetInput(const std::string &s, itk::DataObject *o)
itkTypeMacro(MyPSMProcrustesIterationCommand, Command)
PointType TransformPoint(const PointType &, const TransformType &) const
void SetPSMParticleSystem(PSMParticleSystemType *p)
virtual void Execute(Object *caller, const EventObject &)
void SetTolerance(const std::vector< double > &v)
void RunRegistration(int i)
MyPSMProcrustesIterationCommand Self
vnl_matrix_fixed< double, VDimension+1, VDimension+1 > TransformType