[Insight-users] question ? [ITK]
Luis Ibanez
luis.ibanez at kitware.com
Sun Aug 19 09:36:55 EDT 2012
Hi Agata,
Please find attached the working code.
Here is the hierarchy that you may want to keep in mind:
- One TPS transform has two PointSets
- One PointSet has a PointContainer
- One PointContainer have PointIterators (const and non-const)
You can get points from a PointSet by using GetPointid)
PointType point = pointSet->GetPoint( pointId );
or you can get them from the PointContainer by using GetElement(id)
PointType point = pointContainer->GetElement( pointId );
or you can get them via Iterators from the PointContainer
PointConstIterator pointItr = pointContainer->Begin();
PointType point = pointItr.Value();
++pointItr;
...etc...
Note that you shouldn't use "int" or "short" or "long" for the
PointId type. One must use PointIdType (which is typedefed
from the ElementIdentifier type), this is very important if you
are working on 64 bits.
The attached code shows you both ways of doing this
at the end of the file.
The key lines of code are:
TPSTransformType::InputPointType trackerPointSourcePosition;
TPSTransformType::OutputPointType trackerPointNewPosition;
// Test transforming the input landmark themselves
for( PointIdType i = 0; i < sourceLandMarks->GetNumberOfPoints(); i++ )
{
trackerPointSourcePosition = sourceLandMarks->GetPoint(i);
trackerPointNewPosition =
tps->TransformPoint(trackerPointSourcePosition);
transformed_points_container->InsertElement(i,trackerPointNewPosition);
std::cout << "Input " << trackerPointSourcePosition << " : ";
std::cout << "Output " << trackerPointNewPosition << std::endl;
}
// Do the same using Iterators
std::cout << std::endl << "Now using Iterators " << std::endl;
PointConstIterator inputPointItr = sourceLandMarkContainer->Begin();
PointConstIterator inputPointEnd = sourceLandMarkContainer->End();
PointIdType pointId = itk::NumericTraits< PointIdType >::Zero;
while( inputPointItr != inputPointEnd )
{
trackerPointNewPosition = tps->TransformPoint(inputPointItr.Value());
transformed_points_container->InsertElement(pointId,trackerPointNewPosition);
std::cout << "Input " << inputPointItr.Value() << " : ";
std::cout << "Output " << trackerPointNewPosition << std::endl;
++inputPointItr;
}
Hope this helps,
Luis
------------------------------------------------------------
On Sat, Aug 18, 2012 at 12:17 PM, Agata Krasoń <agatakrason at gmail.com>wrote:
> Dear Luis,
>
> I am trying to use tps transformation in ITK.
> ...
> I created sources and target points.
> I made Tps transformation.
> I need to receive points after registration
>
> ...
>
> I have a problem with iterators.
> I would appreciate for any help please.
>
> int main(int argc, char* argv[])
> {
>
> double tab[30] =
> {
> // Tracker
> -81.29,-31.07, -770.58,
> -83.11, -21.26, -822.64,
> -93.45, -32.44, -858.72,
> -68.08, -126.89,-813.07,
> -61.04, 75.74, -808.36,
>
> // Image
> 140.6, 230.7, -30.5,
> 140.2, 231.7, -71.1,
> 144.8, 235.9, -116.1,
> 45.8, 220.2, -66.7,
> 231.6, 211.3, -66.1
> };
>
> const unsigned int Dimension = 3;
> typedef unsigned char PixelType;
> typedef double CoordinateRepType;
> typedef itk::ThinPlateSplineKernelTransform<
> CoordinateRepType,Dimension> TPSTransformType;
> typedef itk::Point< CoordinateRepType,Dimension > PointType;
> typedef std::vector< PointType > PointArrayType;
> typedef TPSTransformType::PointSetType PointSetType;
> typedef PointSetType::Pointer PointSetPointer;
> typedef PointSetType::PointIdentifier PointIdType;
>
>
> // Landmarks correspondances may be associated with the
> SplineKernelTransforms
> // via Point Set containers. Let us define containers for the landmarks.
> PointSetType::Pointer sourceLandMarks = PointSetType::New();
> PointSetType::Pointer targetLandMarks = PointSetType::New();
> PointType trackerPoint; PointType imagePoint;
> PointSetType::PointsContainer::Pointer sourceLandMarkContainer =
> sourceLandMarks->GetPoints();
> PointSetType::PointsContainer::Pointer targetLandMarkContainer =
> targetLandMarks->GetPoints();
>
>
> // 1 Landmark
> trackerPoint[0] = tab[0];
> trackerPoint[1] = tab[1];
> trackerPoint[2] = tab[2];
> imagePoint[0] = tab[15];
> imagePoint[1] = tab[16];
> imagePoint[2] = tab[17];
> sourceLandMarkContainer->InsertElement( 0,trackerPoint);
> targetLandMarkContainer->InsertElement(0,imagePoint);
>
>
> // 2 Landmark
> trackerPoint[0] = tab[3];
> trackerPoint[1] = tab[4];
> trackerPoint[2] = tab[5];
> imagePoint[0] = tab[18];
> imagePoint[1] = tab[19];
> imagePoint[2] = tab[20];
> sourceLandMarkContainer->InsertElement(1,trackerPoint);
> targetLandMarkContainer->InsertElement(1,imagePoint);
>
> // 3 Landmark
> trackerPoint[0] = tab[6];
> trackerPoint[1] = tab[7];
> trackerPoint[2] = tab[8];
> imagePoint[0] = tab[21];
> imagePoint[1] = tab[22];
> imagePoint[2] = tab[23];
> sourceLandMarkContainer->InsertElement( 2,trackerPoint);
> targetLandMarkContainer->InsertElement(2,imagePoint);
>
> // 4 Landmark
> trackerPoint[0] = tab[9];
> trackerPoint[1] = tab[10];
> trackerPoint[2] = tab[11];
> imagePoint[0] = tab[24];
> imagePoint[1] = tab[25];
> imagePoint[2] = tab[26];
> sourceLandMarkContainer->InsertElement( 3,trackerPoint);
> targetLandMarkContainer->InsertElement(3,imagePoint);
>
> // 5 Landmark
> trackerPoint[0] = tab[12];
> trackerPoint[1] = tab[13];
> trackerPoint[2] = tab[14];
> imagePoint[0] = tab[27];
> imagePoint[1] = tab[28];
> imagePoint[2] = tab[29];
>
> sourceLandMarkContainer->InsertElement(4,trackerPoint);
> targetLandMarkContainer->InsertElement(4,imagePoint);
>
> TPSTransformType::Pointer tps = TPSTransformType::New();
> tps->SetSourceLandmarks(sourceLandMarks);
> tps->SetTargetLandmarks(targetLandMarks);
> tps->ComputeWMatrix();
>
>
> PointSetType::Pointer transformed_points = PointSetType::New();
> PointSetType::PointsContainer::Pointer transformed_points_container =
> transformed_points->GetPoints();
> TPSTransformType::OutputPointType trackerPointNewPosition;
>
> for( int i = 0; i < sourceLandMarks->GetNumberOfPoints(); i++){
> PointSetType src_pnt = sourceLandMarks->GetPoints(); // Here It
> doesn't compile ?
> trackerPointNewPosition = tps->TransformPoint(trackerPoint);
>
> transformed_points_container->InsertElement(i,trackerPointNewPosition);
> std::cout<<" "<<trackerPointNewPosition <<std::endl;
> }
>
>
> return EXIT_SUCCESS;
>
>
> }
>
>
>
>
> I would appreciate for any help please.
>
> Best,
> Agata
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://www.itk.org/pipermail/insight-users/attachments/20120819/a8d757af/attachment.htm>
-------------- next part --------------
# This is the root ITK CMakeLists file.
cmake_minimum_required(VERSION 2.4)
if(COMMAND CMAKE_POLICY)
cmake_policy(SET CMP0003 NEW)
endif(COMMAND CMAKE_POLICY)
# This project is designed to be built outside the Insight source tree.
project(tps)
# Find ITK.
find_package(ITK REQUIRED)
include(${ITK_USE_FILE})
add_executable(tps main.cpp )
target_link_libraries(tps ${ITK_LIBRARIES})
-------------- next part --------------
A non-text attachment was scrubbed...
Name: main.cpp
Type: text/x-c++src
Size: 4693 bytes
Desc: not available
URL: <http://www.itk.org/pipermail/insight-users/attachments/20120819/a8d757af/attachment.cpp>
More information about the Insight-users
mailing list