cmtkAtlasSegmentation.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: 2398 $
00026 //
00027 //  $LastChangedDate: 2010-10-05 14:54:37 -0700 (Tue, 05 Oct 2010) $
00028 //
00029 //  $LastChangedBy: torstenrohlfing $
00030 //
00031 */
00032 
00033 #include "cmtkAtlasSegmentation.h"
00034 
00035 #include <Registration/cmtkRegistrationCallback.h>
00036 #include <Registration/cmtkAffineRegistration.h>
00037 #include <Registration/cmtkElasticRegistration.h>
00038 #include <Registration/cmtkReformatVolume.h>
00039 
00040 cmtk::AtlasSegmentation::AtlasSegmentation
00041 ( UniformVolume::SmartPtr& targetImage, UniformVolume::SmartPtr& atlasImage, UniformVolume::SmartPtr& atlasLabels )  
00042   : m_Verbose( false ),
00043     m_Fast( false ),
00044     m_TargetImage( targetImage ),
00045     m_AtlasImage( atlasImage ),
00046     m_AtlasLabels( atlasLabels ),
00047     m_LabelMap( NULL )
00048 {
00049 }
00050 
00051 void
00052 cmtk::AtlasSegmentation
00053 ::RegisterAffine()
00054 {
00055   AffineRegistration ar;
00056   ar.SetVolume_1( this->m_TargetImage );
00057   ar.SetVolume_2( this->m_AtlasImage );
00058   
00059   // run 6 DOFs first, then 9 at each level.
00060   ar.AddNumberDOFs( 6 );
00061   ar.AddNumberDOFs( 9 );
00062   
00063   ar.SetInitialAlignCenters( true );
00064   ar.SetExploration( 4 * this->m_TargetImage->GetMaxDelta() );
00065   ar.SetAccuracy( .1 * this->m_TargetImage->GetMaxDelta() );
00066   ar.SetSampling( 2 * this->m_TargetImage->GetMaxDelta() );
00067 
00068   ar.SetUseOriginalData( !this->m_Fast );
00069   
00070   if ( this->m_Verbose ) 
00071     {
00072     StdErr << "Affine registration...";
00073     StdErr.flush();
00074     }
00075   ar.Register();
00076   if ( this->m_Verbose ) 
00077     {
00078     StdErr << " done.\n";
00079     }
00080   
00081   this->m_AffineXform = ar.GetTransformation();
00082 }
00083 
00084 void
00085 cmtk::AtlasSegmentation
00086 ::RegisterSpline()
00087 {
00088   ElasticRegistration er;
00089   er.SetVolume_1( this->m_TargetImage );
00090   er.SetVolume_2( this->m_AtlasImage );
00091   er.SetInitialTransformation( this->GetAffineXform() );
00092   
00093   er.SetUseOriginalData( !this->m_Fast );
00094   er.SetFastMode( this->m_Fast );
00095 
00096   const Types::Coordinate minSize = std::min( std::min( this->m_TargetImage->Size[0], this->m_TargetImage->Size[1] ), this->m_TargetImage->Size[2] );
00097   er.SetGridSpacing( minSize / 2 );
00098   er.SetRefineGrid( std::max<int>( 0, static_cast<int>( (log( minSize / this->m_TargetImage->GetMaxDelta() ) / log(2.0)) - 3 ) ) );
00099   er.SetDelayRefineGrid( !this->m_Fast );
00100   
00101   er.SetGridEnergyWeight( 1e-1 );
00102   er.SetAdaptiveFixParameters( true );
00103 
00104   er.SetAlgorithm( 3 );
00105   er.SetExploration( minSize / 8 );
00106   er.SetAccuracy( .1 * this->m_TargetImage->GetMinDelta() );
00107   er.SetSampling( 2 * this->m_TargetImage->GetMaxDelta() );  
00108   
00109   if ( this->m_Verbose ) 
00110     {
00111     StdErr << "Nonrigid registration...";
00112     StdErr.flush();
00113     }
00114   er.Register();
00115   if ( this->m_Verbose ) 
00116     {
00117     StdErr << " done.\n";
00118     }
00119   this->m_WarpXform = er.GetTransformation();
00120 }
00121 
00122 void
00123 cmtk::AtlasSegmentation
00124 ::ReformatLabels()
00125 {
00126   ReformatVolume reformat;
00127   reformat.SetInterpolation( Interpolators::PARTIALVOLUME );
00128   reformat.SetPaddingValue( 0 );
00129   reformat.SetUsePaddingValue( true );
00130   reformat.SetReferenceVolume( this->m_TargetImage );
00131   reformat.SetFloatingVolume( this->m_AtlasLabels );
00132 
00133   WarpXform::SmartPtr warpXform = this->GetWarpXform();
00134   reformat.SetWarpXform( warpXform );
00135 
00136   this->m_LabelMap = UniformVolume::SmartPtr( reformat.PlainReformat() );
00137 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines