CMS 3D CMS Logo

FWPFCandidateWithHitsProxyBuilder.cc
Go to the documentation of this file.
1 #include "TEveBoxSet.h"
2 #include "TEveTrack.h"
3 #include "TEveTrackPropagator.h"
4 #include "TEveCompound.h"
5 #include "TEveStraightLineSet.h"
6 #include "TEveProjectionBases.h"
7 
15 
20 
22 
23 namespace
24 {
25  const std::string cname("particleFlowRecHitHCALUpgrade");
26 
27 void addLineToLineSet(TEveStraightLineSet* ls, const float* p, int i1, int i2)
28 {
29  i1 *= 3;
30  i2 *= 3;
31  ls->AddLine(p[i1], p[i1+1], p[i1+2], p[i2], p[i2+1], p[i2+2]);
32 }
33 
34 
35 void addBoxAsLines(TEveStraightLineSet* lineset, const float* p)
36 {
37  for (int l = 0; l < 5; l+=4)
38  {
39  addLineToLineSet(lineset, p, 0+l, 1+l);
40  addLineToLineSet(lineset, p, 1+l, 2+l);
41  addLineToLineSet(lineset, p, 2+l, 3+l);
42  addLineToLineSet(lineset, p, 3+l, 0+l);
43  }
44  for (int l = 0; l < 4; ++l)
45  addLineToLineSet(lineset, p, 0+l, 4+l);
46 }
47 
48 void editLineInLineSet(TEveChunkManager::iterator& li, const float* p, int i1, int i2)
49 {
50  TEveStraightLineSet::Line_t& line = * (TEveStraightLineSet::Line_t*) li();
51  i1 *= 3;
52  i2 *= 3;
53  for (int i = 0; i < 3 ; ++i) {
54  line.fV1[0+i] = p[i1+i];
55  line.fV2[0+i] = p[i2+i];
56  }
57 
58  li.next();
59 }
60 
61 void editBoxInLineSet(TEveChunkManager::iterator& li, const float* p)
62 {
63 
64  for (int i = 0; i < 5; i+=4)
65  {
66  editLineInLineSet(li, p, 0+i, 1+i);
67 
68  editLineInLineSet(li, p, 1+i, 2+i);
69  editLineInLineSet(li, p, 2+i, 3+i);
70  editLineInLineSet(li, p, 3+i, 0+i);
71  }
72  for (int i = 0; i < 4; ++i)
73  editLineInLineSet(li, p, 0+i, 4+i);
74 }
75 }
76 
77  //______________________________________________________________________________
78 void
79 FWPFCandidateWithHitsProxyBuilder::build(const FWEventItem* iItem, TEveElementList* product, const FWViewContext* vc)
80 {
81  // init PFCandiate collection
83  iItem->get( candidates );
84  if( candidates == 0 ) return;
85 
86 
87  Int_t idx = 0;
89  for( reco::PFCandidateCollection::const_iterator it = candidates->begin(), itEnd = candidates->end(); it != itEnd; ++it, ++idx)
90  {
91  TEveCompound* comp = createCompound();
92  setupAddElement( comp, product );
93  // printf("products size %d/%d \n", (int)iItem->size(), product->NumChildren());
94 
95  const reco::PFCandidate& cand = *it;
96 
97  // track
98  {
99  TEveRecTrack t;
100  t.fBeta = 1.;
101  t.fP = TEveVector( cand.px(), cand.py(), cand.pz() );
102  t.fV = TEveVector( cand.vertex().x(), cand.vertex().y(), cand.vertex().z() );
103  t.fSign = cand.charge();
104  TEveTrack* trk = new TEveTrack(&t, context().getTrackPropagator() );
105  trk->MakeTrack();
106  fireworks::setTrackTypePF( cand, trk );
107  setupAddElement( trk, comp);
108  }
109  // hits
110  {
111  comp->SetMainColor(iItem->defaultDisplayProperties().color());
112  addHitsForCandidate(cand, comp, vc);
113  }
114 
115  }
116 }
117 
118 //______________________________________________________________________________
120 {
121  // ref hcal collections
123 
124 
125  m_collectionHCAL =0;
126  try
127  {
128  // edm::InputTag tag("particleFlowRecHitHCAL");
129  edm::InputTag tag(cname);
130  item()->getEvent()->getByLabel(tag, handle_hits);
131  if (handle_hits.isValid())
132  {
133  m_collectionHCAL = &*handle_hits;
134  }
135  else
136  {
137  fwLog(fwlog::kError) <<"FWPFCandidateWithHitsProxyBuilder, item " << item()->name() <<": Failed to access collection with name " << cname << "." << std::endl;
138  }
139  }
140  catch (...)
141  {
142  fwLog(fwlog::kError) <<"FWPFCandidateWithHitsProxyBuilder, item " << item()->name() <<": Failed to access collection with name " << cname << "." << std::endl;
143  }
144 }
145 
146 //______________________________________________________________________________
147 void FWPFCandidateWithHitsProxyBuilder::viewContextBoxScale( const float* corners, float scale, bool plotEt, std::vector<float>& scaledCorners, const reco::PFRecHit*)
148 {
149  static TEveVector vtmp;
150  vtmp.Set(0.f, 0.f, 0.f);
151  for( unsigned int i = 0; i < 24; i += 3 )
152  {
153  vtmp[0] += corners[i];
154  vtmp[1] += corners[i + 1];
155  vtmp[2] += corners[i + 2];
156  }
157  vtmp *= 1.f/8.f;
158 
159  if (plotEt)
160  {
161  scale *= vtmp.Perp()/vtmp.Mag();
162  }
163 
164  // Coordinates for a scaled version of the original box
165  for( unsigned int i = 0; i < 24; i += 3 )
166  {
167  scaledCorners[i] = vtmp[0] + ( corners[i] - vtmp[0] ) * scale;
168  scaledCorners[i + 1] = vtmp[1] + ( corners[i + 1] - vtmp[1] ) * scale;
169  scaledCorners[i + 2] = vtmp[2] + ( corners[i + 2] - vtmp[2] ) * scale;
170  }
171 }
172 
173 //______________________________________________________________________________
175 {
176 
177  for (reco::PFRecHitCollection::const_iterator it = m_collectionHCAL->begin(); it != m_collectionHCAL->end(); ++it)
178  {
179 
180  if ( it->detId() == candIdx)
181  {
182  return &(*it);
183  }
184  }
185  return 0;
186 }
187 
188 //______________________________________________________________________________
190 {
191  std::vector<float> scaledCorners(24);
192 
193  float scale = vc->getEnergyScale()->getScaleFactor3D()/50;
194  for (TEveElement::List_i i=parent->BeginChildren(); i!=parent->EndChildren(); ++i)
195  {
196  if ((*i)->NumChildren() > 1)
197  {
198  TEveElement::List_i xx = (*i)->BeginChildren(); ++xx;
199  TEveBoxSet* boxset = dynamic_cast<TEveBoxSet*>(*xx);
200  ++xx;
201  TEveStraightLineSet* lineset = dynamic_cast<TEveStraightLineSet*>(*xx);
202  TEveChunkManager::iterator li(lineset->GetLinePlex());
203  li.next();
204 
205 
206  TEveChunkManager* plex = boxset->GetPlex();
207  if (plex->N())
208  {
209  for (int atomIdx=0; atomIdx < plex->Size(); ++atomIdx)
210  {
211 
212  TEveBoxSet::BFreeBox_t* atom = (TEveBoxSet::BFreeBox_t*)boxset->GetPlex()->Atom(atomIdx);
213  reco::PFRecHit* hit = (reco::PFRecHit*)boxset->GetUserData(atomIdx);
214  const float* corners = item()->getGeom()->getCorners(hit->detId());
215  viewContextBoxScale(corners, hit->energy()*scale, vc->getEnergyScale()->getPlotEt(), scaledCorners, hit);
216  memcpy(atom->fVertices, &scaledCorners[0], sizeof(atom->fVertices));
217 
218  editBoxInLineSet(li, &scaledCorners[0]);
219  }
220 
221  for (TEveProjectable::ProjList_i p = lineset->BeginProjecteds(); p != lineset->EndProjecteds(); ++p)
222  {
223  TEveStraightLineSetProjected* projLineSet = (TEveStraightLineSetProjected*)(*p);
224  projLineSet->UpdateProjection();
225  }
226  }
227  }
228  }
229 }
230 //______________________________________________________________________________
231 namespace {
232 TString boxset_tooltip_callback(TEveDigitSet* ds, Int_t idx)
233 {
234  void* ud = ds->GetUserData(idx);
235  if (ud)
236  {
238  // printf("idx %d %p hit data %p\n", idx, (void*)hit, ud);
239  if (hit)
240  return TString::Format("RecHit %d energy '%f'", idx, hit->energy());
241  else
242  return "ERROR";
243  }
244  return "NULL";
245 }
246 }
247 //______________________________________________________________________________
249 {
251 
252  TEveBoxSet* boxset = 0;
253  TEveStraightLineSet* lineset = 0;
254 
255  for(unsigned elIdx=0; elIdx<eleInBlocks.size(); elIdx++)
256  {
257  // unsigned ieTrack = 0;
258  // unsigned ieECAL = 0;
259  unsigned ieHCAL = 0;
260 
261  reco::PFBlockRef blockRef = eleInBlocks[elIdx].first;
262  unsigned indexInBlock = eleInBlocks[elIdx].second;
263  edm::Ptr<reco::PFBlock> myBlock(blockRef.id(),blockRef.get(), blockRef.key());
264  /*
265  if (myBlock->elements()[indexInBlock].type() == 1)
266  ieTrack = indexInBlock;
267  if (myBlock->elements()[indexInBlock].type() == 4)
268  ieECAL = indexInBlock;
269  */
270  if (myBlock->elements()[indexInBlock].type() == 5)
271  ieHCAL = indexInBlock;
272 
273 
274  std::vector<float> scaledCorners(24);
275  float scale = vc->getEnergyScale()->getScaleFactor3D()/50;
276  if (ieHCAL && m_collectionHCAL) {
277  reco::PFClusterRef hcalclusterRef=myBlock->elements()[ieHCAL].clusterRef();
278  edm::Ptr<reco::PFCluster> myCluster(hcalclusterRef.id(),hcalclusterRef.get(), hcalclusterRef.key());
279  if (myCluster.get())
280  {
281  const std::vector< std::pair<DetId, float> > & hitsandfracs = myCluster->hitsAndFractions();
282 
283  if (!boxset)
284  {
285  boxset = new TEveBoxSet();
286  boxset->Reset(TEveBoxSet::kBT_FreeBox, true, hitsandfracs.size());
287  boxset->SetAntiFlick(true);
288  boxset->SetAlwaysSecSelect(1);
289  boxset->SetPickable(1);
290  boxset->SetTooltipCBFoo(boxset_tooltip_callback);
291  }
292 
293  if (!lineset)
294  {
295  lineset = new TEveStraightLineSet();
296  }
297 
298  bool hitsFound = false;
299  for ( int ihandf=0, lastIdx=(int)(hitsandfracs.size()); ihandf<lastIdx; ihandf++)
300  {
301  unsigned int hitDetId = hitsandfracs[ihandf].first;
302  const float* corners = context().getGeom()->getCorners(hitDetId);
303  const reco::PFRecHit* hit = getHitForDetId(hitDetId);
304  if (hit)
305  {
306  viewContextBoxScale( corners, hit->energy()*scale, vc->getEnergyScale()->getPlotEt(), scaledCorners, hit);
307  boxset->AddBox( &scaledCorners[0]);
308  // setup last box
309  boxset->DigitColor(holder->GetMainColor());
310  boxset->DigitUserData((void*)hit);
311  addBoxAsLines(lineset, &scaledCorners[0]);
312  hitsFound = true;
313  }
314  /*
315  // AMT: don't add lines if hit is not found becuse of unconsistency of scaling.
316  else
317  {
318  addBoxAsLines(lineset, corners);
319  }
320  */
321  }
322  if (!hitsFound)
323  fwLog(fwlog::kWarning) << Form("Can't find matching hits with for HCAL block %d in %s collection. Number of hits %d.\n", elIdx, cname.c_str(), (int)hitsandfracs.size());
324 
325 
326  }
327  else
328  {
329  fwLog(fwlog::kInfo) << "empty cluster \n";
330  }
331  }
332  } // endloop cand.elementsInBlocks();
333 
334 
335  if (boxset) {
336  boxset->RefitPlex();
337  setupAddElement(boxset, holder);
338  }
339 
340  if (lineset) {
341  setupAddElement(lineset, holder);
342  }
343 }
344 
type
Definition: HCALResponse.h:21
const fireworks::Context & context() const
float getScaleFactor3D() const
const FWDisplayProperties & defaultDisplayProperties() const
Definition: FWEventItem.cc:453
#define REGISTER_FWPROXYBUILDER(_name_, _type_, _purpose_, _view_)
static const int kAllRPZBits
Definition: FWViewType.h:58
unsigned detId() const
rechit detId
Definition: PFRecHit.h:108
const std::string & name() const
Definition: FWEventItem.cc:502
const FWGeometry * getGeom() const
Definition: Context.h:83
void setupAddElement(TEveElement *el, TEveElement *parent, bool set_color=true) const
void get(const T *&oData) const
Definition: FWEventItem.h:85
FWViewEnergyScale * getEnergyScale() const
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:195
key_type key() const
Accessor for product key.
Definition: Ref.h:265
std::vector< ElementInBlock > ElementsInBlocks
Definition: PFCandidate.h:386
ProductID id() const
Accessor for product ID.
Definition: Ref.h:259
const FWEventItem * item() const
static const int kAll3DBits
Definition: FWViewType.h:59
Color_t color() const
virtual int charge() const final
electric charge
Definition: LeafCandidate.h:91
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
virtual double px() const final
x coordinate of momentum vector
virtual const Point & vertex() const
vertex position (overwritten by PF...)
Definition: PFCandidate.cc:652
double f[11][100]
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:245
float energy() const
rechit energy
Definition: PFRecHit.h:114
bool isValid() const
Definition: HandleBase.h:74
virtual double pz() const final
z coordinate of momentum vector
const edm::EventBase * getEvent() const
Definition: FWEventItem.h:148
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
def ls(path, rec=False)
Definition: eostools.py:348
const float * getCorners(unsigned int id) const
Definition: FWGeometry.cc:286
#define fwLog(_level_)
Definition: fwLog.h:50
void viewContextBoxScale(const float *corners, float scale, bool plotEt, std::vector< float > &scaledCorners, const reco::PFRecHit *)
const reco::PFRecHitCollection * m_collectionHCAL
const reco::PFRecHit * getHitForDetId(unsigned detId)
void setTrackTypePF(const reco::PFCandidate &pfCand, TAttLine *track)
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:94
TEveCompound * createCompound(bool set_color=true, bool propagate_color_to_all_children=false) const
virtual void scaleProduct(TEveElementList *parent, FWViewType::EType type, const FWViewContext *vc)
void addHitsForCandidate(const reco::PFCandidate &c, TEveElement *holder, const FWViewContext *vc)
const ElementsInBlocks & elementsInBlocks() const
Definition: PFCandidate.cc:687
virtual double py() const final
y coordinate of momentum vector
const FWGeometry * getGeom() const
Definition: FWEventItem.cc:683
bool getPlotEt() const