Skip to content
Snippets Groups Projects
Commit 684affeb authored by Grégoire Mercier's avatar Grégoire Mercier
Browse files

COMP Compilation issues

parent 2cd8c193
No related branches found
No related tags found
No related merge requests found
......@@ -22,6 +22,7 @@
#include "itkImageToImageFilter.h"
namespace otb {
/** \class AngularProjectionBinaryImageFilter
* \brief Performs \f$ y_i = \cos \theta_i x_1 + \sin \theta_i x_2\f$
......@@ -36,7 +37,7 @@ class ITK_EXPORT AngularProjectionBinaryImageFilter
{
public:
/** Standard typedefs */
typedef WaveletFilterBank Self;
typedef AngularProjectionBinaryImageFilter Self;
typedef itk::ImageToImageFilter<TInputImage, TOutputImage> Superclass;
typedef itk::SmartPointer<Self> Pointer;
typedef itk::SmartPointer<const Self> ConstPointer;
......@@ -45,7 +46,7 @@ public:
itkNewMacro(Self);
/** Creation through object factory macro */
itkTypeMacro(WaveletFilterBank, ImageToImageFilter);
itkTypeMacro(AngularProjectionBinaryImageFilter, ImageToImageFilter);
/** Template parameters typedefs */
typedef TInputImage InputImageType;
......@@ -69,7 +70,10 @@ public:
itkStaticConstMacro(OutputImageDimension, unsigned int, TOutputImage::ImageDimension);
/** Set/Get Angle set */
itkGetMacro(AngleSet,std::vector<PrecisionType>);
std::vector<PrecisionType> GetAngleSet() const
{
return m_AngleSet;
}
void SetAngleSet ( std::vector<PrecisionType> & angle );
void SetInput1 ( const InputImageType * );
......
......@@ -59,7 +59,7 @@ AngularProjectionBinaryImageFilter< TInputImage, TOutputImage, TPrecision >
return 0;
}
return static_cast<const InputImageType * > (this->GetInput(0) );
return static_cast<const TInputImage * > (this->itk::ProcessObject::GetInput(0) );
}
template < class TInputImage, class TOutputImage, class TPrecision >
......@@ -72,7 +72,7 @@ AngularProjectionBinaryImageFilter< TInputImage, TOutputImage, TPrecision >
return 0;
}
return static_cast<const InputImageType * > (this->GetInput(1) );
return static_cast<const TInputImage * > (this->itk::ProcessObject::GetInput(1));
}
template < class TInputImage, class TOutputImage, class TPrecision >
......@@ -106,7 +106,7 @@ template < class TInputImage, class TOutputImage, class TPrecision >
void
AngularProjectionBinaryImageFilter< TInputImage, TOutputImage, TPrecision >
::ThreadedGenerateData
( const OutputImageRegionType & outputRegionForThread, int threadID )
( const OutputImageRegionType & outputRegionForThread, int threadId )
{
itk::ProgressReporter reporter(this, threadId,
outputRegionForThread.GetNumberOfPixels() );
......@@ -135,8 +135,8 @@ AngularProjectionBinaryImageFilter< TInputImage, TOutputImage, TPrecision >
{
for ( unsigned int i = 0; i < outIter.size(); i++ )
{
outIter[i] = vcl_cos( m_AngleSet[i] ) * *iter1
+ vcl_sin( m_AngleSet[i] ) * *iter2;
outIter[i].Set( vcl_cos( m_AngleSet[i] ) * iter1.Get()
+ vcl_sin( m_AngleSet[i] ) * iter2.Get() );
++outIter[i];
}
......@@ -147,5 +147,8 @@ AngularProjectionBinaryImageFilter< TInputImage, TOutputImage, TPrecision >
}
}
} // end of namespace otb
#endif
......@@ -20,7 +20,9 @@
#include "otbMacro.h"
#include "itkImageToImageFilter.h"
#include "otbWaveletTransformImageFilter.h"
#include "otbWaveletOperator.h"
#include "otbWaveletFilterBank.h"
#include "otbWaveletTransform.h"
#include "otbSparseWvltToAngleMapperListFilter.h"
#include "otbImageList.h"
#include "itkListSample.h"
......@@ -69,14 +71,14 @@ public:
typedef TInputImage InputImageType;
typedef TOutputImage OutputImageType;
typedef TPrecision PrecisionType;
typedef TMotherWaveletOperator MotherWaveletOperatorID;
static const Wavelet::Wavelet MotherWaveletOperatorID = TMotherWaveletOperator;
/** Filter types and related */
typedef Image< PrecisionType, InputImageDimension > InternalImageType;
typedef WaveletOperator< MotherWaveletOperatorID, Wavelet::FORWARD, Precision, InputImageDimension > WaveletOperator;
typedef WaveletFilterBank< InputImageType, InternalImageType, WaveletOperator, Wavelet::FORWARD > FilterBankType;
typedef WaveletTransform< InputImageType, InternalImageType, FilterBankType, Wavelet::FORWARD > WvltFilterType;
typedef WaveletOperator< MotherWaveletOperatorID, Wavelet::FORWARD, PrecisionType, InputImageDimension > WaveletOperatorType;
typedef WaveletFilterBank< InternalImageType, InternalImageType, WaveletOperatorType, Wavelet::FORWARD > FilterBankType;
typedef WaveletTransform< InternalImageType, InternalImageType, FilterBankType, Wavelet::FORWARD > WvltFilterType;
typedef typename WvltFilterType::Pointer WvltFilterPointerType;
typedef typename WvltFilterType::OutputImageListType InternalImageListType;
......@@ -88,7 +90,7 @@ public:
typedef typename itk::Statistics::Histogram< PrecisionType > HistogramType;
typedef typename HistogramType::Pointer HistogramPointerType;
typedef AngularProjectionBinaryImageFilter< InternalImageType, OutputImageType, PrecisionType > TransformFilterType;
typedef AngularProjectionBinaryImageFilter< InputImageType, OutputImageType, PrecisionType > TransformFilterType;
typedef typename TransformFilterType::Pointer TransformFilterPointerType;
/** Sets/Gets */
......@@ -141,7 +143,7 @@ protected:
virtual ~SparseUnmixingImageFilter() { }
virtual void GenerateData();
virtual void GenerateNumberOfComponentsRequired () const;
virtual void GenerateNumberOfComponentsRequired ();
private:
SparseUnmixingImageFilter(const Self &); // not implemented
void operator=(const Self &);
......
......@@ -20,7 +20,7 @@
#include "otbSparseUnmixingImageFilter.h"
#include "otbMath.h"
#include "itkProcessAccumulator.h"
#include "itkProcessObject.h"
namespace otb {
......@@ -34,10 +34,10 @@ SparseUnmixingImageFilter< TInputImage, TOutputImage, TPrecision, TMotherWavelet
m_NumberOfComponentsRequired = 0;
m_WvltFilter1 = WvltFilterType::New();
m_WvltFilter1->SetNumberOfDecomposition(2);
m_WvltFilter1->SetNumberOfDecompositions(2);
m_WvltFilter2 = WvltFilterType::New();
m_WvltFilter2->SetNumberOfDecomposition( m_WvltFilter1->GetNumberOfDecomposition() );
m_WvltFilter2->SetNumberOfDecompositions( m_WvltFilter1->GetNumberOfDecompositions() );
m_ListFilter = ListFilterType::New();
m_ListFilter->GetFunctor().SetLowerThreshold(10.);
......@@ -129,12 +129,10 @@ SparseUnmixingImageFilter< TInputImage, TOutputImage, TPrecision, TMotherWavelet
progress->RegisterInternalFilter(m_ListFilter,0.25f);
m_ListFilter->Update();
typename InternalSampleListType::iterator angleIter = m_ListFilter->GetOutputSampleList()->Begin();
typename InternalSampleListType::Iterator angleIter = m_ListFilter->GetOutputSampleList()->Begin();
while ( angleIter != m_ListFilter->GetOutputSampleList()->End() )
{
typename InternalSampleListType::OutputMeasurementVectorType angleMeasure
= angleIter.GetMeasurementVector();
m_Histogram->IncreaseFrequency( angleMeasure[0], 1 );
m_Histogram->IncreaseFrequency( angleIter.GetMeasurementVector()[0], 1 );
++angleIter;
}
......@@ -143,8 +141,8 @@ SparseUnmixingImageFilter< TInputImage, TOutputImage, TPrecision, TMotherWavelet
m_Transformer->SetInput1( this->GetInput1() );
m_Transformer->SetInput2( this->GetInput2() );
m_Transformer->SetAngularProjectors( m_AngleValue );
progress->RegisterInternalFilter(m_Transformer,0.25f);
m_Transformer->SetAngleSet( m_AngleValue );
progress->RegisterInternalFilter(m_Transformer,0.25f);
this->SetNumberOfOutputs( m_Transformer->GetNumberOfOutputs() );
for ( unsigned int i = 0; i < this->GetNumberOfOutputs(); i++ )
......@@ -163,16 +161,25 @@ SparseUnmixingImageFilter< TInputImage, TOutputImage, TPrecision, TMotherWavelet
template < class TInputImage, class TOutputImage,
class TPrecision, Wavelet::Wavelet TMotherWaveletOperator >
unsigned int
void
SparseUnmixingImageFilter< TInputImage, TOutputImage, TPrecision, TMotherWaveletOperator >
::GenerateNumberOfComponentsRequired () const
::GenerateNumberOfComponentsRequired ()
{
std::vector<PrecisionType> angles;
typename HistogramType::Itertor prevHist = m_Histogram->End();
--prevHist;
typename HistogramType::Itertor curHist = m_Histogram->Begin();
typename HistogramType::Itertor nextHist = m_Histogram->Begin();
typename HistogramType::Iterator prevHist = m_Histogram->Begin();
typename HistogramType::Iterator curHist = m_Histogram->Begin();
/** Since operator-- does not exists in itk::Histgram, we have
* to reach the end step by step
*/
++curHist;
while ( curHist != m_Histogram->End() )
{
++curHist;
++prevHist;
}
curHist = m_Histogram->Begin();
typename HistogramType::Iterator nextHist = m_Histogram->Begin();
++nextHist;
if ( prevHist.GetFrequency() < curHist.GetFrequency()
......
......@@ -68,9 +68,9 @@ int otbAngularProjectionBinaryImageFilterTest ( int argc, char * argv[] )
return EXIT_FAILURE;
}
const char * inputImageName1 = parseResult->GetParameterString("--InputImages",0);
const char * inputImageName2 = parseResult->GetParameterString("--InputImages",1);
const char * outputImageName = parseResult->GetParameterString("--OutputImages");
const char * inputImageName1 = parseResult->GetParameterString("--InputImages",0).c_str();
const char * inputImageName2 = parseResult->GetParameterString("--InputImages",1).c_str();
const char * outputImageName = parseResult->GetParameterString("--OutputImages").c_str();
// Main type definition
const unsigned int Dimension = 2;
......@@ -104,14 +104,14 @@ int otbAngularProjectionBinaryImageFilterTest ( int argc, char * argv[] )
filter->Update();
typedef otb::ImageFileWriter< ImageType > ImageWriterType;
typedef otb::ImageFileWriter< ImageType > WriterType;
std::vector< WriterType::Pointer > writers;
writers.resize( filter->GetNumberOfOutputs() );
for ( unsigned int i = 0; i < filter->GetNumberOfOutputs(); i++ )
{
itk::OStringStream title;
title << outputFileName << "_" << i << ".hdr";
title << outputImageName << "_" << i << ".hdr";
writers[i] = WriterType::New();
writers[i]->SetFileName( title.str().c_str() );
......
......@@ -68,9 +68,9 @@ int otbSparseUnmixingImageFilterTest ( int argc, char * argv[] )
return EXIT_FAILURE;
}
const char * inputImageName1 = parseResult->GetParameterString("--InputImages",0);
const char * inputImageName2 = parseResult->GetParameterString("--InputImages",1);
const char * outputImageName = parseResult->GetParameterString("--OutputImages");
const char * inputImageName1 = parseResult->GetParameterString("--InputImages",0).c_str();
const char * inputImageName2 = parseResult->GetParameterString("--InputImages",1).c_str();
const char * outputImageName = parseResult->GetParameterString("--OutputImages").c_str();
const double threshold = parseResult->IsOptionPresent("--Threshold") ?
parseResult->GetParameterDouble("--Threshold") : 10.;
......@@ -100,14 +100,14 @@ int otbSparseUnmixingImageFilterTest ( int argc, char * argv[] )
filter->Update();
typedef otb::ImageFileWriter< ImageType > ImageWriterType;
typedef otb::ImageFileWriter< ImageType > WriterType;
std::vector< WriterType::Pointer > writers;
writers.resize( filter->GetNumberOfOutputs() );
for ( unsigned int i = 0; i < filter->GetNumberOfOutputs(); i++ )
{
itk::OStringStream title;
title << outputFileName << "_" << i << ".hdr";
title << outputImageName << "_" << i << ".hdr";
writers[i] = WriterType::New();
writers[i]->SetFileName( title.str().c_str() );
......
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