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 "TEveVSDStructs.h"
39 
40 class FWVertexProxyBuilder : public FWSimpleProxyBuilderTemplate<reco::Vertex> {
41 public:
43  ~FWVertexProxyBuilder() override {}
44 
45  void setItem(const FWEventItem* iItem) override {
47  if (iItem) {
48  iItem->getConfig()->assertParam("Draw Tracks", false);
49  iItem->getConfig()->assertParam("Draw Pseudo Track", false);
50  iItem->getConfig()->assertParam("Draw Ellipse", false);
51  iItem->getConfig()->assertParam("Scale Ellipse", 2l, 1l, 10l);
52  iItem->getConfig()->assertParam(
53  "Ellipse Color Index", 6l, 0l, (long)context().colorManager()->numberOfLimitedColors());
54  iItem->getConfig()->assertParam("Event Center Index", -1l, -1l, 500l);
55  }
56  }
57 
59 
60 private:
61  FWVertexProxyBuilder(const FWVertexProxyBuilder&) = delete; // stop default
62  const FWVertexProxyBuilder& operator=(const FWVertexProxyBuilder&) = delete; // stop default
63 
65  void build(const reco::Vertex& iData, unsigned int iIndex, TEveElement& oItemHolder, const FWViewContext*) override;
66 
67  void localModelChanges(const FWModelId& iId,
68  TEveElement* iCompound,
69  FWViewType::EType viewType,
70  const FWViewContext* vc) override;
71 };
72 
74  unsigned int iIndex,
75  TEveElement& oItemHolder,
76  const FWViewContext* vc) {
77  const reco::Vertex& v = iData;
78 
79  // marker
80  TEveGeoManagerHolder gmgr(TEveGeoShape::GetGeoMangeur());
81  TEvePointSet* pointSet = new TEvePointSet();
82  pointSet->SetNextPoint(v.x(), v.y(), v.z());
83  setupAddElement(pointSet, &oItemHolder);
84 
85  // ellipse
86  if (item()->getConfig()->value<bool>("Draw Ellipse")) {
87  TEveEllipsoid* eveEllipsoid = new TEveEllipsoid("Ellipsoid", Form("Ellipsoid %d", iIndex));
88 
89  eveEllipsoid->RefPos().Set(v.x(), v.y(), v.z());
90 
92  TMatrixDSym m(3);
93  for (int i = 0; i < 3; i++)
94  for (int j = 0; j < 3; j++) {
95  m(i, j) = e(i, j);
96  eveEllipsoid->RefEMtx()(i + 1, j + 1) = e(i, j);
97  }
98 
99  // external scaling
100  double ellipseScale = 1.;
101  if (item()->getConfig()->value<long>("Scale Ellipse"))
102  ellipseScale = item()->getConfig()->value<long>("Scale Ellipse");
103 
104  eveEllipsoid->SetScale(ellipseScale);
105 
106  // cache 3D extend used in eval bbox and render 3D
107  TMatrixDEigen eig(m);
108  TVectorD vv(eig.GetEigenValuesRe());
109  eveEllipsoid->RefExtent3D().Set(sqrt(vv(0)) * ellipseScale, sqrt(vv(1)) * ellipseScale, sqrt(vv(2)) * ellipseScale);
110 
111  eveEllipsoid->SetLineWidth(2);
112  setupAddElement(eveEllipsoid, &oItemHolder);
113  eveEllipsoid->SetMainTransparency(TMath::Min(100, 80 + item()->defaultDisplayProperties().transparency() / 5));
114 
115  Color_t color = item()->getConfig()->value<long>("Ellipse Color Index");
116  // eveEllipsoid->SetFillColor(item()->defaultDisplayProperties().color());
117  // eveEllipsoid->SetLineColor(item()->defaultDisplayProperties().color());
118  eveEllipsoid->SetMainColor(color + context().colorManager()->offsetOfLimitedColors());
119  }
120 
121  // tracks
122  if (item()->getConfig()->value<bool>("Draw Tracks")) {
123  for (reco::Vertex::trackRef_iterator it = v.tracks_begin(); it != v.tracks_end(); ++it) {
124  float w = v.trackWeight(*it);
125  if (w < 0.5)
126  continue;
127 
128  const reco::Track& track = *it->get();
129  TEveRecTrack t;
130  t.fBeta = 1.;
131  t.fV = TEveVector(track.vx(), track.vy(), track.vz());
132  t.fP = TEveVector(track.px(), track.py(), track.pz());
133  t.fSign = track.charge();
134  TEveTrack* trk = new TEveTrack(&t, context().getTrackPropagator());
135  trk->SetMainColor(item()->defaultDisplayProperties().color());
136  trk->MakeTrack();
137  setupAddElement(trk, &oItemHolder);
138  }
139  }
140  if (item()->getConfig()->value<bool>("Draw Pseudo Track")) {
141  TEveRecTrack t;
142  t.fBeta = 1.;
143  t.fV = TEveVector(v.x(), v.y(), v.z());
144  t.fP = TEveVector(-v.p4().px(), -v.p4().py(), -v.p4().pz());
145  t.fSign = 1;
146  TEveTrack* trk = new TEveTrack(&t, context().getTrackPropagator());
147  trk->SetLineStyle(7);
148  trk->MakeTrack();
149  setupAddElement(trk, &oItemHolder);
150  }
151 
153  int cIdx = item()->getConfig()->value<long>("Event Center Index");
154  if (cIdx == -1) {
155  context->commonPrefs()->resetEventCenter();
156  }
157  if (cIdx >= 0 && int(iIndex) == int(cIdx)) {
158  fwLog(fwlog::kInfo) << "FWVertexProxyBuilder set event center "
159  << Form("idx [%d], (%f %f %f) \n", cIdx, v.x(), v.y(), v.z());
160  context->commonPrefs()->setEventCenter(v.x(), v.y(), v.z());
161  }
162 }
163 
165  TEveElement* iCompound,
166  FWViewType::EType viewType,
167  const FWViewContext* vc) {
168  increaseComponentTransparency(iId.index(), iCompound, "Ellipsoid", 80);
169  TEveElement* el = iCompound->FindChild("Ellipsoid");
170  if (el)
171  el->SetMainColor(item()->getConfig()->value<long>("Ellipse Color Index") +
172  context().colorManager()->offsetOfLimitedColors());
173 }
174 
175 //
176 // static member functions
177 //
const fireworks::Context & context() const
FWProxyBuilderConfiguration * getConfig() const
Definition: FWEventItem.h:150
TEveVector & RefPos()
Definition: TEveEllipsoid.h:36
void resetEventCenter()
static const int kAllRPZBits
Definition: FWViewType.h:67
const double w
Definition: UKUtility.cc:23
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:73
double y() const
y coordinate
Definition: Vertex.h:117
void setupAddElement(TEveElement *el, TEveElement *parent, bool set_color=true) const
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
void localModelChanges(const FWModelId &iId, TEveElement *iCompound, FWViewType::EType viewType, const FWViewContext *vc) override
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:605
T Min(T a, T b)
Definition: MathUtil.h:39
const FWEventItem * item() const
void SetScale(float x)
Definition: TEveEllipsoid.h:40
const FWVertexProxyBuilder & operator=(const FWVertexProxyBuilder &)=delete
int index() const
Definition: FWModelId.h:41
T sqrt(T t)
Definition: SSEVec.h:19
virtual void setItem(const FWEventItem *iItem)
TEveTrans & RefEMtx()
Definition: TEveEllipsoid.h:38
double z() const
z coordinate
Definition: Vertex.h:119
float trackWeight(const TREF &r) const
returns the weight with which a Track has contributed to the vertex-fit.
Definition: Vertex.h:84
static Context * getInstance()
Definition: Context.cc:193
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:611
FWGenericParameter< T > * assertParam(const std::string &name, T def)
void setEventCenter(float, float, float)
void setItem(const FWEventItem *iItem) override
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:626
double x() const
x coordinate
Definition: Vertex.h:115
TEveVector & RefExtent3D()
Definition: TEveEllipsoid.h:37
#define fwLog(_level_)
Definition: fwLog.h:45
math::XYZTLorentzVectorD p4(float mass=0.13957018, float minWeight=0.5) const
Returns the four momentum of the sum of the tracks, assuming the given mass for the decay products...
Definition: Vertex.cc:118
#define REGISTER_PROXYBUILDER_METHODS()
#define REGISTER_FWPROXYBUILDER(_name_, _type_, _purpose_, _view_)
Error error() const
return SMatrix
Definition: Vertex.h:149
CmsShowCommon * commonPrefs() const
Definition: Context.cc:160
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:623
int charge() const
track electric charge
Definition: TrackBase.h:575
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:71
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:37
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:608
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:620