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
|
/*
//
// Copyright 1997-2009 Torsten Rohlfing
//
// Copyright 2004-2011, 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 <cmtkconfig.h>
#include <System/cmtkCommandLine.h>
#include <System/cmtkExitException.h>
#include <System/cmtkConsole.h>
#include <System/cmtkProgressConsole.h>
#include <IO/cmtkVolumeIO.h>
#include <Segmentation/cmtkAtlasSegmentation.h>
#include <list>
int
doMain( const int argc, const char* argv[] )
{
bool fast = false;
std::string targetImageName;
std::string atlasImageName;
std::string atlasLabelName;
std::string outImageName;
try
{
cmtk::CommandLine cl( cmtk::CommandLine::PROPS_XML );
cl.SetProgramInfo( cmtk::CommandLine::PRG_TITLE, "Atlas-based segmentation" );
cl.SetProgramInfo( cmtk::CommandLine::PRG_DESCR, "Register a target image to an atlas, using affine followed by nonrigid B-spline registration, then reformat the atlas label map to the target image." );
cl.SetProgramInfo( cmtk::CommandLine::PRG_CATEG, "CMTK.Segmentation" );
typedef cmtk::CommandLine::Key Key;
cl.AddSwitch( Key( 'f', "fast" ), &fast, true, "Fast mode." );
cl.AddParameter( &targetImageName, "TargetImage", "Target image path. This is the image to be segmented." )
->SetProperties( cmtk::CommandLine::PROPS_IMAGE );
cl.AddParameter( &atlasImageName, "AtlasImage", "Atlas image path. This is the structural channel (e.g., T1-weighted MRI) of the atlas to be used." )
->SetProperties( cmtk::CommandLine::PROPS_IMAGE );
cl.AddParameter( &atlasLabelName, "AtlasLabels", "Atlas label image path. This is the label map to be reformatted to the target image." )
->SetProperties( cmtk::CommandLine::PROPS_IMAGE | cmtk::CommandLine::PROPS_LABELS );
cl.AddParameter( &outImageName, "OutputImage", "Output image path. This is where the reformatted label map is written." )
->SetProperties( cmtk::CommandLine::PROPS_IMAGE | cmtk::CommandLine::PROPS_LABELS | cmtk::CommandLine::PROPS_OUTPUT );
cl.Parse( argc, argv );
}
catch ( const cmtk::CommandLine::Exception& e )
{
cmtk::StdErr << e << "\n";
throw cmtk::ExitException( 1 );
}
// Instantiate programm progress indicator.
cmtk::ProgressConsole progressIndicator( "Atlas-based Segmentation" );
cmtk::UniformVolume::SmartPtr targetImg( cmtk::VolumeIO::ReadOriented( targetImageName ) );
if ( !targetImg )
{
cmtk::StdErr << "ERROR: could not read target image " << targetImageName << "\n";
throw cmtk::ExitException( 1 );
}
cmtk::UniformVolume::SmartPtr atlasImg( cmtk::VolumeIO::ReadOriented( atlasImageName ) );
if ( !atlasImg )
{
cmtk::StdErr << "ERROR: could not read atlas image " << atlasImageName << "\n";
throw cmtk::ExitException( 1 );
}
cmtk::UniformVolume::SmartPtr atlasLbl( cmtk::VolumeIO::ReadOriented( atlasLabelName ) );
if ( !atlasLbl )
{
cmtk::StdErr << "ERROR: could not read atlas labels " << atlasLabelName << "\n";
throw cmtk::ExitException( 1 );
}
cmtk::AtlasSegmentation segment( targetImg, atlasImg, atlasLbl );
segment.SetFast( fast );
cmtk::VolumeIO::Write( *(segment.GetLabelMap()), outImageName );
return 0;
}
#include "cmtkSafeMain"
|