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
|
/*
//
// Copyright 1997-2009 Torsten Rohlfing
//
// Copyright 2004-2012 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 <cmtkconfig.h>
#include <Base/cmtkUniformVolume.h>
#include <Base/cmtkAffineXform.h>
#include <Base/cmtkSplineWarpXform.h>
namespace
cmtk
{
/** \addtogroup Segmentation */
//@{
/** Class for integrated atlas-based segmentation.
* This class encapsulates the following stages of atlas-based image segmentation: 1) linear image-to-atlas registration, 2) non-linear
* image-to-atlas registration, 3) label map reformatting.
*/
class AtlasSegmentation
{
public:
/// Constructor: compute registrations.
AtlasSegmentation( UniformVolume::SmartPtr& targetImage, UniformVolume::SmartPtr& atlasImage, UniformVolume::SmartPtr& atlasLabels );
/// Get affine transformation.
AffineXform::SmartPtr& GetAffineXform()
{
if ( ! this->m_AffineXform )
this->RegisterAffine();
return this->m_AffineXform;
}
/// Get nonrigid transformation.
WarpXform::SmartPtr GetWarpXform()
{
if ( ! this->m_WarpXform )
this->RegisterSpline();
return this->m_WarpXform;
}
/// Get nonrigid spline transformation.
SplineWarpXform::SmartPtr GetSplineWarpXform()
{
return SplineWarpXform::SmartPtr::DynamicCastFrom( this->GetWarpXform() );
}
/// Get reformatted label map.
UniformVolume::SmartPtr& GetLabelMap()
{
if ( ! this->m_LabelMap )
this->ReformatLabels();
return this->m_LabelMap;
}
/// Set fast flag.
void SetFast( const bool fast )
{
this->m_Fast = fast;
}
private:
/// Flag for "fast" computation.
bool m_Fast;
/// Target image.
UniformVolume::SmartPtr m_TargetImage;
/// Atlas image.
UniformVolume::SmartPtr m_AtlasImage;
/// Atlas labels.
UniformVolume::SmartPtr m_AtlasLabels;
/// Affine registration transformation.
AffineXform::SmartPtr m_AffineXform;
/// Compute affine registration.
void RegisterAffine();
/// Nonrigid, B-spline transformation.
WarpXform::SmartPtr m_WarpXform;
/// Compute spline registration.
void RegisterSpline();
/// Output label map.
UniformVolume::SmartPtr m_LabelMap;
/// Compute label map.
void ReformatLabels();
};
//@}
} // namespace cmtk
|