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