Skip to content
Snippets Groups Projects
Commit ef50e80b authored by Owen Arnold's avatar Owen Arnold
Browse files

refs #6211. Extract and pass through peak radius.

parent abb5c856
No related branches found
No related tags found
No related merge requests found
......@@ -34,17 +34,15 @@
#include "MantidAPI/IMDWorkspace.h"
#include "MantidAPI/IPeaksWorkspace.h"
#include "MantidVatesAPI/vtkDataSetFactory.h"
#include <vtkCellData.h>
#include <vtkFloatArray.h>
#include <vtkHexahedron.h>
#include <vtkUnstructuredGrid.h>
namespace Mantid
{
namespace VATES
{
// Forward dec.
class ProgressAction;
class DLLExport vtkPeakMarkerFactory: public vtkDataSetFactory
class DLLExport vtkPeakMarkerFactory
{
public:
......@@ -79,6 +77,12 @@ public:
return "vtkPeakMarkerFactory";
}
/// Getter for the peak workspace integration radius
double getIntegrationRadius() const;
/// Was the peaks workspace integrated?
bool isPeaksWorkspaceIntegrated() const;
protected:
virtual void validate() const;
......@@ -98,6 +102,8 @@ private:
/// Which peak dimension to use
ePeakDimensions m_dimensionToShow;
/// peak radius value.
double m_peakRadius;
};
......
......@@ -6,6 +6,11 @@
#include "MantidAPI/IPeak.h"
#include "MantidKernel/V3D.h"
#include <vtkVertex.h>
#include <vtkGlyph3D.h>
#include <vtkSphereSource.h>
#include <vtkUnstructuredGrid.h>
#include <vtkFloatArray.h>
#include <vtkCellData.h>
#include "MantidKernel/ReadLock.h"
using Mantid::API::IPeaksWorkspace;
......@@ -18,7 +23,7 @@ namespace VATES
{
vtkPeakMarkerFactory::vtkPeakMarkerFactory(const std::string& scalarName, ePeakDimensions dimensions) :
m_scalarName(scalarName), m_dimensionToShow(dimensions)
m_scalarName(scalarName), m_dimensionToShow(dimensions), m_peakRadius(-1)
{
}
......@@ -53,6 +58,24 @@ namespace VATES
{
m_workspace = boost::dynamic_pointer_cast<IPeaksWorkspace>(workspace);
validateWsNotNull();
try
{
m_peakRadius = atof(m_workspace->run().getProperty("PeakRadius")->value().c_str());
}
catch(Mantid::Kernel::Exception::NotFoundError&)
{
}
}
double vtkPeakMarkerFactory::getIntegrationRadius() const
{
return m_peakRadius;
}
bool vtkPeakMarkerFactory::isPeaksWorkspaceIntegrated() const
{
return (m_peakRadius > 0);
}
void vtkPeakMarkerFactory::validateWsNotNull() const
......@@ -70,7 +93,7 @@ namespace VATES
/**
Create the vtkStructuredGrid from the provided workspace
@param progressUpdating: Reporting object to pass progress information up the stack.
@return vtkUnStructuredGrid with a bunch of points.
@return vtkPolyData glyph.
*/
vtkDataSet* vtkPeakMarkerFactory::create(ProgressAction& progressUpdating) const
{
......@@ -138,6 +161,7 @@ namespace VATES
points->Squeeze();
visualDataSet->Squeeze();
return visualDataSet;
}
......
......@@ -144,18 +144,71 @@ public:
void testCreateWithoutInitializeThrows()
{
FakeProgressAction progressUpdate;
using namespace Mantid::VATES;
vtkPeakMarkerFactory factory("signal");
TS_ASSERT_THROWS(factory.create(progressUpdate), std::runtime_error);
}
void testTypeName()
{
using namespace Mantid::VATES;
vtkPeakMarkerFactory factory ("signal");
TS_ASSERT_EQUALS("vtkPeakMarkerFactory", factory.getFactoryTypeName());
}
void testGetPeakRadiusDefault()
{
vtkPeakMarkerFactory factory("signal");
TS_ASSERT_EQUALS(-1, factory.getIntegrationRadius());
}
void testIsPeaksWorkspaceIntegratedDefault()
{
vtkPeakMarkerFactory factory("signal");
TS_ASSERT_EQUALS(false, factory.isPeaksWorkspaceIntegrated());
}
void testGetPeakRadiusWhenNotIntegrated()
{
IPeaksWorkspace* mockWorkspace = new MockPeaksWorkspace;
const double expectedRadius =-1; // The default
// Note that no PeaksRadius property has been set.
vtkPeakMarkerFactory factory("signal");
factory.initialize(Mantid::API::IPeaksWorkspace_sptr(mockWorkspace));
TS_ASSERT_EQUALS(expectedRadius, factory.getIntegrationRadius());
}
void testIsPeaksWorkspaceIntegratedWhenNotIntegrated()
{
IPeaksWorkspace* mockWorkspace = new MockPeaksWorkspace;
// Note that no PeaksRadius property has been set.
vtkPeakMarkerFactory factory("signal");
factory.initialize(Mantid::API::IPeaksWorkspace_sptr(mockWorkspace));
TS_ASSERT_EQUALS(false, factory.isPeaksWorkspaceIntegrated()); // false is the default
}
void testGetPeakRadiusWhenIntegrated()
{
IPeaksWorkspace* mockWorkspace = new MockPeaksWorkspace;
const double expectedRadius = 4;
mockWorkspace->mutableRun().addProperty("PeakRadius", expectedRadius, true); // Has a PeaksRadius so must have been processed via IntegratePeaksMD
vtkPeakMarkerFactory factory("signal");
factory.initialize(Mantid::API::IPeaksWorkspace_sptr(mockWorkspace));
TS_ASSERT_EQUALS(expectedRadius, factory.getIntegrationRadius());
}
void testIsPeaksWorkspaceIntegratedWhenIntegrated()
{
IPeaksWorkspace* mockWorkspace = new MockPeaksWorkspace;
const double expectedRadius = 4;
mockWorkspace->mutableRun().addProperty("PeakRadius", expectedRadius, true); // Has a PeaksRadius so must have been processed via IntegratePeaksMD
vtkPeakMarkerFactory factory("signal");
factory.initialize(Mantid::API::IPeaksWorkspace_sptr(mockWorkspace));
TS_ASSERT_EQUALS(true, factory.isPeaksWorkspaceIntegrated());
}
};
......
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