cmtkImagePairNonrigidRegistrationCommandLine.cxx

Go to the documentation of this file.
00001 /*
00002 //
00003 //  Copyright 1997-2009 Torsten Rohlfing
00004 //
00005 //  Copyright 2004-2010 SRI International
00006 //
00007 //  This file is part of the Computational Morphometry Toolkit.
00008 //
00009 //  http://www.nitrc.org/projects/cmtk/
00010 //
00011 //  The Computational Morphometry Toolkit is free software: you can
00012 //  redistribute it and/or modify it under the terms of the GNU General Public
00013 //  License as published by the Free Software Foundation, either version 3 of
00014 //  the License, or (at your option) any later version.
00015 //
00016 //  The Computational Morphometry Toolkit is distributed in the hope that it
00017 //  will be useful, but WITHOUT ANY WARRANTY; without even the implied
00018 //  warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00019 //  GNU General Public License for more details.
00020 //
00021 //  You should have received a copy of the GNU General Public License along
00022 //  with the Computational Morphometry Toolkit.  If not, see
00023 //  <http://www.gnu.org/licenses/>.
00024 //
00025 //  $Revision: 2492 $
00026 //
00027 //  $LastChangedDate: 2010-10-21 15:59:13 -0700 (Thu, 21 Oct 2010) $
00028 //
00029 //  $LastChangedBy: torstenrohlfing $
00030 //
00031 */
00032 
00033 #include "cmtkImagePairNonrigidRegistrationCommandLine.h"
00034 
00035 #include <System/cmtkCommandLine.h>
00036 #include <System/cmtkConsole.h>
00037 #include <System/cmtkTimers.h>
00038 #include <System/cmtkThreads.h>
00039 #include <System/cmtkCompressedStream.h>
00040 #include <System/cmtkMountPoints.h>
00041 
00042 #include <IO/cmtkXformIO.h>
00043 #include <IO/cmtkClassStream.h>
00044 #include <IO/cmtkClassStreamAffineXform.h>
00045 #include <IO/cmtkVolumeIO.h>
00046 #include <IO/cmtkSplineWarpXformITKIO.h>
00047 
00048 #include <Base/cmtkAnatomicalOrientation.h>
00049 #include <Base/cmtkTransformChangeFromSpaceAffine.h>
00050 
00051 #include <Registration/cmtkVoxelMatchingElasticFunctional.h>
00052 #include <Registration/cmtkSymmetricElasticFunctional.h>
00053 
00054 #ifdef HAVE_SYS_TYPES_H
00055 #  include <sys/types.h>
00056 #endif
00057 
00058 #ifdef HAVE_SYS_STAT_H
00059 #  include <sys/stat.h>
00060 #endif
00061 
00062 #ifdef HAVE_SYS_UTSNAME_H
00063 #  include <sys/utsname.h>
00064 #endif
00065 
00066 #ifdef _MSC_VER
00067 #  include <direct.h>
00068 #endif // #ifdef _MSC_VER
00069 
00070 #include <stdlib.h>
00071 #include <stdio.h>
00072 #include <string.h>
00073 #include <signal.h>
00074 #include <iostream>
00075 
00076 namespace
00077 cmtk
00078 {
00079 
00082 
00083 ImagePairNonrigidRegistrationCommandLine* ImagePairNonrigidRegistrationCommandLine::StaticThis = NULL;
00084 
00085 ImagePairNonrigidRegistrationCommandLine
00086 ::ImagePairNonrigidRegistrationCommandLine
00087 ( const int argc, const char *argv[] ) :
00088   m_OutputPathITK( NULL ),
00089   m_ReformattedImagePath( NULL )
00090 {
00091   Studylist = Time = NULL;
00092 
00093   this->m_OutputIntermediate = 0;
00094   Verbose = 0;
00095 
00096   InputStudylist = NULL;
00097   const char *initialTransformationFile = NULL;
00098   bool initialTransformationInverse = false;
00099   IntermediateResultIndex = 0;
00100 
00101   bool forceOutsideFlag = false;
00102   Types::DataItem forceOutsideValue = 0;
00103 
00104   const char* clArg1 = NULL; // input studylist or reference image
00105   const char* clArg2 = NULL; // empty or floating image
00106   const char* clArg3 = NULL; // empty or initial transformation
00107 
00108   try
00109     {
00110     CommandLine cl( CommandLine::PROPS_XML );
00111     cl.SetProgramInfo( CommandLine::PRG_TITLE, "B-spline nonrigid registration" );
00112     cl.SetProgramInfo( CommandLine::PRG_DESCR, "This program performs nonrigid image registration using multi-resolution optimization of voxel-based image similarity measures "
00113                        "and a multi-resolution B-spline transformation model." );
00114     cl.SetProgramInfo( CommandLine::PRG_CATEG, "CMTK.Registration" );
00115 
00116     typedef CommandLine::Key Key;
00117     cl.BeginGroup( "Console", "Console output control" )->SetProperties( CommandLine::PROPS_NOXML );
00118     cl.AddSwitch( Key( 'v', "verbose" ), &this->Verbose, true, "Verbose peration" );
00119     cl.AddSwitch( Key( 'q', "quiet" ), &Verbose, false, "Quiet mode" );
00120     cl.EndGroup();
00121 
00122     cl.BeginGroup( "TransformationIO", "Transformation import/export" );
00123     cl.AddOption( Key( "initial" ), &initialTransformationFile, "Initialize transformation from given path" )->SetProperties( CommandLine::PROPS_XFORM );
00124     cl.AddSwitch( Key( "invert-initial" ), &initialTransformationInverse, true, "Invert given (affine) initial transformation." );
00125     cl.AddOption( Key( "write-itk-xform" ), &this->m_OutputPathITK, "Output path for final transformation in ITK format" )
00126       ->SetProperties( CommandLine::PROPS_XFORM | CommandLine::PROPS_OUTPUT )
00127       ->SetAttribute( "type", "bspline" )->SetAttribute( "reference", "FloatingImage" );
00128     cl.AddOption( Key( "write-reformatted" ), &this->m_ReformattedImagePath, "Write reformatted floating image." )
00129       ->SetProperties( CommandLine::PROPS_IMAGE | CommandLine::PROPS_OUTPUT );
00130     cl.EndGroup();
00131 
00132     cl.BeginGroup( "Transformation", "Transformation parameters" );
00133     cl.AddOption( Key( "grid-spacing" ), &this->m_GridSpacing, "Control point grid spacing" );
00134     cl.AddOption( Key( "grid-refine" ), &this->m_RefineGrid, "Number of refinements (control point grid resolution levels)" );
00135     cl.AddSwitch( Key( "delay-refine" ), &this->m_DelayRefineGrid, true, "Delay control point grid refinement; first switch to next higher image resolution" );
00136     cl.AddSwitch( Key( "exact-spacing" ), &this->m_ExactGridSpacing, true, "Use exact control point spacing; do not modify spacing to fit reference image bounding box" );
00137 
00138     cl.AddOption( Key( "ignore-edge" ), &this->IgnoreEdge, "Ignore n control point layers along each image face" );
00139     cl.AddOption( Key( "restrict" ), &this->RestrictToAxes, "Restrict deformation to coordinate dimension(s) [one or more of 'x','y','z']" );
00140 
00141     cl.AddSwitch( Key( "no-adaptive-fix" ), &this->m_AdaptiveFixParameters, false, "Disable adaptive fixing of control points; optimize all deformation parameters" );
00142     cl.AddOption( Key( "adaptive-fix-thresh" ), &this->m_AdaptiveFixThreshFactor, "Threshold factor for entropy criterion to fix local control points" );
00143     cl.AddSwitch( Key( "accurate" ), &this->m_FastMode, false, "Accurate computation mode: may give slightly better results after substantially longer computation" );
00144     cl.AddSwitch( Key( "fast" ), &this->m_FastMode, true, "Fast computation mode: may give slightly worse results than accurate mode, but saves substantial CPU time" );
00145     cl.EndGroup();
00146 
00147     cl.BeginGroup( "Regularization", "Regularization parameters" );
00148     cl.AddOption( Key( "jacobian-constraint-weight" ), &this->m_JacobianConstraintWeight, "Weight for Jacobian-based local volume preservation constraint" );
00149     cl.AddOption( Key( "smoothness-constraint-weight" ), &this->m_GridEnergyWeight, "Weight for smoothness constraint based on second-order grid bending energy." );
00150     cl.AddOption( Key( "landmark-constraint-weight" ), &this->m_LandmarkErrorWeight, "Weight for landmark misregistration registration" );
00151     cl.AddOption( Key( "inverse-consistency-weight" ), &this->m_InverseConsistencyWeight, "Weight for inverse consistency constraint" );
00152     cl.AddOption( Key( "constraint-relaxation-factor" ), &this->m_RelaxWeight, "Weight relaxation factor for alternating under-constrained iterations" );
00153     cl.AddSwitch( Key( "relax-to-unfold" ), &this->m_RelaxToUnfold, true, "Before each resolution level, regularize negative-Jacobian areas of the deformation to unfold them." );
00154     cl.EndGroup();
00155 
00156     cl.BeginGroup( "Optimization", "Optimization parameters" );
00157     cl.AddOption( Key( "max-stepsize" ), &this->m_MaxStepSize, "Maximum optimizer step size, which determines search space exploration." );
00158     cl.AddOption( Key( "min-stepsize" ), &this->m_MinStepSize, "Minimum optimizer step size, which determines precision." );
00159     cl.AddOption( Key( "stepfactor" ), &this->m_OptimizerStepFactor, "Factor for search step size reduction. Must be > 0.0 and < 1.0 [default: 0.5]" );
00160     cl.AddOption( Key( "delta-f-threshold" ), &this->m_DeltaFThreshold, "Optional threshold to terminate optimization (level) if relative change of target function drops below this value." );
00161     cl.AddSwitch( Key( "no-maxnorm" ), &this->m_UseMaxNorm, false, "Use Euclid norm for gradient normalication in optimization, rather than maximum norm" );
00162     cl.EndGroup();
00163 
00164     cl.BeginGroup( "Resolution", "Image resolution parameters" );
00165     cl.AddOption( Key( 's', "sampling" ), &this->m_Sampling, "Image sampling (finest resampled image resolution)" );
00166     cl.AddOption( Key( "coarsest" ), &this->m_CoarsestResolution, "Upper limit for image sampling in multiresolution hierarchy" );
00167 
00168     cl.AddSwitch( Key( "omit-original-data" ), &this->m_UseOriginalData, false, "Do not use original data in full resolution for final registration stage." );
00169     cl.EndGroup();
00170 
00171     cl.BeginGroup( "Images", "Image data" );
00172     CommandLine::EnumGroup<int>::SmartPtr
00173       metricGroup = cl.AddEnum( "registration-metric", &this->m_Metric, "Registration metric for motion estimation by image-to-image registration." );
00174     metricGroup->AddSwitch( Key( "nmi" ), 0, "Normalized Mutual Information metric" );
00175     metricGroup->AddSwitch( Key( "mi" ), 1, "Standard Mutual Information metric" );
00176     metricGroup->AddSwitch( Key( "cr" ), 2, "Correlation Ratio metric" );
00177     metricGroup->AddSwitch( Key( "rms" ), 3, "Root of Mean Squaresa metric (this is the square root of MSD)" );
00178     metricGroup->AddSwitch( Key( "msd" ), 4, "Mean Squared Difference metric" );
00179     metricGroup->AddSwitch( Key( "ncc" ), 5, "Normalized Cross Correlation metric" );
00180 
00181     cl.BeginGroup( "Interpolation", "Floating Image Interpolation Options" );
00182     cmtk::CommandLine::EnumGroup<Interpolators::InterpolationEnum>::SmartPtr kernelGroup = 
00183       cl.AddEnum( "interpolation", &this->m_FloatingImageInterpolation, "Interpolation method for floating image sampling:" );
00184     kernelGroup->AddSwitch( Key( "nearest-neighbor" ), Interpolators::NEAREST_NEIGHBOR, "Nearest neighbor interpolation (for intensity and label data)" );
00185     kernelGroup->AddSwitch( Key( "linear" ), Interpolators::LINEAR, "Trilinear interpolation" );
00186     kernelGroup->AddSwitch( Key( "cubic" ), Interpolators::CUBIC, "Tricubic interpolation" );
00187     kernelGroup->AddSwitch( Key( "cosine-sinc" ), Interpolators::COSINE_SINC, "Cosine-windowed sinc interpolation (most accurate but slowest)" );
00188     kernelGroup->AddSwitch( Key( "partial-volume" ), Interpolators::PARTIALVOLUME, "Partial volume interpolation (for label data)" );
00189     kernelGroup->AddSwitch( Key( "automatic" ), Interpolators::DEFAULT, "Select interpolation automatically based on data type: linear for grey-level data, nearest neighbor for label data." );
00190 
00191     cl.AddSwitch( Key( "match-histograms" ), &this->m_MatchFltToRefHistogram, true, "Match floating image histogram to reference image histogram." );
00192     cl.AddSwitch( Key( "repeat-match-histograms" ), &this->m_RepeatMatchFltToRefHistogram, true, "Repeat histogram matching after every level of the registration to account for volume changes. When registering masked data, it is advisable to also use the --force-outside-value option to prevent poorly matched histograms." );
00193     cl.AddOption( Key( "force-outside-value" ), &forceOutsideValue, "Force values outside field of view to this value rather than drop incomplete pixel pairs", &forceOutsideFlag );
00194 
00195     this->m_PreprocessorRef.AttachToCommandLine( cl );
00196     this->m_PreprocessorFlt.AttachToCommandLine( cl );
00197 
00198     cl.BeginGroup( "Output", "Output parameters" )->SetProperties( CommandLine::PROPS_NOXML );
00199     cl.AddOption( Key( 'o', "outlist" ), &this->Studylist, "Output path for final transformation" );
00200     cl.AddOption( Key( 't', "time" ), &this->Time, "Computation time statistics output file name" );
00201     cl.AddSwitch( Key( "output-intermediate" ), &this->m_OutputIntermediate, true, "Write transformation for each level [default: only write final transformation]" );
00202     cl.EndGroup();
00203 
00204     cl.AddParameter( &clArg1, "ReferenceImage", "Reference (fixed) image path" )
00205       ->SetProperties( CommandLine::PROPS_IMAGE );
00206     cl.AddParameter( &clArg2, "FloatingImage", "Floating (moving) image path" )
00207       ->SetProperties( CommandLine::PROPS_IMAGE | CommandLine::PROPS_OPTIONAL);
00208     cl.AddParameter( &clArg3, "InitialXform", "Initial affine transformation from reference to floating image" )
00209       ->SetProperties( CommandLine::PROPS_NOXML | CommandLine::PROPS_XFORM | CommandLine::PROPS_OPTIONAL );
00210 
00211     cl.Parse( argc, argv );
00212     }
00213   catch ( const CommandLine::Exception& ex )
00214     {
00215     StdErr << ex << "\n";
00216     exit( 1 );
00217     }
00218 
00219   if ( (this->m_OptimizerStepFactor <= 0) || (this->m_OptimizerStepFactor >= 1) ) 
00220     {
00221     StdErr << "ERROR: step factor value " << this->m_OptimizerStepFactor << " is invalid. Must be in range (0..1)\n";
00222     exit( 1 );
00223     }
00224 
00225   if ( clArg2 ) 
00226     {
00227     this->SetInitialTransformation( AffineXform::SmartPtr( new AffineXform() ) );
00228     
00229     Study1 = const_cast<char*>( clArg1 );
00230     Study2 = const_cast<char*>( clArg2 );
00231 
00232     if ( clArg3 )
00233       {
00234       initialTransformationFile = clArg3;
00235       }
00236     }
00237   else 
00238     {
00239     InputStudylist = clArg1;
00240     
00241     if ( Verbose )
00242       fprintf( stderr, "Reading input studylist %s.\n", InputStudylist );
00243     
00244     ClassStream classStream( MountPoints::Translate(InputStudylist),"registration", ClassStream::READ );
00245     if ( ! classStream.IsValid() ) 
00246       {
00247       std::cerr << "Could not open studylist archive " << InputStudylist << ".\n";
00248       exit( 1 );
00249       }
00250     
00251     classStream.Seek ( "registration" );
00252     Study1 = classStream.ReadString( "reference_study" );
00253     Study2 = classStream.ReadString( "floating_study" );
00254     if ( Study2 )
00255       {
00256       AffineXform::SmartPtr affineXform;
00257       classStream >> affineXform;
00258       this->SetInitialTransformation( affineXform );
00259       }
00260     else
00261       {
00262       // legacy studylists have inverse transformation stored in them
00263       Study2 = classStream.ReadString( "model_study" );
00264       AffineXform::SmartPtr affineXform;
00265       classStream >> affineXform;
00266       this->SetInitialTransformation( affineXform->GetInverse() );
00267       }
00268     
00269     classStream.Close();
00270     }
00271 
00273   if ( initialTransformationFile ) 
00274     {
00275     Xform::SmartPtr initialXform( XformIO::Read( initialTransformationFile ) );
00276     AffineXform::SmartPtr affineXform = AffineXform::SmartPtr::DynamicCastFrom( initialXform );
00277     if ( affineXform )
00278       {
00279       if ( initialTransformationInverse )
00280         this->SetInitialTransformation( affineXform->GetInverse() );
00281       else
00282         this->SetInitialTransformation( affineXform );
00283       }
00284     else
00285       {
00286       InitialWarpXform = SplineWarpXform::SmartPtr::DynamicCastFrom( initialXform );
00287       }
00288     }
00289 
00290   UniformVolume::SmartPtr volume( VolumeIO::ReadOriented( Study1, Verbose ) );
00291   if ( !volume ) throw ConstructorFailed();
00292   this->SetVolume_1( UniformVolume::SmartPtr( this->m_PreprocessorRef.GetProcessedImage( volume ) ) );
00293 
00294   volume = UniformVolume::SmartPtr( VolumeIO::ReadOriented( Study2, Verbose ) );
00295   if ( !volume ) throw ConstructorFailed();
00296   this->SetVolume_2( UniformVolume::SmartPtr( this->m_PreprocessorFlt.GetProcessedImage( volume ) ) );
00297 
00298   AffineXform::SmartPtr affineXform( AffineXform::SmartPtr::DynamicCastFrom( this->m_InitialTransformation ) );
00299   if ( affineXform )
00300     {
00301     if ( affineXform->MetaKeyExists( META_SPACE ) && (affineXform->m_MetaInformation[META_SPACE] != AnatomicalOrientation::ORIENTATION_STANDARD ) )
00302       {
00303       TransformChangeFromSpaceAffine toStandardSpace( *affineXform, *(this->m_Volume_1), *(this->m_Volume_2), affineXform->m_MetaInformation[META_SPACE].c_str() );
00304       *affineXform = toStandardSpace.GetTransformation();
00305       affineXform->m_MetaInformation[META_SPACE] = AnatomicalOrientation::ORIENTATION_STANDARD;
00306       this->SetInitialTransformation( affineXform );
00307       }
00308     }
00309 
00310   if ( forceOutsideFlag )
00311     {
00312     this->SetForceOutside( true, forceOutsideValue );
00313     }
00314 }
00315 
00316 ImagePairNonrigidRegistrationCommandLine::~ImagePairNonrigidRegistrationCommandLine ()
00317 {
00318 #ifdef HAVE_SIGRELSE
00319   // release signal handler
00320   sigrelse( SIGUSR1 );
00321 #endif
00322 }
00323 
00324 CallbackResult
00325 ImagePairNonrigidRegistrationCommandLine
00326 ::InitRegistration ()
00327 {
00328   CallbackResult result = this->Superclass::InitRegistration();
00329   if ( result != CALLBACK_OK )
00330     return result;
00331 
00332   if ( this->m_OutputIntermediate )
00333     this->OutputIntermediate();
00334   
00335   // register signal handler for intermediate result output.
00336   Self::StaticThis = this;
00337 #ifndef _MSC_VER
00338   signal( SIGUSR1, cmtkImagePairNonrigidRegistrationCommandLineDispatchSIGUSR1 );
00339 #endif
00340 
00341   return CALLBACK_OK;
00342 }
00343         
00344 void
00345 ImagePairNonrigidRegistrationCommandLine
00346 ::OutputResult
00347 ( const CoordinateVector* )
00348 {
00349   if ( Studylist ) 
00350     this->OutputWarp( Studylist );
00351 
00352   if ( this->m_OutputPathITK ) 
00353     {
00354     SplineWarpXformITKIO::Write( this->m_OutputPathITK, *(this->GetTransformation()), *(this->m_ReferenceVolume), *(this->m_FloatingVolume) );
00355     }
00356 
00357   if ( this->m_ReformattedImagePath )
00358     {
00359     VolumeIO::Write( *(this->GetReformattedFloatingImage() ), this->m_ReformattedImagePath, this->Verbose );
00360     }
00361 }
00362 
00363 void
00364 ImagePairNonrigidRegistrationCommandLine
00365 ::EnterResolution
00366 ( CoordinateVector::SmartPtr& v, Functional::SmartPtr& f, 
00367   const int index, const int total )
00368 {
00369   if ( Verbose )
00370     fprintf( stderr, "\rEntering resolution level %d out of %d...\n", index, total );
00371   
00372   this->Superclass::EnterResolution( v, f, index, total );
00373 }
00374 
00375 int
00376 ImagePairNonrigidRegistrationCommandLine::DoneResolution
00377 ( CoordinateVector::SmartPtr& v, Functional::SmartPtr& f, const int index, const int total )
00378 {
00379   if ( this->m_OutputIntermediate )
00380     this->OutputIntermediate();
00381   return this->Superclass::DoneResolution( v, f, index, total );
00382 }
00383 
00384 void
00385 ImagePairNonrigidRegistrationCommandLine::OutputWarp ( const char* path ) const
00386 {
00387   ClassStream classStream( path, "studylist", ClassStream::WRITE );
00388   if ( ! classStream.IsValid() ) return;
00389 
00390   classStream.Begin( "studylist" );
00391   classStream.WriteInt( "num_sources", 2 );
00392   classStream.End();
00393 
00394   classStream.Begin( "source" );
00395   classStream.WriteString( "studyname", CompressedStream::GetBaseName( Study1 ) );
00396   classStream.End();
00397 
00398   classStream.Begin( "source" );
00399   classStream.WriteString( "studyname", CompressedStream::GetBaseName( Study2 ) );
00400   classStream.End();
00401 
00402   classStream.Close();
00403 
00404   classStream.Open( path, "settings", ClassStream::WRITE );
00405   classStream.WriteInt( "algorithm", this->m_Algorithm );
00406   classStream.WriteBool( "use_maxnorm", this->m_UseMaxNorm );
00407   classStream.WriteDouble( "exploration", this->m_MaxStepSize );
00408   classStream.WriteDouble( "accuracy", this->m_MinStepSize );
00409   classStream.WriteDouble( "min_sampling", this->m_Sampling );
00410   classStream.WriteDouble( "coarsest_resolution", this->m_CoarsestResolution );
00411   classStream.WriteBool( "use_original_data", this->m_UseOriginalData );
00412   classStream.WriteBool( "fast_mode", this->m_FastMode );
00413   classStream.WriteInt( "metric", this->m_Metric );
00414   classStream.WriteDouble( "optimizer_step_factor", this->m_OptimizerStepFactor );
00415   classStream.WriteDouble( "grid_spacing", this->m_GridSpacing );
00416   classStream.WriteInt( "ignore_edge", IgnoreEdge );
00417   classStream.WriteDouble( "jacobian_constraint_weight", this->m_JacobianConstraintWeight );
00418   classStream.WriteDouble( "energy_constraint_weight", this->m_GridEnergyWeight );
00419   classStream.WriteDouble( "inverse_consistency_weight", this->m_InverseConsistencyWeight );
00420   classStream.WriteDouble( "weight_relaxation", this->m_RelaxWeight );
00421   classStream.WriteDouble( "landmark_error_weight", this->m_LandmarkErrorWeight );
00422   classStream.WriteInt( "refine_grid", this->m_RefineGrid );
00423   classStream.WriteBool( "delay_refine_grid", this->m_DelayRefineGrid );
00424   classStream.WriteBool( "adaptive_fix_parameters", this->m_AdaptiveFixParameters );
00425   classStream.WriteDouble( "adaptive_fix_parameters_thresh", this->m_AdaptiveFixThreshFactor );
00426 
00427   this->m_PreprocessorRef.WriteSettings( classStream );  
00428   this->m_PreprocessorFlt.WriteSettings( classStream );  
00429 
00430   classStream.Close();
00431       
00432   classStream.Open( path, "statistics", ClassStream::WRITE );
00433   classStream.WriteDouble( "time_level", this->GetLevelElapsedTime() );
00434   classStream.WriteDouble( "time_total", this->GetTotalElapsedTime() );
00435   classStream.WriteDouble( "walltime_level", this->GetLevelElapsedWalltime() );
00436   classStream.WriteDouble( "walltime_total", this->GetTotalElapsedWalltime() );
00437   classStream.WriteDouble( "thread_time_level", this->GetThreadLevelElapsedTime() );
00438   classStream.WriteDouble( "thread_time_total", this->GetThreadTotalElapsedTime() );
00439   classStream.WriteInt( "number_of_threads", Threads::NumberOfThreads );
00440   classStream.WriteInt( "number_of_cpus", Threads::GetNumberOfProcessors() );
00441 
00442 #ifndef _MSC_VER
00443   struct utsname name;
00444   if ( uname( &name ) >= 0 ) 
00445     {
00446     classStream.WriteString( "host", name.nodename );
00447     classStream.WriteString( "system", name.sysname );
00448     }
00449 #endif
00450   classStream.Close();
00451 
00452   const WarpXform::SmartPtr warp = WarpXform::SmartPtr::DynamicCastFrom( this->m_Xform );
00453   if ( warp ) 
00454     {
00455     classStream.Open( path, "registration", ClassStream::WRITE );
00456     if ( classStream.IsValid() ) 
00457       {
00458       classStream.Begin( "registration" );
00459       classStream.WriteString( "reference_study", CompressedStream::GetBaseName( this->Study1 ) );
00460       classStream.WriteString( "floating_study", CompressedStream::GetBaseName( this->Study2 ) );
00461     
00462       if ( warp->GetInitialAffineXform() ) 
00463         {
00464         classStream << (*warp->GetInitialAffineXform());
00465         } 
00466       else 
00467         {
00468         classStream << *this->m_InitialTransformation;
00469         }
00470       classStream << warp;
00471       classStream.End();
00472       }
00473     classStream.Close();
00474     }
00475 }
00476 
00477 void
00478 ImagePairNonrigidRegistrationCommandLine::OutputIntermediate( const bool incrementCount )
00479 {
00480   char path[PATH_MAX];
00481   if ( Studylist ) 
00482     {
00483     snprintf( path, sizeof( path ), "%s/level-%02d.list", Studylist, IntermediateResultIndex );
00484     } 
00485   else
00486     {
00487     snprintf( path, sizeof( path ), "level-%02d.list", IntermediateResultIndex );
00488     }
00489   this->OutputWarp( path );
00490   
00491   if ( incrementCount )
00492     ++IntermediateResultIndex;
00493 }
00494 
00495 CallbackResult 
00496 ImagePairNonrigidRegistrationCommandLine::Register ()
00497 {
00498   const double baselineTime = Timers::GetTimeProcess();
00499   CallbackResult Result = this->Superclass::Register();
00500   const int elapsed = static_cast<int>( Timers::GetTimeProcess() - baselineTime );
00501 
00502   if ( Time ) 
00503     {
00504     FILE *tfp = fopen( Time, "w" );
00505     
00506     if ( tfp ) 
00507       {
00508       fprintf( tfp, "%d\n", elapsed );
00509       fclose( tfp );
00510       } 
00511     else
00512       {
00513       std::cerr << "Could not open time file " << Time << "\n";
00514       }
00515     }
00516   
00517   return Result;
00518 }
00519 
00520 } // namespace cmtk
00521 
00522 void
00523 cmtkImagePairNonrigidRegistrationCommandLineDispatchSIGUSR1( int sig )
00524 {
00525   fprintf( stderr, "Received USR1 (%d) signal. Writing intermediate result #%d.\nNote that this result is not final.\n", sig, cmtk::ImagePairNonrigidRegistrationCommandLine::StaticThis->IntermediateResultIndex );
00526 
00527 #ifndef _MSC_VER
00528   // set signal handler again.
00529   signal( sig, cmtkImagePairNonrigidRegistrationCommandLineDispatchSIGUSR1 );
00530 #endif
00531   
00532   // write intermediate result. give "false" flag for index increment to 
00533   // preserve to final numbering of levels.
00534   cmtk::ImagePairNonrigidRegistrationCommandLine::StaticThis->OutputIntermediate( true /* Increment count*/ );
00535 }
00536 
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines