CMS 3D CMS Logo

FWTrackingParticleProxyBuilderFullFramework.cc
Go to the documentation of this file.
11 
13 
17 
20 
21 #include "TEveTrack.h"
22 #include "TEveCompound.h"
23 #include "TEveManager.h"
24 #include "TEveBrowser.h"
25 #include "TEveTrackPropagator.h"
26 
27 
29 {
30 public:
33 
34  // virtual void setItem(const FWEventItem* iItem) override;
35 
37 
38 private:
41 
43  void build(const FWEventItem* iItem, TEveElementList* product, const FWViewContext*) override;
44 
45  void build( const TrackingParticle& iData, unsigned int iIndex, TEveElement& oItemHolder, const FWViewContext* ) override;
46 
49 
50 
51 };
52 
53 //______________________________________________________________________________
54 
55 /*
56  void FWTrackingParticleProxyBuilderFullFramework::setItem(const FWEventItem* iItem)
57  {
58  printf("set item\n");
59  FWProxyBuilderBase::setItem(iItem);
60  }
61 */
62 //______________________________________________________________________________
63 void FWTrackingParticleProxyBuilderFullFramework::build(const FWEventItem* iItem, TEveElementList* product, const FWViewContext*)
64 {
65  // setup event handles amd call function from parent class
66 
67  const edm::Event* event = (const edm::Event*)item()->getEvent();
68  if (event) {
69 
70  // get collection handle
71  edm::InputTag coltag(item()->moduleLabel(), item()->productInstanceLabel(), item()->processName());
72  event->getByLabel(coltag, tpch);
73 
74  // AMT todo: check if there is any other way getting the list other than this
75  // ifnot, set proces name as a configurable parameter
77  try {
78  event->getByLabel("xxx", simHitsTPAssoc);
79  m_assocList = &*simHitsTPAssoc;
80  }
81  catch (const std::exception& e) {
82  std::cerr << " FWTrackingParticleProxyBuilderFullFramework::setItem() Can't get hits association list " << e.what() << std::endl;
83  }
84  /*
85  // debug propagator
86  gEve->GetBrowser()->MapWindow();
87  gEve->AddToListTree(context().getTrackPropagator(), true);
88  context().getTrackPropagator()->SetRnrReferences(true);
89  */
90  }
91  FWSimpleProxyBuilder::build(iItem, product, nullptr);
92 }
93 //______________________________________________________________________________
94 void
95 FWTrackingParticleProxyBuilderFullFramework::build( const TrackingParticle& iData, unsigned int tpIdx, TEveElement& comp, const FWViewContext* )
96 {
97  TEveRecTrack t;
98  t.fBeta = 1.0;
99  t.fP = TEveVector( iData.px(), iData.py(), iData.pz() );
100  t.fV = TEveVector( iData.vx(), iData.vy(), iData.vz() );
101  t.fSign = iData.charge();
102 
103  TEveTrack* track = new TEveTrack(&t, context().getTrackPropagator());
104  if( t.fSign == 0 )
105  track->SetLineStyle( 7 );
106 
107  track->MakeTrack();
108  setupAddElement( track, &comp );
109  // printf("add track %d \n", tpIdx);
110 
111 
112  if (m_assocList) {
113  TEvePointSet* pointSet = new TEvePointSet;
114  setupAddElement( pointSet, &comp );
115 
116  const FWGeometry *geom = item()->getGeom();
117  float local[3];
118  float localDir[3];
119  float global[3] = { 0.0, 0.0, 0.0 };
120  float globalDir[3] = { 0.0, 0.0, 0.0 };
121 
122  TrackingParticleRef tpr(tpch, tpIdx);
123  std::pair<TrackingParticleRef, TrackPSimHitRef> clusterTPpairWithDummyTP(tpr,TrackPSimHitRef());
124  auto range = std::equal_range(m_assocList->begin(), m_assocList->end(), clusterTPpairWithDummyTP, SimHitTPAssociationProducer::simHitTPAssociationListGreater);
125  // printf("TrackingParticle[%d] P(%.1f, %.1f, %.1f) matches %d hits\n", tpIdx,iData.px(), iData.py(), iData.pz() ,(int)(range.second-range.first ));
126 
127  std::vector<const PSimHit*> phits;
128  for (auto ri = range.first; ri != range.second; ++ri)
129  phits.push_back(ri->second.get());
130 
131  std::sort(phits.begin(), phits.end(), [](const PSimHit* a, const PSimHit* b){ return a->tof() < b->tof(); });
132  for (auto phi = phits.begin(); phi != phits.end(); ++phi)
133  {
134  const PSimHit* phit = *phi;
135 
136  local[0] = phit->localPosition().x();
137  local[1] = phit->localPosition().y();
138  local[2] = phit->localPosition().z();
139 
140  localDir[0] = phit->momentumAtEntry().x();
141  localDir[1] = phit->momentumAtEntry().y();
142  localDir[2] = phit->momentumAtEntry().z();
143 
144  geom->localToGlobal( phit->detUnitId(), local, global );
145  geom->localToGlobal( phit->detUnitId(), localDir, globalDir, false );
146  pointSet->SetNextPoint( global[0], global[1], global[2] );
147 
148  //printf("localP = (%f, %f, %f) globalP = (%f, %f, %f), loss = %f, tof =%f\n", localDir[0], localDir[1], localDir[2],
149  // globalDir[0], globalDir[1], globalDir[2],
150  // phit->energyLoss(), phit->tof());
151  track->AddPathMark( TEvePathMark( TEvePathMark::kReference, TEveVector( global[0], global[1], global[2] ),
152  TEveVector( globalDir[0], globalDir[1], globalDir[2] )));
153 
154  }
155  }
156 }
157 
const fireworks::Context & context() const
#define REGISTER_PROXYBUILDER_METHODS()
float tof() const
deprecated name for timeOfFlight()
Definition: PSimHit.h:76
#define REGISTER_FWPROXYBUILDER(_name_, _type_, _purpose_, _view_)
static const int kAllRPZBits
Definition: FWViewType.h:58
static bool simHitTPAssociationListGreater(SimHitTPPair i, SimHitTPPair j)
LocalVector momentumAtEntry() const
The momentum of the track that produced the hit, at entry point.
Definition: PSimHit.h:55
void setupAddElement(TEveElement *el, TEveElement *parent, bool set_color=true) const
#define nullptr
T y() const
Definition: PV3DBase.h:63
double py() const
y coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
const FWEventItem * item() const
static const int kAll3DBits
Definition: FWViewType.h:59
double pz() const
z coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
float charge() const
Electric charge. Note this is taken from the first SimTrack only.
Local3DPoint localPosition() const
Definition: PSimHit.h:52
double vy() const
y coordinate of parent vertex position
T z() const
Definition: PV3DBase.h:64
void localToGlobal(unsigned int id, const float *local, float *global, bool translatep=true) const
Definition: FWGeometry.cc:478
const edm::EventBase * getEvent() const
Definition: FWEventItem.h:148
double vx() const
x coordinate of parent vertex position
double b
Definition: hdecay.h:120
const FWTrackingParticleProxyBuilderFullFramework & operator=(const FWTrackingParticleProxyBuilderFullFramework &)=delete
std::vector< SimHitTPPair > SimHitTPAssociationList
double a
Definition: hdecay.h:121
edm::Ref< edm::PSimHitContainer > TrackPSimHitRef
Monte Carlo truth information used for tracking validation.
double px() const
x coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
T x() const
Definition: PV3DBase.h:62
const FWGeometry * getGeom() const
Definition: FWEventItem.cc:686
unsigned int detUnitId() const
Definition: PSimHit.h:97
Definition: event.py:1
double vz() const
const SimHitTPAssociationProducer::SimHitTPAssociationList * m_assocList