CMS 3D CMS Logo

FWVertexProxyBuilder.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Vertexs
4 // Class : FWVertexProxyBuilder
5 //
6 // Implementation:
7 // <Notes on implementation>
8 //
9 // Original Author: Chris Jones
10 // Created: Tue Dec 2 14:17:03 EST 2008
11 //
12 // user include files// user include files
22 
25 
26 #include "TEvePointSet.h"
27 #include "TMatrixDEigen.h"
28 #include "TMatrixDSym.h"
29 #include "TDecompSVD.h"
30 #include "TVectorD.h"
31 #include "TEveTrans.h"
32 #include "TEveTrack.h"
33 #include "TEveTrackPropagator.h"
34 #include "TEveStraightLineSet.h"
35 #include "TEveBoxSet.h"
36 #include "TGeoSphere.h"
37 #include "TEveGeoNode.h"
38 #include "TEveGeoShape.h"
39 #include "TEveVSDStructs.h"
40 
41 class FWVertexProxyBuilder : public FWSimpleProxyBuilderTemplate<reco::Vertex> {
42 public:
44  ~FWVertexProxyBuilder() override {}
45 
46  void setItem(const FWEventItem* iItem) override {
48  if (iItem) {
49  iItem->getConfig()->assertParam("Draw Tracks", false);
50  iItem->getConfig()->assertParam("Draw Pseudo Track", false);
51  iItem->getConfig()->assertParam("Draw Ellipse", false);
52  iItem->getConfig()->assertParam("Scale Ellipse", 2l, 1l, 10l);
53  iItem->getConfig()->assertParam(
54  "Ellipse Color Index", 6l, 0l, (long)context().colorManager()->numberOfLimitedColors());
55  iItem->getConfig()->assertParam("Event Center Index", -1l, -1l, 500l);
56  }
57  }
58 
60 
61  FWVertexProxyBuilder(const FWVertexProxyBuilder&) = delete; // stop default
62  const FWVertexProxyBuilder& operator=(const FWVertexProxyBuilder&) = delete; // stop default
63 
64 private:
66  void build(const reco::Vertex& iData, unsigned int iIndex, TEveElement& oItemHolder, const FWViewContext*) override;
67 
68  void localModelChanges(const FWModelId& iId,
69  TEveElement* iCompound,
70  FWViewType::EType viewType,
71  const FWViewContext* vc) override;
72 };
73 
75  unsigned int iIndex,
76  TEveElement& oItemHolder,
77  const FWViewContext* vc) {
78  const reco::Vertex& v = iData;
79 
80  // marker
81  TEveGeoManagerHolder gmgr(TEveGeoShape::GetGeoMangeur());
82  TEvePointSet* pointSet = new TEvePointSet();
83  pointSet->SetNextPoint(v.x(), v.y(), v.z());
84  setupAddElement(pointSet, &oItemHolder);
85 
86  // ellipse
87  if (item()->getConfig()->value<bool>("Draw Ellipse")) {
88  TEveEllipsoid* eveEllipsoid = new TEveEllipsoid("Ellipsoid", Form("Ellipsoid %d", iIndex));
89 
90  eveEllipsoid->RefPos().Set(v.x(), v.y(), v.z());
91 
92  reco::Vertex::Error e = v.error();
93  TMatrixDSym m(3);
94  for (int i = 0; i < 3; i++)
95  for (int j = 0; j < 3; j++) {
96  m(i, j) = e(i, j);
97  eveEllipsoid->RefEMtx()(i + 1, j + 1) = e(i, j);
98  }
99 
100  // external scaling
101  double ellipseScale = 1.;
102  if (item()->getConfig()->value<long>("Scale Ellipse"))
103  ellipseScale = item()->getConfig()->value<long>("Scale Ellipse");
104 
105  eveEllipsoid->SetScale(ellipseScale);
106 
107  // cache 3D extend used in eval bbox and render 3D
108  TMatrixDEigen eig(m);
109  TVectorD vv(eig.GetEigenValuesRe());
110  eveEllipsoid->RefExtent3D().Set(sqrt(vv(0)) * ellipseScale, sqrt(vv(1)) * ellipseScale, sqrt(vv(2)) * ellipseScale);
111 
112  eveEllipsoid->SetLineWidth(2);
113  setupAddElement(eveEllipsoid, &oItemHolder);
114  eveEllipsoid->SetMainTransparency(TMath::Min(100, 80 + item()->defaultDisplayProperties().transparency() / 5));
115 
116  Color_t color = item()->getConfig()->value<long>("Ellipse Color Index");
117  // eveEllipsoid->SetFillColor(item()->defaultDisplayProperties().color());
118  // eveEllipsoid->SetLineColor(item()->defaultDisplayProperties().color());
119  eveEllipsoid->SetMainColor(color + context().colorManager()->offsetOfLimitedColors());
120  }
121 
122  // tracks
123  if (item()->getConfig()->value<bool>("Draw Tracks")) {
124  for (reco::Vertex::trackRef_iterator it = v.tracks_begin(); it != v.tracks_end(); ++it) {
125  float w = v.trackWeight(*it);
126  if (w < 0.5)
127  continue;
128 
129  const reco::Track& track = *it->get();
130  TEveRecTrack t;
131  t.fBeta = 1.;
132  t.fV = TEveVector(track.vx(), track.vy(), track.vz());
133  t.fP = TEveVector(track.px(), track.py(), track.pz());
134  t.fSign = track.charge();
135  TEveTrack* trk = new TEveTrack(&t, context().getTrackPropagator());
136  trk->SetMainColor(item()->defaultDisplayProperties().color());
137  trk->MakeTrack();
138  setupAddElement(trk, &oItemHolder);
139  }
140  }
141  if (item()->getConfig()->value<bool>("Draw Pseudo Track")) {
142  TEveRecTrack t;
143  t.fBeta = 1.;
144  t.fV = TEveVector(v.x(), v.y(), v.z());
145  t.fP = TEveVector(-v.p4().px(), -v.p4().py(), -v.p4().pz());
146  t.fSign = 1;
147  TEveTrack* trk = new TEveTrack(&t, context().getTrackPropagator());
148  trk->SetLineStyle(7);
149  trk->MakeTrack();
150  setupAddElement(trk, &oItemHolder);
151  }
152 
154  int cIdx = item()->getConfig()->value<long>("Event Center Index");
155  if (cIdx == -1) {
157  }
158  if (cIdx >= 0 && int(iIndex) == int(cIdx)) {
159  fwLog(fwlog::kInfo) << "FWVertexProxyBuilder set event center "
160  << Form("idx [%d], (%f %f %f) \n", cIdx, v.x(), v.y(), v.z());
161  context->commonPrefs()->setEventCenter(v.x(), v.y(), v.z());
162  }
163 }
164 
166  TEveElement* iCompound,
167  FWViewType::EType viewType,
168  const FWViewContext* vc) {
169  increaseComponentTransparency(iId.index(), iCompound, "Ellipsoid", 80);
170  TEveElement* el = iCompound->FindChild("Ellipsoid");
171  if (el)
172  el->SetMainColor(item()->getConfig()->value<long>("Ellipse Color Index") +
173  context().colorManager()->offsetOfLimitedColors());
174 }
175 
176 //
177 // static member functions
178 //
TEveVector & RefPos()
Definition: TEveEllipsoid.h:36
#define REGISTER_PROXYBUILDER_METHODS()
void resetEventCenter()
#define REGISTER_FWPROXYBUILDER(_name_, _type_, _purpose_, _view_)
static const int kAllRPZBits
Definition: FWViewType.h:57
void setupAddElement(TEveElement *el, TEveElement *parent, bool set_color=true) const
CmsShowCommon * commonPrefs() const
Definition: Context.cc:167
T w() const
FWProxyBuilderConfiguration * getConfig() const
Definition: FWEventItem.h:150
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:45
void localModelChanges(const FWModelId &iId, TEveElement *iCompound, FWViewType::EType viewType, const FWViewContext *vc) override
int index() const
Definition: FWModelId.h:41
const fireworks::Context & context() const
void SetScale(float x)
Definition: TEveEllipsoid.h:40
const FWVertexProxyBuilder & operator=(const FWVertexProxyBuilder &)=delete
T sqrt(T t)
Definition: SSEVec.h:23
virtual void setItem(const FWEventItem *iItem)
TEveTrans & RefEMtx()
Definition: TEveEllipsoid.h:38
static Context * getInstance()
Definition: Context.cc:209
FWGenericParameter< T > * assertParam(const std::string &name, T def)
void setEventCenter(float, float, float)
void setItem(const FWEventItem *iItem) override
TEveVector & RefExtent3D()
Definition: TEveEllipsoid.h:37
#define fwLog(_level_)
Definition: fwLog.h:45
static constexpr unsigned int k3DBit
Definition: FWViewType.h:49
const FWEventItem * item() const
void increaseComponentTransparency(unsigned int index, TEveElement *holder, const std::string &name, Char_t transpOffset)
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:38