Files
uLib/src/Vtk/HEP/Detectors/vtkDetectorChamber.cxx

129 lines
4.5 KiB
C++

/*//////////////////////////////////////////////////////////////////////////////
// CMT Cosmic Muon Tomography project //////////////////////////////////////////
////////////////////////////////////////////////////////////////////////////////
Copyright (c) 2014, Universita' degli Studi di Padova, INFN sez. di Padova
All rights reserved
Authors: Andrea Rigoni Garola < andrea.rigoni@pd.infn.it >
------------------------------------------------------------------
This library is free software; you can redistribute it and/or
modify it under the terms of the GNU Lesser General Public
License as published by the Free Software Foundation; either
version 3.0 of the License, or (at your option) any later version.
This library is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public
License along with this library.
//////////////////////////////////////////////////////////////////////////////*/
#include <vtkAbstractTransform.h>
#include <vtkAxes.h>
#include <vtkCubeSource.h>
#include <vtkLineSource.h>
#include <vtkMatrix4x4.h>
#include <vtkPolyDataMapper.h>
#include <vtkPropPicker.h>
#include <vtkProperty.h>
#include <vtkRenderWindow.h>
#include <vtkRenderWindowInteractor.h>
#include <vtkRenderer.h>
#include <vtkRendererCollection.h>
#include <vtkSmartPointer.h>
#include <vtkTransform.h>
#include "Vtk/HEP/Detectors/vtkDetectorChamber.h"
#include <vtkBoxWidget.h>
#include <vtkTransformPolyDataFilter.h>
#include <vtkPlaneSource.h>
#include <vtkProperty.h>
#include <vtkNew.h>
namespace uLib {
namespace Vtk {
DetectorChamber::DetectorChamber(uLib::DetectorChamber *content)
: ContainerBox(content) {
m_PlaneSource = vtkPlaneSource::New();
vtkNew<vtkPolyDataMapper> mapper;
mapper->SetInputConnection(m_PlaneSource->GetOutputPort());
m_PlaneActor = vtkActor::New();
m_PlaneActor->SetMapper(mapper);
m_PlaneActor->GetProperty()->SetColor(0.2, 0.8, 0.2); // Light green
m_PlaneActor->GetProperty()->SetOpacity(0.3);
m_PlaneActor->GetProperty()->SetAmbient(1.0);
m_PlaneActor->GetProperty()->SetDiffuse(0.0);
this->SetProp(m_PlaneActor);
this->Update();
}
DetectorChamber::~DetectorChamber() {
m_PlaneSource->Delete();
m_PlaneActor->Delete();
}
DetectorChamber::Content *DetectorChamber::GetContent() const {
return static_cast<Content *>(this->m_model.get());
}
void DetectorChamber::Update() {
this->BaseClass::Update();
if (!this->m_model) return;
Content *c = this->GetContent();
Vector3f size = c->GetSize();
HLine3f plane = c->GetProjectionPlane();
// Normalized local space of the chamber
float Lx = plane.origin.x() / size.x();
float Ly = plane.origin.y() / size.y();
float Lz = plane.origin.z() / size.z();
float Dx = plane.direction.x() * size.x();
float Dy = plane.direction.y() * size.y();
float Dz = plane.direction.z() * size.z();
Vector3f normal(Dx, Dy, Dz);
normal.normalize();
// Find spans
Vector3f v1;
if (std::abs(normal.z()) < 0.9) v1 = Vector3f(0,0,1);
else v1 = Vector3f(1,0,0);
Vector3f p1 = normal.cross(v1).normalized();
Vector3f p2 = normal.cross(p1).normalized();
// Center the visual representation on the box extents (unit box [0,1]^3)
// instead of the plane origin (which might be a corner).
Vector3f boxCenter(0.5f, 0.5f, 0.5f);
Vector3f planePt(Lx, Ly, Lz);
Vector3f visualCenter = boxCenter - (boxCenter - planePt).dot(normal) * normal;
float scale = 1.5; // Slightly larger than the diagonal of a face (1.41)
m_PlaneSource->SetOrigin(visualCenter.x() - 0.5f*scale*p1.x() - 0.5f*scale*p2.x(),
visualCenter.y() - 0.5f*scale*p1.y() - 0.5f*scale*p2.y(),
visualCenter.z() - 0.5f*scale*p1.z() - 0.5f*scale*p2.z());
m_PlaneSource->SetPoint1(visualCenter.x() + 0.5f*scale*p1.x() - 0.5f*scale*p2.x(),
visualCenter.y() + 0.5f*scale*p1.y() - 0.5f*scale*p2.y(),
visualCenter.z() + 0.5f*scale*p1.z() - 0.5f*scale*p2.z());
m_PlaneSource->SetPoint2(visualCenter.x() - 0.5f*scale*p1.x() + 0.5f*scale*p2.x(),
visualCenter.y() - 0.5f*scale*p1.y() + 0.5f*scale*p2.y(),
visualCenter.z() - 0.5f*scale*p1.z() + 0.5f*scale*p2.z());
}
} // namespace Vtk
} // namespace uLib