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

ENH: (Visu Refactoring) Updating classes to use GlComponent and ImageToScreenTransform

parent e8d0d07a
No related branches found
No related tags found
No related merge requests found
......@@ -69,11 +69,11 @@ public:
&& event == FL_PUSH)
{
otbMsgDevMacro(<<"ChangeExtractRegionActionHandler::HandleWidgetEvent(): handling ("<<widgetId<<", "<<event<<")");
// Get the clicked index
typename ViewType::ImageWidgetType::PointType screenPoint, imagePoint;
screenPoint[0] = Fl::event_x();
screenPoint[1] = Fl::event_y();
screenPoint = m_View->GetScrollWidget()->GetMousePosition();
// Transform to image point
imagePoint = m_View->GetScrollWidget()->GetScreenToImageTransform()->TransformPoint(screenPoint);
......
......@@ -71,8 +71,7 @@ public:
otbMsgDevMacro(<<"ChangeScaledExtractRegionActionHandler::HandleWidgetEvent(): handling ("<<widgetId<<", "<<event<<")");
// Get the clicked index
typename ViewType::ImageWidgetType::PointType screenPoint, imagePoint;
screenPoint[0] = Fl::event_x();
screenPoint[1] = Fl::event_y();
screenPoint = m_View->GetFullWidget()->GetMousePosition();
// Transform to image point
imagePoint = m_View->GetFullWidget()->GetScreenToImageTransform()->TransformPoint(screenPoint);
......
......@@ -131,5 +131,16 @@ int GlWidget::handle(int event)
return 0;
}
}
GlWidget::PointType GlWidget::GetMousePosition()
{
// Get the cursor position
PointType index;
index[0] = Fl::event_x();
index[1] = Fl::event_y();
// Flip the y axis
index[1]= this->h()-index[1];
return index;
}
}
#endif
......@@ -20,12 +20,14 @@
// FLTK includes
#include <FL/gl.h>
#include <FL/Fl.h>
#include "FL/Fl_Gl_Window.H"
// This include is needed to get the OTB_GL_USE_ACCEL definition
#include "otbConfigure.h"
#include "itkFixedArray.h"
#include "itkPoint.h"
#include "otbImageWidgetController.h"
......@@ -62,6 +64,9 @@ public:
/** Controller typedef */
typedef otb::ImageWidgetController ControllerType;
typedef ControllerType::Pointer ControllerPointerType;
/** Index typedef */
typedef itk::Point<double,2> PointType;
/** Color typedef (used to draw the rectangle, 4th channel is alpha) */
typedef itk::FixedArray<float,4> ColorType;
......@@ -89,6 +94,11 @@ public:
itkSetMacro(BackgroundColor,ColorType);
itkGetMacro(BackgroundColor,ColorType);
/** Fltk y axis is flipped, therefore we use this function to get
* the cursor position using gl axis */
PointType GetMousePosition();
protected:
/** Constructor */
GlWidget();
......
......@@ -21,6 +21,7 @@
#include "otbImageWidget.h"
#include "otbImageLayerRenderingModelListener.h"
#include "otbImageWidgetController.h"
#include "otbRegionGlComponent.h"
namespace otb
{
......@@ -67,6 +68,10 @@ public:
typedef otb::ImageWidget<ImageType> ImageWidgetType;
typedef typename ImageWidgetType::Pointer ImageWidgetPointerType;
/** Region gl component typedef */
typedef RegionGlComponent RegionGlComponentType;
typedef typename RegionGlComponentType::Pointer RegionGlComponentPointerType;
/**
* This method unregister with previous model if any, and
* register with the new one.
......@@ -123,6 +128,11 @@ private:
/** Controller pointer */
ControllerPointerType m_Controller;
/** Viewed region gl components */
RegionGlComponentPointerType m_ExtractRegionGlComponent;
RegionGlComponentPointerType m_ScaledExtractRegionGlComponent;
}; // end class
} // end namespace otb
......
......@@ -25,13 +25,21 @@ namespace otb
template < class TInputImage >
ImageView<TInputImage>
::ImageView() : m_ScrollWidget(), m_FullWidget(), m_ZoomWidget(),
m_Model(), m_Controller()
m_Model(), m_Controller(), m_ExtractRegionGlComponent(), m_ScaledExtractRegionGlComponent()
{
// Initializing the widgets
m_ScrollWidget = ImageWidgetType::New();
m_FullWidget = ImageWidgetType::New();
m_ZoomWidget = ImageWidgetType::New();
// Extract regions gl components
m_ExtractRegionGlComponent = RegionGlComponentType::New();
m_ScaledExtractRegionGlComponent = RegionGlComponentType::New();
m_ExtractRegionGlComponent->SetVisible(false);
m_ScaledExtractRegionGlComponent->SetVisible(false);
m_ScrollWidget->AddGlComponent(m_ExtractRegionGlComponent);
m_FullWidget->AddGlComponent(m_ScaledExtractRegionGlComponent);
// Set the widget identifiers
m_ScrollWidget->SetIdentifier("Scroll");
m_FullWidget->SetIdentifier("Full");
......@@ -131,12 +139,12 @@ ImageView<TInputImage>
// display the zoom rectangle if necessary
if(m_Model->GetHasExtract())
{
m_ScrollWidget->SetDisplayRectangle(true);
m_ScrollWidget->SetRectangle(m_Model->GetExtractRegion());
m_ExtractRegionGlComponent->SetVisible(true);
m_ExtractRegionGlComponent->SetRegion(m_Model->GetExtractRegion());
}
else
{
m_ScrollWidget->SetDisplayRectangle(false);
m_ExtractRegionGlComponent->SetVisible(false);
}
}
......@@ -165,12 +173,12 @@ ImageView<TInputImage>
// display the zoom rectangle if necessary
if(m_Model->GetHasScaledExtract())
{
m_FullWidget->SetDisplayRectangle(true);
m_FullWidget->SetRectangle(m_Model->GetScaledExtractRegion());
m_ScaledExtractRegionGlComponent->SetVisible(true);
m_ScaledExtractRegionGlComponent->SetRegion(m_Model->GetScaledExtractRegion());
}
else
{
m_FullWidget->SetDisplayRectangle(false);
m_ScaledExtractRegionGlComponent->SetVisible(false);
}
// redraw the widget
m_FullWidget->redraw();
......
......@@ -88,19 +88,6 @@ public:
itkSetMacro(IsotropicZoom,double);
itkGetMacro(IsotropicZoom,double);
/** Enable/disable rectangle drawing */
itkSetMacro(DisplayRectangle,bool);
itkGetMacro(DisplayRectangle,bool);
itkBooleanMacro(DisplayRectangle);
/** Set/Get the rectangle to display */
itkSetMacro(Rectangle,RegionType);
itkGetConstReferenceMacro(Rectangle,RegionType);
/** Set/Get the color of the rectangle */
itkSetMacro(RectangleColor,ColorType);
itkGetConstReferenceMacro(RectangleColor,ColorType);
/** Set/Get the subsampling rate */
itkSetMacro(SubsamplingRate,unsigned int);
itkGetMacro(SubsamplingRate,unsigned int);
......@@ -111,7 +98,7 @@ public:
/** Add a GlComponent */
unsigned int AddGlComponent(const GlComponent * glComponent)
unsigned int AddGlComponent(GlComponent * glComponent)
{
m_GlComponents->PushBack(glComponent);
return m_GlComponents->Size()-1;
......@@ -185,11 +172,6 @@ private:
/** OpenGl buffered region */
RegionType m_OpenGlBufferedRegion;
/** Rectangle region */
RegionType m_Rectangle;
bool m_DisplayRectangle;
ColorType m_RectangleColor;
/** The display extent */
RegionType m_Extent;
......
......@@ -26,15 +26,9 @@ namespace otb
template <class TInputImage>
ImageWidget<TInputImage>
::ImageWidget() : m_IsotropicZoom(1.0), m_OpenGlBuffer(NULL), m_OpenGlBufferedRegion(),
m_Rectangle(),m_DisplayRectangle(false),m_RectangleColor(), m_Extent(),
m_SubsamplingRate(1), m_ImageToScreenTransform(),m_ScreenToImageTransform(),
m_GlComponents()
m_Extent(), m_SubsamplingRate(1), m_ImageToScreenTransform(),
m_ScreenToImageTransform(), m_GlComponents()
{
// Default color for rectangle
m_RectangleColor.Fill(0.);
m_RectangleColor[0]=1.0;
m_RectangleColor[3]=1.0;
// Initialize space to screen transform and inverse
m_ImageToScreenTransform = AffineTransformType::New();
m_ScreenToImageTransform = AffineTransformType::New();
......@@ -150,17 +144,18 @@ ImageWidget<TInputImage>
s2iMatrix.Fill(0);
const double s2iSpacing =(m_SubsamplingRate)/m_IsotropicZoom;
s2iMatrix(0,0)=s2iSpacing;
s2iMatrix(1,1)=s2iSpacing;
s2iMatrix(1,1)=-s2iSpacing;
m_ScreenToImageTransform->SetMatrix(s2iMatrix);
// Image to screen translation
typename AffineTransformType::OutputVectorType translation;
translation[0]= m_SubsamplingRate * (m_OpenGlBufferedRegion.GetIndex()[0]-m_Extent.GetIndex()[0]/m_IsotropicZoom);
translation[1]= m_SubsamplingRate * (m_OpenGlBufferedRegion.GetIndex()[0]-m_Extent.GetIndex()[0]/m_IsotropicZoom);
translation[1]= m_SubsamplingRate * (((m_Extent.GetIndex()[1]+m_Extent.GetSize()[1])/m_IsotropicZoom) + m_OpenGlBufferedRegion.GetIndex()[1]);
m_ScreenToImageTransform->SetTranslation(translation);
// Compute the inverse transform
m_ScreenToImageTransform->GetInverse(m_ImageToScreenTransform);
bool couldInvert = m_ScreenToImageTransform->GetInverse(m_ImageToScreenTransform);
assert(couldInvert && "Could not invert ScreenToImageTransform");
}
template <class TInputImage>
......@@ -181,7 +176,6 @@ ImageWidget<TInputImage>
// Update transforms
this->UpdateTransforms();
if(!this->GetUseGlAcceleration())
{
// Set the pixel Zoom
......@@ -228,31 +222,6 @@ ImageWidget<TInputImage>
it.Get()->Render(m_Extent,m_ImageToScreenTransform);
}
}
// // Draw the rectangle if necessary
// if(m_DisplayRectangle)
// {
// typename RegionType::IndexType index;
// typename RegionType::SizeType size;
// index = m_Rectangle.GetIndex();
// // UL left in image space is LR in opengl space, so we need to get
// // the real upper left
// index[1]+=m_Rectangle.GetSize()[1];
// index = ImageIndexToScreenIndex(index);
// size[0]= static_cast<unsigned int>(static_cast<double>(m_Rectangle.GetSize()[0]/m_SubsamplingRate)*m_IsotropicZoom);
// size[1] = static_cast<unsigned int>(static_cast<double>(m_Rectangle.GetSize()[1]/m_SubsamplingRate)*m_IsotropicZoom);
// glEnable(GL_BLEND);
// glBlendFunc(GL_SRC_ALPHA, GL_ONE_MINUS_SRC_ALPHA);
// glColor4f(m_RectangleColor[0],m_RectangleColor[1],m_RectangleColor[2],m_RectangleColor[3]);
// glBegin(GL_LINE_LOOP);
// gl_rect(index[0],index[1],size[0],size[1]);
// glEnd();
// glDisable(GL_BLEND);
// }
}
}
#endif
......@@ -98,8 +98,7 @@ public:
{
// Get the clicked index
typename ViewType::ImageWidgetType::PointType screenPoint, imagePoint;
screenPoint[0] = Fl::event_x();
screenPoint[1] = Fl::event_y();
screenPoint = sourceWidget->GetMousePosition();
// Transform to image point
imagePoint = sourceWidget->GetScreenToImageTransform()->TransformPoint(screenPoint);
......
/*=========================================================================
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.
=========================================================================*/
#ifndef __otbRegionGlComponent_h
#define __otbRegionGlComponent_h
#include "otbGlComponent.h"
#include "itkImageRegion.h"
namespace otb
{
/** \class RegionGlComponent
* \brief This Gl Component represents a region.
*/
class RegionGlComponent
: public GlComponent
{
public:
/** Standard class typedefs */
typedef RegionGlComponent Self;
typedef GlComponent Superclass;
typedef itk::SmartPointer<Self> Pointer;
typedef itk::SmartPointer<const Self> ConstPointer;
typedef itk::ImageRegion<2> RegionType;
// affine transform
typedef Superclass::AffineTransformType AffineTransformType;
typedef AffineTransformType::InputPointType PointType;
typedef AffineTransformType::InputVectorType VectorType;
typedef Superclass::ColorType ColorType;
/** Runtime information */
itkTypeMacro(RegionGlComponent,GlComponent);
/** New macro */
itkNewMacro(Self);
/// Render the curve according to display extent and axis characteristics
virtual void Render(const RegionType& extent,const AffineTransformType * space2ScreenTransform)
{
PointType ip1, ip2, sp1,sp2;
// Get the input region points
ip1[0]=m_Region.GetIndex()[0]+1;
ip1[1]=m_Region.GetIndex()[1]+1;
ip2[0]=m_Region.GetIndex()[0] + m_Region.GetSize()[0]-1;
ip2[1]=m_Region.GetIndex()[1] + m_Region.GetSize()[1]-1;
// Convert to screen points
sp1 = space2ScreenTransform->TransformPoint(ip1);
sp2 = space2ScreenTransform->TransformPoint(ip2);
// draw the box
glEnable(GL_BLEND);
glBlendFunc(GL_SRC_ALPHA, GL_ONE_MINUS_SRC_ALPHA);
glColor4d(m_Color[0],m_Color[1],m_Color[2],m_Color[3]);
glBegin(GL_LINE_LOOP);
glVertex2d(sp1[0],sp1[1]);
glVertex2d(sp1[0],sp2[1]);
glVertex2d(sp2[0],sp2[1]);
glVertex2d(sp2[0],sp1[1]);
glEnd();
glDisable(GL_BLEND);
}
/** Set/Get the region to render */
itkSetMacro(Region,RegionType);
itkGetConstReferenceMacro(Region,RegionType);
/** Set/Get the color of the region to render */
itkSetMacro(Color,ColorType);
itkGetConstReferenceMacro(Color,ColorType);
protected:
/** Constructor */
RegionGlComponent() : m_Region(), m_Color()
{
m_Color.Fill(0);
m_Color[0] = 1.;
m_Color[3] = 1.;
}
/** Destructor */
virtual ~RegionGlComponent(){}
/** Printself method */
void PrintSelf(std::ostream& os, itk::Indent indent) const
{
Superclass::PrintSelf(os,indent);
}
private:
RegionGlComponent(const Self&); // purposely not implemented
void operator=(const Self&); // purposely not implemented
// The region to render
RegionType m_Region;
// The color of the region to render
ColorType m_Color;
}; // end class
} // end namespace otb
#endif
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