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  void build(const FWEventItem* iItem, TEveElementList* product, const FWViewContext*) override;
42 
43  void build( const TrackingParticle& iData, unsigned int iIndex, TEveElement& oItemHolder, const FWViewContext* ) override;
44 
47 
48 
49 };
50 
51 //______________________________________________________________________________
52 
53 /*
54  void FWTrackingParticleProxyBuilderFullFramework::setItem(const FWEventItem* iItem)
55  {
56  printf("set item\n");
57  FWProxyBuilderBase::setItem(iItem);
58  }
59 */
60 //______________________________________________________________________________
61 void FWTrackingParticleProxyBuilderFullFramework::build(const FWEventItem* iItem, TEveElementList* product, const FWViewContext*)
62 {
63  // setup event handles amd call function from parent class
64 
65  const edm::Event* event = (const edm::Event*)item()->getEvent();
66  if (event) {
67 
68  // get collection handle
69  edm::InputTag coltag(item()->moduleLabel(), item()->productInstanceLabel(), item()->processName());
70  event->getByLabel(coltag, tpch);
71 
72  // AMT todo: check if there is any other way getting the list other than this
73  // ifnot, set proces name as a configurable parameter
75  try {
76  event->getByLabel("xxx", simHitsTPAssoc);
77  m_assocList = &*simHitsTPAssoc;
78  }
79  catch (const std::exception& e) {
80  std::cerr << " FWTrackingParticleProxyBuilderFullFramework::setItem() Can't get hits association list " << e.what() << std::endl;
81  }
82  /*
83  // debug propagator
84  gEve->GetBrowser()->MapWindow();
85  gEve->AddToListTree(context().getTrackPropagator(), true);
86  context().getTrackPropagator()->SetRnrReferences(true);
87  */
88  }
89  FWSimpleProxyBuilder::build(iItem, product, nullptr);
90 }
91 //______________________________________________________________________________
92 void
93 FWTrackingParticleProxyBuilderFullFramework::build( const TrackingParticle& iData, unsigned int tpIdx, TEveElement& comp, const FWViewContext* )
94 {
95  TEveRecTrack t;
96  t.fBeta = 1.0;
97  t.fP = TEveVector( iData.px(), iData.py(), iData.pz() );
98  t.fV = TEveVector( iData.vx(), iData.vy(), iData.vz() );
99  t.fSign = iData.charge();
100 
101  TEveTrack* track = new TEveTrack(&t, context().getTrackPropagator());
102  if( t.fSign == 0 )
103  track->SetLineStyle( 7 );
104 
105  track->MakeTrack();
106  setupAddElement( track, &comp );
107  // printf("add track %d \n", tpIdx);
108 
109 
110  if (m_assocList) {
111  TEvePointSet* pointSet = new TEvePointSet;
112  setupAddElement( pointSet, &comp );
113 
114  const FWGeometry *geom = item()->getGeom();
115  float local[3];
116  float localDir[3];
117  float global[3] = { 0.0, 0.0, 0.0 };
118  float globalDir[3] = { 0.0, 0.0, 0.0 };
119 
120  TrackingParticleRef tpr(tpch, tpIdx);
121  std::pair<TrackingParticleRef, TrackPSimHitRef> clusterTPpairWithDummyTP(tpr,TrackPSimHitRef());
122  auto range = std::equal_range(m_assocList->begin(), m_assocList->end(), clusterTPpairWithDummyTP, SimHitTPAssociationProducer::simHitTPAssociationListGreater);
123  // printf("TrackingParticle[%d] P(%.1f, %.1f, %.1f) matches %d hits\n", tpIdx,iData.px(), iData.py(), iData.pz() ,(int)(range.second-range.first ));
124 
125  std::vector<const PSimHit*> phits;
126  for (auto ri = range.first; ri != range.second; ++ri)
127  phits.push_back(ri->second.get());
128 
129  std::sort(phits.begin(), phits.end(), [](const PSimHit* a, const PSimHit* b){ return a->tof() < b->tof(); });
130  for (auto phi = phits.begin(); phi != phits.end(); ++phi)
131  {
132  const PSimHit* phit = *phi;
133 
134  local[0] = phit->localPosition().x();
135  local[1] = phit->localPosition().y();
136  local[2] = phit->localPosition().z();
137 
138  localDir[0] = phit->momentumAtEntry().x();
139  localDir[1] = phit->momentumAtEntry().y();
140  localDir[2] = phit->momentumAtEntry().z();
141 
142  geom->localToGlobal( phit->detUnitId(), local, global );
143  geom->localToGlobal( phit->detUnitId(), localDir, globalDir, false );
144  pointSet->SetNextPoint( global[0], global[1], global[2] );
145 
146  //printf("localP = (%f, %f, %f) globalP = (%f, %f, %f), loss = %f, tof =%f\n", localDir[0], localDir[1], localDir[2],
147  // globalDir[0], globalDir[1], globalDir[2],
148  // phit->energyLoss(), phit->tof());
149  track->AddPathMark( TEvePathMark( TEvePathMark::kReference, TEveVector( global[0], global[1], global[2] ),
150  TEveVector( globalDir[0], globalDir[1], globalDir[2] )));
151 
152  }
153  }
154 }
155 
const fireworks::Context & context() const
#define REGISTER_PROXYBUILDER_METHODS()
float tof() const
deprecated name for timeOfFlight()
Definition: PSimHit.h:72
#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:47
void setupAddElement(TEveElement *el, TEveElement *parent, bool set_color=true) const
T y() const
Definition: PV3DBase.h:63
double py() const
y coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
#define nullptr
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:44
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:683
unsigned int detUnitId() const
Definition: PSimHit.h:93
Definition: event.py:1
double vz() const
const SimHitTPAssociationProducer::SimHitTPAssociationList * m_assocList