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 //
12 // user include files// user include files
21 
23 
24 #include "TEvePointSet.h"
25 #include "TMatrixDEigen.h"
26 #include "TMatrixDSym.h"
27 #include "TDecompSVD.h"
28 #include "TVectorD.h"
29 #include "TEveTrans.h"
30 #include "TEveTrack.h"
31 #include "TEveTrackPropagator.h"
32 #include "TEveStraightLineSet.h"
33 #include "TEveBoxSet.h"
34 #include "TGeoSphere.h"
35 #include "TEveGeoNode.h"
36 #include "TEveVSDStructs.h"
37 
39 {
40 public:
42  virtual ~FWVertexProxyBuilder() {}
43 
44  virtual void setItem(const FWEventItem* iItem) override
45  {
47  if (iItem)
48  {
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("Ellipse Color Index", 6l, 0l, (long)context().colorManager()->numberOfLimitedColors());
54  }
55  }
56 
58 
59 private:
60  FWVertexProxyBuilder(const FWVertexProxyBuilder&); // stop default
61  const FWVertexProxyBuilder& operator=(const FWVertexProxyBuilder&); // stop default
62 
63  virtual void build(const reco::Vertex& iData, unsigned int iIndex,TEveElement& oItemHolder, const FWViewContext*) override;
64 
65  virtual void localModelChanges(const FWModelId& iId, TEveElement* iCompound,
66  FWViewType::EType viewType, const FWViewContext* vc) override;
67 
68 };
69 
70 
71 void
72 FWVertexProxyBuilder::build(const reco::Vertex& iData, unsigned int iIndex, TEveElement& oItemHolder , const FWViewContext*)
73 {
74  const reco::Vertex & v = iData;
75 
76  // marker
77  TEveGeoManagerHolder gmgr(TEveGeoShape::GetGeoMangeur());
78  TEvePointSet* pointSet = new TEvePointSet();
79  pointSet->SetNextPoint( v.x(), v.y(), v.z() );
80  setupAddElement(pointSet, &oItemHolder);
81 
82 
83  // ellipse
84  if ( item()->getConfig()->value<bool>("Draw Ellipse"))
85  {
86 
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  {
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 
117 
118  Color_t color = item()->getConfig()->value<long>("Ellipse Color Index");
119  // eveEllipsoid->SetFillColor(item()->defaultDisplayProperties().color());
120  // eveEllipsoid->SetLineColor(item()->defaultDisplayProperties().color());
121  eveEllipsoid->SetMainColor(color + context().colorManager()->offsetOfLimitedColors());
122  }
123 
124  // tracks
125  if ( item()->getConfig()->value<bool>("Draw Tracks"))
126  {
128  it != v.tracks_end() ; ++it)
129  {
130  float w = v.trackWeight(*it);
131  if (w < 0.5) continue;
132 
133  const reco::Track & track = *it->get();
134  TEveRecTrack t;
135  t.fBeta = 1.;
136  t.fV = TEveVector(track.vx(), track.vy(), track.vz());
137  t.fP = TEveVector(track.px(), track.py(), track.pz());
138  t.fSign = track.charge();
139  TEveTrack* trk = new TEveTrack(&t, context().getTrackPropagator());
140  trk->SetMainColor(item()->defaultDisplayProperties().color());
141  trk->MakeTrack();
142  setupAddElement(trk, &oItemHolder);
143  }
144  }
145  if ( item()->getConfig()->value<bool>("Draw Pseudo Track"))
146  {
147  TEveRecTrack t;
148  t.fBeta = 1.;
149  t.fV = TEveVector(v.x(),v.y(),v.z());
150  t.fP = TEveVector(-v.p4().px(), -v.p4().py(), -v.p4().pz());
151  t.fSign = 1;
152  TEveTrack* trk = new TEveTrack(&t, context().getTrackPropagator());
153  trk->SetLineStyle(7);
154  trk->MakeTrack();
155  setupAddElement(trk, &oItemHolder);
156 
157  }
158 }
159 
160 void
161 FWVertexProxyBuilder::localModelChanges(const FWModelId& iId, TEveElement* iCompound,
162  FWViewType::EType viewType, const FWViewContext* vc)
163 {
164  increaseComponentTransparency(iId.index(), iCompound, "Ellipsoid", 80);
165  TEveElement* el = iCompound->FindChild("Ellipsoid");
166  if (el)
167  el->SetMainColor(item()->getConfig()->value<long>("Ellipse Color Index") + context().colorManager()->offsetOfLimitedColors());
168 }
169 
170 //
171 // static member functions
172 //
const FWVertexProxyBuilder & operator=(const FWVertexProxyBuilder &)
const fireworks::Context & context() const
FWProxyBuilderConfiguration * getConfig() const
Definition: FWEventItem.h:165
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:58
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:44
double y() const
y coordinate
Definition: Vertex.h:96
void setupAddElement(TEveElement *el, TEveElement *parent, bool set_color=true) const
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
virtual 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:131
const FWEventItem * item() const
void SetScale(float x)
Definition: TEveEllipsoid.h:42
virtual void setItem(const FWEventItem *iItem) override
int index() const
Definition: FWModelId.h:49
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
int j
Definition: DBlmapReader.cc:9
double z() const
y coordinate
Definition: Vertex.h:98
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:135
FWGenericParameter< T > * assertParam(const std::string &name, T def)
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:145
double x() const
x coordinate
Definition: Vertex.h:94
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:112
Error error() const
return SMatrix
Definition: Vertex.h:115
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:143
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:37
int charge() const
track electric charge
Definition: TrackBase.h:111
T w() const
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:39
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:133
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:141