CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 // $Id: FWVertexProxyBuilder.cc,v 1.16 2012/03/13 04:22:06 amraktad Exp $
12 //
13 // user include files// user include files
22 
24 
25 #include "TEvePointSet.h"
26 #include "TMatrixDEigen.h"
27 #include "TMatrixDSym.h"
28 #include "TDecompSVD.h"
29 #include "TVectorD.h"
30 #include "TEveTrans.h"
31 #include "TEveTrack.h"
32 #include "TEveTrackPropagator.h"
33 #include "TEveStraightLineSet.h"
34 #include "TEveBoxSet.h"
35 #include "TGeoSphere.h"
36 #include "TEveGeoNode.h"
37 #include "TEveVSDStructs.h"
38 
40 {
41 public:
43  virtual ~FWVertexProxyBuilder() {}
44 
45  virtual void setItem(const FWEventItem* iItem)
46  {
48  if (iItem)
49  {
50  iItem->getConfig()->assertParam("Draw Tracks", false);
51  iItem->getConfig()->assertParam("Draw Pseudo Track", false);
52  iItem->getConfig()->assertParam("Draw Ellipse", false);
53  iItem->getConfig()->assertParam("Scale Ellipse",2l, 1l, 10l);
54  iItem->getConfig()->assertParam("Ellipse Color Index", 6l, 0l, (long)context().colorManager()->numberOfLimitedColors());
55  }
56  }
57 
59 
60 private:
61  FWVertexProxyBuilder(const FWVertexProxyBuilder&); // stop default
62  const FWVertexProxyBuilder& operator=(const FWVertexProxyBuilder&); // stop default
63 
64  virtual void build(const reco::Vertex& iData, unsigned int iIndex,TEveElement& oItemHolder, const FWViewContext*);
65 
66  virtual void localModelChanges(const FWModelId& iId, TEveElement* iCompound,
67  FWViewType::EType viewType, const FWViewContext* vc);
68 
69 };
70 
71 
72 void
73 FWVertexProxyBuilder::build(const reco::Vertex& iData, unsigned int iIndex, TEveElement& oItemHolder , const FWViewContext*)
74 {
75  const reco::Vertex & v = iData;
76 
77  // marker
78  TEveGeoManagerHolder gmgr(TEveGeoShape::GetGeoMangeur());
79  TEvePointSet* pointSet = new TEvePointSet();
80  pointSet->SetNextPoint( v.x(), v.y(), v.z() );
81  setupAddElement(pointSet, &oItemHolder);
82 
83 
84  // ellipse
85  if ( item()->getConfig()->value<bool>("Draw Ellipse"))
86  {
87 
88  TEveEllipsoid* eveEllipsoid = new TEveEllipsoid("Ellipsoid", Form("Ellipsoid %d", iIndex));
89 
90  eveEllipsoid->RefPos().Set(v.x(),v.y(),v.z());
91 
93  TMatrixDSym m(3);
94  for(int i=0;i<3;i++)
95  for(int j=0;j<3;j++)
96  {
97  m(i,j) = e(i,j);
98  eveEllipsoid->RefEMtx()(i+1, j+1) = e(i,j);
99  }
100 
101  // external scaling
102  double ellipseScale = 1.;
103  if ( item()->getConfig()->value<long>("Scale Ellipse"))
104  ellipseScale = item()->getConfig()->value<long>("Scale Ellipse");
105 
106  eveEllipsoid->SetScale(ellipseScale);
107 
108  // cache 3D extend used in eval bbox and render 3D
109  TMatrixDEigen eig(m);
110  TVectorD vv ( eig.GetEigenValuesRe());
111  eveEllipsoid->RefExtent3D().Set(sqrt(vv(0))*ellipseScale,sqrt(vv(1))*ellipseScale,sqrt(vv(2))*ellipseScale);
112 
113  eveEllipsoid->SetLineWidth(2);
114  setupAddElement(eveEllipsoid, &oItemHolder);
115  eveEllipsoid->SetMainTransparency(TMath::Min(100, 80 + item()->defaultDisplayProperties().transparency() / 5));
116 
117 
118 
119  Color_t color = item()->getConfig()->value<long>("Ellipse Color Index");
120  // eveEllipsoid->SetFillColor(item()->defaultDisplayProperties().color());
121  // eveEllipsoid->SetLineColor(item()->defaultDisplayProperties().color());
122  eveEllipsoid->SetMainColor(color + context().colorManager()->offsetOfLimitedColors());
123  }
124 
125  // tracks
126  if ( item()->getConfig()->value<bool>("Draw Tracks"))
127  {
129  it != v.tracks_end() ; ++it)
130  {
131  float w = v.trackWeight(*it);
132  if (w < 0.5) continue;
133 
134  const reco::Track & track = *it->get();
135  TEveRecTrack t;
136  t.fBeta = 1.;
137  t.fV = TEveVector(track.vx(), track.vy(), track.vz());
138  t.fP = TEveVector(track.px(), track.py(), track.pz());
139  t.fSign = track.charge();
140  TEveTrack* trk = new TEveTrack(&t, context().getTrackPropagator());
141  trk->SetMainColor(item()->defaultDisplayProperties().color());
142  trk->MakeTrack();
143  setupAddElement(trk, &oItemHolder);
144  }
145  }
146  if ( item()->getConfig()->value<bool>("Draw Pseudo Track"))
147  {
148  TEveRecTrack t;
149  t.fBeta = 1.;
150  t.fV = TEveVector(v.x(),v.y(),v.z());
151  t.fP = TEveVector(-v.p4().px(), -v.p4().py(), -v.p4().pz());
152  t.fSign = 1;
153  TEveTrack* trk = new TEveTrack(&t, context().getTrackPropagator());
154  trk->SetLineStyle(7);
155  trk->MakeTrack();
156  setupAddElement(trk, &oItemHolder);
157 
158  }
159 }
160 
161 void
162 FWVertexProxyBuilder::localModelChanges(const FWModelId& iId, TEveElement* iCompound,
163  FWViewType::EType viewType, const FWViewContext* vc)
164 {
165  increaseComponentTransparency(iId.index(), iCompound, "Ellipsoid", 80);
166  TEveElement* el = iCompound->FindChild("Ellipsoid");
167  if (el)
168  el->SetMainColor(item()->getConfig()->value<long>("Ellipse Color Index") + context().colorManager()->offsetOfLimitedColors());
169 }
170 
171 //
172 // static member functions
173 //
const FWVertexProxyBuilder & operator=(const FWVertexProxyBuilder &)
const fireworks::Context & context() const
FWProxyBuilderConfiguration * getConfig() const
Definition: FWEventItem.h:166
int i
Definition: DBlmapReader.cc:9
TEveVector & RefPos()
Definition: TEveEllipsoid.h:38
#define REGISTER_FWPROXYBUILDER(_name_, _type_, _purpose_, _view_)
static const int kAllRPZBits
Definition: FWViewType.h:59
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:45
double y() const
y coordinate
Definition: Vertex.h:97
void setupAddElement(TEveElement *el, TEveElement *parent, bool set_color=true) const
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:133
const FWEventItem * item() const
void SetScale(float x)
Definition: TEveEllipsoid.h:42
int index() const
Definition: FWModelId.h:50
float trackWeight(const TrackBaseRef &r) const
returns the weight with which a Track has contributed to the vertex-fit.
T sqrt(T t)
Definition: SSEVec.h:48
virtual void setItem(const FWEventItem *iItem)
TEveTrans & RefEMtx()
Definition: TEveEllipsoid.h:40
virtual void localModelChanges(const FWModelId &iId, TEveElement *iCompound, FWViewType::EType viewType, const FWViewContext *vc)
int j
Definition: DBlmapReader.cc:9
double z() const
y coordinate
Definition: Vertex.h:99
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:137
FWGenericParameter< T > * assertParam(const std::string &name, T def)
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:147
virtual void setItem(const FWEventItem *iItem)
double x() const
x coordinate
Definition: Vertex.h:95
TEveVector & RefExtent3D()
Definition: TEveEllipsoid.h:39
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:113
Error error() const
return SMatrix
Definition: Vertex.h:116
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:145
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:38
int charge() const
track electric charge
Definition: TrackBase.h:113
T w() const
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:40
void increaseComponentTransparency(unsigned int index, TEveElement *holder, const std::string &name, Char_t transpOffset)
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:135
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:143