ITK  5.4.0
Insight Toolkit
Examples/SpatialObjects/SpatialObjectToImageStatisticsCalculator.cxx
/*=========================================================================
*
* Copyright NumFOCUS
*
* 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
*
* https://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::Spatial\-Object\-To\-Image\-Statistics\-Calculator}
// This example describes how to use the
// \doxygen{SpatialObjectToImageStatisticsCalculator} to compute statistics
// of an \doxygen{Image} only in a region defined inside a given
// \doxygen{SpatialObject}.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
// Software Guide : EndCodeSnippet
#include "itkImage.h"
int
main(int, char *[])
{
// Software Guide : BeginLatex
//
// We first create a test image using the \doxygen{RandomImageSource}
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
using ImageType = itk::Image<unsigned char, 2>;
using RandomImageSourceType = itk::RandomImageSource<ImageType>;
auto randomImageSource = RandomImageSourceType::New();
size[0] = 10;
size[1] = 10;
randomImageSource->SetSize(size);
randomImageSource->Update();
ImageType::Pointer image = randomImageSource->GetOutput();
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// Next we create an \doxygen{EllipseSpatialObject} with a radius of 2.
// We also move the ellipse to the center of the image.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
using EllipseType = itk::EllipseSpatialObject<2>;
auto ellipse = EllipseType::New();
ellipse->SetRadiusInObjectSpace(2);
offset.Fill(5);
ellipse->SetCenterInObjectSpace(offset);
ellipse->Update();
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// Then we can create the
// \doxygen{SpatialObjectToImageStatisticsCalculator}.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
using CalculatorType =
auto calculator = CalculatorType::New();
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// We pass a pointer to the image to the calculator.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
calculator->SetImage(image);
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// We also pass the SpatialObject. The statistics will be computed inside
// the SpatialObject (Internally the calculator is using the
// \code{IsInside()} function).
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
calculator->SetSpatialObject(ellipse);
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// At the end we trigger the computation via the \code{Update()} function
// and we can retrieve the mean and the covariance matrix using
// \code{GetMean()} and \code{GetCovarianceMatrix()} respectively.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
calculator->Update();
std::cout << "Sample mean = " << calculator->GetMean() << std::endl;
std::cout << "Sample covariance = " << calculator->GetCovarianceMatrix();
// Software Guide : EndCodeSnippet
return EXIT_SUCCESS;
}
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itkEllipseSpatialObject.h
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itkSpatialObjectToImageStatisticsCalculator.h
itkImage.h
itk::EllipseSpatialObject
Definition: itkEllipseSpatialObject.h:38
itk::SpatialObjectToImageStatisticsCalculator
Definition: itkSpatialObjectToImageStatisticsCalculator.h:37
itk::RandomImageSource
Generate an n-dimensional image of random pixel values.
Definition: itkRandomImageSource.h:55
itkRandomImageSource.h
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
New
static Pointer New()
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:83
itk::FixedArray::Fill
void Fill(const ValueType &)