File: cmtkAtlasSegmentation.cxx

package info (click to toggle)
cmtk 3.3.1p2%2Bdfsg-4
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 10,524 kB
  • sloc: cpp: 87,098; ansic: 23,347; sh: 3,896; xml: 1,551; perl: 707; makefile: 334
file content (134 lines) | stat: -rw-r--r-- 4,176 bytes parent folder | download | duplicates (5)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
/*
//
//  Copyright 1997-2009 Torsten Rohlfing
//
//  Copyright 2004-2013 SRI International
//
//  This file is part of the Computational Morphometry Toolkit.
//
//  http://www.nitrc.org/projects/cmtk/
//
//  The Computational Morphometry Toolkit is free software: you can
//  redistribute it and/or modify it under the terms of the GNU General Public
//  License as published by the Free Software Foundation, either version 3 of
//  the License, or (at your option) any later version.
//
//  The Computational Morphometry Toolkit is distributed in the hope that it
//  will be useful, but WITHOUT ANY WARRANTY; without even the implied
//  warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
//  GNU General Public License for more details.
//
//  You should have received a copy of the GNU General Public License along
//  with the Computational Morphometry Toolkit.  If not, see
//  <http://www.gnu.org/licenses/>.
//
//  $Revision: 5436 $
//
//  $LastChangedDate: 2018-12-10 19:01:20 -0800 (Mon, 10 Dec 2018) $
//
//  $LastChangedBy: torstenrohlfing $
//
*/

#include "cmtkAtlasSegmentation.h"

#include <Registration/cmtkRegistrationCallback.h>
#include <Registration/cmtkAffineRegistration.h>
#include <Registration/cmtkElasticRegistration.h>
#include <Registration/cmtkReformatVolume.h>

#include <System/cmtkDebugOutput.h>
#include <System/cmtkExitException.h>

cmtk::AtlasSegmentation::AtlasSegmentation
( UniformVolume::SmartPtr& targetImage, UniformVolume::SmartPtr& atlasImage, UniformVolume::SmartPtr& atlasLabels ) :
    m_Fast( false ),
    m_TargetImage( targetImage ),
    m_AtlasImage( atlasImage ),
    m_AtlasLabels( atlasLabels ),
    m_LabelMap( NULL )
{
}

void
cmtk::AtlasSegmentation
::RegisterAffine()
{
  AffineRegistration ar;
  ar.SetVolume_1( this->m_TargetImage );
  ar.SetVolume_2( this->m_AtlasImage );
  
  // run 6 DOFs first, then 9 at each level.
  ar.AddNumberDOFs( 6 );
  ar.AddNumberDOFs( 9 );
  
  ar.SetInitialAlignCenters( true );
  ar.SetExploration( 4 * this->m_TargetImage->GetMaxDelta() );
  ar.SetAccuracy( .1 * this->m_TargetImage->GetMaxDelta() );
  ar.SetSampling( 2 * this->m_TargetImage->GetMaxDelta() );

  ar.SetUseOriginalData( !this->m_Fast );
  
  (DebugOutput( 1 ) << "Affine registration...").flush();
  ar.Register();
  DebugOutput( 1 ) << " done.\n";

  try
    {
    this->m_AffineXform = ar.GetTransformation();
    }
  catch ( const AffineXform::MatrixType::SingularMatrixException& )
    {
    StdErr << "ERROR: singular matrix in AffineRegistration::GetTransformation()\n";
    throw ExitException( 1 );
    }
}

void
cmtk::AtlasSegmentation
::RegisterSpline()
{
  ElasticRegistration er;
  er.SetVolume_1( this->m_TargetImage );
  er.SetVolume_2( this->m_AtlasImage );
  er.SetInitialTransformation( this->GetAffineXform() );
  
  er.SetUseOriginalData( !this->m_Fast );
  er.SetFastMode( this->m_Fast );

  const Types::Coordinate minSize = std::min( std::min( this->m_TargetImage->m_Size[0], this->m_TargetImage->m_Size[1] ), this->m_TargetImage->m_Size[2] );
  er.SetGridSpacing( minSize / 2 );
  er.SetRefineGrid( std::max<int>( 0, static_cast<int>( (log( minSize / this->m_TargetImage->GetMaxDelta() ) / log(2.0)) - 3 ) ) );
  er.SetDelayRefineGrid( !this->m_Fast );
  
  er.SetGridEnergyWeight( 1e-1f );
  er.SetAdaptiveFixParameters( true );

  er.SetAlgorithm( 3 );
  er.SetExploration( minSize / 8 );
  er.SetAccuracy( .1 * this->m_TargetImage->GetMinDelta() );
  er.SetSampling( 2 * this->m_TargetImage->GetMaxDelta() );  
  
  (DebugOutput( 1 ) << "Nonrigid registration...").flush();
  er.Register();
  DebugOutput( 1 ) << " done.\n";

  this->m_WarpXform = er.GetTransformation();
}

void
cmtk::AtlasSegmentation
::ReformatLabels()
{
  ReformatVolume reformat;
  reformat.SetInterpolation( Interpolators::PARTIALVOLUME );
  reformat.SetPaddingValue( 0 );
  reformat.SetUsePaddingValue( true );
  reformat.SetReferenceVolume( this->m_TargetImage );
  reformat.SetFloatingVolume( this->m_AtlasLabels );

  WarpXform::SmartPtr warpXform = this->GetWarpXform();
  reformat.SetWarpXform( warpXform );

  this->m_LabelMap = UniformVolume::SmartPtr( reformat.PlainReformat() );
}