00001
00002
00003
00004
00005
00006
00007
00008
00009 #include "Fireworks/Core/interface/FWSimpleProxyBuilderTemplate.h"
00010 #include "Fireworks/Core/interface/Context.h"
00011 #include "Fireworks/Core/interface/FWEventItem.h"
00012 #include "Fireworks/Core/interface/FWGeometry.h"
00013 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
00014
00015 #include "Fireworks/Core/interface/FWProxyBuilderConfiguration.h"
00016 #include "Fireworks/Core/interface/FWParameters.h"
00017
00018 #include "TEveTrack.h"
00019
00020 class FWTrackingParticleProxyBuilder : public FWSimpleProxyBuilderTemplate<TrackingParticle>
00021 {
00022 public:
00023 FWTrackingParticleProxyBuilder( void ) {}
00024 virtual ~FWTrackingParticleProxyBuilder( void ) {}
00025
00026 virtual void setItem(const FWEventItem* iItem) {
00027 FWProxyBuilderBase::setItem(iItem);
00028 iItem->getConfig()->assertParam("Point Size", 1l, 3l, 1l);
00029 }
00030
00031 REGISTER_PROXYBUILDER_METHODS();
00032
00033 private:
00034
00035 FWTrackingParticleProxyBuilder( const FWTrackingParticleProxyBuilder& );
00036
00037 const FWTrackingParticleProxyBuilder& operator=( const FWTrackingParticleProxyBuilder& );
00038
00039 void build( const TrackingParticle& iData, unsigned int iIndex, TEveElement& oItemHolder, const FWViewContext* );
00040 };
00041
00042 void
00043 FWTrackingParticleProxyBuilder::build( const TrackingParticle& iData, unsigned int iIndex, TEveElement& oItemHolder, const FWViewContext* )
00044 {
00045 TEveRecTrack t;
00046 t.fBeta = 1.0;
00047 t.fP = TEveVector( iData.px(), iData.py(), iData.pz() );
00048 t.fV = TEveVector( iData.vx(), iData.vy(), iData.vz() );
00049 t.fSign = iData.charge();
00050
00051 TEveTrack* track = new TEveTrack(&t, context().getTrackPropagator());
00052 if( t.fSign == 0 )
00053 track->SetLineStyle( 7 );
00054
00055 TEvePointSet* pointSet = new TEvePointSet;
00056 setupAddElement( pointSet, track );
00057 pointSet->SetMarkerSize(item()->getConfig()->value<long>("Point Size"));
00058 const FWGeometry *geom = item()->getGeom();
00059 const std::vector<PSimHit>& hits = iData.trackPSimHit();
00060
00061 float local[3];
00062 float localDir[3];
00063 float global[3] = { 0.0, 0.0, 0.0 };
00064 float globalDir[3] = { 0.0, 0.0, 0.0 };
00065 std::vector<PSimHit>::const_iterator it = hits.begin();
00066 std::vector<PSimHit>::const_iterator end = hits.end();
00067 if( it != end )
00068 {
00069 unsigned int trackid = hits.begin()->trackId();
00070
00071 for( ; it != end; ++it )
00072 {
00073 const PSimHit& phit = (*it);
00074 if( phit.trackId() != trackid )
00075 {
00076 trackid = phit.trackId();
00077 track->AddPathMark( TEvePathMark( TEvePathMark::kDecay, TEveVector( global[0], global[1], global[2] ),
00078 TEveVector( globalDir[0], globalDir[1], globalDir[2] )));
00079 }
00080 local[0] = phit.localPosition().x();
00081 local[1] = phit.localPosition().y();
00082 local[2] = phit.localPosition().z();
00083 localDir[0] = phit.momentumAtEntry().x();
00084 localDir[1] = phit.momentumAtEntry().y();
00085 localDir[2] = phit.momentumAtEntry().z();
00086 geom->localToGlobal( phit.detUnitId(), local, global );
00087 geom->localToGlobal( phit.detUnitId(), localDir, globalDir );
00088 pointSet->SetNextPoint( global[0], global[1], global[2] );
00089 track->AddPathMark( TEvePathMark( TEvePathMark::kReference, TEveVector( global[0], global[1], global[2] ),
00090 TEveVector( globalDir[0], globalDir[1], globalDir[2] )));
00091 }
00092 if( hits.size() > 1 )
00093 track->AddPathMark( TEvePathMark( TEvePathMark::kDecay, TEveVector( global[0], global[1], global[2] ),
00094 TEveVector( globalDir[0], globalDir[1], globalDir[2] )));
00095 }
00096
00097 track->MakeTrack();
00098 setupAddElement( track, &oItemHolder );
00099 }
00100
00101 REGISTER_FWPROXYBUILDER( FWTrackingParticleProxyBuilder, TrackingParticle, "TrackingParticles", FWViewType::kAll3DBits | FWViewType::kAllRPZBits );