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 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259
|
/*=========================================================================
*
* Copyright Insight Software Consortium
*
* Licensed under the Apache License, Version 2.0 (the "License");
* you may not use this file except in compliance with the License.
* You may obtain a copy of the License at
*
* http://www.apache.org/licenses/LICENSE-2.0.txt
*
* Unless required by applicable law or agreed to in writing, software
* distributed under the License is distributed on an "AS IS" BASIS,
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
* See the License for the specific language governing permissions and
* limitations under the License.
*
*=========================================================================*/
// Software Guide : BeginLatex
//
// \index{itk::bio::CellularAggregate}
//
// The following example illustrates the use of Cellular Algorithms for performing image segmentation.
// Cellular algorithms are implemented by combining the following classes
//
// \subdoxygen{bio}{CellularAggregate}
// \subdoxygen{bio}{Cell}
//
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
#include "itkBioCellularAggregate.h"
// Software Guide : EndCodeSnippet
#include "itkImageFileReader.h"
#include "itkVTKPolyDataWriter.h"
int main( int argc, char *argv[] )
{
if( argc < 8 )
{
std::cerr << "Missing Parameters " << std::endl;
std::cerr << "Usage: " << argv[0];
std::cerr << " inputImage seedX seedY lowThreshold highThreshold iterations outputMesh.vtk" << std::endl;
return EXIT_FAILURE;
}
// Software Guide : BeginLatex
//
// We now define the image type using a pixel type and a particular
// dimension. In this case the \code{float} type is used for the pixels due
// to the requirements of the smoothing filter.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
typedef float InternalPixelType;
const unsigned int Dimension = 2;
typedef itk::Image< InternalPixelType, Dimension > ImageType;
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// The \subdoxygen{bio}{CellularAggregate} class must be instantiated using
// the dimension of the image to be segmented.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
typedef itk::bio::CellularAggregate< Dimension > CellularAggregateType;
typedef CellularAggregateType::BioCellType CellType;
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// Then an object of this class can be constructed by invoking the
// \code{New} operator and receiving the result in a \code{SmartPointer},
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
CellularAggregateType::Pointer cellularAggregate
= CellularAggregateType::New();
// Software Guide : EndCodeSnippet
// We instantiate reader and writer types
typedef itk::ImageFileReader< ImageType > ReaderType;
ReaderType::Pointer reader = ReaderType::New();
reader->SetFileName( argv[1] );
std::cout << "Filename = " << argv[1] << std::endl;
try
{
reader->Update();
}
catch( itk::ExceptionObject & excep )
{
std::cerr << "Exception caught !" << std::endl;
std::cerr << excep << std::endl;
return EXIT_FAILURE;
}
// Software Guide : BeginLatex
//
// The CellularAggregate considers the image as a chemical substrate in
// which the Cells are going to develop. The intensity values of the image
// will influence the behavior of the Cells, in particular they will
// intervine to regulate the Cell Cycle. A Cellular Aggregate could be
// gathering information from several images simultaneously, in this context
// each image can bee seen as a map of concentration of a particular
// chemical compound. The set of images will describe the chemical
// composition of the extra cellular matrix.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
cellularAggregate->AddSubstrate( reader->GetOutput() );
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// The initialization of the algorithm requires the user to provide a seed
// point. It is convenient to select this point to be placed in a
// \emph{typical} region of the anatomical structure to be segmented. A
// small neighborhood around the seed point will be used to compute the
// initial mean and standard deviation for the inclusion criterion. The
// seed is passed in the form of a \doxygen{Index} to the \code{SetSeed()}
// method.
//
// \index{itk::ConfidenceConnectedImageFilter!SetSeed()}
// \index{itk::ConfidenceConnectedImageFilter!SetInitialNeighborhoodRadius()}
//
// Software Guide : EndLatex
ImageType::IndexType index;
index[0] = atoi( argv[2] );
index[1] = atoi( argv[3] );
CellType::PointType position;
reader->GetOutput()->TransformIndexToPhysicalPoint( index, position );
std::cout << "Egg position index = " << index << std::endl;
std::cout << "Egg position point = " << position << std::endl;
// Software Guide : BeginLatex
//
// Individual Cell do not derive from the \doxygen{Object} class in order to
// avoid the penalties of Mutex operations when passing pointers to them.
// The Creation of a new cell is done by invoking the normal \code{new}
// operator.
//
// \index{itk::bio::Cell!Creation}
// \index{itk::bio::Cell!pointer}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
CellType * egg = CellType::CreateEgg();
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// In this particular example, the Cell cycle is going to be controled
// mostly by the intensity values of the image. These values are asimilated
// to concentrations of a particular chemical compound. Cell will feel
// compfortable at when the concentration of this chemical is inside a
// particular range. In this circumstances cells will be able to
// proliferate. When the chemical concentration is out of the range, cell
// will not enter their division stage and will anchor to the cellular
// matrix. The values defining this range can be set by invoking the methods
// \code{SetChemoAttractantHighThreshold} and
// \code{SetChemoAttractantLowThreshold). These to methods are static and
// set the values to be used by all the cells.
//
// \index{itk::bio::Cell!SetChemoAttractantLowThreshold}
// \index{itk::bio::Cell!SetChemoAttractantHighThreshold}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
CellType::SetChemoAttractantLowThreshold( atof( argv[4] ) );
CellType::SetChemoAttractantHighThreshold( atof( argv[5] ) );
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// The newly created Cell is passed to the \code{CellularAggregate} object
// that will take care of controling the development of the cells.
//
// \index{itk::bio::CellularAggregate!SetEgg}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
cellularAggregate->SetEgg( egg, position );
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// The CellularAggregate will update the life cycle of all the cells in an
// iterative way. The User must select how many iterations to run.
// CellularAlgorithms can in principle run forever. It is up to the User to
// define an stopping criterion. One of the simplest options is to set a
// limit to the number of iterations, by invoking the AdvanceTimeStep()
// method inside a for loop.
//
// \index{itk::bio::CellularAggregate!SetEgg}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
unsigned int numberOfIterations = atoi( argv[6] );
std::cout << "numberOfIterations " << numberOfIterations << std::endl;
for(unsigned int i=0; i<numberOfIterations; i++)
{
cellularAggregate->AdvanceTimeStep();
}
// Software Guide : EndCodeSnippet
std::cout << " Final number of Cells = " << cellularAggregate->GetNumberOfCells() << std::endl;
// Write the mesh to a file
//
typedef itk::VTKPolyDataWriter< CellularAggregateType::MeshType > WriterType;
WriterType::Pointer writer = WriterType::New();
writer->SetInput( cellularAggregate->GetMesh() );
writer->SetFileName( argv[7] );
try
{
writer->Update();
}
catch( itk::ExceptionObject & excep )
{
std::cerr << "Exception caught !" << std::endl;
std::cerr << excep << std::endl;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
}
|