ITK  6.0.0
Insight Toolkit
SphinxExamples/src/Filtering/ImageGrid/FitSplineIntoPointSet/Code.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.
*
*=========================================================================*/
#include "itkPointSet.h"
#include "itkImage.h"
#include "itkVectorImage.h"
int
main()
{
constexpr unsigned int ParametricDimension = 1;
constexpr unsigned int DataDimension = 2;
auto pointSet = PointSetType::New();
PointSetType::PointType param0, param1, param2;
param0[0] = 0.0;
DataType p0;
p0[0] = 10.0;
p0[1] = 10.0;
pointSet->SetPoint(0, param0);
pointSet->SetPointData(0, p0);
param1[0] = 1.0;
DataType p1;
p1[0] = 80.0;
p1[1] = 50.0;
pointSet->SetPoint(1, param1);
pointSet->SetPointData(1, p1);
param2[0] = 2.0;
DataType p2;
p2[0] = 180.0;
p2[1] = 180.0;
pointSet->SetPoint(2, param2);
pointSet->SetPointData(2, p2);
auto splineFilter = SplineFilterType::New();
int splineorder = 2; // complexity of the spline
SplineFilterType::ArrayType ncontrol;
ncontrol[0] = splineorder + 1;
SplineFilterType::ArrayType closedim;
closedim[0] = 0;
ImageType::PointType parametricDomainOrigin;
parametricDomainOrigin[0] = 0.0;
ImageType::SpacingType parametricDomainSpacing;
parametricDomainSpacing[0] = 0.0001; // this determines the sampling of the continuous B-spline object.
ImageType::SizeType parametricDomainSize;
parametricDomainSize[0] = 2.0 / parametricDomainSpacing[0] + 1;
splineFilter->SetGenerateOutputImage(true); // the only reason to turn this off is if one only wants to use the
// control point lattice for further processing
splineFilter->SetInput(pointSet);
splineFilter->SetSplineOrder(splineorder);
splineFilter->SetNumberOfControlPoints(ncontrol);
splineFilter->SetNumberOfLevels(3);
splineFilter->SetCloseDimension(closedim);
splineFilter->SetSize(parametricDomainSize);
splineFilter->SetSpacing(parametricDomainSpacing);
splineFilter->SetOrigin(parametricDomainOrigin);
splineFilter->Update();
// The output will consist of a 1-D image where each voxel contains the
// (x,y,z) locations of the points
using OutputImageType = itk::Image<unsigned char, 2>;
auto outputImage = OutputImageType::New();
size.Fill(200);
start.Fill(0);
OutputImageType::RegionType region(start, size);
outputImage->SetRegions(region);
outputImage->Allocate();
outputImage->FillBuffer(0);
for (unsigned int i = 0; i < splineFilter->GetOutput()->GetLargestPossibleRegion().GetSize()[0]; ++i)
{
ImageType::IndexType splineIndex;
splineIndex[0] = i;
DataType outputPixel = splineFilter->GetOutput()->GetPixel(splineIndex);
index[0] = outputPixel[0];
index[1] = outputPixel[1];
outputImage->SetPixel(index, 255);
}
itk::WriteImage(outputImage, "spline.png");
return EXIT_SUCCESS;
};
itk::PointSet
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:81
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::BSplineScatteredDataPointSetToImageFilter
Image filter which provides a B-spline output approximation.
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:132
itk::Vector
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itk::Index::Fill
void Fill(IndexValueType value)
Definition: itkIndex.h:272
itk::Size::Fill
void Fill(SizeValueType value)
Definition: itkSize.h:211
itkVectorImage.h
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itkBSplineScatteredDataPointSetToImageFilter.h
itkImageFileWriter.h
itk::Size::SetSize
void SetSize(const SizeValueType val[VDimension])
Definition: itkSize.h:179
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itkPointSet.h
New
static Pointer New()
itk::WriteImage
ITK_TEMPLATE_EXPORT void WriteImage(TImagePointer &&image, const std::string &filename, bool compress=false)
Definition: itkImageFileWriter.h:256