Skip to content
Snippets Groups Projects
Commit 104a3bcd authored by Julien Michel's avatar Julien Michel
Browse files

MRG

parents 79e80945 57558ec7
No related branches found
No related tags found
No related merge requests found
......@@ -138,8 +138,8 @@ namespace otb
/** For each octave, we fill the imageList for the extremum search*/
m_ImageList->PushBack(m_determinantImage);
}
}
/*----------------------------------------------------*/
/* extremum Search over octave's scales */
/*----------------------------------------------------*/
......@@ -351,7 +351,7 @@ namespace otb
{
col = i%Largeur - rayon ;
raw = i/Largeur - rayon ;
dist = vcl_sqrt(col *col + raw * raw );
dist = vcl_sqrt(static_cast<double>(col *col + raw * raw) );
col +=rayon;
raw +=rayon; // Backup to the image coordinate axes
......
......@@ -59,6 +59,14 @@ template<class TOutputPointSet>
itkSetStringMacro(FileName);
itkGetStringMacro(FileName);
/** Get Macro*/
itkGetMacro( NumberOfPoints, long int );
itkGetMacro( MinX, double );
itkGetMacro( MaxX, double );
itkGetMacro( MinY, double );
itkGetMacro( MaxY, double );
virtual void GenerateOutputInformation(void);
/** Does the real work. */
......@@ -74,6 +82,10 @@ template<class TOutputPointSet>
std::string m_FileName; // The file to be read
long int m_NumberOfPoints;
double m_MinX;
double m_MaxX;
double m_MinY;
double m_MaxY;
private:
PointSetFileReader(const Self&); //purposely not implemented
......
......@@ -34,6 +34,10 @@ namespace otb
::PointSetFileReader() : otb::PointSetSource<TOutputPointSet>()
{
m_NumberOfPoints=-1;
m_MinX=0;
m_MaxX=0;
m_MinY=0;
m_MaxY=0;
}
template <class TOutputPointSet>
......@@ -85,7 +89,10 @@ namespace otb
otbDebugMacro(<< "Points count: " << header.GetPointRecordsCount());
m_NumberOfPoints = header.GetPointRecordsCount();
m_MinX = header.GetMinX();
m_MaxX = header.GetMaxX();
m_MinY = header.GetMinY();
m_MaxY = header.GetMaxY();
ifs.close();
}
......@@ -176,6 +183,11 @@ namespace otb
{
Superclass::PrintSelf(os, indent);
os << indent << "Number of points: " << this->m_NumberOfPoints << std::endl;
os << indent << std::setprecision(15);
os << indent << "Min X: " << this->m_MinX << std::endl;
os << indent << "Max X: " << this->m_MaxX << std::endl;
os << indent << "Min Y: " << this->m_MinY << std::endl;
os << indent << "Max Y: " << this->m_MaxY << std::endl;
os << indent << "m_FileName: " << this->m_FileName << "\n";
}
......
......@@ -61,6 +61,14 @@ TARGET_LINK_LIBRARIES(DEMToImageGenerator OTBProjections OTBIO OTBCommon ITKComm
ADD_EXECUTABLE(DEMToOrthoImageGenerator DEMToOrthoImageGenerator.cxx )
TARGET_LINK_LIBRARIES(DEMToOrthoImageGenerator OTBProjections OTBIO OTBCommon ITKCommon ITKIO otbossim )
ADD_EXECUTABLE(LidarReaderExample LidarReaderExample.cxx )
TARGET_LINK_LIBRARIES(LidarReaderExample OTBIO OTBCommon ITKCommon ITKIO otbossim )
IF(ITK_USE_REVIEW)
ADD_EXECUTABLE(LidarToImageExample LidarToImageExample.cxx )
TARGET_LINK_LIBRARIES(LidarToImageExample OTBIO OTBCommon ITKCommon ITKIO otbossim )
ENDIF(ITK_USE_REVIEW)
#ADD_EXECUTABLE(ImageReadExportVTK ImageReadExportVTK.cxx )
#TARGET_LINK_LIBRARIES(ImageReadExportVTK ITKCommon ITKIO)
......@@ -104,6 +112,7 @@ IF( NOT OTB_DISABLE_CXX_TESTING AND BUILD_TESTING )
SET(BASELINE ${OTB_DATA_ROOT}/Baseline/Examples/IO)
SET(DATA ${OTB_DATA_ROOT}/Input)
SET(INPUTDATA ${OTB_SOURCE_DIR}/Examples/Data)
SET(TEMP ${OTB_BINARY_DIR}/Testing/Temporary)
......@@ -142,9 +151,27 @@ ADD_TEST(DEMToImageGeneratorTest ${EXE_TESTS}
-0.002
${INPUTDATA}/DEM_srtm
)
IF(ITK_USE_REVIEW)
ADD_TEST(LidarToImageExampleTest ${EXE_TESTS}
--compare-image ${TOL} ${BASELINE}/lidar-image-4.hdr
${TEMP}/lidar-image-4.hdr
LidarToImageExampleTest
${DATA}/TO_core_last_zoom.las
${TEMP}/lidar-image-4.hdr
${TEMP}/lidar-image-4.png
1.0
5
4
)
ENDIF(ITK_USE_REVIEW)
ADD_EXECUTABLE(otbIOExamplesTests otbIOExamplesTests.cxx)
TARGET_LINK_LIBRARIES(otbIOExamplesTests otbossim OTBBasicFilters OTBCommon OTBDisparityMap OTBIO OTBSpatialReasoning OTBChangeDetection OTBFeatureExtraction OTBLearning OTBMultiScale OTBProjections ITKIO ITKAlgorithms ITKStatistics ITKCommon ${OTB_IO_UTILITIES_DEPENDENT_LIBRARIES})
IF(ITK_USE_REVIEW)
ADD_EXECUTABLE(otbIOExamplesTests2 otbIOExamplesTests2.cxx)
TARGET_LINK_LIBRARIES(otbIOExamplesTests2 otbossim OTBBasicFilters OTBCommon OTBDisparityMap OTBIO OTBSpatialReasoning OTBChangeDetection OTBFeatureExtraction OTBLearning OTBMultiScale OTBProjections ITKIO ITKAlgorithms ITKStatistics ITKCommon ${OTB_IO_UTILITIES_DEPENDENT_LIBRARIES})
ENDIF(ITK_USE_REVIEW)
ENDIF( NOT OTB_DISABLE_CXX_TESTING AND BUILD_TESTING )
/*=========================================================================
Program: ORFEO Toolbox
Language: C++
Date: $Date$
Version: $Revision$
Copyright (c) Centre National d'Etudes Spatiales. All rights reserved.
See OTBCopyright.txt for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notices for more information.
=========================================================================*/
// Software Guide : BeginCommandLineArgs
// INPUTS: {srs.las}
// Software Guide : EndCommandLineArgs
// Software Guide : BeginLatex
//
// This example describes how to read a lidar data file and to display the basic
// information about it. Here, OTB make use of libLAS.
//
// The first step toward the use of these filters is to include the proper header files.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
#include "otbPointSetFileReader.h"
#include "itkPointSet.h"
#include <fstream>
int main(int argc, char * argv[])
{
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// We need now to declare the data types that we will be using and instanciate the
// reader (which is a \doxygen{otb}{PointSetFileReader}).
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
typedef itk::PointSet <double, 2> PointSetType;
typedef otb::PointSetFileReader<PointSetType> PointSetFileReaderType;
PointSetFileReaderType::Pointer reader = PointSetFileReaderType::New();
reader->SetFileName(argv[1]);
reader->Update();
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// We can already display some interesting information such as the
// data extension and the number of points:
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
std::cout << "*** Data area *** " << std::endl;
std::cout << std::setprecision(15);
std::cout << " - Easting min: " << reader->GetMinX() << std::endl;
std::cout << " - Easting max: " << reader->GetMaxX() << std::endl;
std::cout << " - Northing min: " << reader->GetMinY() << std::endl;
std::cout << " - Northing max: " << reader->GetMaxY() << std::endl;
std::cout << "Data points: " << reader->GetNumberOfPoints() << std::endl;
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// We can also loop on the point to see each point with its coordinates and
// value. Be careful here, data set can have hundred of thousands of points:
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
PointSetType::Pointer data = reader->GetOutput();
unsigned long nPoints = data->GetNumberOfPoints();
for(unsigned long i=0; i < nPoints; ++i)
{
PointSetType::PointType point;
data->GetPoint(i,&point);
std::cout << point << " : ";
PointSetType::PixelType value;
data->GetPointData(i,&value);
std::cout << value << std::endl;
}
return EXIT_SUCCESS;
}
// Software Guide : EndCodeSnippet
/*=========================================================================
Program: ORFEO Toolbox
Language: C++
Date: $Date$
Version: $Revision$
Copyright (c) Centre National d'Etudes Spatiales. All rights reserved.
See OTBCopyright.txt for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notices for more information.
=========================================================================*/
#if defined(_MSC_VER)
#pragma warning ( disable : 4786 )
#endif
// Software Guide : BeginCommandLineArgs
// INPUTS: {TO_core_last_zoom.las}
// OUTPUTS: {lidar-image-4.hdr}, {lidar-image-4.png}
// 1.0 5 4
// Software Guide : EndCommandLineArgs
// Software Guide : BeginCommandLineArgs
// INPUTS: {TO_core_last_zoom.las}
// OUTPUTS: {lidar-image-8.hdr}, {lidar-image-8.png}
// 1.0 5 8
// Software Guide : EndCommandLineArgs
// Software Guide : BeginLatex
//
// This example describes how to convert a point set obtained from some lidar data
// to an image file. Lidar produce a point set which is irregular in terms of spatial
// sampling. To be able to generate an image, an interpolation is required.
//
// The interpolation is done using the
// \doxygen{itk}{BSplineScatteredDataPointSetToImageFilter} which uses BSplines. The
// method is fully describes in \cite{Tustison2005} and \cite{Lee1997}.
//
// The first step toward the use of these filters is to include the proper header files.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
#include "itkPointSet.h"
#include "otbImage.h"
#include "otbPointSetFileReader.h"
#include "itkBSplineScatteredDataPointSetToImageFilter.h"
#include "otbImageFileWriter.h"
// Software Guide : EndCodeSnippet
#include "itkRescaleIntensityImageFilter.h"
int main( int argc, char* argv[] )
{
if(argc!=7)
{
std::cout << argv[0] <<" <input_lidar_filename> <output_image_filename(double)>"
<< " <output_image_filename(unsigned char)>"
<< " <output_resolution> <spline_order>"
<< " <number_of_level>"
<< std::endl;
return EXIT_FAILURE;
}
// Software Guide : BeginLatex
//
// Then, we declare the type of the data that we are going to use. As
// lidar decribes the altitude of the point (often to centimeter accuracy),
// it is required to use real type to represent it.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
typedef double RealType;
typedef itk::Vector<RealType, 1> VectorType;
typedef itk::PointSet <VectorType, 2> PointSetType;
typedef otb::Image<RealType, 2> ImageType;
typedef itk::Image<VectorType, 2> VectorImageType;
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// Lidar data are read into a point set using the
// \doxygen{otb}{PointSetFileReader}. Its usage is very similar to
// the \doxygen{otb}{ImagePointSetFileReader}:
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
typedef otb::PointSetFileReader<PointSetType> PointSetFileReaderType;
PointSetFileReaderType::Pointer reader = PointSetFileReaderType::New();
reader->SetFileName(argv[1]);
reader->Update();
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// We can now prepare the parameters to pass to the interpolation filter:
// you have to be aware that the origin of the image is on the upper left
// corner and thus correspond to the minimun easting but the maximum northing.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
double resolution = atof(argv[4]);
int splineOrder = atoi(argv[5]);
int level = atoi(argv[6]);
// Software Guide : EndCodeSnippet
std::cout << "*** Data area *** " << std::endl;
std::cout << std::setprecision(15);
std::cout << " - Easting min: " << reader->GetMinX() << std::endl;
std::cout << " - Easting max: " << reader->GetMaxX() << std::endl;
std::cout << " - Northing min: " << reader->GetMinY() << std::endl;
std::cout << " - Northing max: " << reader->GetMaxY() << std::endl;
std::cout << "Data points: " << reader->GetNumberOfPoints() << std::endl;
// Software Guide : BeginCodeSnippet
ImageType::IndexType start;
start[0] = 0;
start[1] = 0;
ImageType::SizeType size;
size[0] = static_cast<long int >(
ceil((reader->GetMaxX()-reader->GetMinX()+1)/resolution))+1;
size[1] = static_cast<long int >(
ceil((reader->GetMaxY()-reader->GetMinY()+1)/resolution))+1;
ImageType::PointType origin;
origin[0] = reader->GetMinX();
origin[1] = reader->GetMaxY();
ImageType::SpacingType spacing;
spacing[0] = resolution;
spacing[1] = -resolution;
// Software Guide : EndCodeSnippet
// std::cout << "Size: " << size << std::endl;
// std::cout << "Origin: " << origin << std::endl;
// std::cout << "Spacing: " << spacing << std::endl;
// Software Guide : BeginLatex
//
// All these parameters are passed to the interpolation filter:
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
typedef itk::BSplineScatteredDataPointSetToImageFilter
<PointSetType, VectorImageType> FilterType;
FilterType::Pointer filter = FilterType::New();
filter->SetSplineOrder( splineOrder );
FilterType::ArrayType ncps;
ncps.Fill( 6 );
filter->SetNumberOfControlPoints( ncps );
filter->SetNumberOfLevels( level );
filter->SetOrigin( origin );
filter->SetSpacing( spacing );
filter->SetSize( size );
filter->SetInput( reader->GetOutput() );
filter->Update();
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// The result of this filter is an image in which every pixel
// is a vector (with only one element). For now, the otb writer
// does not know how to process that (hopefully soon!). So we
// have to manually copy the element in a standard image:
//
// Software Guide : EndLatex
//TODO ImageOfVectorsToVectorImage needed !
// Write the output to an image.
// Software Guide : BeginCodeSnippet
typedef otb::Image<RealType, 2> RealImageType;
RealImageType::Pointer image = RealImageType::New();
ImageType::RegionType region;
region.SetSize( size );
region.SetIndex( start );
image->SetRegions( region );
image->Allocate();
itk::ImageRegionIteratorWithIndex<RealImageType>
Itt( image, image->GetLargestPossibleRegion() );
for ( Itt.GoToBegin(); !Itt.IsAtEnd(); ++Itt )
{
Itt.Set( filter->GetOutput()->GetPixel( Itt.GetIndex() )[0] );
}
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// Everything is ready so we can just write the image:
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
typedef otb::ImageFileWriter<ImageType> WriterType;
WriterType::Pointer writer = WriterType::New();
writer->SetInput( image );
writer->SetFileName( argv[2] );
writer->Update();
// Software Guide : EndCodeSnippet
typedef otb::Image<unsigned char, 2> UCharImageType;
typedef itk::RescaleIntensityImageFilter<ImageType,
UCharImageType> RescalerType;
RescalerType::Pointer rescaler = RescalerType::New();
rescaler->SetInput(image);
rescaler->SetOutputMaximum(255);
rescaler->SetOutputMinimum(0);
typedef otb::ImageFileWriter<UCharImageType> WriterUCharType;
WriterUCharType::Pointer writerUChar = WriterUCharType::New();
writerUChar->SetInput(rescaler->GetOutput());
writerUChar->SetFileName(argv[3]);
writerUChar->Update();
// Software Guide : BeginLatex
//
// Figure~\ref{fig:LIDARTOIMAGEEXAMPLE} shows the output images with two sets of parameters
//
//
// \begin{figure}
// \center
// \includegraphics[width=0.40\textwidth]{lidar-image-4.png}
// \includegraphics[width=0.40\textwidth]{lidar-image-8.png}
// \itkcaption[Image from lidar data]{Image obtained with 4 level spline interpolation (left) and 8 levels (right)}
// \label{fig:LIDARTOIMAGEEXAMPLE}
// \end{figure}
//
// Software Guide : EndLatex
return EXIT_SUCCESS;
}
......@@ -36,4 +36,3 @@ REGISTER_TEST(DEMToImageGeneratorTest);
#undef main
#define main DEMToImageGeneratorTest
#include "DEMToImageGenerator.cxx"
/*=========================================================================
Program: ORFEO Toolbox
Language: C++
Date: $Date$
Version: $Revision$
Copyright (c) Centre National d'Etudes Spatiales. All rights reserved.
See OTBCopyright.txt for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notices for more information.
=========================================================================*/
// this file defines the otbMultiScaleTest for the test driver
// and all it expects is that you have a function called RegisterTests
#if defined(_MSC_VER)
#pragma warning ( disable : 4786 )
#endif
#include <iostream>
#include "otbTestMain.h"
void RegisterTests()
{
REGISTER_TEST(LidarToImageExampleTest);
}
#undef main
#define main LidarToImageExampleTest
#include "LidarToImageExample.cxx"
# do not do coverage in this directory
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment