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