Skip to content
Snippets Groups Projects
Commit 411d5735 authored by Rashad Kanavath's avatar Rashad Kanavath
Browse files

Merge remote-tracking branch 'origin/develop' into rfc-9-sarcalibration

parents 2624e213 9369ccf6
No related branches found
No related tags found
No related merge requests found
Showing
with 1400 additions and 22 deletions
......@@ -68,7 +68,7 @@ otb_add_test(NAME prTeEstimateRPCSensorModelExampleTest COMMAND ${OTB_TEST_DRIVE
--compare-ascii ${EPSILON_4}
${BASELINE}/otbGCPsToRPCSensorModelImageFilterWithoutDEMOutput.txt
${TEMP}/otbGCPsToRPCSensorModelImageFilterWithoutDEMOutput.txt
--ignore-lines-with 5 PipelineMTime ImportImageContaine Source: Image Time:
--ignore-lines-with 6 PipelineMTime ImportImageContaine Source: Image Time: Pointer:
Execute $<TARGET_FILE:EstimateRPCSensorModelExample>
LARGEINPUT{SPOT4/RIO_DE_JANEIRO/IMAG_01.DAT}
${TEMP}/otbGCPsToRPCSensorModelImageFilterWithoutDEMOutput.txt
......
project(OTBAppSARDecompositions)
otb_module_impl()
set(OTBAppFiltering_LINK_LIBS
${OTBPolarimetry_LIBRARIES}
${OTBImageManipulation_LIBRARIES}
${OTBApplicationEngine_LIBRARIES}
${OTBImageBase_LIBRARIES}
)
otb_create_application(
NAME SARDecompositions
SOURCES otbSARDecompositions.cxx
LINK_LIBRARIES ${${otb-module}_LIBRARIES})
/*=========================================================================
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.
=========================================================================*/
#include "otbWrapperApplication.h"
#include "otbWrapperApplicationFactory.h"
#include "otbReciprocalHAlphaImageFilter.h"
#include "otbSinclairReciprocalImageFilter.h"
#include "otbSinclairToReciprocalCoherencyMatrixFunctor.h"
#include "otbPerBandVectorImageFilter.h"
#include "itkMeanImageFilter.h"
namespace otb
{
namespace Wrapper
{
class SARDecompositions : public Application
{
public:
/** Standard class typedefs. */
typedef SARDecompositions Self;
typedef Application Superclass;
typedef itk::SmartPointer<Self> Pointer;
typedef itk::SmartPointer<const Self> ConstPointer;
typedef otb::Functor::SinclairToReciprocalCoherencyMatrixFunctor<ComplexDoubleImageType::PixelType,
ComplexDoubleImageType::PixelType,
ComplexDoubleImageType::PixelType,
ComplexDoubleVectorImageType::PixelType> FunctorType;
typedef SinclairReciprocalImageFilter<ComplexDoubleImageType,
ComplexDoubleImageType,
ComplexDoubleImageType,
ComplexDoubleVectorImageType,
FunctorType > SRFilterType;
typedef otb::ReciprocalHAlphaImageFilter<ComplexDoubleVectorImageType, DoubleVectorImageType> HAFilterType;
typedef itk::MeanImageFilter<ComplexDoubleImageType, ComplexDoubleImageType> MeanFilterType;
typedef otb::PerBandVectorImageFilter<ComplexDoubleVectorImageType, ComplexDoubleVectorImageType, MeanFilterType> PerBandMeanFilterType;
//FloatImageType
/** Standard macro */
itkNewMacro(Self);
itkTypeMacro(SARDecompositions, otb::Application);
private:
void DoInit()
{
SetName("SARDecompositions");
SetDescription("From one-band complex images (each one related to an element of the Sinclair matrix), returns the selected decomposition.");
// Documentation
SetDocName("SARDecompositions");
SetDocLongDescription("From one-band complex images (HH, HV, VH, VV), returns the selected decomposition.\n \n"
"The H-alpha-A decomposition is currently the only one available; it is implemented for the monostatic case (transmitter and receiver are co-located).\n"
"User must provide three one-band complex images HH, HV or VH, and VV (monostatic case <=> HV = VH).\n"
"The H-alpha-A decomposition consists in averaging 3x3 complex coherency matrices (incoherent analysis); the user must provide the size of the averaging window, thanks to the parameter inco.kernelsize.\n "
"The applications returns a float vector image, made up of three channels : H (entropy), Alpha, A (Anisotropy)." );
SetDocLimitations("None");
SetDocAuthors("OTB-Team");
SetDocSeeAlso("SARPolarMatrixConvert, SARPolarSynth");
AddDocTag(Tags::SAR);
AddParameter(ParameterType_ComplexInputImage, "inhh", "Input Image");
SetParameterDescription("inhh", "Input image (HH)");
AddParameter(ParameterType_ComplexInputImage, "inhv", "Input Image");
SetParameterDescription("inhv", "Input image (HV)");
MandatoryOff("inhv");
AddParameter(ParameterType_ComplexInputImage, "invh", "Input Image");
SetParameterDescription("invh", "Input image (VH)");
MandatoryOff("invh");
AddParameter(ParameterType_ComplexInputImage, "invv", "Input Image");
SetParameterDescription("invv", "Input image (VV)");
AddParameter(ParameterType_OutputImage, "out", "Output Image");
SetParameterDescription("out", "Output image");
AddParameter(ParameterType_Choice, "decomp", "Decompositions");
AddChoice("decomp.haa","H-alpha-A decomposition");
SetParameterDescription("decomp.haa","H-alpha-A decomposition");
AddParameter(ParameterType_Group,"inco","Incoherent decompositions");
SetParameterDescription("inco","This group allows to set parameters related to the incoherent decompositions.");
AddParameter(ParameterType_Int, "inco.kernelsize", "Kernel size for spatial incoherent averaging.");
SetParameterDescription("inco.kernelsize", "Minute (0-59)");
SetMinimumParameterIntValue("inco.kernelsize", 1);
SetDefaultParameterInt("inco.kernelsize", 3);
MandatoryOff("inco.kernelsize");
AddRAMParameter();
// Default values
SetDefaultParameterInt("decomp", 0); // H-alpha-A
// Doc example parameter settings
SetDocExampleParameterValue("inhh", "HH.tif");
SetDocExampleParameterValue("invh", "VH.tif");
SetDocExampleParameterValue("invv", "VV.tif");
SetDocExampleParameterValue("decomp", "haa");
SetDocExampleParameterValue("out", "HaA.tif");
}
void DoUpdateParameters()
{
// Nothing to do here : all parameters are independent
}
void DoExecute()
{
bool inhv = HasUserValue("inhv");
bool invh = HasUserValue("invh");
if ( (!inhv) && (!invh) )
otbAppLogFATAL( << "Parameter inhv or invh not set. Please provide a HV or a VH complex image.");
switch (GetParameterInt("decomp"))
{
case 0: // H-alpha-A
m_SRFilter = SRFilterType::New();
m_HAFilter = HAFilterType::New();
m_MeanFilter = PerBandMeanFilterType::New();
if (inhv)
m_SRFilter->SetInputHV_VH(GetParameterComplexDoubleImage("inhv"));
else if (invh)
m_SRFilter->SetInputHV_VH(GetParameterComplexDoubleImage("invh"));
m_SRFilter->SetInputHH(GetParameterComplexDoubleImage("inhh"));
m_SRFilter->SetInputVV(GetParameterComplexDoubleImage("invv"));
MeanFilterType::InputSizeType radius;
radius.Fill( GetParameterInt("inco.kernelsize") );
m_MeanFilter->GetFilter()->SetRadius(radius);
m_MeanFilter->SetInput(m_SRFilter->GetOutput());
m_HAFilter->SetInput(m_MeanFilter->GetOutput());
SetParameterOutputImage("out", m_HAFilter->GetOutput() );
break;
}
}
//MCPSFilterType::Pointer m_MCPSFilter;
SRFilterType::Pointer m_SRFilter;
HAFilterType::Pointer m_HAFilter;
PerBandMeanFilterType::Pointer m_MeanFilter;
};
} //end namespace Wrapper
} //end namespace otb
OTB_APPLICATION_EXPORT(otb::Wrapper::SARDecompositions)
set(DOCUMENTATION "Basic filters application.")
otb_module(OTBAppSARDecompositions
DEPENDS
OTBPolarimetry
OTBImageManipulation
OTBITK
OTBApplicationEngine
OTBImageBase
TEST_DEPENDS
OTBTestKernel
OTBCommandLine
DESCRIPTION
"${DOCUMENTATION}"
)
otb_module_test()
#----------- SARDecompositions TESTS ----------------
otb_test_application(NAME apTvSARDecompositions
APP SARDecompositions
OPTIONS
-inhh ${INPUTDATA}/RSAT_imageryC_HH.tif
-inhv ${INPUTDATA}/RSAT_imageryC_HV.tif
-invv ${INPUTDATA}/RSAT_imageryC_VV.tif
-decomp haa
-out ${TEMP}/apTvReciprocalHAlpha.tif
)
project(OTBAppSARPolarMatrixConvert)
otb_module_impl()
set(OTBAppFiltering_LINK_LIBS
${OTBPolarimetry_LIBRARIES}
${OTBImageManipulation_LIBRARIES}
${OTBApplicationEngine_LIBRARIES}
${OTBImageBase_LIBRARIES}
)
otb_create_application(
NAME SARPolarMatrixConvert
SOURCES otbSARPolarMatrixConvert.cxx
LINK_LIBRARIES ${${otb-module}_LIBRARIES})
set(DOCUMENTATION "Basic filters application.")
otb_module(OTBAppSARPolarMatrixConvert
DEPENDS
OTBPolarimetry
OTBImageManipulation
OTBITK
OTBApplicationEngine
OTBImageBase
TEST_DEPENDS
OTBTestKernel
OTBCommandLine
DESCRIPTION
"${DOCUMENTATION}"
)
otb_module_test()
#----------- SARPolarMatrixConvert TESTS ----------------
#1
otb_test_application(NAME apTvSARPolarMatrixConvertRecCoherency
APP SARPolarMatrixConvert
OPTIONS
-inhh ${INPUTDATA}/RSAT_imageryC_HH.tif
-inhv ${INPUTDATA}/RSAT_imageryC_HV.tif
-invv ${INPUTDATA}/RSAT_imageryC_VV.tif
-conv msinclairtocoherency
-outc ${TEMP}/apTvSARPolarMatrixConvertRecCoherency.tif
)
#2
otb_test_application(NAME apTvSARPolarMatrixConvertRecCovariance
APP SARPolarMatrixConvert
OPTIONS
-inhh ${INPUTDATA}/RSAT_imageryC_HH.tif
-inhv ${INPUTDATA}/RSAT_imageryC_HV.tif
-invv ${INPUTDATA}/RSAT_imageryC_VV.tif
-conv msinclairtocovariance
-outc ${TEMP}/apTvSARPolarMatrixConvertRecCovariance.tif
)
#3
otb_test_application(NAME apTvSARPolarMatrixConvertRecCirCovariance
APP SARPolarMatrixConvert
OPTIONS
-inhh ${INPUTDATA}/RSAT_imageryC_HH.tif
-inhv ${INPUTDATA}/RSAT_imageryC_HV.tif
-invv ${INPUTDATA}/RSAT_imageryC_VV.tif
-conv msinclairtocircovariance
-outc ${TEMP}/apTvSARPolarMatrixConvertRecCirCovariance.tif
)
#4
otb_test_application(NAME apTvSARPolarMatrixConvertRecCohToMueller
APP SARPolarMatrixConvert
OPTIONS
-inc ${BASELINE}/saTvSinclairImageFilter_SinclairToReciprocalCovariance.tif
-conv mcoherencytomueller
-outf ${TEMP}/apTvSARPolarMatrixConvertRecCohToMueller.tif
)
#5
otb_test_application(NAME apTvSARPolarMatrixConvertRecCovToCohDeg
APP SARPolarMatrixConvert
OPTIONS
-inc ${BASELINE}/saTvSinclairImageFilter_SinclairToReciprocalCovariance.tif
-conv mcovariancetocoherencydegree
-outc ${TEMP}/apTvSARPolarMatrixConvertRecCovToCohDeg.tif
)
#6
otb_test_application(NAME apTvSARPolarMatrixConvertRecCovToRecCoh
APP SARPolarMatrixConvert
OPTIONS
-inc ${BASELINE}/saTvSinclairImageFilter_SinclairToReciprocalCovariance.tif
-conv mcovariancetocoherency
-outc ${TEMP}/apTvSARPolarMatrixConvertRecCovToRecCoh.tif
)
#7
otb_test_application(NAME apTvSARPolarMatrixConvertRecLinCovToRecCirCov
APP SARPolarMatrixConvert
OPTIONS
-inc ${BASELINE}/saTvSinclairImageFilter_SinclairToReciprocalCovariance.tif
-conv mlinearcovariancetocircularcovariance
-outc ${TEMP}/apTvSARPolarMatrixConvertRecCovToRecCoh.tif
)
#8
otb_test_application(NAME apTvSARPolarMatrixConvertMuellerToRecCov
APP SARPolarMatrixConvert
OPTIONS
-inf ${BASELINE}/saTvSinclairImageFilter_SinclairToMueller.tif
-conv muellertomcovariance
-outc ${TEMP}/apTvSARPolarMatrixConvertMuellerToRecCov.tif
)
#9
otb_test_application(NAME apTvSARPolarMatrixConvertBiSincToCoherency
APP SARPolarMatrixConvert
OPTIONS
-inhh ${INPUTDATA}/RSAT_imageryC_HH.tif
-inhv ${INPUTDATA}/RSAT_imageryC_HV.tif
-invh ${INPUTDATA}/RSAT_imageryC_VH.tif
-invv ${INPUTDATA}/RSAT_imageryC_VV.tif
-conv bsinclairtocoherency
-outc ${TEMP}/apTvSARPolarMatrixConvertBiSincToCoherency.tif
)
#10
otb_test_application(NAME apTvSARPolarMatrixConvertBiSincToCovariance
APP SARPolarMatrixConvert
OPTIONS
-inhh ${INPUTDATA}/RSAT_imageryC_HH.tif
-inhv ${INPUTDATA}/RSAT_imageryC_HV.tif
-invh ${INPUTDATA}/RSAT_imageryC_VH.tif
-invv ${INPUTDATA}/RSAT_imageryC_VV.tif
-conv bsinclairtocovariance
-outc ${TEMP}/apTvSARPolarMatrixConvertBiSincToCovariance.tif
)
#11
otb_test_application(NAME apTvSARPolarMatrixConvertBiSincToCirCovariance
APP SARPolarMatrixConvert
OPTIONS
-inhh ${INPUTDATA}/RSAT_imageryC_HH.tif
-inhv ${INPUTDATA}/RSAT_imageryC_HV.tif
-invh ${INPUTDATA}/RSAT_imageryC_VH.tif
-invv ${INPUTDATA}/RSAT_imageryC_VV.tif
-conv bsinclairtocircovariance
-outc ${TEMP}/apTvSARPolarMatrixConvertBiSincToCirCovariance.tif
)
#12
otb_test_application(NAME apTvSARPolarMatrixConvertSincToMueller
APP SARPolarMatrixConvert
OPTIONS
-inhh ${INPUTDATA}/RSAT_imageryC_HH.tif
-inhv ${INPUTDATA}/RSAT_imageryC_HV.tif
-invh ${INPUTDATA}/RSAT_imageryC_HV.tif #monostatic
-invv ${INPUTDATA}/RSAT_imageryC_VV.tif
-conv sinclairtomueller
-outf ${TEMP}/apTvSARPolarMatrixConvertSincToMueller.tif
)
#13
otb_test_application(NAME apTvSARPolarMatrixConvertMuellerToPolDeGPow
APP SARPolarMatrixConvert
OPTIONS
-inf ${BASELINE}/saTvSinclairImageFilter_SinclairToMueller.tif
-conv muellertopoldegandpower
-outf ${TEMP}/apTvSARPolarMatrixConvertMuellerToPolDeGPow.tif
)
project(OTBAppSARPolarSynth)
otb_module_impl()
set(OTBAppFiltering_LINK_LIBS
${OTBPolarimetry_LIBRARIES}
${OTBImageManipulation_LIBRARIES}
${OTBApplicationEngine_LIBRARIES}
${OTBImageBase_LIBRARIES}
)
otb_create_application(
NAME SARPolarSynth
SOURCES otbSARPolarSynth.cxx
LINK_LIBRARIES ${${otb-module}_LIBRARIES})
/*=========================================================================
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.
=========================================================================*/
#include "otbWrapperApplication.h"
#include "otbWrapperApplicationFactory.h"
#include "otbMultiChannelsPolarimetricSynthesisFilter.h"
namespace otb
{
namespace Wrapper
{
class SARPolarSynth : public Application
{
public:
/** Standard class typedefs. */
typedef SARPolarSynth Self;
typedef Application Superclass;
typedef itk::SmartPointer<Self> Pointer;
typedef itk::SmartPointer<const Self> ConstPointer;
typedef MultiChannelsPolarimetricSynthesisFilter<ComplexDoubleVectorImageType, FloatImageType> MCPSFilterType;
/** Standard macro */
itkNewMacro(Self);
itkTypeMacro(SARPolarSynth, otb::Application);
private:
void DoInit()
{
SetName("SARPolarSynth");
SetDescription("Gives, for each pixel, the power that would have been received by a SAR system with a basis different from the classical (H,V) one (polarimetric synthetis).");
// Documentation
SetDocName("SARPolarSynth");
SetDocLongDescription("This application gives, for each pixel, the power that would have been received by a SAR system with a basis different from the classical (H,V) one (polarimetric synthetis).\n"
"The new basis A and B are indicated through two Jones vectors, defined by the user thanks to orientation (psi) and ellipticity (khi) parameters.\n"
"These parameters are namely psii, khii, psir and khir. The suffixes (i) and (r) refer to the transmiting antenna and the receiving antenna respectively.\n"
"Orientations and ellipticities are given in degrees, and are between -90°/90° and -45°/45° respectively.\n "
"\n"
"Four polarization architectures can be processed : \n"
"1) HH_HV_VH_VV : full polarization, general bistatic case.\n"
"2) HH_HV_VV or HH_VH_VV : full polarization, monostatic case (transmitter and receiver are co-located).\n"
"3) HH_HV : dual polarization.\n"
"4) VH_VV : dual polarization.\n"
"The application takes a complex vector image as input, where each band correspond to a particular emission/reception polarization scheme.\n"
"User must comply with the band order given above, since the bands are used to build the Sinclair matrix.\n"
"\n"
"In order to determine the architecture, the application first relies on the number of bands of the input image.\n"
"1) Architecture HH_HV_VH_VV is the only one with four bands, there is no possible confusion.\n"
"2) Concerning HH_HV_VV and HH_VH_VV architectures, both correspond to a three channels image. But they are processed in the same way, as the Sinclair matrix is symetric in the monostatic case.\n"
"3) Finally, the two last architectures (dual polarizations), can't be distinguished only by the number of bands of the input image.\n"
" User must then use the parameters emissionh and emissionv to indicate the architecture of the system : emissionh=1 and emissionv=0 --> HH_HV, emissionh=0 and emissionv=1 --> VH_VV.\n"
"Note : if the architecture is HH_HV, khii and psii are automatically set to 0°/0°; if the architecture is VH_VV, khii and psii are automatically set to 0°/90°.\n"
"\n"
"It is also possible to force the calculation to co-polar or cross-polar modes.\n"
"In the co-polar case, values for psir and khir will be ignored and forced to psii and khii; same as the cross-polar mode, where khir and psir will be forced to psii+90° and -khii.\n"
"\n"
"Finally, the result of the polarimetric synthetis is expressed in the power domain, through a one-band scalar image.\n"
"Note: this application doesn't take into account the terms which do not depend on the polarization of the antennas. \n"
"The parameter gain can be used for this purpose.\n"
"\n"
"The final formula is thus : P = | B^T . [S] . A |², where A ans B are two Jones vectors and S is a Sinclair matrix.");
SetDocLimitations("None");
SetDocAuthors("OTB-Team");
SetDocSeeAlso("SARDecompositions, SARPolarMatrixConvert");
AddDocTag(Tags::SAR);
AddParameter(ParameterType_ComplexInputImage, "in", "Input Image");
SetParameterDescription("in", "Input image.");
AddParameter(ParameterType_OutputImage, "out", "Output Image");
SetParameterDescription("out", "Output image.");
AddParameter(ParameterType_Float,"psii","psii");
SetParameterDescription("psii","Orientation (transmitting antenna)");
SetMinimumParameterFloatValue("psii",-90.0);
SetMaximumParameterFloatValue("psii",90.0);
AddParameter(ParameterType_Float,"khii","khii");
SetParameterDescription("khii","Ellipticity (transmitting antenna)");
SetMinimumParameterFloatValue("khii",-45.0);
SetMaximumParameterFloatValue("khii",45.0);
AddParameter(ParameterType_Float,"psir","psir");
SetParameterDescription("psir","Orientation (receiving antenna)");
SetMinimumParameterFloatValue("psir",-90.0);
SetMaximumParameterFloatValue("psir",90.0);
AddParameter(ParameterType_Float,"khir","khir");
SetParameterDescription("khir","Ellipticity (receiving antenna)");
SetMinimumParameterFloatValue("khir",-45.0);
SetMaximumParameterFloatValue("khir",45.0);
AddParameter(ParameterType_Int,"emissionh","Emission H");
SetParameterDescription("emissionh","This parameter is useful in determining the polarization architecture (dual polarization case).");
SetMinimumParameterIntValue("emissionh",0);
SetMaximumParameterIntValue("emissionh",1);
MandatoryOff("emissionh");
AddParameter(ParameterType_Int,"emissionv","Emission V");
SetParameterDescription("emissionv","This parameter is useful in determining the polarization architecture (dual polarization case).");
SetMinimumParameterIntValue("emissionv",0);
SetMaximumParameterIntValue("emissionv",1);
MandatoryOff("emissionv");
AddParameter(ParameterType_Choice, "mode", "Forced mode");
AddChoice("mode.none","None");
SetParameterDescription("mode.none","None");
AddChoice("mode.co","Copolarization");
SetParameterDescription("mode.none","Copolarization");
AddChoice("mode.cross","Crosspolarization");
SetParameterDescription("mode.cross","Crosspolarization");
AddRAMParameter();
// Default values
SetDefaultParameterFloat("psii", 0.);
SetDefaultParameterFloat("khii", 0.);
SetDefaultParameterFloat("psir", 0.);
SetDefaultParameterFloat("khir", 0.);
SetDefaultParameterInt("emissionh", 0);
SetDefaultParameterInt("emissionv", 0);
SetDefaultParameterFloat("mode", 0);
// Doc example parameter settings
SetDocExampleParameterValue("in", "sar.tif");
SetDocExampleParameterValue("psii","15.");
SetDocExampleParameterValue("khii", "5.");
SetDocExampleParameterValue("psir","-25.");
SetDocExampleParameterValue("khir", "10.");
SetDocExampleParameterValue("out", "newbasis.tif");
}
void DoUpdateParameters()
{
// Nothing to do here : all parameters are independent
}
void DoExecute()
{
m_MCPSFilter = MCPSFilterType::New();
m_MCPSFilter->SetPsiI(GetParameterFloat("psii"));
m_MCPSFilter->SetKhiI(GetParameterFloat("khii"));
m_MCPSFilter->SetPsiR(GetParameterFloat("psir"));
m_MCPSFilter->SetKhiR(GetParameterFloat("khir"));
m_MCPSFilter->SetEmissionH(GetParameterInt("emissionh"));
m_MCPSFilter->SetEmissionV(GetParameterInt("emissionv"));
m_MCPSFilter->SetMode(GetParameterInt("mode"));
ComplexDoubleVectorImageType* inVImage = GetParameterComplexDoubleVectorImage("in");
inVImage->UpdateOutputInformation();
int nbBands = inVImage->GetNumberOfComponentsPerPixel();
otbAppLogINFO( << "nbBands = " << nbBands);
m_MCPSFilter->SetInput(inVImage);
SetParameterOutputImage("out", m_MCPSFilter->GetOutput());
}
//std::vector<itk::ProcessObject::Pointer> m_Ref;
MCPSFilterType::Pointer m_MCPSFilter;
};
} //end namespace Wrapper
} //end namespace otb
OTB_APPLICATION_EXPORT(otb::Wrapper::SARPolarSynth)
set(DOCUMENTATION "Basic filters application.")
otb_module(OTBAppSARPolarSynth
DEPENDS
OTBPolarimetry
OTBImageManipulation
OTBITK
OTBApplicationEngine
OTBImageBase
TEST_DEPENDS
OTBTestKernel
OTBCommandLine
DESCRIPTION
"${DOCUMENTATION}"
)
otb_module_test()
#----------- PolarSynth TESTS ----------------
otb_test_application(NAME apTvSARPolarSynth
APP SARPolarSynth
OPTIONS -in ${INPUTDATA}/RSAT2_AltonaExtract_1000_1000_100_100.hdr
-out ${TEMP}/resApMultiPolarimetricSynthesis1.tif
-psii 10.0
-khii 0.0
-psir 0.0
-khir 0.0
)
......@@ -59,7 +59,7 @@ public:
typedef double RealType;
typedef PointSetType::PointType PointType;
typedef SarCalibrationLookupData LookupDataType;
typedef typename LookupDataType::Pointer LookupDataPointerType;
typedef LookupDataType::Pointer LookupDataPointerType;
virtual void CreateCalibrationLookupData(const short t);
......
......@@ -137,6 +137,7 @@ public:
plog[k]=0.0;
else
plog[k]=-p[k]*log(p[k])/log(3.0);
}
entropy = 0.0;
......
......@@ -18,8 +18,10 @@
#ifndef __otbSinclairToReciprocalCoherencyMatrixFunctor_h
#define __otbSinclairToReciprocalCoherencyMatrixFunctor_h
#include "vcl_complex.h"
#include "itkMacro.h"
#include "vcl_complex.h"
#include "otbMath.h"
#include "vnl/vnl_matrix.h"
namespace otb
{
......@@ -62,6 +64,7 @@ class SinclairToReciprocalCoherencyMatrixFunctor
public:
/** Some typedefs. */
typedef typename std::complex <double> ComplexType;
typedef vnl_matrix<ComplexType> VNLMatrixType;
typedef typename TOutput::ValueType OutputValueType;
itkStaticConstMacro(NumberOfComponentsPerPixel, unsigned int, 6);
......@@ -75,19 +78,23 @@ public:
const ComplexType S_hh = static_cast<ComplexType>(Shh);
const ComplexType S_hv = static_cast<ComplexType>(Shv);
const ComplexType S_vv = static_cast<ComplexType>(Svv);
const ComplexType HHPlusVV = S_hh + S_vv;
const ComplexType HHMinusVV = S_hh - S_vv;
const ComplexType twoHV = ComplexType( 2.0 ) * S_hv;
result[0] = static_cast<OutputValueType>( std::norm(HHPlusVV) );
result[1] = static_cast<OutputValueType>( HHPlusVV * vcl_conj(HHMinusVV) );
result[2] = static_cast<OutputValueType>( HHPlusVV * vcl_conj(twoHV) );
result[3] = static_cast<OutputValueType>( std::norm(HHMinusVV) );
result[4] = static_cast<OutputValueType>( HHMinusVV *vcl_conj(twoHV) );
result[5] = static_cast<OutputValueType>( std::norm(twoHV) );
result /= 2.0;
VNLMatrixType f3p(3, 1, 0.);
f3p[0][0]= (S_hh + S_vv) / ComplexType( std::sqrt(2.0) , 0.0);
f3p[1][0]= (S_hh - S_vv) / ComplexType( std::sqrt(2.0) , 0.0);
f3p[2][0]= ComplexType( std::sqrt(2.0) , 0.0) * S_hv;
VNLMatrixType res = f3p*f3p.conjugate_transpose();
result[0] = static_cast<OutputValueType>( res[0][0] );
result[1] = static_cast<OutputValueType>( res[0][1] );
result[2] = static_cast<OutputValueType>( res[0][2] );
result[3] = static_cast<OutputValueType>( res[1][1] );
result[4] = static_cast<OutputValueType>( res[1][2] );
result[5] = static_cast<OutputValueType>( res[2][2] );
return (result);
}
......
......@@ -19,6 +19,8 @@
#define __otbSinclairToReciprocalCovarianceMatrixFunctor_h
#include "vcl_complex.h"
#include "otbMath.h"
#include "vnl/vnl_matrix.h"
namespace otb
{
......@@ -61,6 +63,7 @@ class SinclairToReciprocalCovarianceMatrixFunctor
public:
/** Some typedefs. */
typedef typename std::complex <double> ComplexType;
typedef vnl_matrix<ComplexType> VNLMatrixType;
typedef typename TOutput::ValueType OutputValueType;
inline TOutput operator ()(const TInput1& Shh, const TInput2& Shv, const TInput3& Svv)
{
......@@ -72,12 +75,19 @@ public:
const ComplexType S_hv = static_cast<ComplexType>(Shv);
const ComplexType S_vv = static_cast<ComplexType>(Svv);
result[0] = static_cast<OutputValueType>( std::norm( S_hh ) );
result[1] = static_cast<OutputValueType>( vcl_sqrt(2.0)*S_hh*vcl_conj(S_hv) );
result[2] = static_cast<OutputValueType>( S_hh*vcl_conj(S_vv) );
result[3] = static_cast<OutputValueType>( 2.0*std::norm( S_hv ) );
result[4] = static_cast<OutputValueType>( vcl_sqrt(2.0)*S_hv*vcl_conj(S_vv) );
result[5] = static_cast<OutputValueType>( std::norm( S_vv ) );
VNLMatrixType f3l(3, 1, 0.);
f3l[0][0]=S_hh;
f3l[1][0]=ComplexType(std::sqrt(2.0),0.0)*S_hv;
f3l[2][0]=S_vv;
VNLMatrixType res = f3l*f3l.conjugate_transpose();
result[0] = static_cast<OutputValueType>( res[0][0] );
result[1] = static_cast<OutputValueType>( res[0][1] );
result[2] = static_cast<OutputValueType>( res[0][2] );
result[3] = static_cast<OutputValueType>( res[1][1] );
result[4] = static_cast<OutputValueType>( res[1][2] );
result[5] = static_cast<OutputValueType>( res[2][2] );
return (result);
}
......
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