Belle II Software  release-08-01-10
EVEVisualization.cc
1 //For GFTrack visualisation:
2 /* Copyright 2011, Technische Universitaet Muenchen,
3  Author: Karl Bicker
4 
5  This file is part of GENFIT.
6 
7  GENFIT is free software: you can redistribute it and/or modify
8  it under the terms of the GNU Lesser General Public License as published
9  by the Free Software Foundation, either version 3 of the License, or
10  (at your option) any later version.
11 
12  GENFIT is distributed in the hope that it will be useful,
13  but WITHOUT ANY WARRANTY; without even the implied warranty of
14  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  GNU Lesser General Public License for more details.
16 
17  You should have received a copy of the GNU Lesser General Public License
18  along with GENFIT. If not, see <http://www.gnu.org/licenses/>.
19 */
20 #include <display/EVEVisualization.h>
21 
22 #include <display/VisualRepMap.h>
23 #include <display/EveGeometry.h>
24 #include <display/EveVisBField.h>
25 #include <display/ObjectInfo.h>
26 
27 #include <vxd/geometry/GeoCache.h>
28 #include <klm/bklm/geometry/GeometryPar.h>
29 #include <cdc/geometry/CDCGeometryPar.h>
30 #include <cdc/dataobjects/CDCRecoHit.h>
31 #include <cdc/translators/RealisticTDCCountTranslator.h>
32 #include <arich/dbobjects/ARICHGeometryConfig.h>
33 #include <simulation/dataobjects/MCParticleTrajectory.h>
34 #include <svd/reconstruction/SVDRecoHit.h>
35 #include <top/geometry/TOPGeometryPar.h>
36 
37 #include <genfit/AbsMeasurement.h>
38 #include <genfit/PlanarMeasurement.h>
39 #include <genfit/SpacepointMeasurement.h>
40 #include <genfit/WireMeasurement.h>
41 #include <genfit/WireMeasurementNew.h>
42 #include <genfit/WirePointMeasurement.h>
43 #include <genfit/DetPlane.h>
44 #include <genfit/Exception.h>
45 #include <genfit/KalmanFitterInfo.h>
46 #include <genfit/GblFitterInfo.h>
47 
48 #include <framework/dataobjects/DisplayData.h>
49 #include <framework/logging/Logger.h>
50 #include <framework/utilities/ColorPalette.h>
51 #include <framework/geometry/VectorUtil.h>
52 
53 #include <Math/VectorUtil.h>
54 #include <TEveArrow.h>
55 #include <TEveBox.h>
56 #include <TEveCalo.h>
57 #include <TEveManager.h>
58 #include <TEveGeoShape.h>
59 #include <TEveLine.h>
60 #include <TEvePointSet.h>
61 #include <TEveSelection.h>
62 #include <TEveStraightLineSet.h>
63 #include <TEveTriangleSet.h>
64 #include <TEveText.h>
65 #include <TEveTrack.h>
66 #include <TEveTrackPropagator.h>
67 #include <TGeoEltu.h>
68 #include <TGeoMatrix.h>
69 #include <TGeoManager.h>
70 #include <TGeoSphere.h>
71 #include <TGeoTube.h>
72 #include <TGLLogicalShape.h>
73 #include <TParticle.h>
74 #include <TMath.h>
75 #include <TVectorD.h>
76 #include <TMatrixD.h>
77 #include <TMatrixDSymEigen.h>
78 
79 #include <boost/math/special_functions/fpclassify.hpp>
80 #include <boost/scoped_ptr.hpp>
81 #include <boost/algorithm/string/find.hpp>
82 #include <boost/algorithm/string/trim.hpp>
83 #include <boost/algorithm/string/classification.hpp>
84 
85 #include <cassert>
86 #include <cmath>
87 
88 using namespace Belle2;
89 using namespace Belle2::TangoPalette;
90 
91 namespace {
93  template <class T> void destroyEveElement(T*& el)
94  {
95  if (!el) return;
96  if (el->GetDenyDestroy() > 1)
97  B2WARNING("destroyEveElement(): Element " << el->GetName() << " has unexpected refcount " << el->GetDenyDestroy());
98 
99  el->DecDenyDestroy(); //also deletes el when refcount <= 0
100  el = nullptr;
101  }
102 
109  void fixGeoShapeRefCount(TEveGeoShape* eveshape)
110  {
111  TGeoShape* s = eveshape->GetShape();
112  //Remove(obj) usually takes about 20us, but since usually s is at the end, we can shorten this
113  if (gGeoManager->GetListOfShapes()->Last() == s)
114  gGeoManager->GetListOfShapes()->RemoveAt(gGeoManager->GetListOfShapes()->GetLast());
115  else
116  gGeoManager->GetListOfShapes()->Remove(s);
117  }
118 }
119 
120 const int EVEVisualization::c_recoHitColor = getTColorID("Orange", 1);
122 const int EVEVisualization::c_trackColor = getTColorID("Sky Blue", 2);
123 const int EVEVisualization::c_recoTrackColor = getTColorID("Sky Blue", 1);
124 const int EVEVisualization::c_trackMarkerColor = getTColorID("Chameleon", 3);
125 const int EVEVisualization::c_klmClusterColor = getTColorID("Chameleon", 1);
126 
128  m_assignToPrimaries(false),
129  m_eclData(0),
130  m_bfield(new EveVisBField())
131 {
132  setErrScale();
133 
134  TGLLogicalShape::SetIgnoreSizeForCameraInterest(kTRUE); // Allows the visualization of the "small" error ellipsoid.
135 
136  //create new containers
137  m_trackpropagator = new TEveTrackPropagator();
138  m_trackpropagator->IncDenyDestroy();
139  m_trackpropagator->SetMagFieldObj(m_bfield, false);
140  m_trackpropagator->SetFitDaughters(false); //most secondaries are no longer immediate daughters since we might discard those!
141  m_trackpropagator->SetMaxR(EveGeometry::getMaxR()); //don't draw tracks outside detector
142  //TODO is this actually needed?
143  m_trackpropagator->SetMaxStep(1.0); //make sure to reeval magnetic field often enough
144 
145  m_tracklist = new TEveTrackList(m_trackpropagator);
146  m_tracklist->IncDenyDestroy();
147  m_tracklist->SetName("MCParticles");
148  m_tracklist->SelectByP(c_minPCut, FLT_MAX); //don't show too many particles by default...
149 
150  m_gftrackpropagator = new TEveTrackPropagator();
151  m_gftrackpropagator->IncDenyDestroy();
152  m_gftrackpropagator->SetMagFieldObj(m_bfield, false);
153  m_gftrackpropagator->SetMaxOrbs(0.5); //stop after track markers
154 
155  m_consttrackpropagator = new TEveTrackPropagator();
156  m_consttrackpropagator->IncDenyDestroy();
157  m_consttrackpropagator->SetMagField(0, 0, -1.5);
159 
160  m_calo3d = new TEveCalo3D(NULL, "ECLClusters");
161  m_calo3d->SetBarrelRadius(125.80); //inner radius of ECL barrel
162  m_calo3d->SetForwardEndCapPos(196.5); //inner edge of forward endcap
163  m_calo3d->SetBackwardEndCapPos(-102.0); //inner edge of backward endcap
164  m_calo3d->SetMaxValAbs(2.1);
165  m_calo3d->SetRnrFrame(false, false); //don't show crystal grid
166  m_calo3d->IncDenyDestroy();
167 
168  //Stop eve from deleting contents... (which might already be deleted)
169  gEve->GetSelection()->IncDenyDestroy();
170  gEve->GetHighlight()->IncDenyDestroy();
171 
172  clearEvent();
173 }
174 
175 void EVEVisualization::setOptions(const std::string& opts) { m_options = opts; }
176 
177 void EVEVisualization::setErrScale(double errScale) { m_errorScale = errScale; }
178 
180 {
181  if (!gEve)
182  return; //objects are probably already freed by Eve
183 
184  //Eve objects
185  destroyEveElement(m_eclData);
186  destroyEveElement(m_unassignedRecoHits);
187  destroyEveElement(m_tracklist);
188  destroyEveElement(m_trackpropagator);
189  destroyEveElement(m_gftrackpropagator);
190  destroyEveElement(m_consttrackpropagator);
191  destroyEveElement(m_calo3d);
192  delete m_bfield;
193 }
194 
195 void EVEVisualization::addTrackCandidate(const std::string& collectionName,
196  const RecoTrack& recoTrack)
197 {
198  const TString label = ObjectInfo::getIdentifier(&recoTrack);
199  // parse the option string ------------------------------------------------------------------------
200  bool drawHits = false;
201 
202  if (m_options != "") {
203  for (size_t i = 0; i < m_options.length(); i++) {
204  if (m_options.at(i) == 'H') drawHits = true;
205  }
206  }
207  // finished parsing the option string -------------------------------------------------------------
208 
209 
210  //track seeds
211  const B2Vector3D& track_pos = recoTrack.getPositionSeed();
212  const B2Vector3D& track_mom = recoTrack.getMomentumSeed();
213 
214  TEveStraightLineSet* lines = new TEveStraightLineSet("RecoHits for " + label);
215  lines->SetMainColor(c_recoTrackColor);
216  lines->SetMarkerColor(c_recoTrackColor);
217  lines->SetMarkerStyle(6);
218  lines->SetMainTransparency(60);
219 
220  if (drawHits) {
221  // Loop over all hits in the track (three different types)
222  for (const RecoHitInformation::UsedPXDHit* pxdHit : recoTrack.getPXDHitList()) {
223  addRecoHit(pxdHit, lines);
224  }
225 
226  for (const RecoHitInformation::UsedSVDHit* svdHit : recoTrack.getSVDHitList()) {
227  addRecoHit(svdHit, lines);
228  }
229 
230  for (const RecoHitInformation::UsedCDCHit* cdcHit : recoTrack.getCDCHitList()) {
231  addRecoHit(cdcHit, lines);
232  }
233  }
234 
235  TEveRecTrack rectrack;
236  rectrack.fP.Set(track_mom);
237  rectrack.fV.Set(track_pos);
238 
239  TEveTrack* track_lines = new TEveTrack(&rectrack, m_gftrackpropagator);
240  track_lines->SetName(label); //popup label set at end of function
241  track_lines->SetPropagator(m_gftrackpropagator);
242  track_lines->SetLineColor(c_recoTrackColor);
243  track_lines->SetLineWidth(1);
244  track_lines->SetTitle(ObjectInfo::getTitle(&recoTrack));
245 
246  track_lines->SetCharge((int)recoTrack.getChargeSeed());
247 
248 
249  track_lines->AddElement(lines);
250  addToGroup(collectionName, track_lines);
251  addObject(&recoTrack, track_lines);
252 }
253 
254 void EVEVisualization::addTrackCandidateImproved(const std::string& collectionName,
255  const RecoTrack& recoTrack)
256 {
257  const TString label = ObjectInfo::getIdentifier(&recoTrack);
258  // parse the option string ------------------------------------------------------------------------
259  bool drawHits = false;
260 
261  if (m_options != "") {
262  for (size_t i = 0; i < m_options.length(); i++) {
263  if (m_options.at(i) == 'H') drawHits = true;
264  }
265  }
266  // finished parsing the option string -------------------------------------------------------------
267 
268  // Create a track as a polyline through reconstructed points
269  // FIXME this is snatched from PrimitivePlotter, need to add extrapolation out of CDC
270  TEveLine* track = new TEveLine(); // We are going to just add points with SetNextPoint
271  std::vector<TVector3> posPoints; // But first we'll have to sort them as in RecoHits axial and stereo come in blocks
272  track->SetName(label); //popup label set at end of function
273  track->SetLineColor(c_recoTrackColor);
274  track->SetLineWidth(3);
275  track->SetTitle(ObjectInfo::getTitle(&recoTrack));
276  track->SetSmooth(true);
277 
278  for (auto recoHit : recoTrack.getRecoHitInformations()) {
279  // skip for reco hits which have not been used in the fit (and therefore have no fitted information on the plane
280  if (!recoHit->useInFit())
281  continue;
282 
283  TVector3 pos;
284  TVector3 mom;
285  TMatrixDSym cov;
286 
287  try {
288  const auto* trackPoint = recoTrack.getCreatedTrackPoint(recoHit);
289  const auto* fittedResult = trackPoint->getFitterInfo();
290  if (not fittedResult) {
291  B2WARNING("Skipping unfitted track point");
292  continue;
293  }
294  const genfit::MeasuredStateOnPlane& state = fittedResult->getFittedState();
295  state.getPosMomCov(pos, mom, cov);
296  } catch (const genfit::Exception&) {
297  B2WARNING("Skipping state with strange pos, mom or cov");
298  continue;
299  }
300 
301  posPoints.push_back(TVector3(pos.X(), pos.Y(), pos.Z()));
302  }
303 
304  sort(posPoints.begin(), posPoints.end(),
305  [](const TVector3 & a, const TVector3 & b) -> bool {
306  return a.X() * a.X() + a.Y() * a.Y() > b.X() * b.X() + b.Y() * b.Y();
307  });
308  for (auto vec : posPoints) {
309  track->SetNextPoint(vec.X(), vec.Y(), vec.Z());
310  }
311  // add corresponding hits ---------------------------------------------------------------------
312  TEveStraightLineSet* lines = new TEveStraightLineSet("RecoHits for " + label);
313  lines->SetMainColor(c_recoTrackColor);
314  lines->SetMarkerColor(c_recoTrackColor);
315  lines->SetMarkerStyle(6);
316  lines->SetMainTransparency(60);
317 
318  if (drawHits) {
319  // Loop over all hits in the track (three different types)
320  for (const RecoHitInformation::UsedPXDHit* pxdHit : recoTrack.getPXDHitList()) {
321  addRecoHit(pxdHit, lines);
322  }
323 
324  for (const RecoHitInformation::UsedSVDHit* svdHit : recoTrack.getSVDHitList()) {
325  addRecoHit(svdHit, lines);
326  }
327 
328  for (const RecoHitInformation::UsedCDCHit* cdcHit : recoTrack.getCDCHitList()) {
329  addRecoHit(cdcHit, lines);
330  }
331  }
332 
333  track->AddElement(lines);
334  addToGroup(collectionName, track);
335  addObject(&recoTrack, track);
336 }
337 
338 void EVEVisualization::addCDCTriggerTrack(const std::string& collectionName,
339  const CDCTriggerTrack& trgTrack)
340 {
341  const TString label = ObjectInfo::getIdentifier(&trgTrack);
342 
343  TVector3 track_pos = TVector3(0, 0, trgTrack.getZ0());
344  B2Vector3D track_mom = (trgTrack.getChargeSign() == 0) ?
345  trgTrack.getDirection() * 1000 :
346  trgTrack.getMomentum(1.5);
347 
348  TEveRecTrack rectrack;
349  rectrack.fP.Set(track_mom);
350  rectrack.fV.Set(track_pos);
351 
352  TEveTrack* track_lines = new TEveTrack(&rectrack, m_consttrackpropagator);
353  track_lines->SetName(label); //popup label set at end of function
354  track_lines->SetPropagator(m_consttrackpropagator);
355  track_lines->SetLineColor(kOrange + 2);
356  track_lines->SetLineWidth(1);
357  track_lines->SetTitle(ObjectInfo::getTitle(&trgTrack) +
358  TString::Format("\ncharge: %d, phi: %.2fdeg, pt: %.2fGeV, theta: %.2fdeg, z: %.2fcm",
359  trgTrack.getChargeSign(),
360  trgTrack.getPhi0() * 180 / M_PI,
361  trgTrack.getTransverseMomentum(1.5),
362  trgTrack.getDirection().Theta() * 180 / M_PI,
363  trgTrack.getZ0()));
364 
365 
366  track_lines->SetCharge(trgTrack.getChargeSign());
367 
368  // show 2D tracks with dashed lines
369  if (trgTrack.getZ0() == 0 && trgTrack.getCotTheta() == 0)
370  track_lines->SetLineStyle(2);
371 
372  addToGroup(collectionName, track_lines);
373  addObject(&trgTrack, track_lines);
374 }
375 
377 {
378  // load the pion fit hypothesis or the hypothesis which is the closest in mass to a pion
379  // the tracking will not always successfully fit with a pion hypothesis
380  const TrackFitResult* fitResult = belle2Track->getTrackFitResultWithClosestMass(Const::pion);
381  if (!fitResult) {
382  B2ERROR("Track without TrackFitResult skipped.");
383  return;
384  }
385  const RecoTrack* track = belle2Track->getRelated<RecoTrack>();
386 
387  const TString label = ObjectInfo::getIdentifier(belle2Track) + " (" + ObjectInfo::getIdentifier(fitResult) + ")";
388 
389  // parse the option string ------------------------------------------------------------------------
390  bool drawDetectors = false;
391  bool drawHits = false;
392  bool drawPlanes = false;
393 
394  if (m_options != "") {
395  for (size_t i = 0; i < m_options.length(); i++) {
396  if (m_options.at(i) == 'D') drawDetectors = true;
397  if (m_options.at(i) == 'H') drawHits = true;
398  if (m_options.at(i) == 'P') drawPlanes = true;
399  }
400  }
401  // finished parsing the option string -------------------------------------------------------------
402 
403  // We loop over all points (scattering non-measurement points for GBL)
404  // and for Kalman we skip those with no measurements, which should
405  // not be there
406  bool isPruned = (track == nullptr);
407 
408 
409  TEveRecTrackD recTrack;
410  const B2Vector3D& poca = fitResult->getPosition();
411  recTrack.fV.Set(poca);
412 
413  const B2Vector3D& poca_momentum = fitResult->getMomentum();
414  if (std::isfinite(poca_momentum.Mag()))
415  recTrack.fP.Set(poca_momentum);
416  else //use 1TeV momentum for tracks without curvature
417  recTrack.fP.Set(B2Vector3D(fitResult->getHelix().getDirection() * 1000));
418 
419  recTrack.fSign = fitResult->getChargeSign();
420  TEveTrack* eveTrack = new TEveTrack(&recTrack, m_gftrackpropagator);
421  eveTrack->SetName(label);
422 
423 
424  if (track) {
425  const genfit::AbsTrackRep* representation;
426 
427  if (m_drawCardinalRep) {
428  representation = track->getCardinalRepresentation();
429  B2DEBUG(100, "Draw cardinal rep");
430  } else {
431  const auto& representations = track->getRepresentations();
432  if (representations.empty()) {
433  B2ERROR("No representations found in the reco track!");
434  return;
435  }
436  B2DEBUG(100, "Draw representation number 0.");
437  representation = representations.front();
438  }
439 
440  if (!track->hasTrackFitStatus(representation)) {
441  B2ERROR("RecoTrack without FitStatus: will be skipped!");
442  return;
443  }
444 
445  const genfit::FitStatus* fitStatus = track->getTrackFitStatus(representation);
446 
447  isPruned = fitStatus->isTrackPruned();
448 
449  // GBL and Kalman cannot mix in a track.
450  // What is 0 after first loop will stay 0:
451  genfit::KalmanFitterInfo* fi = 0;
452  genfit::KalmanFitterInfo* prevFi = 0;
453  genfit::GblFitterInfo* gfi = 0;
454  genfit::GblFitterInfo* prevGFi = 0;
455  const genfit::MeasuredStateOnPlane* fittedState(NULL);
456  const genfit::MeasuredStateOnPlane* prevFittedState(NULL);
457 
458 
459  const auto& hitPoints = track->getHitPointsWithMeasurement();
460  const unsigned int numpoints = hitPoints.size();
461 
462  int hitCounter = -1;
463  for (genfit::TrackPoint* tp : hitPoints) { // loop over all points in the track
464  hitCounter++;
465 
466  // get the fitter infos ------------------------------------------------------------------
467  if (! tp->hasFitterInfo(representation)) {
468  B2ERROR("trackPoint has no fitterInfo for rep");
469  continue;
470  }
471 
472  genfit::AbsFitterInfo* fitterInfo = tp->getFitterInfo(representation);
473 
474  fi = dynamic_cast<genfit::KalmanFitterInfo*>(fitterInfo);
475  gfi = dynamic_cast<genfit::GblFitterInfo*>(fitterInfo);
476 
477  if (!gfi && !fi) {
478  B2ERROR("Can only display KalmanFitterInfo or GblFitterInfo");
479  continue;
480  }
481 
482  if (gfi && fi)
483  B2FATAL("AbsFitterInfo dynamic-casted to both KalmanFitterInfo and GblFitterInfo!");
484 
485 
486  if (fi && ! tp->hasRawMeasurements()) {
487  B2ERROR("trackPoint has no raw measurements");
488  continue;
489  }
490 
491  if (fi && ! fi->hasPredictionsAndUpdates()) {
492  B2ERROR("KalmanFitterInfo does not have all predictions and updates");
493  continue;
494  }
495 
496  try {
497  if (gfi)
498  fittedState = &(gfi->getFittedState(true));
499  if (fi)
500  fittedState = &(fi->getFittedState(true));
501  } catch (genfit::Exception& e) {
502  B2ERROR(e.what() << " - can not get fitted state");
503  continue;
504  }
505 
506  TVector3 track_pos = representation->getPos(*fittedState);
507 
508  // draw track if corresponding option is set ------------------------------------------
509  if (prevFittedState != NULL) {
510 
511  TEvePathMark::EType_e markType = TEvePathMark::kReference;
512  if (hitCounter + 1 == static_cast<int>(numpoints)) //track should stop here.
513  markType = TEvePathMark::kDecay;
514 
515  // Kalman: non-null prevFi ensures that the previous fitter info was also KalmanFitterInfo
516  if (fi && prevFi) {
517  makeLines(eveTrack, prevFittedState, fittedState, representation, markType, m_drawErrors);
518  if (m_drawErrors) { // make sure to draw errors in both directions
519  makeLines(eveTrack, prevFittedState, fittedState, representation, markType, m_drawErrors, 0);
520  }
521  //these are currently disabled.
522  //TODO: if activated, I want to have a separate TEveStraightLineSet instead of eveTrack (different colors/options)
523  if (m_drawForward)
524  makeLines(eveTrack, prevFi->getForwardUpdate(), fi->getForwardPrediction(), representation, markType, m_drawErrors, 0);
525  if (m_drawBackward)
526  makeLines(eveTrack, prevFi->getBackwardPrediction(), fi->getBackwardUpdate(), representation, markType, m_drawErrors);
527  // draw reference track if corresponding option is set ------------------------------------------
528  if (m_drawRefTrack && fi->hasReferenceState() && prevFi->hasReferenceState())
529  makeLines(eveTrack, prevFi->getReferenceState(), fi->getReferenceState(), representation, markType, false);
530  }
531 
532  // GBL: non-null prevGFi ensures that the previous fitter info was also GblFitterInfo
533  if (gfi && prevGFi) {
534  makeLines(eveTrack, prevFittedState, fittedState, representation, markType, m_drawErrors);
535  if (m_drawErrors) {
536  makeLines(eveTrack, prevFittedState, fittedState, representation, markType, m_drawErrors, 0);
537  }
538 
539  if (m_drawRefTrack && gfi->hasReferenceState() && prevGFi->hasReferenceState()) {
540  genfit::StateOnPlane prevSop = prevGFi->getReferenceState();
542  makeLines(eveTrack, &prevSop, &sop, representation, markType, false);
543  }
544  }
545 
546  }
547  prevFi = fi; // Kalman
548  prevGFi = gfi; // GBL
549  prevFittedState = fittedState;
550 
551 
552  //loop over all measurements for this point (e.g. u and v-type strips)
553  const int numMeasurements = tp->getNumRawMeasurements();
554  for (int iMeasurement = 0; iMeasurement < numMeasurements; iMeasurement++) {
555  const genfit::AbsMeasurement* m = tp->getRawMeasurement(iMeasurement);
556 
557  TVectorT<double> hit_coords;
558  TMatrixTSym<double> hit_cov;
559 
560  if (fi) {
561  // Kalman
562  genfit::MeasurementOnPlane* mop = fi->getMeasurementOnPlane(iMeasurement);
563  hit_coords.ResizeTo(mop->getState());
564  hit_cov.ResizeTo(mop->getCov());
565  hit_coords = mop->getState();
566  hit_cov = mop->getCov();
567  }
568  if (gfi) {
569  // GBL - has only one measurement (other should be already merged before the track is constructed)
570  // That means tp->getNumRawMeasurements() returns 1 for tracks fitted by GBL, because GBLfit Module
571  // merges all corresponding SVD strips to 2D combined clusters.
573  hit_coords.ResizeTo(gblMeas.getState());
574  hit_cov.ResizeTo(gblMeas.getCov());
575  hit_coords = gblMeas.getState();
576  hit_cov = gblMeas.getCov();
577  }
578 
579  // finished getting the hit infos -----------------------------------------------------
580 
581  // sort hit infos into variables ------------------------------------------------------
582  TVector3 o = fittedState->getPlane()->getO();
583  TVector3 u = fittedState->getPlane()->getU();
584  TVector3 v = fittedState->getPlane()->getV();
585 
586  bool planar_hit = false;
587  bool planar_pixel_hit = false;
588  bool space_hit = false;
589  bool wire_hit = false;
590  bool wirepoint_hit = false;
591  double_t hit_u = 0;
592  double_t hit_v = 0;
593  double_t plane_size = 4;
594 
595  int hit_coords_dim = m->getDim();
596 
597  if (dynamic_cast<const genfit::PlanarMeasurement*>(m) != NULL) {
598  planar_hit = true;
599  if (hit_coords_dim == 1) {
600  hit_u = hit_coords(0);
601  } else if (hit_coords_dim == 2) {
602  planar_pixel_hit = true;
603  hit_u = hit_coords(0);
604  hit_v = hit_coords(1);
605  }
606  } else if (dynamic_cast<const genfit::SpacepointMeasurement*>(m) != NULL) {
607  space_hit = true;
608  } else if (dynamic_cast<const CDCRecoHit*>(m)
609  || dynamic_cast<const genfit::WireMeasurement*>(m)
610  || dynamic_cast<const genfit::WireMeasurementNew*>(m)) {
611  wire_hit = true;
612  hit_u = fabs(hit_coords(0));
613  hit_v = v * (track_pos - o); // move the covariance tube so that the track goes through it
614  if (dynamic_cast<const genfit::WirePointMeasurement*>(m) != NULL) {
615  wirepoint_hit = true;
616  hit_v = hit_coords(1);
617  }
618  } else {
619  B2ERROR("Hit " << hitCounter << ", Measurement " << iMeasurement << ": Unknown measurement type: skipping hit!");
620  break;
621  }
622 
623  // finished setting variables ---------------------------------------------------------
624 
625  // draw planes if corresponding option is set -----------------------------------------
626  if (drawPlanes || (drawDetectors && planar_hit)) {
627  TVector3 move(0, 0, 0);
628  if (wire_hit) move = v * (v * (track_pos - o)); // move the plane along the wire until the track goes through it
629  TEveBox* box = boxCreator(o + move, u, v, plane_size, plane_size, 0.01);
630  if (drawDetectors && planar_hit) {
631  box->SetMainColor(kCyan);
632  } else {
633  box->SetMainColor(kGray);
634  }
635  box->SetMainTransparency(50);
636  eveTrack->AddElement(box);
637  }
638  // finished drawing planes ------------------------------------------------------------
639 
640  // draw detectors if option is set, only important for wire hits ----------------------
641  if (drawDetectors) {
642 
643  if (wire_hit) {
644  TEveGeoShape* det_shape = new TEveGeoShape("det_shape");
645  det_shape->SetShape(new TGeoTube(std::max(0., (double)(hit_u - 0.0105 / 2.)), hit_u + 0.0105 / 2., plane_size));
646  fixGeoShapeRefCount(det_shape);
647 
648  TVector3 norm = u.Cross(v);
649  TGeoRotation det_rot("det_rot", (u.Theta() * 180) / TMath::Pi(), (u.Phi() * 180) / TMath::Pi(),
650  (norm.Theta() * 180) / TMath::Pi(), (norm.Phi() * 180) / TMath::Pi(),
651  (v.Theta() * 180) / TMath::Pi(), (v.Phi() * 180) / TMath::Pi()); // move the tube to the right place and rotate it correctly
652  TVector3 move = v * (v * (track_pos - o)); // move the tube along the wire until the track goes through it
653  TGeoCombiTrans det_trans(o(0) + move.X(),
654  o(1) + move.Y(),
655  o(2) + move.Z(),
656  &det_rot);
657  det_shape->SetTransMatrix(det_trans);
658  det_shape->SetMainColor(kCyan);
659  det_shape->SetMainTransparency(25);
660  if ((drawHits && (hit_u + 0.0105 / 2 > 0)) || !drawHits) {
661  eveTrack->AddElement(det_shape);
662  }
663  }
664 
665  }
666  // finished drawing detectors ---------------------------------------------------------
667 
668  if (drawHits) {
669 
670  // draw planar hits, with distinction between strip and pixel hits ----------------
671  if (planar_hit) {
672  if (!planar_pixel_hit) {
673  //strip hit
675  const SVDRecoHit* recoHit = dynamic_cast<const SVDRecoHit*>(m);
676  if (!recoHit) {
677  B2WARNING("SVD recohit couldn't be converted... ");
678  continue;
679  }
680 
681  const VXD::SensorInfoBase& sensor = geo.get(recoHit->getSensorID());
682  double du, dv;
683  TVector3 a = o; //defines position of sensor plane
684  double hit_res_u = hit_cov(0, 0);
685  if (recoHit->isU()) {
686  du = std::sqrt(hit_res_u);
687  dv = sensor.getLength();
688  a += hit_u * u;
689  } else {
690  du = sensor.getWidth();
691  dv = std::sqrt(hit_res_u);
692  a += hit_u * v;
693  }
694  double depth = sensor.getThickness();
695  TEveBox* hit_box = boxCreator(a, u, v, du, dv, depth);
696  hit_box->SetName("SVDRecoHit");
697  hit_box->SetMainColor(c_recoHitColor);
698  hit_box->SetMainTransparency(0);
699  eveTrack->AddElement(hit_box);
700  } else {
701  //pixel hit
702  // calculate eigenvalues to draw error-ellipse ----------------------------
703  TMatrixDSymEigen eigen_values(hit_cov);
704  TEveGeoShape* cov_shape = new TEveGeoShape("PXDRecoHit");
705  const TVectorD& ev = eigen_values.GetEigenValues();
706  const TMatrixD& eVec = eigen_values.GetEigenVectors();
707  double pseudo_res_0 = m_errorScale * std::sqrt(ev(0));
708  double pseudo_res_1 = m_errorScale * std::sqrt(ev(1));
709  // finished calcluating, got the values -----------------------------------
710 
711  // calculate the semiaxis of the error ellipse ----------------------------
712  cov_shape->SetShape(new TGeoEltu(pseudo_res_0, pseudo_res_1, 0.0105));
713  fixGeoShapeRefCount(cov_shape);
714  TVector3 pix_pos = o + hit_u * u + hit_v * v;
715  TVector3 u_semiaxis = (pix_pos + eVec(0, 0) * u + eVec(1, 0) * v) - pix_pos;
716  TVector3 v_semiaxis = (pix_pos + eVec(0, 1) * u + eVec(1, 1) * v) - pix_pos;
717  TVector3 norm = u.Cross(v);
718  // finished calculating ---------------------------------------------------
719 
720  // rotate and translate everything correctly ------------------------------
721  TGeoRotation det_rot("det_rot", (u_semiaxis.Theta() * 180) / TMath::Pi(), (u_semiaxis.Phi() * 180) / TMath::Pi(),
722  (v_semiaxis.Theta() * 180) / TMath::Pi(), (v_semiaxis.Phi() * 180) / TMath::Pi(),
723  (norm.Theta() * 180) / TMath::Pi(), (norm.Phi() * 180) / TMath::Pi());
724  TGeoCombiTrans det_trans(pix_pos(0), pix_pos(1), pix_pos(2), &det_rot);
725  cov_shape->SetTransMatrix(det_trans);
726  // finished rotating and translating --------------------------------------
727 
728  cov_shape->SetMainColor(c_recoHitColor);
729  cov_shape->SetMainTransparency(0);
730  eveTrack->AddElement(cov_shape);
731  }
732  }
733  // finished drawing planar hits ---------------------------------------------------
734 
735  // draw spacepoint hits -----------------------------------------------------------
736  if (space_hit) {
737 
738  // get eigenvalues of covariance to know how to draw the ellipsoid ------------
739  TMatrixDSymEigen eigen_values(m->getRawHitCov());
740  TEveGeoShape* cov_shape = new TEveGeoShape("SpacePoint Hit");
741  cov_shape->SetShape(new TGeoSphere(0., 1.));
742  fixGeoShapeRefCount(cov_shape);
743  const TVectorD& ev = eigen_values.GetEigenValues();
744  const TMatrixD& eVec = eigen_values.GetEigenVectors();
745  TVector3 eVec1(eVec(0, 0), eVec(1, 0), eVec(2, 0));
746  TVector3 eVec2(eVec(0, 1), eVec(1, 1), eVec(2, 1));
747  TVector3 eVec3(eVec(0, 2), eVec(1, 2), eVec(2, 2));
748  TVector3 norm = u.Cross(v);
749  // got everything we need -----------------------------------------------------
750 
751 
752  TGeoRotation det_rot("det_rot", (eVec1.Theta() * 180) / TMath::Pi(), (eVec1.Phi() * 180) / TMath::Pi(),
753  (eVec2.Theta() * 180) / TMath::Pi(), (eVec2.Phi() * 180) / TMath::Pi(),
754  (eVec3.Theta() * 180) / TMath::Pi(), (eVec3.Phi() * 180) / TMath::Pi()); // the rotation is already clear
755 
756  // set the scaled eigenvalues -------------------------------------------------
757  double pseudo_res_0 = m_errorScale * std::sqrt(ev(0));
758  double pseudo_res_1 = m_errorScale * std::sqrt(ev(1));
759  double pseudo_res_2 = m_errorScale * std::sqrt(ev(2));
760  // finished scaling -----------------------------------------------------------
761 
762  // rotate and translate -------------------------------------------------------
763  TGeoGenTrans det_trans(o(0), o(1), o(2),
764  //std::sqrt(pseudo_res_0/pseudo_res_1/pseudo_res_2), std::sqrt(pseudo_res_1/pseudo_res_0/pseudo_res_2), std::sqrt(pseudo_res_2/pseudo_res_0/pseudo_res_1), // this workaround is necessary due to the "normalization" performed in TGeoGenTrans::SetScale
765  //1/(pseudo_res_0),1/(pseudo_res_1),1/(pseudo_res_2),
766  pseudo_res_0, pseudo_res_1, pseudo_res_2,
767  &det_rot);
768  cov_shape->SetTransMatrix(det_trans);
769  // finished rotating and translating ------------------------------------------
770 
771  cov_shape->SetMainColor(c_recoHitColor);
772  cov_shape->SetMainTransparency(10);
773  eveTrack->AddElement(cov_shape);
774  }
775  // finished drawing spacepoint hits -----------------------------------------------
776 
777  // draw wire hits -----------------------------------------------------------------
778  if (wire_hit) {
779  const double cdcErrorScale = 1.0;
780  TEveGeoShape* cov_shape = new TEveGeoShape("CDCRecoHit");
781  double pseudo_res_0 = cdcErrorScale * std::sqrt(hit_cov(0, 0));
782  double pseudo_res_1 = plane_size;
783  if (wirepoint_hit) pseudo_res_1 = cdcErrorScale * std::sqrt(hit_cov(1, 1));
784 
785  cov_shape->SetShape(new TGeoTube(std::max(0., (double)(hit_u - pseudo_res_0)), hit_u + pseudo_res_0, pseudo_res_1));
786  fixGeoShapeRefCount(cov_shape);
787  TVector3 norm = u.Cross(v);
788 
789  // rotate and translate -------------------------------------------------------
790  TGeoRotation det_rot("det_rot", (u.Theta() * 180) / TMath::Pi(), (u.Phi() * 180) / TMath::Pi(),
791  (norm.Theta() * 180) / TMath::Pi(), (norm.Phi() * 180) / TMath::Pi(),
792  (v.Theta() * 180) / TMath::Pi(), (v.Phi() * 180) / TMath::Pi());
793  TGeoCombiTrans det_trans(o(0) + hit_v * v.X(),
794  o(1) + hit_v * v.Y(),
795  o(2) + hit_v * v.Z(),
796  &det_rot);
797  cov_shape->SetTransMatrix(det_trans);
798  // finished rotating and translating ------------------------------------------
799 
800  cov_shape->SetMainColor(c_recoHitColor);
801  cov_shape->SetMainTransparency(50);
802  eveTrack->AddElement(cov_shape);
803  }
804  // finished drawing wire hits -----------------------------------------------------
805  }
806 
807  }
808  }
809 
810  auto& firstref = eveTrack->RefPathMarks().front();
811  auto& lastref = eveTrack->RefPathMarks().back();
812  double f = firstref.fV.Distance(recTrack.fV);
813  double b = lastref.fV.Distance(recTrack.fV);
814  if (f > 100 and f > b) {
815  B2WARNING("Decay vertex is much closer to POCA than first vertex, reversing order of track points... (this is intended for cosmic tracks, if you see this message in other context it might indicate a problem)");
816  //last ref is better than first...
817  lastref.fType = TEvePathMarkD::kReference;
818  firstref.fType = TEvePathMarkD::kDecay;
819  std::reverse(eveTrack->RefPathMarks().begin(), eveTrack->RefPathMarks().end());
820  }
821  }
822  eveTrack->SetTitle(TString::Format("%s\n"
823  "pruned: %s\n"
824  "pT=%.3f, pZ=%.3f\n"
825  "pVal: %e",
826  label.Data(),
827  isPruned ? " yes" : "no",
828  poca_momentum.Pt(), poca_momentum.Pz(),
829  fitResult->getPValue()
830  ));
831  eveTrack->SetLineColor(c_trackColor);
832  eveTrack->SetLineStyle(1);
833  eveTrack->SetLineWidth(3.0);
834 
835 
836  addToGroup("Fitted Tracks", eveTrack);
837  if (track)
838  addObject(track, eveTrack);
839  addObject(belle2Track, eveTrack);
840 }
841 
842 TEveBox* EVEVisualization::boxCreator(const TVector3& o, TVector3 u, TVector3 v, float ud, float vd, float depth)
843 {
844  //force minimum width of polygon to deal with Eve limits
845  float min = 0.04;
846  if (vd < min)
847  vd = min;
848  if (ud < min)
849  ud = min;
850  if (depth < min)
851  depth = min;
852 
853  TEveBox* box = new TEveBox;
854  box->SetPickable(true);
855 
856  TVector3 norm = u.Cross(v);
857  u *= (0.5 * ud);
858  v *= (0.5 * vd);
859  norm *= (0.5 * depth);
860 
861 
862  for (int k = 0; k < 8; ++k) {
863  // Coordinates for all eight corners of the box
864  // as two parallel rectangles, with vertices specified in clockwise direction
865  int signU = ((k + 1) & 2) ? -1 : 1;
866  int signV = (k & 4) ? -1 : 1;
867  int signN = (k & 2) ? -1 : 1;
868  float vertex[3];
869  for (int i = 0; i < 3; ++i) {
870  vertex[i] = o(i) + signU * u(i) + signV * v(i) + signN * norm(i);
871  }
872  box->SetVertex(k, vertex);
873  }
874 
875  return box;
876 }
877 
878 void EVEVisualization::makeLines(TEveTrack* eveTrack, const genfit::StateOnPlane* prevState, const genfit::StateOnPlane* state,
879  const genfit::AbsTrackRep* rep,
880  TEvePathMark::EType_e markType, bool drawErrors, int markerPos)
881 {
882  using namespace genfit;
883 
884  TVector3 pos, dir, oldPos, oldDir;
885  rep->getPosDir(*state, pos, dir);
886  rep->getPosDir(*prevState, oldPos, oldDir);
887 
888  double distA = (pos - oldPos).Mag();
889  double distB = distA;
890  if ((pos - oldPos)*oldDir < 0)
891  distA *= -1.;
892  if ((pos - oldPos)*dir < 0)
893  distB *= -1.;
894 
895  TEvePathMark mark(
896  markType,
897  TEveVector(pos(0), pos(1), pos(2)),
898  TEveVector(dir(0), dir(1), dir(2))
899  );
900  eveTrack->AddPathMark(mark);
901 
902 
903  if (drawErrors) {
904  const MeasuredStateOnPlane* measuredState;
905  if (markerPos == 0)
906  measuredState = dynamic_cast<const MeasuredStateOnPlane*>(prevState);
907  else
908  measuredState = dynamic_cast<const MeasuredStateOnPlane*>(state);
909 
910  if (measuredState != NULL) {
911 
912  // step for evaluate at a distance from the original plane
913  TVector3 eval;
914  if (markerPos == 0)
915  eval = 0.2 * distA * oldDir;
916  else
917  eval = -0.2 * distB * dir;
918 
919 
920  // get cov at first plane
921  TMatrixDSym cov;
922  TVector3 position, direction;
923  rep->getPosMomCov(*measuredState, position, direction, cov);
924 
925  // get eigenvalues & -vectors
926  TMatrixDSymEigen eigen_values(cov.GetSub(0, 2, 0, 2));
927  const TVectorD& ev = eigen_values.GetEigenValues();
928  const TMatrixD& eVec = eigen_values.GetEigenVectors();
929  TVector3 eVec1, eVec2;
930  // limit
931  static const double maxErr = 1000.;
932  double ev0 = std::min(ev(0), maxErr);
933  double ev1 = std::min(ev(1), maxErr);
934  double ev2 = std::min(ev(2), maxErr);
935 
936  // get two largest eigenvalues/-vectors
937  if (ev0 < ev1 && ev0 < ev2) {
938  eVec1.SetXYZ(eVec(0, 1), eVec(1, 1), eVec(2, 1));
939  eVec1 *= sqrt(ev1);
940  eVec2.SetXYZ(eVec(0, 2), eVec(1, 2), eVec(2, 2));
941  eVec2 *= sqrt(ev2);
942  } else if (ev1 < ev0 && ev1 < ev2) {
943  eVec1.SetXYZ(eVec(0, 0), eVec(1, 0), eVec(2, 0));
944  eVec1 *= sqrt(ev0);
945  eVec2.SetXYZ(eVec(0, 2), eVec(1, 2), eVec(2, 2));
946  eVec2 *= sqrt(ev2);
947  } else {
948  eVec1.SetXYZ(eVec(0, 0), eVec(1, 0), eVec(2, 0));
949  eVec1 *= sqrt(ev0);
950  eVec2.SetXYZ(eVec(0, 1), eVec(1, 1), eVec(2, 1));
951  eVec2 *= sqrt(ev1);
952  }
953 
954  if (eVec1.Cross(eVec2)*eval < 0)
955  eVec2 *= -1;
956  //assert(eVec1.Cross(eVec2)*eval > 0);
957 
958  TVector3 oldEVec1(eVec1);
959 
960  const int nEdges = 24;
961  std::vector<TVector3> vertices;
962 
963  vertices.push_back(position);
964 
965  // vertices at plane
966  for (int i = 0; i < nEdges; ++i) {
967  const double angle = 2 * TMath::Pi() / nEdges * i;
968  vertices.push_back(position + cos(angle)*eVec1 + sin(angle)*eVec2);
969  }
970 
971 
972 
973  SharedPlanePtr newPlane(new DetPlane(*(measuredState->getPlane())));
974  newPlane->setO(position + eval);
975 
976  MeasuredStateOnPlane stateCopy(*measuredState);
977  try {
978  rep->extrapolateToPlane(stateCopy, newPlane);
979  } catch (Exception& e) {
980  B2ERROR(e.what());
981  return;
982  }
983 
984  // get cov at 2nd plane
985  rep->getPosMomCov(stateCopy, position, direction, cov);
986 
987  // get eigenvalues & -vectors
988  {
989  TMatrixDSymEigen eigen_values2(cov.GetSub(0, 2, 0, 2));
990  const TVectorD& eVal = eigen_values2.GetEigenValues();
991  const TMatrixD& eVect = eigen_values2.GetEigenVectors();
992  // limit
993  ev0 = std::min(eVal(0), maxErr);
994  ev1 = std::min(eVal(1), maxErr);
995  ev2 = std::min(eVal(2), maxErr);
996 
997  // get two largest eigenvalues/-vectors
998  if (ev0 < ev1 && ev0 < ev2) {
999  eVec1.SetXYZ(eVect(0, 1), eVect(1, 1), eVect(2, 1));
1000  eVec1 *= sqrt(ev1);
1001  eVec2.SetXYZ(eVect(0, 2), eVect(1, 2), eVect(2, 2));
1002  eVec2 *= sqrt(ev2);
1003  } else if (ev1 < ev0 && ev1 < ev2) {
1004  eVec1.SetXYZ(eVect(0, 0), eVect(1, 0), eVect(2, 0));
1005  eVec1 *= sqrt(ev0);
1006  eVec2.SetXYZ(eVect(0, 2), eVect(1, 2), eVect(2, 2));
1007  eVec2 *= sqrt(ev2);
1008  } else {
1009  eVec1.SetXYZ(eVect(0, 0), eVect(1, 0), eVect(2, 0));
1010  eVec1 *= sqrt(ev0);
1011  eVec2.SetXYZ(eVect(0, 1), eVect(1, 1), eVect(2, 1));
1012  } eVec2 *= sqrt(ev1);
1013  }
1014 
1015  if (eVec1.Cross(eVec2)*eval < 0)
1016  eVec2 *= -1;
1017  //assert(eVec1.Cross(eVec2)*eval > 0);
1018 
1019  if (oldEVec1 * eVec1 < 0) {
1020  eVec1 *= -1;
1021  eVec2 *= -1;
1022  }
1023 
1024  // vertices at 2nd plane
1025  double angle0 = eVec1.Angle(oldEVec1);
1026  if (eVec1 * (eval.Cross(oldEVec1)) < 0)
1027  angle0 *= -1;
1028  for (int i = 0; i < nEdges; ++i) {
1029  const double angle = 2 * TMath::Pi() / nEdges * i - angle0;
1030  vertices.push_back(position + cos(angle)*eVec1 + sin(angle)*eVec2);
1031  }
1032 
1033  vertices.push_back(position);
1034 
1035 
1036  TEveTriangleSet* error_shape = new TEveTriangleSet(vertices.size(), nEdges * 2);
1037  for (unsigned int k = 0; k < vertices.size(); ++k) {
1038  error_shape->SetVertex(k, vertices[k].X(), vertices[k].Y(), vertices[k].Z());
1039  }
1040 
1041  assert(vertices.size() == 2 * nEdges + 2);
1042 
1043  int iTri(0);
1044  for (int i = 0; i < nEdges; ++i) {
1045  //error_shape->SetTriangle(iTri++, 0, i+1, (i+1)%nEdges+1);
1046  error_shape->SetTriangle(iTri++, i + 1, i + 1 + nEdges, (i + 1) % nEdges + 1);
1047  error_shape->SetTriangle(iTri++, (i + 1) % nEdges + 1, i + 1 + nEdges, (i + 1) % nEdges + 1 + nEdges);
1048  //error_shape->SetTriangle(iTri++, 2*nEdges+1, i+1+nEdges, (i+1)%nEdges+1+nEdges);
1049  }
1050 
1051  //assert(iTri == nEdges*4);
1052 
1053  error_shape->SetMainColor(c_trackColor);
1054  error_shape->SetMainTransparency(25);
1055  eveTrack->AddElement(error_shape);
1056  }
1057  }
1058 }
1059 
1060 void EVEVisualization::addSimHit(const CDCSimHit* hit, const MCParticle* particle)
1061 {
1062  addSimHit(ROOT::Math::XYZVector(hit->getPosWire()), particle);
1063 }
1064 void EVEVisualization::addSimHit(const PXDSimHit* hit, const MCParticle* particle)
1065 {
1066  static VXD::GeoCache& geo = VXD::GeoCache::getInstance();
1067  const ROOT::Math::XYZVector& global_pos = geo.get(hit->getSensorID()).pointToGlobal(hit->getPosIn());
1068  addSimHit(global_pos, particle);
1069 }
1070 void EVEVisualization::addSimHit(const SVDSimHit* hit, const MCParticle* particle)
1071 {
1072  static VXD::GeoCache& geo = VXD::GeoCache::getInstance();
1073  const ROOT::Math::XYZVector& global_pos = geo.get(hit->getSensorID()).pointToGlobal(hit->getPosIn());
1074  addSimHit(global_pos, particle);
1075 }
1076 void EVEVisualization::addSimHit(const KLMSimHit* hit, const MCParticle* particle)
1077 {
1078  const TVector3& global_pos = hit->getPosition();
1079  addSimHit(global_pos, particle);
1080 }
1081 void EVEVisualization::addSimHit(const TVector3& v, const MCParticle* particle)
1082 {
1083  MCTrack* track = addMCParticle(particle);
1084  if (!track)
1085  return; //hide hits from this particle
1086  track->simhits->SetNextPoint(v.X(), v.Y(), v.Z());
1087 }
1088 
1090 {
1091  if (!particle) {
1092  if (!m_mcparticleTracks[nullptr].simhits) {
1093  const TString pointsTitle("Unassigned SimHits");
1094  m_mcparticleTracks[nullptr].simhits = new TEvePointSet(pointsTitle);
1095  m_mcparticleTracks[nullptr].simhits->SetTitle(pointsTitle);
1096  m_mcparticleTracks[nullptr].simhits->SetMarkerStyle(6);
1097  m_mcparticleTracks[nullptr].simhits->SetMainColor(c_unassignedHitColor);
1098  //m_mcparticleTracks[nullptr].simhits->SetMainTransparency(50);
1099  m_mcparticleTracks[nullptr].track = NULL;
1100  }
1101  return &m_mcparticleTracks[nullptr];
1102  }
1103 
1104  if (m_hideSecondaries and !particle->hasStatus(MCParticle::c_PrimaryParticle)) {
1105  return NULL;
1106  }
1107  if (m_assignToPrimaries) {
1108  while (!particle->hasStatus(MCParticle::c_PrimaryParticle) and particle->getMother())
1109  particle = particle->getMother();
1110  }
1111 
1112  if (!m_mcparticleTracks[particle].track) {
1113  const ROOT::Math::XYZVector& p = particle->getMomentum();
1114  const ROOT::Math::XYZVector& vertex = particle->getProductionVertex();
1115  const int pdg = particle->getPDG();
1116  TParticle tparticle(pdg, particle->getStatus(),
1117  (particle->getMother() ? particle->getMother()->getIndex() : 0), 0, particle->getFirstDaughter(), particle->getLastDaughter(),
1118  p.X(), p.Y(), p.Z(), particle->getEnergy(),
1119  vertex.X(), vertex.Y(), vertex.Z(), particle->getProductionTime());
1120  TEveMCTrack mctrack;
1121  mctrack = tparticle;
1122  mctrack.fTDecay = particle->getDecayTime();
1123  mctrack.fVDecay.Set(B2Vector3D(particle->getDecayVertex()));
1124  mctrack.fDecayed = !boost::math::isinf(mctrack.fTDecay);
1125  mctrack.fIndex = particle->getIndex();
1126  m_mcparticleTracks[particle].track = new TEveTrack(&mctrack, m_trackpropagator);
1127 
1128  //Check if there is a trajectory stored for this particle
1129  const auto mcTrajectories = particle->getRelationsTo<MCParticleTrajectory>();
1130  bool hasTrajectory(false);
1131  for (auto rel : mcTrajectories.relations()) {
1132  //Trajectories with negative weight are from secondary daughters which
1133  //were ignored so we don't use them.
1134  if (rel.weight <= 0) continue;
1135  //Found one, let's add tose point as reference points to the TEveTrack.
1136  //This will force the track propagation to visit all points in order but
1137  //provide smooth helix interpolation between the points
1138  const MCParticleTrajectory& trajectory = dynamic_cast<const MCParticleTrajectory&>(*rel.object);
1139  for (const MCTrajectoryPoint& pt : trajectory) {
1140  m_mcparticleTracks[particle].track->AddPathMark(
1141  TEvePathMark(
1142  //Add the last trajectory point as decay point to prevent TEve to
1143  //propagate beyond the end of the track. So lets compare the adress
1144  //to the address of last point and choose the pathmark accordingly
1145  (&pt == &trajectory.back()) ? TEvePathMark::kDecay : TEvePathMark::kReference,
1146  TEveVector(pt.x, pt.y, pt.z),
1147  TEveVector(pt.px, pt.py, pt.pz)
1148  ));
1149  }
1150  //"There can only be One" -> found a trajectory, stop the loop
1151  hasTrajectory = true;
1152  break;
1153  }
1154 
1155  //If we have the full trajectory there is no need to add additional path marks
1156  if (!hasTrajectory) {
1157  //add daughter vertices - improves track rendering as lost momentum is taken into account
1158  for (int iDaughter = particle->getFirstDaughter(); iDaughter <= particle->getLastDaughter(); iDaughter++) {
1159  if (iDaughter == 0)
1160  continue; //no actual daughter
1161 
1162  const MCParticle* daughter = StoreArray<MCParticle>()[iDaughter - 1];
1163 
1164  TEvePathMarkD refMark(TEvePathMarkD::kDaughter);
1165  refMark.fV.Set(B2Vector3D(daughter->getProductionVertex()));
1166  refMark.fP.Set(B2Vector3D(daughter->getMomentum()));
1167  refMark.fTime = daughter->getProductionTime();
1168  m_mcparticleTracks[particle].track->AddPathMark(refMark);
1169  }
1170 
1171  //neutrals and very short-lived particles should stop somewhere
1172  //(can result in wrong shapes for particles stopped in the detector, so not used there)
1173  if ((TMath::Nint(particle->getCharge()) == 0 or !particle->hasStatus(MCParticle::c_StoppedInDetector))
1174  and mctrack.fDecayed) {
1175  TEvePathMarkD decayMark(TEvePathMarkD::kDecay);
1176  decayMark.fV.Set(B2Vector3D(particle->getDecayVertex()));
1177  m_mcparticleTracks[particle].track->AddPathMark(decayMark);
1178  }
1179  }
1180  TString particle_name(mctrack.GetName());
1181 
1182  //set track title (for popup)
1183  const MCParticle* mom = particle->getMother();
1184  if (mom) {
1185  m_mcparticleTracks[particle].parentParticle = mom;
1186  addMCParticle(mom);
1187  }
1188 
1189  TString title = ObjectInfo::getTitle(particle);
1190  if (!hasTrajectory) {
1191  //Hijack the mother label to show that the track position is only
1192  //extrapolated, not known from simulation
1193  title += "\n(track estimated from initial momentum)";
1194  //Also, show those tracks with dashed lines
1195  m_mcparticleTracks[particle].track->SetLineStyle(2);
1196  }
1197 
1198  m_mcparticleTracks[particle].track->SetTitle(title);
1199 
1200  //add some color (avoid black & white)
1201  switch (abs(pdg)) {
1202  case 11:
1203  m_mcparticleTracks[particle].track->SetLineColor(kAzure);
1204  break;
1205  case 13:
1206  m_mcparticleTracks[particle].track->SetLineColor(kCyan + 1);
1207  break;
1208  case 22:
1209  m_mcparticleTracks[particle].track->SetLineColor(kSpring);
1210  break;
1211  case 211:
1212  m_mcparticleTracks[particle].track->SetLineColor(kGray + 1);
1213  break;
1214  case 321:
1215  m_mcparticleTracks[particle].track->SetLineColor(kRed + 1);
1216  break;
1217  case 2212:
1218  m_mcparticleTracks[particle].track->SetLineColor(kOrange - 2);
1219  break;
1220  default:
1221  m_mcparticleTracks[particle].track->SetLineColor(kMagenta);
1222  }
1223 
1224  //create point set for hits
1225  const TString pointsTitle = "SimHits for " + ObjectInfo::getIdentifier(particle) + " - " + particle_name;
1226  m_mcparticleTracks[particle].simhits = new TEvePointSet(pointsTitle);
1227  m_mcparticleTracks[particle].simhits->SetTitle(pointsTitle);
1228  m_mcparticleTracks[particle].simhits->SetMarkerStyle(6);
1229  m_mcparticleTracks[particle].simhits->SetMainColor(m_mcparticleTracks[particle].track->GetLineColor());
1230  //m_mcparticleTracks[particle].simhits->SetMainTransparency(50);
1231  addObject(particle, m_mcparticleTracks[particle].track);
1232  }
1233  return &m_mcparticleTracks[particle];
1234 }
1235 
1237 {
1238  for (auto& mcTrackPair : m_mcparticleTracks) {
1239  MCTrack& mcTrack = mcTrackPair.second;
1240  if (mcTrack.track) {
1241  if (mcTrack.simhits->Size() > 0) {
1242  mcTrack.track->AddElement(mcTrack.simhits);
1243  } else {
1244  //if we don't add it, remove empty collection
1245  destroyEveElement(mcTrack.simhits);
1246  }
1247 
1248  TEveElement* parent = m_tracklist;
1249  if (mcTrack.parentParticle) {
1250  const auto& parentIt = m_mcparticleTracks.find(mcTrack.parentParticle);
1251  if (parentIt != m_mcparticleTracks.end()) {
1252  parent = parentIt->second.track;
1253  }
1254  }
1255  parent->AddElement(mcTrack.track);
1256  } else { //add simhits directly
1257  gEve->AddElement(mcTrack.simhits);
1258  }
1259  }
1260  gEve->AddElement(m_tracklist);
1261  m_tracklist->MakeTracks();
1262  m_tracklist->SelectByP(c_minPCut, FLT_MAX); //don't show too many particles by default...
1263 
1264  for (size_t i = 0; i < m_options.length(); i++) {
1265  if (m_options.at(i) == 'M') {
1266  m_gftrackpropagator->SetRnrDaughters(true);
1267  m_gftrackpropagator->SetRnrReferences(true);
1268  //m_gftrackpropagator->SetRnrFV(true); //TODO: this crashes?
1269  TMarker m;
1270  m.SetMarkerColor(c_trackMarkerColor);
1271  m.SetMarkerStyle(1); //a single pixel
1272  m.SetMarkerSize(1); //ignored.
1273  m_gftrackpropagator->RefPMAtt() = m;
1274  m_gftrackpropagator->RefFVAtt() = m;
1275  }
1276  }
1277 
1278  m_consttrackpropagator->SetMagField(0, 0, -1.5);
1279 
1280  m_eclData->DataChanged(); //update limits (Empty() won't work otherwise)
1281  if (!m_eclData->Empty()) {
1282  m_eclData->SetAxisFromBins(0.0,
1283  0.0); //epsilon_x/y = 0 so we don't merge neighboring bins. This avoids some rendering issues with projections of small clusters.
1284  m_calo3d->SetData(m_eclData);
1285  }
1286  gEve->AddElement(m_calo3d);
1287 
1288  if (m_unassignedRecoHits) {
1290  gEve->AddElement(m_unassignedRecoHits);
1291  }
1292 
1293 }
1294 
1296 {
1297  if (!gEve)
1298  return;
1299 
1301  for (auto& groupPair : m_groups) {
1302  //store visibility, invalidate pointers
1303  if (groupPair.second.group)
1304  groupPair.second.visible = groupPair.second.group->GetRnrState();
1305  groupPair.second.group = nullptr;
1306  }
1307 
1308  m_mcparticleTracks.clear();
1309  m_shownRecohits.clear();
1310  m_tracklist->DestroyElements();
1311 
1312  //remove ECL data from event
1313  m_calo3d->SetData(NULL);
1314  m_calo3d->DestroyElements();
1315 
1316  //lower energy threshold for ECL
1317  float ecl_threshold = 0.01;
1318  if (m_eclData)
1319  ecl_threshold = m_eclData->GetSliceThreshold(0);
1320 
1321  destroyEveElement(m_eclData);
1322  m_eclData = new TEveCaloDataVec(1); //#slices
1323  m_eclData->IncDenyDestroy();
1324  m_eclData->RefSliceInfo(0).Setup("ECL", ecl_threshold, kRed);
1325 
1328  destroyEveElement(m_unassignedRecoHits);
1329 
1330  gEve->GetSelection()->RemoveElements();
1331  gEve->GetHighlight()->RemoveElements();
1332  //other things are cleaned up by TEve...
1333 }
1334 
1335 
1336 
1337 
1339 {
1340  TVector3 v = vertex->getPos();
1341  TEvePointSet* vertexPoint = new TEvePointSet(ObjectInfo::getInfo(vertex));
1342  //sadly, setting a title for a TEveGeoShape doesn't result in a popup...
1343  vertexPoint->SetTitle(ObjectInfo::getTitle(vertex));
1344  vertexPoint->SetMainColor(c_recoHitColor);
1345  vertexPoint->SetNextPoint(v.X(), v.Y(), v.Z());
1346 
1347  TMatrixDSymEigen eigen_values(vertex->getCov());
1348  TEveGeoShape* det_shape = new TEveGeoShape(ObjectInfo::getInfo(vertex) + " Error");
1349  det_shape->SetShape(new TGeoSphere(0., 1.)); //Initially created as a sphere, then "scaled" into an ellipsoid.
1350  fixGeoShapeRefCount(det_shape);
1351  const TVectorD& ev = eigen_values.GetEigenValues(); //Assigns the eigenvalues into the "ev" matrix.
1352  const TMatrixD& eVec = eigen_values.GetEigenVectors(); //Assigns the eigenvalues into the "eVec" matrix.
1353  TVector3 eVec1(eVec(0, 0), eVec(1, 0), eVec(2,
1354  0)); //Define the 3 eigenvectors of the covariance matrix as objects of the TVector3 class using constructor.
1355  TVector3 eVec2(eVec(0, 1), eVec(1, 1), eVec(2,
1356  1)); //eVec(i,j) uses the method/overloaded operator ( . ) of the TMatrixT class to return the matrix entry.
1357  TVector3 eVec3(eVec(0, 2), eVec(1, 2), eVec(2, 2));
1358  // got everything we need ----------------------------------------------------- //Eigenvalues(semi axis) of the covariance matrix accquired!
1359 
1360 
1361  TGeoRotation det_rot("det_rot", (eVec1.Theta() * 180) / TMath::Pi(), (eVec1.Phi() * 180) / TMath::Pi(),
1362  (eVec2.Theta() * 180) / TMath::Pi(), (eVec2.Phi() * 180) / TMath::Pi(),
1363  (eVec3.Theta() * 180) / TMath::Pi(), (eVec3.Phi() * 180) / TMath::Pi()); // the rotation is already clear
1364 
1365  // set the scaled eigenvalues -------------------------------------------------
1366  //"Scaled" eigenvalues pseudo_res (lengths of the semi axis) are the sqrt of the eigenvalues.
1367  double pseudo_res_0 = std::sqrt(ev(0));
1368  double pseudo_res_1 = std::sqrt(ev(1));
1369  double pseudo_res_2 = std::sqrt(ev(2));
1370 
1371  //B2INFO("The pseudo_res_0/1/2 are " << pseudo_res_0 << "," << pseudo_res_1 << "," << pseudo_res_2); //shows the scaled eigenvalues
1372 
1373 
1374 
1375  // rotate and translate -------------------------------------------------------
1376  TGeoGenTrans det_trans(v(0), v(1), v(2), pseudo_res_0, pseudo_res_1, pseudo_res_2,
1377  &det_rot); //Puts the ellipsoid at the position of the vertex, v(0)=v.X(), operator () overloaded.
1378  det_shape->SetTransMatrix(det_trans);
1379  // finished rotating and translating ------------------------------------------
1380 
1381  det_shape->SetMainColor(kOrange); //The color of the error ellipsoid.
1382  det_shape->SetMainTransparency(0);
1383 
1384  vertexPoint->AddElement(det_shape);
1385  addToGroup("Vertices", vertexPoint);
1386  addObject(vertex, vertexPoint);
1387 }
1388 
1389 
1391 {
1392 
1393  // only display c_nPhotons hypothesis clusters
1394  if (cluster->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) {
1395 
1396  const float phi = cluster->getPhi();
1397  float dPhi = cluster->getUncertaintyPhi();
1398  float dTheta = cluster->getUncertaintyTheta();
1399  if (dPhi >= M_PI / 4 or dTheta >= M_PI / 4 or cluster->getUncertaintyEnergy() == 1.0) {
1400  B2WARNING("Found ECL cluster with broken errors (unit matrix or too large). Using 0.05 as error in phi/theta. The 3x3 error matrix previously was:");
1401  cluster->getCovarianceMatrix3x3().Print();
1402  dPhi = dTheta = 0.05;
1403  }
1404 
1405  if (!std::isfinite(dPhi) or !std::isfinite(dTheta)) {
1406  B2ERROR("ECLCluster phi or theta error is NaN or infinite, skipping cluster!");
1407  return;
1408  }
1409 
1410  //convert theta +- dTheta into eta +- dEta
1411  TVector3 thetaLow;
1412  thetaLow.SetPtThetaPhi(1.0, cluster->getTheta() - dTheta, phi);
1413  TVector3 thetaHigh;
1414  thetaHigh.SetPtThetaPhi(1.0, cluster->getTheta() + dTheta, phi);
1415  float etaLow = thetaLow.Eta();
1416  float etaHigh = thetaHigh.Eta();
1417  if (etaLow > etaHigh) {
1418  std::swap(etaLow, etaHigh);
1419  }
1420 
1421  int id = m_eclData->AddTower(etaLow, etaHigh, phi - dPhi, phi + dPhi);
1422  m_eclData->FillSlice(0, cluster->getEnergy(ECLCluster::EHypothesisBit::c_nPhotons));
1424  }
1425 }
1426 
1428 {
1429  const double layerThicknessCm = 3.16; //TDR, Fig 10.2
1430  const double layerDistanceCm = 9.1 - layerThicknessCm;
1431 
1432  // Pposition of first RPC plane.
1433  ROOT::Math::XYZVector position = cluster->getClusterPosition();
1434  TVector3 startPos(position.X(), position.Y(), position.Z());
1435  TVector3 dir; //direction of cluster stack, Mag() == distance between planes
1436  TVector3 a, b; //defines RPC plane
1437  bool isBarrel = (startPos.Z() > -175.0 and startPos.Z() < 270.0);
1438  if (isBarrel) {
1439  //barrel
1440  b = TVector3(0, 0, 1);
1441  a = startPos.Cross(b).Unit();
1442  double c = M_PI / 4.0;
1443  double offset = c / 2.0 + M_PI;
1444  a.SetPhi(int((a.Phi() + offset) / (c))*c - M_PI);
1445  TVector3 perp = b.Cross(a);
1446 
1447  const double barrelRadiusCm = 204.0;
1448  startPos.SetMag(barrelRadiusCm / perp.Dot(startPos.Unit()));
1449 
1450  dir = startPos.Unit();
1451  dir.SetMag((layerDistanceCm + layerThicknessCm) / perp.Dot(dir));
1452  } else {
1453  //endcap
1454  b = TVector3(startPos.X(), startPos.Y(), 0).Unit();
1455  a = startPos.Cross(b).Unit();
1456  double endcapStartZ = 284;
1457  if (startPos.Z() < 0)
1458  endcapStartZ = -189.5;
1459 
1460  double scaleFac = endcapStartZ / startPos.Z();
1461  startPos.SetMag(startPos.Mag() * scaleFac);
1462 
1463  dir = startPos.Unit();
1464  dir.SetMag((layerDistanceCm + layerThicknessCm) / fabs(dir.Z()));
1465  }
1466 
1467  for (int i = 0; i < cluster->getLayers(); i++) {
1468  TVector3 layerPos = startPos;
1469  layerPos += (cluster->getInnermostLayer() + i) * dir;
1470  auto* layer = boxCreator(layerPos, a, b, 20.0, 20.0, layerThicknessCm / 2);
1471  layer->SetMainColor(c_klmClusterColor);
1472  layer->SetMainTransparency(70);
1473  layer->SetName(ObjectInfo::getIdentifier(cluster));
1474  layer->SetTitle(ObjectInfo::getTitle(cluster));
1475 
1476  addToGroup(std::string("KLMClusters/") + ObjectInfo::getIdentifier(cluster).Data(), layer);
1477  addObject(cluster, layer);
1478  }
1479 }
1480 
1482 {
1483  //TVector3 globalPosition= bklm2dhit->getGlobalPosition();
1485  const bklm::Module* module = m_GeoPar->findModule(bklm2dhit->getSection(), bklm2dhit->getSector(), bklm2dhit->getLayer());
1486 
1487  CLHEP::Hep3Vector global;
1488  //+++ global coordinates of the hit
1489  global[0] = bklm2dhit->getPositionX();
1490  global[1] = bklm2dhit->getPositionY();
1491  global[2] = bklm2dhit->getPositionZ();
1492 
1493  //+++ local coordinates of the hit
1494  CLHEP::Hep3Vector local = module->globalToLocal(global);
1495  //double localU = local[1]; //phi
1496  //double localV = local[2]; //z
1497  int Nphistrip = bklm2dhit->getPhiStripMax() - bklm2dhit->getPhiStripMin() + 1;
1498  int Nztrip = bklm2dhit->getZStripMax() - bklm2dhit->getZStripMin() + 1;
1499  double du = module->getPhiStripWidth() * Nphistrip;
1500  double dv = module->getZStripWidth() * Nztrip;
1501 
1502  //Let's do some simple thing
1503  CLHEP::Hep3Vector localU(local[0], local[1] + 1.0, local[2]);
1504  CLHEP::Hep3Vector localV(local[0], local[1], local[2] + 1.0);
1505 
1506  CLHEP::Hep3Vector globalU = module->localToGlobal(localU);
1507  CLHEP::Hep3Vector globalV = module->localToGlobal(localV);
1508 
1509  TVector3 o(global[0], global[1], global[2]);
1510  TVector3 u(globalU[0], globalU[1], globalU[2]);
1511  TVector3 v(globalV[0], globalV[1], globalV[2]);
1512 
1513  //Lest's just assign the depth is 1.0 cm (thickness of a layer), better to update
1514  TEveBox* bklmbox = boxCreator(o, u - o, v - o, du, dv, 1.0);
1515 
1516  bklmbox->SetMainColor(kGreen);
1517  //bklmbox->SetName((std::to_string(hitModule)).c_str());
1518  bklmbox->SetName("BKLMHit2d");
1519 
1520  addToGroup("BKLM2dHits", bklmbox);
1521  addObject(bklm2dhit, bklmbox);
1522 }
1523 
1525 {
1526  const double du = 2.0;
1527  const double dv = 2.0;
1528  ROOT::Math::XYZVector hitPosition = eklm2dhit->getPosition();
1529  TVector3 o(hitPosition.X(), hitPosition.Y(), hitPosition.Z());
1530  TVector3 u(1.0, 0.0, 0.0);
1531  TVector3 v(0.0, 1.0, 0.0);
1532  TEveBox* eklmbox = boxCreator(o, u, v, du, dv, 4.0);
1533  eklmbox->SetMainColor(kGreen);
1534  eklmbox->SetName("EKLMHit2d");
1535  addToGroup("EKLM2dHits", eklmbox);
1536  addObject(eklm2dhit, eklmbox);
1537 }
1538 
1540 {
1541  const VXD::GeoCache& aGeometry = VXD::GeoCache::getInstance();
1542 
1543  VxdID sensorID = roi->getSensorID();
1544  const VXD::SensorInfoBase& aSensorInfo = aGeometry.getSensorInfo(sensorID);
1545 
1546  double minU = aSensorInfo.getUCellPosition(roi->getMinUid(), roi->getMinVid());
1547  double minV = aSensorInfo.getVCellPosition(roi->getMinVid());
1548  double maxU = aSensorInfo.getUCellPosition(roi->getMaxUid(), roi->getMaxVid());
1549  double maxV = aSensorInfo.getVCellPosition(roi->getMaxVid());
1550 
1551 
1552  ROOT::Math::XYZVector localA(minU, minV, 0);
1553  ROOT::Math::XYZVector localB(minU, maxV, 0);
1554  ROOT::Math::XYZVector localC(maxU, minV, 0);
1555 
1556  ROOT::Math::XYZVector globalA = aSensorInfo.pointToGlobal(localA);
1557  ROOT::Math::XYZVector globalB = aSensorInfo.pointToGlobal(localB);
1558  ROOT::Math::XYZVector globalC = aSensorInfo.pointToGlobal(localC);
1559 
1560  TEveBox* ROIbox = boxCreator(XYZToTVector(globalB + globalC) * 0.5, XYZToTVector(globalB - globalA),
1561  XYZToTVector(globalC - globalA), 1, 1, 0.01);
1562 
1563  ROIbox->SetName(ObjectInfo::getIdentifier(roi));
1564  ROIbox->SetMainColor(kSpring - 9);
1565  ROIbox->SetMainTransparency(50);
1566 
1567  addToGroup("ROIs", ROIbox);
1568  addObject(roi, ROIbox);
1569 }
1570 
1571 void EVEVisualization::addRecoHit(const SVDCluster* hit, TEveStraightLineSet* lines)
1572 {
1573  static VXD::GeoCache& geo = VXD::GeoCache::getInstance();
1574  const VXD::SensorInfoBase& sensor = geo.get(hit->getSensorID());
1575 
1576  ROOT::Math::XYZVector a, b;
1577  if (hit->isUCluster()) {
1578  const float u = hit->getPosition();
1579  a = sensor.pointToGlobal(ROOT::Math::XYZVector(sensor.getBackwardWidth() / sensor.getWidth(0) * u, -0.5 * sensor.getLength(), 0.0));
1580  b = sensor.pointToGlobal(ROOT::Math::XYZVector(sensor.getForwardWidth() / sensor.getWidth(0) * u, +0.5 * sensor.getLength(), 0.0));
1581  } else {
1582  const float v = hit->getPosition();
1583  a = sensor.pointToGlobal(ROOT::Math::XYZVector(-0.5 * sensor.getWidth(v), v, 0.0));
1584  b = sensor.pointToGlobal(ROOT::Math::XYZVector(+0.5 * sensor.getWidth(v), v, 0.0));
1585  }
1586 
1587  lines->AddLine(a.X(), a.Y(), a.Z(), b.X(), b.Y(), b.Z());
1588  m_shownRecohits.insert(hit);
1589 }
1590 
1591 void EVEVisualization::addRecoHit(const CDCHit* hit, TEveStraightLineSet* lines)
1592 {
1594  const TVector3& wire_pos_f = cdcgeo.wireForwardPosition(WireID(hit->getID()));
1595  const TVector3& wire_pos_b = cdcgeo.wireBackwardPosition(WireID(hit->getID()));
1596 
1597  lines->AddLine(wire_pos_f.X(), wire_pos_f.Y(), wire_pos_f.Z(), wire_pos_b.X(), wire_pos_b.Y(), wire_pos_b.Z());
1598  m_shownRecohits.insert(hit);
1599 
1600 }
1601 
1602 void EVEVisualization::addCDCHit(const CDCHit* hit, bool showTriggerHits)
1603 {
1605  const TVector3& wire_pos_f = cdcgeo.wireForwardPosition(WireID(hit->getID()));
1606  const TVector3& wire_pos_b = cdcgeo.wireBackwardPosition(WireID(hit->getID()));
1607  static CDC::RealisticTDCCountTranslator tdcTranslator;
1608  TEveGeoShape* cov_shape = new TEveGeoShape("cov_shape");
1609  //TODO: leftrightflag not set! (same for other parameters, unsure which ones should be set)
1610  double driftLength = tdcTranslator.getDriftLength(hit->getTDCCount(), WireID(hit->getID()));
1611  double driftLengthRes = tdcTranslator.getDriftLengthResolution(driftLength, WireID(hit->getID()));
1612  driftLengthRes = std::max(driftLengthRes, 0.005);
1613  const double lengthOfWireSection = 3.0;
1614 
1615  //z in wire direction, x,y orthogonal
1616  const TVector3 zaxis = wire_pos_b - wire_pos_f;
1617  const TVector3 xaxis = zaxis.Orthogonal();
1618  const TVector3 yaxis = xaxis.Cross(zaxis);
1619 
1620  // move to z=0
1621  const TVector3 midPoint = wire_pos_f - zaxis * (wire_pos_f.Z() / zaxis.Z());
1622 
1623  cov_shape->SetShape(new TGeoTube(std::max(0., (double)(driftLength - driftLengthRes)), driftLength + driftLengthRes,
1624  lengthOfWireSection));
1625  fixGeoShapeRefCount(cov_shape);
1626 
1627  TGeoRotation det_rot("det_rot",
1628  xaxis.Theta() * 180 / TMath::Pi(), xaxis.Phi() * 180 / TMath::Pi(),
1629  yaxis.Theta() * 180 / TMath::Pi(), yaxis.Phi() * 180 / TMath::Pi(),
1630  zaxis.Theta() * 180 / TMath::Pi(), zaxis.Phi() * 180 / TMath::Pi()
1631  );
1632 
1633  TGeoCombiTrans det_trans(midPoint(0), midPoint(1), midPoint(2), &det_rot);
1634  cov_shape->SetTransMatrix(det_trans);
1635 
1636  // get relation to trigger track segments
1637  bool isPartOfTS = false;
1638  const auto segments = hit->getRelationsFrom<CDCTriggerSegmentHit>();
1639  if (showTriggerHits && segments.size() > 0) {
1640  isPartOfTS = true;
1641  }
1642 
1643  if (hit->getISuperLayer() % 2 == 0) {
1644  if (isPartOfTS)
1645  cov_shape->SetMainColor(kCyan + 3);
1646  else
1647  cov_shape->SetMainColor(kCyan);
1648  } else {
1649  if (isPartOfTS)
1650  cov_shape->SetMainColor(kPink + 6);
1651  else
1652  cov_shape->SetMainColor(kPink + 7);
1653  }
1654 
1655  cov_shape->SetMainTransparency(50);
1656  cov_shape->SetName(ObjectInfo::getIdentifier(hit));
1657  cov_shape->SetTitle(ObjectInfo::getInfo(hit) + TString::Format("\nWire ID: %d\nADC: %d\nTDC: %d",
1658  hit->getID(), hit->getADCCount(), hit->getTDCCount()));
1659 
1660  addToGroup("CDCHits", cov_shape);
1661  addObject(hit, cov_shape);
1662  if (isPartOfTS) {
1663  addToGroup("CDCTriggerSegmentHits", cov_shape);
1664  for (auto rel : segments.relations()) {
1665  addObject(rel.object, cov_shape);
1666  }
1667  }
1668 }
1669 
1670 void EVEVisualization::addCDCTriggerSegmentHit(const std::string& collectionName, const CDCTriggerSegmentHit* hit)
1671 {
1673  TEveStraightLineSet* shape = new TEveStraightLineSet();
1674 
1675  // get center wire
1676  unsigned iL = WireID(hit->getID()).getICLayer();
1677  if (hit->getPriorityPosition() < 3) iL -= 1;
1678  unsigned nWires = cdcgeo.nWiresInLayer(iL);
1679  unsigned iCenter = hit->getIWire();
1680  if (hit->getPriorityPosition() == 1) iCenter += 1;
1681 
1682  // a track segment consists of 11 wires (15 in SL0) in a special configuration
1683  // -> get the shift with respect to the center wire (*) for all wires
1684  // SL 1-8:
1685  // _ _ _
1686  // |_|_|_|
1687  // |_|_|
1688  // |*|
1689  // |_|_|
1690  // |_|_|_|
1691  std::vector<int> layershift = { -2, -1, 0, 1, 2};
1692  std::vector<std::vector<float>> cellshift = {
1693  { -1, 0, 1},
1694  { -0.5, 0.5},
1695  { 0},
1696  { -0.5, 0.5},
1697  { -1, 0, 1}
1698  };
1699  // SL 0:
1700  // _ _ _ _ _
1701  // |_|_|_|_|_|
1702  // |_|_|_|_|
1703  // |_|_|_|
1704  // |_|_|
1705  // |*|
1706  if (hit->getISuperLayer() == 0) {
1707  layershift = { 0, 1, 2, 3, 4};
1708  cellshift = {
1709  { 0},
1710  { -0.5, 0.5},
1711  { -1, 0, 1},
1712  { -1.5, -0.5, 0.5, 1.5},
1713  { -2, -1, 0, 1, 2}
1714  };
1715  }
1716 
1717  // draw all cells in segment
1718  for (unsigned il = 0; il < layershift.size(); ++il) {
1719  for (unsigned ic = 0; ic < cellshift[il].size(); ++ic) {
1720  TVector3 corners[2][2];
1721  for (unsigned ir = 0; ir < 2; ++ir) {
1722  double r = cdcgeo.fieldWireR(iL + layershift[il] - ir);
1723  double fz = cdcgeo.fieldWireFZ(iL + layershift[il] - ir);
1724  double bz = cdcgeo.fieldWireBZ(iL + layershift[il] - ir);
1725  for (unsigned iphi = 0; iphi < 2; ++iphi) {
1726  double phib = (iCenter + cellshift[il][ic] + iphi - 0.5) * 2 * M_PI / nWires;
1727  double phif = phib + cdcgeo.nShifts(iL + layershift[il]) * M_PI / nWires;
1728 
1729  TVector3 pos_f = TVector3(cos(phif) * r, sin(phif) * r, fz);
1730  TVector3 pos_b = TVector3(cos(phib) * r, sin(phib) * r, bz);
1731  TVector3 zaxis = pos_b - pos_f;
1732  corners[ir][iphi] = pos_f - zaxis * (pos_f.Z() / zaxis.Z());
1733  }
1734  }
1735 
1736  shape->AddLine(corners[0][0].X(), corners[0][0].Y(), 0,
1737  corners[0][1].X(), corners[0][1].Y(), 0);
1738  shape->AddLine(corners[0][1].X(), corners[0][1].Y(), 0,
1739  corners[1][1].X(), corners[1][1].Y(), 0);
1740  shape->AddLine(corners[1][1].X(), corners[1][1].Y(), 0,
1741  corners[1][0].X(), corners[1][0].Y(), 0);
1742  shape->AddLine(corners[1][0].X(), corners[1][0].Y(), 0,
1743  corners[0][0].X(), corners[0][0].Y(), 0);
1744  }
1745  }
1746 
1747  if (hit->getISuperLayer() % 2 == 0) {
1748  shape->SetMainColor(kCyan + 3);
1749  } else {
1750  shape->SetMainColor(kPink + 6);
1751  }
1752 
1753  shape->SetName(ObjectInfo::getIdentifier(hit));
1754  shape->SetTitle(ObjectInfo::getTitle(hit) +
1755  TString::Format("\nPriority: %d\nLeft/Right: %d",
1756  hit->getPriorityPosition(), hit->getLeftRight()));
1757  addToGroup(collectionName, shape);
1758  addObject(hit, shape);
1759 }
1760 
1762 {
1764 
1765  int hitModule = hit->getModule();
1766  float fi = arichGeo->getDetectorPlane().getSlotPhi(hitModule);
1767 
1768  TVector3 centerPos3D = hit->getPosition();
1769 
1770  TVector3 channelX(1, 0, 0); channelX.RotateZ(fi);
1771  TVector3 channelY(0, 1, 0); channelY.RotateZ(fi);
1772 
1773  auto* arichbox = boxCreator(centerPos3D, arichGeo->getMasterVolume().momentumToGlobal(channelX),
1774  arichGeo->getMasterVolume().momentumToGlobal(channelY), 0.49, 0.49, 0.05);
1775  arichbox->SetMainColor(kOrange + 10);
1776  arichbox->SetName((std::to_string(hitModule)).c_str());
1777 
1778  addToGroup("ARICHHits", arichbox);
1779  addObject(hit, arichbox);
1780 }
1781 
1783 {
1784  /* TOP module ID -> #digits */
1785  std::map<int, int> m_topSummary;
1786  for (const TOPDigit& hit : digits) {
1787  int mod = hit.getModuleID();
1788  ++m_topSummary[mod];
1789  }
1790  int maxcount = 0;
1791  for (auto modCountPair : m_topSummary) {
1792  if (modCountPair.second > maxcount)
1793  maxcount = modCountPair.second;
1794  }
1795  for (auto modCountPair : m_topSummary) {
1796  const auto& topmod = TOP::TOPGeometryPar::Instance()->getGeometry()->getModule(modCountPair.first);
1797  double phi = topmod.getPhi();
1798  double r_center = topmod.getRadius();
1799  double z = topmod.getZc();
1800 
1801  TVector3 centerPos3D;
1802  centerPos3D.SetMagThetaPhi(r_center, M_PI / 2, phi);
1803  centerPos3D.SetZ(z);
1804 
1805  TVector3 channelX(1, 0, 0); channelX.RotateZ(phi);
1806  TVector3 channelY(0, 1, 0); channelY.RotateZ(phi);
1807 
1808  //bar is a bit thicker so we can mouse over without getting the geometry
1809  auto* moduleBox = boxCreator(centerPos3D, channelX, channelY,
1810  3.0 * topmod.getBarThickness(), topmod.getBarWidth(), topmod.getBarLength());
1811  moduleBox->SetMainColor(kAzure + 10);
1812  double weight = double(modCountPair.second) / maxcount;
1813  moduleBox->SetMainTransparency(90 - weight * 50);
1814  moduleBox->SetName(("TOP module " + std::to_string(modCountPair.first)).c_str());
1815  moduleBox->SetTitle(TString::Format("#TOPDigits: %d ", modCountPair.second));
1816 
1817  addToGroup("TOP Modules", moduleBox);
1818  //associate all TOPDigits with this module.
1819  for (const TOPDigit& hit : digits) {
1820  if (modCountPair.first == hit.getModuleID())
1821  addObject(&hit, moduleBox);
1822  }
1823  }
1824 }
1825 
1826 
1828 {
1829  for (const auto& labelPair : displayData.m_labels) {
1830  TEveText* text = new TEveText(labelPair.first.c_str());
1831  text->SetName(labelPair.first.c_str());
1832  text->SetTitle(labelPair.first.c_str());
1833  text->SetMainColor(kGray + 1);
1834  const ROOT::Math::XYZVector& p = labelPair.second;
1835  text->PtrMainTrans()->SetPos(p.X(), p.Y(), p.Z());
1836  addToGroup("DisplayData", text);
1837  }
1838 
1839  for (const auto& pointPair : displayData.m_pointSets) {
1840  TEvePointSet* points = new TEvePointSet(pointPair.first.c_str());
1841  points->SetTitle(pointPair.first.c_str());
1842  points->SetMarkerStyle(7);
1843  points->SetMainColor(kGreen);
1844  for (const auto& p : pointPair.second) {
1845  points->SetNextPoint(p.X(), p.Y(), p.Z());
1846  }
1847  addToGroup("DisplayData", points);
1848  }
1849 
1850  int randomColor = 2; //primary colours, changing rapidly with index
1851  for (const auto& arrow : displayData.m_arrows) {
1852  const ROOT::Math::XYZVector pos = arrow.start;
1853  const ROOT::Math::XYZVector dir = arrow.end - pos;
1854  TEveArrow* eveArrow = new TEveArrow(dir.X(), dir.Y(), dir.Z(), pos.X(), pos.Y(), pos.Z());
1855  eveArrow->SetName(arrow.name.c_str());
1856  eveArrow->SetTitle(arrow.name.c_str());
1857  int arrowColor = arrow.color;
1858  if (arrowColor == -1) {
1859  arrowColor = randomColor;
1860  randomColor++;
1861  }
1862  eveArrow->SetMainColor(arrowColor);
1863 
1864  //add label
1865  TEveText* text = new TEveText(arrow.name.c_str());
1866  text->SetMainColor(arrowColor);
1867  //place label in middle of arrow, with some slight offset
1868  // orthogonal direction is arbitrary, set smallest component zero
1869  ROOT::Math::XYZVector orthogonalDir;
1870  if (std::abs(dir.X()) < std::abs(dir.Y())) {
1871  if (std::abs(dir.X()) < std::abs(dir.Z())) {
1872  orthogonalDir.SetCoordinates(0, dir.Z(), -dir.Y());
1873  } else {
1874  orthogonalDir.SetCoordinates(dir.Y(), -dir.X(), 0);
1875  }
1876  } else {
1877  if (std::abs(dir.Y()) < std::abs(dir.Z())) {
1878  orthogonalDir.SetCoordinates(-dir.Z(), 0, dir.X());
1879  } else {
1880  orthogonalDir.SetCoordinates(dir.Y(), -dir.X(), 0);
1881  }
1882  }
1883  const ROOT::Math::XYZVector& labelPos = pos + 0.5 * dir + 0.1 * orthogonalDir;
1884  text->PtrMainTrans()->SetPos(labelPos.X(), labelPos.Y(), labelPos.Z());
1885  eveArrow->AddElement(text);
1886  addToGroup("DisplayData", eveArrow);
1887  }
1888 
1889 }
1890 void EVEVisualization::addObject(const TObject* dataStoreObject, TEveElement* visualRepresentation)
1891 {
1892  VisualRepMap::getInstance()->add(dataStoreObject, visualRepresentation);
1893 }
1894 
1895 void EVEVisualization::addToGroup(const std::string& name, TEveElement* elem)
1896 {
1897  //slashes at beginning and end are ignored
1898  const std::string& groupName = boost::algorithm::trim_copy_if(name, boost::algorithm::is_any_of("/"));
1899 
1900  TEveElementList* group = m_groups[groupName].group;
1901  if (!group) {
1902  group = new TEveElementList(groupName.c_str(), groupName.c_str());
1903  group->SetRnrState(m_groups[groupName].visible);
1904  m_groups[groupName].group = group;
1905 
1906  //if groupName contains '/', remove last bit and add to parent group
1907  //e.g. if groupName=A/B/C, call addToGroup("A/B", groupC)
1908  auto lastSlash = boost::algorithm::find_last(groupName, "/");
1909  if (lastSlash) {
1910  const std::string parentGroup(groupName.begin(), lastSlash.begin());
1911  const std::string thisGroup(lastSlash.end(), groupName.end());
1912  group->SetElementName(thisGroup.c_str());
1913  addToGroup(parentGroup, group);
1914  } else {
1915  gEve->AddElement(group);
1916  }
1917  }
1918  group->AddElement(elem);
1919 }
1920 
Datastore class that holds photon hits. Input to the reconstruction.
Definition: ARICHHit.h:23
DataType Pz() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:441
DataType Mag() const
The magnitude (rho in spherical coordinate system).
Definition: B2Vector3.h:159
DataType Pt() const
The transverse component (R in cylindrical coordinate system).
Definition: B2Vector3.h:198
Class containing the result of the unpacker in raw data and the result of the digitizer in simulation...
Definition: CDCHit.h:40
This class is used to transfer CDC information to the track fit.
Definition: CDCRecoHit.h:32
Example Detector.
Definition: CDCSimHit.h:21
Combination of several CDCHits to a track segment hit for the trigger.
Track created by the CDC trigger.
The Class for CDC Geometry Parameters.
int nShifts(int layerId) const
Returns number shift.
double fieldWireR(int layerId) const
Returns radius of field wire in each layer.
const B2Vector3D wireForwardPosition(uint layerId, int cellId, EWirePosition set=c_Base) const
Returns the forward position of the input sense wire.
const B2Vector3D wireBackwardPosition(uint layerId, int cellId, EWirePosition set=c_Base) const
Returns the backward position of the input sense wire.
unsigned nWiresInLayer(int layerId) const
Returns wire numbers in a layer.
double fieldWireBZ(int layerId) const
Returns backward z position of field wire in each layer.
double fieldWireFZ(int layerId) const
Returns forward z position of field wire in each layer.
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
Translator mirroring the realistic Digitization.
double getDriftLength(unsigned short tdcCount, const WireID &wireID=WireID(), double timeOfFlightEstimator=0, bool leftRight=false, double z=0, double alpha=0, double theta=static_cast< double >(TMath::Pi()/2.), unsigned short adcCount=0) override
Get Drift length.
double getDriftLengthResolution(double driftLength, const WireID &wireID=WireID(), bool leftRight=false, double z=0, double alpha=0, double=static_cast< double >(TMath::Pi()/2.)) override
Get position resolution^2 corresponding to the drift length from getDriftLength of this class.
static const ChargedStable pion
charged pion particle
Definition: Const.h:652
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
Add custom information to the display.
Definition: DisplayData.h:55
std::vector< Arrow > m_arrows
List of arrows.
Definition: DisplayData.h:114
std::map< std::string, std::vector< ROOT::Math::XYZVector > > m_pointSets
name -> points map
Definition: DisplayData.h:108
std::vector< std::pair< std::string, ROOT::Math::XYZVector > > m_labels
text labels (to be shown at a given position).
Definition: DisplayData.h:109
ECL cluster data.
Definition: ECLCluster.h:27
@ c_nPhotons
CR is split into n photons (N1)
void clearEvent()
clear event data.
void setOptions(const std::string &opts)
Set the display options.
static constexpr double c_minPCut
don't show MCParticles with momentum below this cutoff.
void addTrackCandidateImproved(const std::string &collectionName, const RecoTrack &recoTrack)
Add a RecoTrack, but use stored genfit track representation to make visualisation objects.
EveVisBField * m_bfield
The global magnetic field.
TEveStraightLineSet * m_unassignedRecoHits
Unassigned recohits.
static const int c_recoHitColor
Color for reco hits.
void addCDCHit(const CDCHit *hit, bool showTriggerHits=false)
show CDCHits directly.
std::map< const MCParticle *, MCTrack > m_mcparticleTracks
map MCParticles to MCTrack (so hits can be added to the correct track).
void addBKLMHit2d(const KLMHit2d *bklm2dhit)
Add a reconstructed 2d hit in the BKLM.
void setErrScale(double errScale=1.)
Set the scaling factor for the visualization of track hit errors.
bool m_drawRefTrack
Draw reference track in addTrack.
void addVertex(const genfit::GFRaveVertex *vertex)
Add a vertex point and its covariance matrix.
bool m_drawForward
draw forward in addTrack
TEveCalo3D * m_calo3d
Object for the energy bar visualisation.
bool m_drawCardinalRep
Draw cardinal representation in addTrack.
void addCDCTriggerTrack(const std::string &collectionName, const CDCTriggerTrack &track)
Add a CDCTriggerTrack.
void addCDCTriggerSegmentHit(const std::string &collectionName, const CDCTriggerSegmentHit *hit)
show outline of track segments.
TEveTrackPropagator * m_trackpropagator
Track propagator for MCParticles.
std::map< std::string, ElementGroup > m_groups
name -> grouping element.
void addSimHit(const CDCSimHit *hit, const MCParticle *particle)
Add a CDCSimHit.
void addECLCluster(const ECLCluster *cluster)
Add a reconstructed cluster in the ECL.
TEveTrackList * m_tracklist
parent object for MC tracks.
void addToGroup(const std::string &name, TEveElement *elem)
Add 'elem' to the element group 'name' (created if necessary).
void makeTracks()
Create visual representation of all tracks.
void addRecoHit(const SomeVXDHit *hit, TEveStraightLineSet *lines)
adds given VXD hit to lines.
void addTrack(const Belle2::Track *belle2Track)
Add this genfit::Track to event data.
bool m_unassignedRecoHitsVisibility
is m_unassignedRecoHits visible?
static const int c_trackMarkerColor
Color for track markers.
TEveTrackPropagator * m_gftrackpropagator
Track propagator for genfit::Tracks (different mainly because of drawing options)
std::string m_options
Option string for genfit::Track visualisation.
static const int c_trackColor
Color for tracks.
static const int c_unassignedHitColor
Color for unassigned (reco)hits.
bool m_hideSecondaries
If true, secondary MCParticles (and hits created by them) will not be shown.
void addEKLMHit2d(const KLMHit2d *eklm2dhit)
Add a reconstructed 2d hit in the EKLM.
TEveBox * boxCreator(const TVector3 &o, TVector3 u, TVector3 v, float ud, float vd, float depth)
Create a box around o, oriented along u and v with widths ud, vd and depth and return a pointer to th...
void addARICHHit(const ARICHHit *hit)
Add recontructed hit in ARICH.
void showUserData(const DisplayData &displayData)
Add user-defined data (labels, points, etc.)
std::set< const TObject * > m_shownRecohits
List of shown recohits (PXDCluster, SVDCluster, CDCHit).
bool m_drawBackward
draw backward in addTrack
TEveTrackPropagator * m_consttrackpropagator
Track propagator for CDCTriggerTracks (uses constant B field)
bool m_assignToPrimaries
If true, hits created by secondary particles (e.g.
bool m_drawErrors
Draw errors in addTrack.
double m_errorScale
Rescale PXD/SVD errors with this factor to ensure visibility.
void addTOPDigits(const StoreArray< TOPDigit > &digits)
Add TOPDigits (shown aggregated per module).
void addROI(const ROIid *roi)
Add a Region Of Interest, computed by the PXDDataReduction module.
void addKLMCluster(const KLMCluster *cluster)
Add a reconstructed cluster in the KLM.
static const int c_klmClusterColor
Color for KLMCluster objects.
MCTrack * addMCParticle(const MCParticle *particle)
Return MCTrack for given particle, add it if it doesn't exist yet.
TEveCaloDataVec * m_eclData
ECL cluster data.
static const int c_recoTrackColor
Color for TrackCandidates.
void addTrackCandidate(const std::string &collectionName, const RecoTrack &recoTrack)
Add a RecoTrack, to evaluate track finding.
void addObject(const TObject *dataStoreObject, TEveElement *visualRepresentation)
Generic function to keep track of which objects have which visual representation.
void makeLines(TEveTrack *eveTrack, const genfit::StateOnPlane *prevState, const genfit::StateOnPlane *state, const genfit::AbsTrackRep *rep, TEvePathMark::EType_e markType, bool drawErrors, int markerPos=1)
Create hit visualisation for the given options, and add them to 'eveTrack'.
Provide magnetic field values for TEveTrackPropagator.
Definition: EveVisBField.h:19
KLM cluster data.
Definition: KLMCluster.h:28
KLM 2d hit.
Definition: KLMHit2d.h:33
int getLayer() const
Get layer number.
Definition: KLMHit2d.h:132
int getZStripMax() const
Get last strip number for z plane.
Definition: KLMHit2d.h:202
int getSection() const
Get section number.
Definition: KLMHit2d.h:96
float getPositionZ() const
Get hit global position z coordinate.
Definition: KLMHit2d.h:312
int getSector() const
Get sector number.
Definition: KLMHit2d.h:114
float getPositionX() const
Set hit global position.
Definition: KLMHit2d.h:294
int getPhiStripMin() const
Get strip number for phi plane.
Definition: KLMHit2d.h:218
int getZStripMin() const
Get strip number for z plane.
Definition: KLMHit2d.h:194
ROOT::Math::XYZVector getPosition() const
Get hit global position.
Definition: KLMHit2d.h:321
int getPhiStripMax() const
Get last strip number for phi plane.
Definition: KLMHit2d.h:226
float getPositionY() const
Get hit global position y coordinate.
Definition: KLMHit2d.h:303
KLM simulation hit.
Definition: KLMSimHit.h:31
Class to save the full simulated trajectory of a particle.
const MCTrajectoryPoint & back() const
return reference to the last point
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
@ c_PrimaryParticle
bit 0: Particle is primary particle.
Definition: MCParticle.h:47
@ c_StoppedInDetector
bit 3: Particle was stopped in the detector (the simulation volume).
Definition: MCParticle.h:53
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
Definition: PXDCluster.h:30
Class PXDSimHit - Geant4 simulated hit for the PXD.
Definition: PXDSimHit.h:24
ROIid stores the U and V ids and the sensor id of the Region Of Interest.
Definition: ROIid.h:25
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
std::vector< Belle2::RecoTrack::UsedSVDHit * > getSVDHitList() const
Return an unsorted list of svd hits.
Definition: RecoTrack.h:452
ROOT::Math::XYZVector getPositionSeed() const
Return the position seed stored in the reco track. ATTENTION: This is not the fitted position.
Definition: RecoTrack.h:480
std::vector< Belle2::RecoTrack::UsedPXDHit * > getPXDHitList() const
Return an unsorted list of pxd hits.
Definition: RecoTrack.h:449
std::vector< Belle2::RecoTrack::UsedCDCHit * > getCDCHitList() const
Return an unsorted list of cdc hits.
Definition: RecoTrack.h:455
const genfit::TrackPoint * getCreatedTrackPoint(const RecoHitInformation *recoHitInformation) const
Get a pointer to the TrackPoint that was created from this hit.
Definition: RecoTrack.cc:230
std::vector< RecoHitInformation * > getRecoHitInformations(bool getSorted=false) const
Return a list of all RecoHitInformations associated with the RecoTrack.
Definition: RecoTrack.cc:557
short int getChargeSeed() const
Return the charge seed stored in the reco track. ATTENTION: This is not the fitted charge.
Definition: RecoTrack.h:508
ROOT::Math::XYZVector getMomentumSeed() const
Return the momentum seed stored in the reco track. ATTENTION: This is not the fitted momentum.
Definition: RecoTrack.h:487
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
The SVD Cluster class This class stores all information about reconstructed SVD clusters.
Definition: SVDCluster.h:29
SVDRecoHit - an extended form of SVDHit containing geometry information.
Definition: SVDRecoHit.h:47
bool isU() const
Is the coordinate u or v?
Definition: SVDRecoHit.h:91
VxdID getSensorID() const
Get the compact ID.
Definition: SVDRecoHit.h:82
Class SVDSimHit - Geant4 simulated hit for the SVD.
Definition: SVDSimHit.h:26
Class to store TOP digitized hits (output of TOPDigitizer or raw data unpacker) relations to TOPSimHi...
Definition: TOPDigit.h:24
const TOPGeometry * getGeometry() const
Returns pointer to geometry object using basf2 units.
static TOPGeometryPar * Instance()
Static method to obtain the pointer to its instance.
Values of the result of a track fit with a given particle hypothesis.
Helix getHelix() const
Conversion to framework Helix (without covariance).
short getChargeSign() const
Return track charge (1 or -1).
double getPValue() const
Getter for Chi2 Probability of the track fit.
ROOT::Math::XYZVector getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
ROOT::Math::XYZVector getPosition() const
Getter for vector of position at closest approach of track in r/phi projection.
Class that bundles various TrackFitResults.
Definition: Track.h:25
const TrackFitResult * getTrackFitResultWithClosestMass(const Const::ChargedStable &requestedType) const
Return the track fit for a fit hypothesis with the closest mass.
Definition: Track.cc:104
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:67
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
static const SensorInfoBase & get(Belle2::VxdID id)
Return a reference to the SensorInfo of a given SensorID.
Definition: GeoCache.h:139
Base class to provide Sensor Information for PXD and SVD.
double getVCellPosition(int vID) const
Return the position of a specific strip/pixel in v direction.
double getUCellPosition(int uID, int vID=-1) const
Return the position of a specific strip/pixel in u direction.
ROOT::Math::XYZVector pointToGlobal(const ROOT::Math::XYZVector &local, bool reco=false) const
Convert a point from local to global coordinates.
void addCluster(const TObject *dataStoreObject, TEveCaloData *caloData, int towerID)
Selection inside TEveCalo* is complicated, use this to keep track of ECL clusters.
static VisualRepMap * getInstance()
get instance pointer.
Definition: VisualRepMap.cc:36
void clear()
Remove all contents in map.
Definition: VisualRepMap.cc:46
void add(const TObject *dataStoreObject, TEveElement *visualRepresentation)
Generic function to keep track of which objects have which visual representation.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
Class to identify a wire inside the CDC.
Definition: WireID.h:34
unsigned short getICLayer() const
Getter for continuous layer numbering.
Definition: WireID.cc:24
Provides BKLM geometry parameters for simulation, reconstruction etc (from Gearbox or DataBase)
Definition: GeometryPar.h:37
const Module * findModule(int section, int sector, int layer) const
Get the pointer to the definition of a module.
Definition: GeometryPar.cc:721
static GeometryPar * instance(void)
Static method to get a reference to the singleton GeometryPar instance.
Definition: GeometryPar.cc:27
Define the geometry of a BKLM module Each sector [octant] contains Modules.
Definition: Module.h:76
This class collects all information needed and produced by a specific AbsFitter and is specific to on...
Definition: AbsFitterInfo.h:42
Contains the measurement and covariance in raw detector coordinates.
Abstract base class for a track representation.
Definition: AbsTrackRep.h:66
void getPosDir(const StateOnPlane &state, TVector3 &pos, TVector3 &dir) const
Get cartesian position and direction vector of a state.
Definition: AbsTrackRep.h:252
virtual double extrapolateToPlane(StateOnPlane &state, const genfit::SharedPlanePtr &plane, bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0
Extrapolates the state to plane, and returns the extrapolation length and, via reference,...
virtual TVector3 getPos(const StateOnPlane &state) const =0
Get the cartesian position of a state.
virtual void getPosMomCov(const MeasuredStateOnPlane &state, TVector3 &pos, TVector3 &mom, TMatrixDSym &cov) const =0
Translates MeasuredStateOnPlane into 3D position, momentum and 6x6 covariance.
Detector plane.
Definition: DetPlane.h:59
Exception class for error handling in GENFIT (provides storage for diagnostic information)
Definition: Exception.h:48
Class where important numbers and properties of a fit can be stored.
Definition: FitStatus.h:80
bool isTrackPruned() const
Has the track been pruned after the fit?
Definition: FitStatus.h:116
GFRaveVertex class.
Definition: GFRaveVertex.h:48
Collects information needed and produced by a GblFitter/GBL and is specific to one AbsTrackRep of the...
Definition: GblFitterInfo.h:52
const MeasuredStateOnPlane & getFittedState(bool afterKink=true) const override
Get the prediction at this point Always biased in GBL (global fit) There are 2 states,...
MeasurementOnPlane getMeasurement() const
Get the measurement on plane from stored measurement data (from last construction/update)
StateOnPlane getReferenceState() const
Returns (copy of) the stored reference 5D state at current plane with internal rep.
Collects information needed and produced by a AbsKalmanFitter implementations and is specific to one ...
#StateOnPlane with additional covariance matrix.
Measured coordinates on a plane.
Measurement class implementing a planar hit geometry (1 or 2D).
Class for measurements implementing a space point hit geometry.
A state with arbitrary dimension defined in a DetPlane.
Definition: StateOnPlane.h:47
Object containing AbsMeasurement and AbsFitterInfo objects.
Definition: TrackPoint.h:46
AbsFitterInfo * getFitterInfo(const AbsTrackRep *rep=nullptr) const
Get fitterInfo for rep. Per default, use cardinal rep.
Definition: TrackPoint.cc:170
Class for measurements in wire detectors (Straw tubes and drift chambers) which do not measure the co...
Class for measurements in wire detectors (Straw tubes and drift chambers) which do not measure the co...
Class for measurements in wire detectors (Straw tubes and drift chambers) which can measure the coord...
static constexpr auto XYZToTVector
Helper function to convert XYZVector to TVector3.
Definition: VectorUtil.h:24
B2Vector3< double > B2Vector3D
typedef for common usage with double
Definition: B2Vector3.h:516
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
const TOPGeoModule & getModule(int moduleID) const
Returns module.
Definition: TOPGeometry.cc:42
double eval(const std::vector< double > &spl, const std::vector< double > &vals, double x)
Evaluate spline (zero order or first order) in point x.
Definition: tools.h:115
double getMaxR()
find a point that is inside the top node.
Definition: EveGeometry.cc:184
TString getIdentifier(const TObject *obj)
Where is this object in the datastore?
Definition: ObjectInfo.cc:105
TString getInfo(const TObject *obj)
Get object info HTML (e.g.
Definition: ObjectInfo.cc:55
TString getTitle(const TObject *obj)
Get plain text for TEve object titles (shown on mouse-over).
Definition: ObjectInfo.cc:68
Implements a colour palette, see http://sobac.com/sobac/tangocolors.htm.
Definition: ColorPalette.h:23
int getTColorID(const std::string &tangoName, int tangoId=1)
Get TColor ID for given name in tango colour palette.
Definition: ColorPalette.cc:39
Abstract base class for different kinds of events.
Defines for I/O streams used for error and debug printing.
std::shared_ptr< genfit::DetPlane > SharedPlanePtr
Shared Pointer to a DetPlane.
Hold MC tracks and associated visualisation objects.
TEvePointSet * simhits
simhit positions.
const MCParticle * parentParticle
parent particle, or nullptr.
TEveTrack * track
the actual MC track.
Small struct to encode a position/momentum without additional overhead.