2 #include "EventDisplay.h"
11 #include "AbsMeasurement.h"
12 #include "FullMeasurement.h"
13 #include "PlanarMeasurement.h"
14 #include "ProlateSpacepointMeasurement.h"
15 #include "SpacepointMeasurement.h"
16 #include "WireMeasurement.h"
17 #include "WirePointMeasurement.h"
18 #include "AbsTrackRep.h"
19 #include "ConstField.h"
21 #include "Exception.h"
22 #include "FieldManager.h"
24 #include "KalmanFitterInfo.h"
25 #include "KalmanFitter.h"
27 #include "KalmanFitterRefTrack.h"
28 #include "RKTrackRep.h"
30 #include <TApplication.h>
31 #include <TEveBrowser.h>
32 #include <TEveManager.h>
33 #include <TEveEventManager.h>
34 #include <TEveGeoNode.h>
35 #include <TEveGeoShape.h>
36 #include <TEveStraightLineSet.h>
37 #include <TEveTriangleSet.h>
38 #include <TDecompSVD.h>
41 #include <TGNumberEntry.h>
43 #include <TGeoManager.h>
44 #include <TGeoMatrix.h>
46 #include <TGeoSphere.h>
50 #include <TMatrixTSym.h>
51 #include <TMatrixDSymEigen.h>
62 EventDisplay* EventDisplay::eventDisplay_ =
nullptr;
64 EventDisplay::EventDisplay() :
71 drawTrackMarkers_(true),
79 drawCardinalRep_(true),
85 fitterId_(SimpleKalman),
87 squareRootFormalism_(false),
97 if((!gApplication) || (gApplication && gApplication->TestBit(TApplication::kDefaultApplication))) {
98 std::cout <<
"In EventDisplay ctor: gApplication not found, creating..." << std::flush;
99 new TApplication(
"ROOT_application", 0, 0);
100 std::cout <<
"done!" << std::endl;
103 std::cout <<
"In EventDisplay ctor: gEve not found, creating..." << std::flush;
104 TEveManager::Create();
105 std::cout <<
"done!" << std::endl;
112 void EventDisplay::setOptions(std::string opts) {
115 for(
size_t i = 0; i < opts.length(); ++i) {
116 if(opts[i] ==
'A') drawAutoScale_ =
true;
117 if(opts[i] ==
'B') drawBackward_ =
true;
118 if(opts[i] ==
'D') drawDetectors_ =
true;
119 if(opts[i] ==
'E') drawErrors_ =
true;
120 if(opts[i] ==
'F') drawForward_ =
true;
121 if(opts[i] ==
'H') drawHits_ =
true;
122 if(opts[i] ==
'M') drawTrackMarkers_ =
true;
123 if(opts[i] ==
'P') drawPlanes_ =
true;
124 if(opts[i] ==
'S') drawScaleMan_ =
true;
125 if(opts[i] ==
'T') drawTrack_ =
true;
126 if(opts[i] ==
'X') drawSilent_ =
true;
127 if(opts[i] ==
'G') drawGeometry_ =
true;
133 void EventDisplay::setErrScale(
double errScale) { errorScale_ = errScale; }
135 double EventDisplay::getErrScale() {
return errorScale_; }
139 if(eventDisplay_ ==
nullptr) {
142 return eventDisplay_;
146 EventDisplay::~EventDisplay() { reset(); }
148 void EventDisplay::reset() {
150 for(
unsigned int i = 0; i < events_.size(); i++) {
152 for(
unsigned int j = 0; j < events_[i]->size(); j++) {
154 delete events_[i]->at(j);
164 void EventDisplay::addEvent(std::vector<Track*>& tracks) {
166 std::vector<Track*>* vec =
new std::vector<Track*>;
168 for(
unsigned int i = 0; i < tracks.size(); i++) {
169 vec->push_back(
new Track(*(tracks[i])));
172 events_.push_back(vec);
176 void EventDisplay::addEvent(std::vector<const Track*>& tracks) {
178 std::vector<Track*>* vec =
new std::vector<Track*>;
180 for(
unsigned int i = 0; i < tracks.size(); i++) {
181 vec->push_back(
new Track(*(tracks[i])));
184 events_.push_back(vec);
188 void EventDisplay::addEvent(
const Track* tr) {
190 std::vector<Track*>* vec =
new std::vector<Track*>;
191 vec->push_back(
new Track(*tr));
192 events_.push_back(vec);
196 void EventDisplay::next(
unsigned int stp) {
198 gotoEvent(eventId_ + stp);
202 void EventDisplay::prev(
unsigned int stp) {
204 if(events_.size() == 0)
return;
208 gotoEvent(eventId_ - stp);
213 int EventDisplay::getNEvents() {
return events_.size(); }
216 void EventDisplay::gotoEvent(
unsigned int id) {
218 if (events_.size() == 0)
220 else if(
id >= events_.size())
221 id = events_.size() - 1;
223 bool resetCam =
true;
230 std::cout <<
"At event " <<
id << std::endl;
231 if (gEve->GetCurrentEvent()) {
232 gEve->GetCurrentEvent()->DestroyElements();
234 double old_error_scale = errorScale_;
235 drawEvent(eventId_, resetCam);
236 if(old_error_scale != errorScale_) {
237 if (gEve->GetCurrentEvent()) {
238 gEve->GetCurrentEvent()->DestroyElements();
240 drawEvent(eventId_, resetCam);
242 errorScale_ = old_error_scale;
246 void EventDisplay::open() {
248 std::cout <<
"EventDisplay::open(); " << getNEvents() <<
" events loaded" << std::endl;
250 if(getNEvents() > 0) {
251 double old_error_scale = errorScale_;
253 if(old_error_scale != errorScale_) {
254 std::cout <<
"autoscaling changed the error, draw again." << std::endl;
257 errorScale_ = old_error_scale;
263 gApplication->Run(kTRUE);
266 std::cout <<
"opened" << std::endl;
271 void EventDisplay::drawEvent(
unsigned int id,
bool resetCam) {
273 std::cout <<
"EventDisplay::drawEvent(" <<
id <<
")" << std::endl;
278 TGeoNode* top_node = gGeoManager->GetTopNode();
279 assert(top_node !=
nullptr);
282 TObjArray* volumes = gGeoManager->GetListOfVolumes();
283 for(
int i = 0; i < volumes->GetEntriesFast(); i++) {
284 TGeoVolume* volume =
dynamic_cast<TGeoVolume*
>(volumes->At(i));
285 assert(volume !=
nullptr);
286 volume->SetLineColor(12);
287 volume->SetTransparency(50);
290 TEveGeoTopNode* eve_top_node =
new TEveGeoTopNode(gGeoManager, top_node);
291 eve_top_node->IncDenyDestroy();
292 gEve->AddElement(eve_top_node);
296 for(
unsigned int i = 0; i < events_.at(
id)->size(); i++) {
298 if (!drawAllTracks_ && trackId_ != i)
301 Track* track = events_[id]->at(i);
303 track->checkConsistency();
305 std::cerr<< e.getExcString() <<std::endl;
309 std::unique_ptr<Track> refittedTrack(
nullptr);
312 std::cout <<
"Refit track:" << std::endl;
314 std::unique_ptr<AbsKalmanFitter> fitter;
318 fitter->setMultipleMeasurementHandling(mmHandling_);
319 (
static_cast<KalmanFitter*
>(fitter.get()))->useSquareRootFormalism(squareRootFormalism_);
324 fitter->setMultipleMeasurementHandling(mmHandling_);
329 fitter.reset(
new DAF(
false));
330 (
static_cast<KalmanFitter*
>( (
static_cast<DAF*
>(fitter.get()))->getKalman() ) )->useSquareRootFormalism(squareRootFormalism_);
333 fitter.reset(
new DAF());
334 (
static_cast<KalmanFitterRefTrack*
>( (
static_cast<DAF*
>(fitter.get()))->getKalman() ) )->setDeltaChi2Ref(dChi2Ref_);
338 fitter->setDebugLvl(std::max(0, (
int)debugLvl_-1));
339 fitter->setMinIterations(nMinIter_);
340 fitter->setMaxIterations(nMaxIter_);
341 fitter->setRelChi2Change(dRelChi2_);
342 fitter->setMaxFailedHits(nMaxFailed_);
345 refittedTrack.reset(
new Track(*track));
346 refittedTrack->deleteFitterInfo();
349 refittedTrack->Print(
"C");
351 timeval startcputime, endcputime;
354 gettimeofday(&startcputime,
nullptr);
355 fitter->processTrack(refittedTrack.get(), resort_);
356 gettimeofday(&endcputime,
nullptr);
359 std::cerr << e.what();
360 std::cerr <<
"Exception, could not refit track" << std::endl;
364 int microseconds = 1000000*(endcputime.tv_sec - startcputime.tv_sec) + (endcputime.tv_usec - startcputime.tv_usec);
365 std::cout <<
"it took " << double(microseconds) / 1000 <<
" ms of CPU to fit the track\n";
368 refittedTrack->checkConsistency();
370 std::cerr<< e.getExcString() <<std::endl;
374 track = refittedTrack.get();
382 if (drawCardinalRep_) {
383 rep = track->getCardinalRep();
384 std::cout <<
"Draw cardinal rep" << std::endl;
387 if (repId_ >= track->getNumReps())
388 repId_ = track->getNumReps() - 1;
389 rep = track->getTrackRep(repId_);
390 std::cout <<
"Draw rep" << repId_ << std::endl;
394 std::cout <<
"track " << i << std::endl;
397 track->getFitStatus(rep)->Print();
399 if (track->getFitStatus(rep)->isFitted()) {
401 std::cout <<
"fitted state: \n";
402 track->getFittedState().Print();
405 std::cerr << e.what();
414 unsigned int numhits = track->getNumPointsWithMeasurement();
421 for(
unsigned int j = 0; j < numhits; j++) {
423 fittedState =
nullptr;
425 TrackPoint* tp = track->getPointWithMeasurement(j);
426 if (! tp->hasRawMeasurements()) {
427 std::cerr<<
"trackPoint has no raw measurements"<<std::endl;
432 int hit_coords_dim = m->getDim();
435 if (tp->getNumRawMeasurements() > 1) {
436 bool sameTypes(
true);
437 for (
unsigned int iM=1; iM<tp->getNumRawMeasurements(); ++iM) {
438 auto& rawMeasurement = *(tp->getRawMeasurement(iM));
439 if (
typeid(rawMeasurement) !=
typeid(*m))
443 std::cerr<<
"cannot draw trackpoint containing multiple Measurements of differend types"<<std::endl;
451 if (! tp->hasFitterInfo(rep)) {
452 std::cerr<<
"trackPoint has no fitterInfo for rep"<<std::endl;
460 std::cerr<<
"can only display KalmanFitterInfo"<<std::endl;
463 if (! fi->hasPredictionsAndUpdates()) {
464 std::cerr<<
"KalmanFitterInfo does not have all predictions and updates"<<std::endl;
469 fittedState = &(fi->getFittedState(
true));
472 std::cerr << e.what();
473 std::cerr<<
"can not get fitted state"<<std::endl;
474 fittedState =
nullptr;
476 prevFittedState = fittedState;
481 if (fittedState ==
nullptr) {
482 if (fi->hasForwardUpdate()) {
483 fittedState = fi->getForwardUpdate();
485 else if (fi->hasBackwardUpdate()) {
486 fittedState = fi->getBackwardUpdate();
488 else if (fi->hasForwardPrediction()) {
489 fittedState = fi->getForwardPrediction();
491 else if (fi->hasBackwardPrediction()) {
492 fittedState = fi->getBackwardPrediction();
496 if (fittedState ==
nullptr) {
497 std::cout <<
"cannot get any state from fitterInfo, continue.\n";
499 prevFittedState = fittedState;
503 TVector3 track_pos = fittedState->getPos();
504 double charge = fittedState->getCharge();
512 bool planar_pixel_hit = planar_hit && hit_coords_dim == 2;
514 bool wire_hit = m && m->isLeftRightMeasurement();
516 if (!full_hit && !planar_hit && !planar_pixel_hit && !space_hit && !wire_hit && !wirepoint_hit) {
517 std::cout <<
"Track " << i <<
", Hit " << j <<
": Unknown measurement type: skipping hit!" << std::endl;
523 unsigned int nMeas = fi->getNumMeasurements();
524 for (
unsigned int iMeas = 0; iMeas < nMeas; ++iMeas) {
526 if (iMeas > 0 && wire_hit)
530 const TVectorT<double>& hit_coords = mop->getState();
531 const TMatrixTSym<double>& hit_cov = mop->getCov();
536 TVector3 o = fittedState->getPlane()->getO();
537 TVector3 u = fittedState->getPlane()->getU();
538 TVector3 v = fittedState->getPlane()->getV();
542 double_t plane_size = 4;
543 TVector2 stripDir(1,0);
546 if(!planar_pixel_hit) {
547 if (
dynamic_cast<RKTrackRep*
>(rep) !=
nullptr) {
548 const TMatrixD& H = mop->getHMatrix()->
getMatrix();
549 stripDir.Set(H(0,3), H(0,4));
551 hit_u = hit_coords(0);
553 hit_u = hit_coords(0);
554 hit_v = hit_coords(1);
556 }
else if (wire_hit) {
557 hit_u = fabs(hit_coords(0));
558 hit_v = v*(track_pos-o);
560 hit_v = hit_coords(1);
564 if(plane_size < 4) plane_size = 4;
569 (drawPlanes_ || (drawDetectors_ && planar_hit))) {
570 TVector3 move(0,0,0);
571 if (planar_hit) move = track_pos-o;
572 if (wire_hit) move = v*(v*(track_pos-o));
573 TEveBox* box = boxCreator(o + move, u, v, plane_size, plane_size, 0.01);
574 if (drawDetectors_ && planar_hit) {
575 box->SetMainColor(kCyan);
577 box->SetMainColor(kGray);
579 box->SetMainTransparency(50);
580 gEve->AddElement(box);
589 update.extrapolateBy(-3.);
590 makeLines(&update, fi->getBackwardUpdate(), rep, kMagenta, 1, drawTrackMarkers_, drawErrors_, 1);
593 if (j > 0 && prevFi !=
nullptr) {
595 makeLines(prevFittedState, fittedState, rep, charge > 0 ? kRed : kBlue, 1, drawTrackMarkers_, drawErrors_, 3);
597 makeLines(prevFittedState, fittedState, rep, charge > 0 ? kRed : kBlue, 1,
false, drawErrors_, 0, 0);
601 makeLines(prevFi->getForwardUpdate(), fi->getForwardPrediction(), rep, kCyan, 1, drawTrackMarkers_, drawErrors_, 1, 0);
602 if (j == numhits-1) {
604 update.extrapolateBy(3.);
605 makeLines(fi->getForwardUpdate(), &update, rep, kCyan, 1, drawTrackMarkers_, drawErrors_, 1, 0);
609 makeLines(prevFi->getBackwardPrediction(), fi->getBackwardUpdate(), rep, kMagenta, 1, drawTrackMarkers_, drawErrors_, 1);
612 if(drawRefTrack_ && fi->hasReferenceState() && prevFi->hasReferenceState())
613 makeLines(prevFi->getReferenceState(), fi->getReferenceState(), rep, charge > 0 ? kRed + 2 : kBlue + 2, 2, drawTrackMarkers_,
false, 3);
615 else if (j > 0 && prevFi ==
nullptr) {
616 std::cout <<
"previous FitterInfo == nullptr \n";
620 std::cerr <<
"extrapolation failed, cannot draw track" << std::endl;
621 std::cerr << e.what();
628 TEveGeoShape* det_shape =
new TEveGeoShape(
"det_shape");
629 det_shape->IncDenyDestroy();
630 det_shape->SetShape(
new TGeoTube(std::max(0., (
double)(hit_u-0.0105/2.)), hit_u+0.0105/2., plane_size));
632 TVector3 norm = u.Cross(v);
633 TGeoRotation det_rot(
"det_rot", (u.Theta()*180)/TMath::Pi(), (u.Phi()*180)/TMath::Pi(),
634 (norm.Theta()*180)/TMath::Pi(), (norm.Phi()*180)/TMath::Pi(),
635 (v.Theta()*180)/TMath::Pi(), (v.Phi()*180)/TMath::Pi());
636 TVector3 move = v*(v*(track_pos-o));
637 TGeoCombiTrans det_trans(o(0) + move.X(),
641 det_shape->SetTransMatrix(det_trans);
642 det_shape->SetMainColor(kCyan);
643 det_shape->SetMainTransparency(25);
644 if((drawHits_ && (hit_u+0.0105/2 > 0)) || !drawHits_) {
645 gEve->AddElement(det_shape);
660 sop.getCov()*=errorScale_;
663 prevSop.extrapolateBy(-3);
664 makeLines(&sop, &prevSop, rep, kYellow, 1,
false,
true, 0, 0);
667 prevSop.extrapolateBy(3);
668 makeLines(&sop, &prevSop, rep, kYellow, 1,
false,
true, 0, 0);
672 if(!planar_pixel_hit) {
674 TVector3 stripDir3 = stripDir.X()*u + stripDir.Y()*v;
675 TVector3 stripDir3perp = stripDir.Y()*u - stripDir.X()*v;
676 TVector3 move = stripDir3perp*(stripDir3perp*(track_pos-o));
677 hit_box = boxCreator((o + move + hit_u*stripDir3), stripDir3, stripDir3perp, errorScale_*
std::sqrt(hit_cov(0,0)), plane_size, 0.0105);
678 hit_box->SetMainColor(kYellow);
679 hit_box->SetMainTransparency(0);
680 gEve->AddElement(hit_box);
683 TMatrixDSymEigen eigen_values(hit_cov);
684 TEveGeoShape* cov_shape =
new TEveGeoShape(
"cov_shape");
685 cov_shape->IncDenyDestroy();
686 TVectorT<double> ev = eigen_values.GetEigenValues();
687 TMatrixT<double> eVec = eigen_values.GetEigenVectors();
688 double pseudo_res_0 = errorScale_*
std::sqrt(ev(0));
689 double pseudo_res_1 = errorScale_*
std::sqrt(ev(1));
694 double min_cov = std::min(pseudo_res_0,pseudo_res_1);
696 std::cout <<
"Track " << i <<
", Hit " << j <<
": Invalid covariance matrix (Eigenvalue < 1e-5), autoscaling not possible!" << std::endl;
698 if(min_cov < 0.049) {
699 double cor = 0.05 / min_cov;
700 std::cout <<
"Track " << i <<
", Hit " << j <<
": Pixel covariance too small, rescaling by " << cor;
704 std::cout <<
" to " << errorScale_ << std::endl;
711 cov_shape->SetShape(
new TGeoEltu(pseudo_res_0, pseudo_res_1, 0.0105));
712 TVector3 pix_pos = o + hit_u*u + hit_v*v;
713 TVector3 u_semiaxis = (pix_pos + eVec(0,0)*u + eVec(1,0)*v)-pix_pos;
714 TVector3 v_semiaxis = (pix_pos + eVec(0,1)*u + eVec(1,1)*v)-pix_pos;
715 TVector3 norm = u.Cross(v);
719 TGeoRotation det_rot(
"det_rot", (u_semiaxis.Theta()*180)/TMath::Pi(), (u_semiaxis.Phi()*180)/TMath::Pi(),
720 (v_semiaxis.Theta()*180)/TMath::Pi(), (v_semiaxis.Phi()*180)/TMath::Pi(),
721 (norm.Theta()*180)/TMath::Pi(), (norm.Phi()*180)/TMath::Pi());
722 TGeoCombiTrans det_trans(pix_pos(0),pix_pos(1),pix_pos(2), &det_rot);
723 cov_shape->SetTransMatrix(det_trans);
726 cov_shape->SetMainColor(kYellow);
727 cov_shape->SetMainTransparency(0);
728 gEve->AddElement(cov_shape);
737 TMatrixDSymEigen eigen_values(m->getRawHitCov());
738 TEveGeoShape* cov_shape =
new TEveGeoShape(
"cov_shape");
739 cov_shape->IncDenyDestroy();
740 cov_shape->SetShape(
new TGeoSphere(0.,1.));
741 TVectorT<double> ev = eigen_values.GetEigenValues();
742 TMatrixT<double> eVec = eigen_values.GetEigenVectors();
743 TVector3 eVec1(eVec(0,0),eVec(1,0),eVec(2,0));
744 TVector3 eVec2(eVec(0,1),eVec(1,1),eVec(2,1));
745 TVector3 eVec3(eVec(0,2),eVec(1,2),eVec(2,2));
746 const TVector3 norm = u.Cross(v);
749 static const double radDeg(180./TMath::Pi());
750 TGeoRotation det_rot(
"det_rot", eVec1.Theta()*radDeg, eVec1.Phi()*radDeg,
751 eVec2.Theta()*radDeg, eVec2.Phi()*radDeg,
752 eVec3.Theta()*radDeg, eVec3.Phi()*radDeg);
754 if (! det_rot.IsValid()){
756 if (fabs(eVec2*eVec3) > 1.e-10)
757 eVec3 = eVec1.Cross(eVec2);
759 det_rot.SetAngles(eVec1.Theta()*radDeg, eVec1.Phi()*radDeg,
760 eVec2.Theta()*radDeg, eVec2.Phi()*radDeg,
761 eVec3.Theta()*radDeg, eVec3.Phi()*radDeg);
765 double pseudo_res_0 = errorScale_*
std::sqrt(ev(0));
766 double pseudo_res_1 = errorScale_*
std::sqrt(ev(1));
767 double pseudo_res_2 = errorScale_*
std::sqrt(ev(2));
769 pseudo_res_0 = errorScale_*0.5;
770 pseudo_res_1 = errorScale_*0.5;
771 pseudo_res_2 = errorScale_*0.5;
777 double min_cov = std::min(pseudo_res_0,std::min(pseudo_res_1,pseudo_res_2));
779 std::cout <<
"Track " << i <<
", Hit " << j <<
": Invalid covariance matrix (Eigenvalue < 1e-5), autoscaling not possible!" << std::endl;
781 if(min_cov <= 0.149) {
782 double cor = 0.15 / min_cov;
783 std::cout <<
"Track " << i <<
", Hit " << j <<
": Space hit covariance too small, rescaling by " << cor;
788 std::cout <<
" to " << errorScale_ << std::endl;
796 TGeoGenTrans det_trans(o(0),o(1),o(2),
799 pseudo_res_0, pseudo_res_1, pseudo_res_2,
801 cov_shape->SetTransMatrix(det_trans);
804 cov_shape->SetMainColor(kYellow);
805 cov_shape->SetMainTransparency(10);
806 gEve->AddElement(cov_shape);
812 TMatrixDSymEigen eigen_values(hit_cov);
813 TEveGeoShape* cov_shape =
new TEveGeoShape(
"cov_shape");
814 cov_shape->IncDenyDestroy();
815 TVectorT<double> ev = eigen_values.GetEigenValues();
816 TMatrixT<double> eVec = eigen_values.GetEigenVectors();
817 double pseudo_res_0 = errorScale_*
std::sqrt(ev(0));
818 double pseudo_res_1 = errorScale_*
std::sqrt(ev(1));
823 double min_cov = std::min(pseudo_res_0,pseudo_res_1);
825 std::cout <<
"Track " << i <<
", Hit " << j <<
": Invalid covariance matrix (Eigenvalue < 1e-5), autoscaling not possible!" << std::endl;
827 if(min_cov < 0.049) {
828 double cor = 0.05 / min_cov;
829 std::cout <<
"Track " << i <<
", Hit " << j <<
": Pixel covariance too small, rescaling by " << cor;
833 std::cout <<
" to " << errorScale_ << std::endl;
840 cov_shape->SetShape(
new TGeoEltu(pseudo_res_0, pseudo_res_1, 0.0105));
841 TVector3 pix_pos = o + hit_u*u + hit_v*v;
842 TVector3 u_semiaxis = (pix_pos + eVec(0,0)*u + eVec(1,0)*v)-pix_pos;
843 TVector3 v_semiaxis = (pix_pos + eVec(0,1)*u + eVec(1,1)*v)-pix_pos;
844 TVector3 norm = u.Cross(v);
848 static const double radDeg(180./TMath::Pi());
849 TGeoRotation det_rot(
"det_rot", u_semiaxis.Theta()*radDeg, u_semiaxis.Phi()*radDeg,
850 v_semiaxis.Theta()*radDeg, v_semiaxis.Phi()*radDeg,
851 norm.Theta()*radDeg, norm.Phi()*radDeg);
857 TGeoCombiTrans det_trans(pix_pos(0),pix_pos(1),pix_pos(2), &det_rot);
858 cov_shape->SetTransMatrix(det_trans);
861 cov_shape->SetMainColor(kYellow);
862 cov_shape->SetMainTransparency(0);
863 gEve->AddElement(cov_shape);
870 TEveGeoShape* cov_shape =
new TEveGeoShape(
"cov_shape");
871 cov_shape->IncDenyDestroy();
872 double pseudo_res_0 = errorScale_*
std::sqrt(hit_cov(0,0));
873 double pseudo_res_1 = plane_size;
874 if (wirepoint_hit) pseudo_res_1 = errorScale_*
std::sqrt(hit_cov(1,1));
878 if(pseudo_res_0 < 1e-5) {
879 std::cout <<
"Track " << i <<
", Hit " << j <<
": Invalid wire resolution (< 1e-5), autoscaling not possible!" << std::endl;
881 if(pseudo_res_0 < 0.0049) {
882 double cor = 0.005 / pseudo_res_0;
883 std::cout <<
"Track " << i <<
", Hit " << j <<
": Wire covariance too small, rescaling by " << cor;
886 std::cout <<
" to " << errorScale_ << std::endl;
890 if(wirepoint_hit && pseudo_res_1 < 1e-5) {
891 std::cout <<
"Track " << i <<
", Hit " << j <<
": Invalid wire resolution (< 1e-5), autoscaling not possible!" << std::endl;
893 if(pseudo_res_1 < 0.0049) {
894 double cor = 0.005 / pseudo_res_1;
895 std::cout <<
"Track " << i <<
", Hit " << j <<
": Wire covariance too small, rescaling by " << cor;
898 std::cout <<
" to " << errorScale_ << std::endl;
905 TVector3 move = v*(v*(track_pos-o));
906 hit_box = boxCreator((o + move + hit_u*u), u, v, errorScale_*
std::sqrt(hit_cov(0,0)), pseudo_res_1, 0.0105);
907 hit_box->SetMainColor(kYellow);
908 hit_box->SetMainTransparency(0);
909 gEve->AddElement(hit_box);
911 hit_box = boxCreator((o + move - hit_u*u), u, v, errorScale_*
std::sqrt(hit_cov(0,0)), pseudo_res_1, 0.0105);
912 hit_box->SetMainColor(kYellow);
913 hit_box->SetMainTransparency(0);
914 gEve->AddElement(hit_box);
924 prevFittedState = fittedState;
930 gEve->Redraw3D(resetCam);
937 TEveBox* EventDisplay::boxCreator(TVector3 o, TVector3 u, TVector3 v,
float ud,
float vd,
float depth) {
939 TEveBox* box =
new TEveBox(
"detPlane_shape");
942 TVector3 norm = u.Cross(v);
947 vertices[0] = o(0) - u(0) - v(0) - norm(0);
948 vertices[1] = o(1) - u(1) - v(1) - norm(1);
949 vertices[2] = o(2) - u(2) - v(2) - norm(2);
951 vertices[3] = o(0) + u(0) - v(0) - norm(0);
952 vertices[4] = o(1) + u(1) - v(1) - norm(1);
953 vertices[5] = o(2) + u(2) - v(2) - norm(2);
955 vertices[6] = o(0) + u(0) - v(0) + norm(0);
956 vertices[7] = o(1) + u(1) - v(1) + norm(1);
957 vertices[8] = o(2) + u(2) - v(2) + norm(2);
959 vertices[9] = o(0) - u(0) - v(0) + norm(0);
960 vertices[10] = o(1) - u(1) - v(1) + norm(1);
961 vertices[11] = o(2) - u(2) - v(2) + norm(2);
963 vertices[12] = o(0) - u(0) + v(0) - norm(0);
964 vertices[13] = o(1) - u(1) + v(1) - norm(1);
965 vertices[14] = o(2) - u(2) + v(2) - norm(2);
967 vertices[15] = o(0) + u(0) + v(0) - norm(0);
968 vertices[16] = o(1) + u(1) + v(1) - norm(1);
969 vertices[17] = o(2) + u(2) + v(2) - norm(2);
971 vertices[18] = o(0) + u(0) + v(0) + norm(0);
972 vertices[19] = o(1) + u(1) + v(1) + norm(1);
973 vertices[20] = o(2) + u(2) + v(2) + norm(2);
975 vertices[21] = o(0) - u(0) + v(0) + norm(0);
976 vertices[22] = o(1) - u(1) + v(1) + norm(1);
977 vertices[23] = o(2) - u(2) + v(2) + norm(2);
980 for(
int k = 0; k < 24; k += 3) box->SetVertex((k/3), vertices[k], vertices[k+1], vertices[k+2]);
988 const Color_t& color,
const Style_t& style,
bool drawMarkers,
bool drawErrors,
double lineWidth,
int markerPos)
990 if (prevState ==
nullptr || state ==
nullptr) {
991 std::cerr <<
"prevState == nullptr || state == nullptr\n";
995 TVector3 pos, dir, oldPos, oldDir;
997 rep->
getPosDir(*prevState, oldPos, oldDir);
999 double distA = (pos-oldPos).Mag();
1000 double distB = distA;
1001 if ((pos-oldPos)*oldDir < 0)
1003 if ((pos-oldPos)*dir < 0)
1005 TVector3 intermediate1 = oldPos + 0.3 * distA * oldDir;
1006 TVector3 intermediate2 = pos - 0.3 * distB * dir;
1007 TEveStraightLineSet* lineSet =
new TEveStraightLineSet;
1008 lineSet->AddLine(oldPos(0), oldPos(1), oldPos(2), intermediate1(0), intermediate1(1), intermediate1(2));
1009 lineSet->AddLine(intermediate1(0), intermediate1(1), intermediate1(2), intermediate2(0), intermediate2(1), intermediate2(2));
1010 lineSet->AddLine(intermediate2(0), intermediate2(1), intermediate2(2), pos(0), pos(1), pos(2));
1011 lineSet->SetLineColor(color);
1012 lineSet->SetLineStyle(style);
1013 lineSet->SetLineWidth(lineWidth);
1016 lineSet->AddMarker(oldPos(0), oldPos(1), oldPos(2));
1018 lineSet->AddMarker(pos(0), pos(1), pos(2));
1022 gEve->AddElement(lineSet);
1026 const MeasuredStateOnPlane* measuredState;
1028 measuredState =
dynamic_cast<const MeasuredStateOnPlane*
>(prevState);
1030 measuredState =
dynamic_cast<const MeasuredStateOnPlane*
>(state);
1032 if (measuredState !=
nullptr) {
1036 if (markerPos == 0) {
1037 if (fabs(distA) < 1.) {
1038 distA < 0 ? distA = -1 : distA = 1;
1040 eval = 0.2 * distA * oldDir;
1043 if (fabs(distB) < 1.) {
1044 distB < 0 ? distB = -1 : distB = 1;
1046 eval = -0.2 * distB * dir;
1052 TVector3 position, direction;
1053 rep->
getPosMomCov(*measuredState, position, direction, cov);
1056 TMatrixDSymEigen eigen_values(cov.GetSub(0,2, 0,2));
1057 TVectorT<double> ev = eigen_values.GetEigenValues();
1058 TMatrixT<double> eVec = eigen_values.GetEigenVectors();
1059 TVector3 eVec1, eVec2;
1061 static const double maxErr = 1000.;
1062 double ev0 = std::min(ev(0), maxErr);
1063 double ev1 = std::min(ev(1), maxErr);
1064 double ev2 = std::min(ev(2), maxErr);
1067 if (ev0 < ev1 && ev0 < ev2) {
1068 eVec1.SetXYZ(eVec(0,1),eVec(1,1),eVec(2,1));
1070 eVec2.SetXYZ(eVec(0,2),eVec(1,2),eVec(2,2));
1073 else if (ev1 < ev0 && ev1 < ev2) {
1074 eVec1.SetXYZ(eVec(0,0),eVec(1,0),eVec(2,0));
1076 eVec2.SetXYZ(eVec(0,2),eVec(1,2),eVec(2,2));
1080 eVec1.SetXYZ(eVec(0,0),eVec(1,0),eVec(2,0));
1082 eVec2.SetXYZ(eVec(0,1),eVec(1,1),eVec(2,1));
1086 if (eVec1.Cross(eVec2)*
eval < 0)
1090 const TVector3 oldEVec1(eVec1);
1091 const TVector3 oldEVec2(eVec2);
1093 const int nEdges = 24;
1094 std::vector<TVector3> vertices;
1096 vertices.push_back(position);
1099 for (
int i=0; i<nEdges; ++i) {
1100 const double angle = 2*TMath::Pi()/nEdges * i;
1101 vertices.push_back(position + cos(angle)*eVec1 + sin(angle)*eVec2);
1106 DetPlane* newPlane =
new DetPlane(*(measuredState->getPlane()));
1107 newPlane->setO(position +
eval);
1109 MeasuredStateOnPlane stateCopy(*measuredState);
1113 catch(Exception& e){
1114 std::cerr<<e.what();
1119 rep->
getPosMomCov(stateCopy, position, direction, cov);
1122 TMatrixDSymEigen eigen_values2(cov.GetSub(0,2, 0,2));
1123 ev = eigen_values2.GetEigenValues();
1124 eVec = eigen_values2.GetEigenVectors();
1126 ev0 = std::min(ev(0), maxErr);
1127 ev1 = std::min(ev(1), maxErr);
1128 ev2 = std::min(ev(2), maxErr);
1131 if (ev0 < ev1 && ev0 < ev2) {
1132 eVec1.SetXYZ(eVec(0,1),eVec(1,1),eVec(2,1));
1134 eVec2.SetXYZ(eVec(0,2),eVec(1,2),eVec(2,2));
1137 else if (ev1 < ev0 && ev1 < ev2) {
1138 eVec1.SetXYZ(eVec(0,0),eVec(1,0),eVec(2,0));
1140 eVec2.SetXYZ(eVec(0,2),eVec(1,2),eVec(2,2));
1144 eVec1.SetXYZ(eVec(0,0),eVec(1,0),eVec(2,0));
1146 eVec2.SetXYZ(eVec(0,1),eVec(1,1),eVec(2,1));
1150 if (eVec1.Cross(eVec2)*
eval < 0)
1154 if (oldEVec1*eVec1 < 0) {
1160 double angle0 = eVec1.Angle(oldEVec1);
1161 if (eVec1*(
eval.Cross(oldEVec1)) < 0)
1163 for (
int i=0; i<nEdges; ++i) {
1164 const double angle = 2*TMath::Pi()/nEdges * i - angle0;
1165 vertices.push_back(position + cos(angle)*eVec1 + sin(angle)*eVec2);
1168 vertices.push_back(position);
1171 TEveTriangleSet* error_shape =
new TEveTriangleSet(vertices.size(), nEdges*2);
1172 for(
unsigned int k = 0; k < vertices.size(); ++k) {
1173 error_shape->SetVertex(k, vertices[k].X(), vertices[k].Y(), vertices[k].Z());
1176 assert(vertices.size() == 2*nEdges+2);
1179 for (
int i=0; i<nEdges; ++i) {
1181 error_shape->SetTriangle(iTri++, i+1, i+1+nEdges, (i+1)%nEdges+1);
1182 error_shape->SetTriangle(iTri++, (i+1)%nEdges+1, i+1+nEdges, (i+1)%nEdges+1+nEdges);
1188 error_shape->SetMainColor(color);
1189 error_shape->SetMainTransparency(25);
1190 gEve->AddElement(error_shape);
1196 void EventDisplay::makeGui() {
1198 TEveBrowser* browser = gEve->GetBrowser();
1199 browser->StartEmbedding(TRootBrowser::kLeft);
1201 TGMainFrame* frmMain =
new TGMainFrame(gClient->GetRoot(), 1000, 600);
1202 frmMain->SetWindowName(
"XX GUI");
1203 frmMain->SetCleanup(kDeepCleanup);
1206 TGTextButton* tb = 0;
1209 TGHorizontalFrame* hf =
new TGHorizontalFrame(frmMain); {
1211 lbl =
new TGLabel(hf,
"Go to event: ");
1213 guiEvent =
new TGNumberEntry(hf, 0, 9,999, TGNumberFormat::kNESInteger,
1214 TGNumberFormat::kNEANonNegative,
1215 TGNumberFormat::kNELLimitMinMax,
1217 hf->AddFrame(guiEvent);
1218 guiEvent->Connect(
"ValueSet(Long_t)",
"genfit::EventDisplay", fh,
"guiGoto()");
1221 tb =
new TGTextButton(hf,
"Redraw Event");
1223 tb->Connect(
"Clicked()",
"genfit::EventDisplay", fh,
"guiGoto()");
1225 frmMain->AddFrame(hf);
1228 hf =
new TGHorizontalFrame(frmMain); {
1229 lbl =
new TGLabel(hf,
"\n Draw Options");
1232 frmMain->AddFrame(hf);
1234 hf =
new TGHorizontalFrame(frmMain); {
1235 guiDrawGeometry_ =
new TGCheckButton(hf,
"Draw geometry");
1236 if(drawGeometry_) guiDrawGeometry_->Toggle();
1237 hf->AddFrame(guiDrawGeometry_);
1238 guiDrawGeometry_->Connect(
"Toggled(Bool_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1240 frmMain->AddFrame(hf);
1242 hf =
new TGHorizontalFrame(frmMain); {
1243 guiDrawDetectors_ =
new TGCheckButton(hf,
"Draw detectors");
1244 if(drawDetectors_) guiDrawDetectors_->Toggle();
1245 hf->AddFrame(guiDrawDetectors_);
1246 guiDrawDetectors_->Connect(
"Toggled(Bool_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1248 frmMain->AddFrame(hf);
1250 hf =
new TGHorizontalFrame(frmMain); {
1251 guiDrawHits_ =
new TGCheckButton(hf,
"Draw hits");
1252 if(drawHits_) guiDrawHits_->Toggle();
1253 hf->AddFrame(guiDrawHits_);
1254 guiDrawHits_->Connect(
"Toggled(Bool_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1256 frmMain->AddFrame(hf);
1260 hf =
new TGHorizontalFrame(frmMain); {
1261 guiDrawPlanes_ =
new TGCheckButton(hf,
"Draw planes");
1262 if(drawPlanes_) guiDrawPlanes_->Toggle();
1263 hf->AddFrame(guiDrawPlanes_);
1264 guiDrawPlanes_->Connect(
"Toggled(Bool_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1266 frmMain->AddFrame(hf);
1268 hf =
new TGHorizontalFrame(frmMain); {
1269 guiDrawTrackMarkers_ =
new TGCheckButton(hf,
"Draw track markers");
1270 if(drawTrackMarkers_) guiDrawTrackMarkers_->Toggle();
1271 hf->AddFrame(guiDrawTrackMarkers_);
1272 guiDrawTrackMarkers_->Connect(
"Toggled(Bool_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1274 frmMain->AddFrame(hf);
1277 hf =
new TGHorizontalFrame(frmMain); {
1278 guiDrawTrack_ =
new TGCheckButton(hf,
"Draw track");
1279 if(drawTrack_) guiDrawTrack_->Toggle();
1280 hf->AddFrame(guiDrawTrack_);
1281 guiDrawTrack_->Connect(
"Toggled(Bool_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1283 frmMain->AddFrame(hf);
1285 hf =
new TGHorizontalFrame(frmMain); {
1286 guiDrawRefTrack_ =
new TGCheckButton(hf,
"Draw reference track");
1287 if(drawRefTrack_) guiDrawRefTrack_->Toggle();
1288 hf->AddFrame(guiDrawRefTrack_);
1289 guiDrawRefTrack_->Connect(
"Toggled(Bool_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1291 frmMain->AddFrame(hf);
1293 hf =
new TGHorizontalFrame(frmMain); {
1294 guiDrawErrors_ =
new TGCheckButton(hf,
"Draw track errors");
1295 if(drawErrors_) guiDrawErrors_->Toggle();
1296 hf->AddFrame(guiDrawErrors_);
1297 guiDrawErrors_->Connect(
"Toggled(Bool_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1299 frmMain->AddFrame(hf);
1301 hf =
new TGHorizontalFrame(frmMain); {
1302 guiDrawForward_ =
new TGCheckButton(hf,
"Draw forward fit");
1303 if(drawForward_) guiDrawForward_->Toggle();
1304 hf->AddFrame(guiDrawForward_);
1305 guiDrawForward_->Connect(
"Toggled(Bool_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1307 frmMain->AddFrame(hf);
1309 hf =
new TGHorizontalFrame(frmMain); {
1310 guiDrawBackward_ =
new TGCheckButton(hf,
"Draw backward fit");
1311 if(drawBackward_) guiDrawBackward_->Toggle();
1312 hf->AddFrame(guiDrawBackward_);
1313 guiDrawBackward_->Connect(
"Toggled(Bool_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1315 frmMain->AddFrame(hf);
1318 hf =
new TGHorizontalFrame(frmMain); {
1319 guiDrawAutoScale_ =
new TGCheckButton(hf,
"Auto-scale errors");
1320 if(drawAutoScale_) guiDrawAutoScale_->Toggle();
1321 hf->AddFrame(guiDrawAutoScale_);
1322 guiDrawAutoScale_->Connect(
"Toggled(Bool_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1324 frmMain->AddFrame(hf);
1326 hf =
new TGHorizontalFrame(frmMain); {
1327 guiDrawScaleMan_ =
new TGCheckButton(hf,
"Manually scale errors");
1328 if(drawScaleMan_) guiDrawScaleMan_->Toggle();
1329 hf->AddFrame(guiDrawScaleMan_);
1330 guiDrawScaleMan_->Connect(
"Toggled(Bool_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1332 frmMain->AddFrame(hf);
1334 hf =
new TGHorizontalFrame(frmMain); {
1335 guiErrorScale_ =
new TGNumberEntry(hf, errorScale_, 6,999, TGNumberFormat::kNESReal,
1336 TGNumberFormat::kNEANonNegative,
1337 TGNumberFormat::kNELLimitMinMax,
1339 hf->AddFrame(guiErrorScale_);
1340 guiErrorScale_->Connect(
"ValueSet(Long_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1341 lbl =
new TGLabel(hf,
"Error scale");
1344 frmMain->AddFrame(hf);
1348 hf =
new TGHorizontalFrame(frmMain); {
1349 lbl =
new TGLabel(hf,
"\n TrackRep options");
1352 frmMain->AddFrame(hf);
1354 hf =
new TGHorizontalFrame(frmMain); {
1355 guiDrawCardinalRep_ =
new TGCheckButton(hf,
"Draw cardinal rep");
1356 if(drawCardinalRep_) guiDrawCardinalRep_->Toggle();
1357 hf->AddFrame(guiDrawCardinalRep_);
1358 guiDrawCardinalRep_->Connect(
"Toggled(Bool_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1360 frmMain->AddFrame(hf);
1362 hf =
new TGHorizontalFrame(frmMain); {
1363 guiRepId_ =
new TGNumberEntry(hf, repId_, 6,999, TGNumberFormat::kNESInteger,
1364 TGNumberFormat::kNEANonNegative,
1365 TGNumberFormat::kNELLimitMinMax,
1367 hf->AddFrame(guiRepId_);
1368 guiRepId_->Connect(
"ValueSet(Long_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1369 lbl =
new TGLabel(hf,
"Else draw rep with id");
1372 frmMain->AddFrame(hf);
1374 hf =
new TGHorizontalFrame(frmMain); {
1375 guiDrawAllTracks_ =
new TGCheckButton(hf,
"Draw all tracks");
1376 if(drawAllTracks_) guiDrawAllTracks_->Toggle();
1377 hf->AddFrame(guiDrawAllTracks_);
1378 guiDrawAllTracks_->Connect(
"Toggled(Bool_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1380 frmMain->AddFrame(hf);
1382 hf =
new TGHorizontalFrame(frmMain); {
1383 guiTrackId_ =
new TGNumberEntry(hf, trackId_, 6,999, TGNumberFormat::kNESInteger,
1384 TGNumberFormat::kNEANonNegative,
1385 TGNumberFormat::kNELLimitMinMax,
1387 hf->AddFrame(guiTrackId_);
1388 guiTrackId_->Connect(
"ValueSet(Long_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1389 lbl =
new TGLabel(hf,
"Else draw track nr. ");
1392 frmMain->AddFrame(hf);
1396 frmMain->MapSubwindows();
1398 frmMain->MapWindow();
1400 browser->StopEmbedding();
1401 browser->SetTabTitle(
"Draw Control", 0);
1404 browser->StartEmbedding(TRootBrowser::kLeft);
1405 TGMainFrame* frmMain2 =
new TGMainFrame(gClient->GetRoot(), 1000, 600);
1406 frmMain2->SetWindowName(
"XX GUI");
1407 frmMain2->SetCleanup(kDeepCleanup);
1409 hf =
new TGHorizontalFrame(frmMain2); {
1411 lbl =
new TGLabel(hf,
"Go to event: ");
1413 guiEvent2 =
new TGNumberEntry(hf, 0, 9,999, TGNumberFormat::kNESInteger,
1414 TGNumberFormat::kNEANonNegative,
1415 TGNumberFormat::kNELLimitMinMax,
1417 hf->AddFrame(guiEvent2);
1418 guiEvent2->Connect(
"ValueSet(Long_t)",
"genfit::EventDisplay", fh,
"guiGoto2()");
1421 tb =
new TGTextButton(hf,
"Redraw Event");
1423 tb->Connect(
"Clicked()",
"genfit::EventDisplay", fh,
"guiGoto()");
1425 frmMain2->AddFrame(hf);
1427 hf =
new TGHorizontalFrame(frmMain2); {
1428 lbl =
new TGLabel(hf,
"\n Fitting options");
1431 frmMain2->AddFrame(hf);
1433 hf =
new TGHorizontalFrame(frmMain2); {
1434 guiRefit_ =
new TGCheckButton(hf,
"Refit");
1435 if(refit_) guiRefit_->Toggle();
1436 hf->AddFrame(guiRefit_);
1437 guiRefit_->Connect(
"Toggled(Bool_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1439 frmMain2->AddFrame(hf);
1441 hf =
new TGHorizontalFrame(frmMain2); {
1442 guiDebugLvl_ =
new TGNumberEntry(hf, debugLvl_, 6,999, TGNumberFormat::kNESInteger,
1443 TGNumberFormat::kNEANonNegative,
1444 TGNumberFormat::kNELLimitMinMax,
1446 hf->AddFrame(guiDebugLvl_);
1447 guiDebugLvl_->Connect(
"ValueSet(Long_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1448 lbl =
new TGLabel(hf,
"debug level");
1451 frmMain2->AddFrame(hf);
1453 hf =
new TGHorizontalFrame(frmMain2); {
1454 guiFitterId_ =
new TGButtonGroup(hf,
"Fitter type:");
1455 guiFitterId_->Connect(
"Clicked(Int_t)",
"genfit::EventDisplay", fh,
"guiSelectFitterId(int)");
1456 hf->AddFrame(guiFitterId_,
new TGLayoutHints(kLHintsTop));
1457 TGRadioButton* fitterId_button =
new TGRadioButton(guiFitterId_,
"Simple Kalman");
1458 new TGRadioButton(guiFitterId_,
"Reference Kalman");
1459 new TGRadioButton(guiFitterId_,
"DAF w/ simple Kalman");
1460 new TGRadioButton(guiFitterId_,
"DAF w/ reference Kalman");
1461 fitterId_button->SetDown(
true,
false);
1462 guiFitterId_->Show();
1464 frmMain2->AddFrame(hf);
1466 hf =
new TGHorizontalFrame(frmMain2); {
1467 guiMmHandling_ =
new TGButtonGroup(hf,
"Multiple measurement handling in Kalman:");
1468 guiMmHandling_->Connect(
"Clicked(Int_t)",
"genfit::EventDisplay", fh,
"guiSelectMmHandling(int)");
1469 hf->AddFrame(guiMmHandling_,
new TGLayoutHints(kLHintsTop));
1470 TGRadioButton* mmHandling_button =
new TGRadioButton(guiMmHandling_,
"weighted average");
1471 new TGRadioButton(guiMmHandling_,
"unweighted average");
1472 new TGRadioButton(guiMmHandling_,
"weighted, closest to reference");
1473 new TGRadioButton(guiMmHandling_,
"unweighted, closest to reference");
1474 new TGRadioButton(guiMmHandling_,
"weighted, closest to prediction");
1475 new TGRadioButton(guiMmHandling_,
"unweighted, closest to prediction");
1476 new TGRadioButton(guiMmHandling_,
"weighted, closest to reference for WireMeasurements, weighted average else");
1477 new TGRadioButton(guiMmHandling_,
"unweighted, closest to reference for WireMeasurements, unweighted average else");
1478 new TGRadioButton(guiMmHandling_,
"weighted, closest to prediction for WireMeasurements, weighted average else");
1479 new TGRadioButton(guiMmHandling_,
"unweighted, closest to prediction for WireMeasurements, unweighted average else");
1480 mmHandling_button->SetDown(
true,
false);
1481 guiMmHandling_->Show();
1483 frmMain2->AddFrame(hf);
1485 hf =
new TGHorizontalFrame(frmMain2); {
1486 guiSquareRootFormalism_ =
new TGCheckButton(hf,
"Use square root formalism (simple Kalman/simple DAF)");
1487 if(squareRootFormalism_) guiSquareRootFormalism_->Toggle();
1488 hf->AddFrame(guiSquareRootFormalism_);
1489 guiSquareRootFormalism_->Connect(
"Toggled(Bool_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1491 frmMain2->AddFrame(hf);
1493 hf =
new TGHorizontalFrame(frmMain2); {
1494 guiDPVal_ =
new TGNumberEntry(hf, dPVal_, 6,9999, TGNumberFormat::kNESReal,
1495 TGNumberFormat::kNEANonNegative,
1496 TGNumberFormat::kNELLimitMinMax,
1498 hf->AddFrame(guiDPVal_);
1499 guiDPVal_->Connect(
"ValueSet(Long_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1500 lbl =
new TGLabel(hf,
"delta pVal (convergence criterium)");
1503 frmMain2->AddFrame(hf);
1505 hf =
new TGHorizontalFrame(frmMain2); {
1506 guiRelChi2_ =
new TGNumberEntry(hf, dRelChi2_, 6,9999, TGNumberFormat::kNESReal,
1507 TGNumberFormat::kNEANonNegative,
1508 TGNumberFormat::kNELLimitMinMax,
1510 hf->AddFrame(guiRelChi2_);
1511 guiRelChi2_->Connect(
"ValueSet(Long_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1512 lbl =
new TGLabel(hf,
"rel chi^2 change (non-convergence criterium)");
1515 frmMain2->AddFrame(hf);
1517 hf =
new TGHorizontalFrame(frmMain2); {
1518 guiDChi2Ref_ =
new TGNumberEntry(hf, dChi2Ref_, 6,9999, TGNumberFormat::kNESReal,
1519 TGNumberFormat::kNEANonNegative,
1520 TGNumberFormat::kNELLimitMinMax,
1522 hf->AddFrame(guiDChi2Ref_);
1523 guiDChi2Ref_->Connect(
"ValueSet(Long_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1524 lbl =
new TGLabel(hf,
"min chi^2 change for re-calculating reference track (Ref Kalman)");
1527 frmMain2->AddFrame(hf);
1529 hf =
new TGHorizontalFrame(frmMain2); {
1530 guiNMinIter_ =
new TGNumberEntry(hf, nMinIter_, 6,999, TGNumberFormat::kNESInteger,
1531 TGNumberFormat::kNEANonNegative,
1532 TGNumberFormat::kNELLimitMinMax,
1534 hf->AddFrame(guiNMinIter_);
1535 guiNMinIter_->Connect(
"ValueSet(Long_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1536 lbl =
new TGLabel(hf,
"Minimum nr of iterations");
1539 frmMain2->AddFrame(hf);
1541 hf =
new TGHorizontalFrame(frmMain2); {
1542 guiNMaxIter_ =
new TGNumberEntry(hf, nMaxIter_, 6,999, TGNumberFormat::kNESInteger,
1543 TGNumberFormat::kNEANonNegative,
1544 TGNumberFormat::kNELLimitMinMax,
1546 hf->AddFrame(guiNMaxIter_);
1547 guiNMaxIter_->Connect(
"ValueSet(Long_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1548 lbl =
new TGLabel(hf,
"Maximum nr of iterations");
1551 frmMain2->AddFrame(hf);
1553 hf =
new TGHorizontalFrame(frmMain2); {
1554 guiNMaxFailed_ =
new TGNumberEntry(hf, nMaxFailed_, 6,999, TGNumberFormat::kNESInteger,
1555 TGNumberFormat::kNEAAnyNumber,
1556 TGNumberFormat::kNELLimitMinMax,
1558 hf->AddFrame(guiNMaxFailed_);
1559 guiNMaxFailed_->Connect(
"ValueSet(Long_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1560 lbl =
new TGLabel(hf,
"Maximum nr of failed hits");
1563 frmMain2->AddFrame(hf);
1566 hf =
new TGHorizontalFrame(frmMain2); {
1567 guiResort_ =
new TGCheckButton(hf,
"Resort track");
1568 if(resort_) guiResort_->Toggle();
1569 hf->AddFrame(guiResort_);
1570 guiResort_->Connect(
"Toggled(Bool_t)",
"genfit::EventDisplay", fh,
"guiSetDrawParams()");
1572 frmMain2->AddFrame(hf);
1577 frmMain2->MapSubwindows();
1579 frmMain2->MapWindow();
1581 browser->StopEmbedding();
1582 browser->SetTabTitle(
"Refit Control", 0);
1586 void EventDisplay::guiGoto(){
1587 Long_t n = guiEvent->GetNumberEntry()->GetIntNumber();
1588 guiEvent2->SetIntNumber(n);
1592 void EventDisplay::guiGoto2(){
1593 Long_t n = guiEvent2->GetNumberEntry()->GetIntNumber();
1594 guiEvent->SetIntNumber(n);
1599 void EventDisplay::guiSetDrawParams(){
1601 drawGeometry_ = guiDrawGeometry_->IsOn();
1602 drawDetectors_ = guiDrawDetectors_->IsOn();
1603 drawHits_ = guiDrawHits_->IsOn();
1604 drawErrors_ = guiDrawErrors_->IsOn();
1606 drawPlanes_ = guiDrawPlanes_->IsOn();
1607 drawTrackMarkers_ = guiDrawTrackMarkers_->IsOn();
1608 drawTrack_ = guiDrawTrack_->IsOn();
1609 drawRefTrack_ = guiDrawRefTrack_->IsOn();
1610 drawForward_ = guiDrawForward_->IsOn();
1611 drawBackward_ = guiDrawBackward_->IsOn();
1613 drawAutoScale_ = guiDrawAutoScale_->IsOn();
1614 drawScaleMan_ = guiDrawScaleMan_->IsOn();
1616 errorScale_ = guiErrorScale_->GetNumberEntry()->GetNumber();
1618 drawCardinalRep_ = guiDrawCardinalRep_->IsOn();
1619 repId_ = guiRepId_->GetNumberEntry()->GetNumber();
1621 drawAllTracks_ = guiDrawAllTracks_->IsOn();
1622 trackId_ = guiTrackId_->GetNumberEntry()->GetNumber();
1625 refit_ = guiRefit_->IsOn();
1626 debugLvl_ = guiDebugLvl_->GetNumberEntry()->GetNumber();
1628 squareRootFormalism_ = guiSquareRootFormalism_->IsOn();
1629 dPVal_ = guiDPVal_->GetNumberEntry()->GetNumber();
1630 dRelChi2_ = guiRelChi2_->GetNumberEntry()->GetNumber();
1631 dChi2Ref_ = guiDChi2Ref_->GetNumberEntry()->GetNumber();
1632 nMinIter_ = guiNMinIter_->GetNumberEntry()->GetNumber();
1633 nMaxIter_ = guiNMaxIter_->GetNumberEntry()->GetNumber();
1634 nMaxFailed_ = guiNMaxFailed_->GetNumberEntry()->GetNumber();
1635 resort_ = guiResort_->IsOn();
1637 gotoEvent(eventId_);
1641 void EventDisplay::guiSelectFitterId(
int val){
1642 fitterId_ = eFitterType(val-1);
1643 gotoEvent(eventId_);
1646 void EventDisplay::guiSelectMmHandling(
int val){
1648 gotoEvent(eventId_);
This class collects all information needed and produced by a specific AbsFitter and is specific to on...
virtual const TMatrixD & getMatrix() const =0
Get the actual matrix representation.
Contains the measurement and covariance in raw detector coordinates.
Abstract base class for a track representation.
virtual unsigned int getDim() const =0
Get the dimension of the state vector used by the track representation.
void getPosDir(const StateOnPlane &state, TVector3 &pos, TVector3 &dir) const
Get cartesian position and direction vector of a state.
virtual double extrapolateToPlane(StateOnPlane &state, const genfit::SharedPlanePtr &plane, bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0
Extrapolates the state to plane, and returns the extrapolation length and, via reference,...
void setPropDir(int dir)
Set propagation direction. (-1, 0, 1) -> (backward, auto, forward).
virtual void getPosMomCov(const MeasuredStateOnPlane &state, TVector3 &pos, TVector3 &mom, TMatrixDSym &cov) const =0
Translates MeasuredStateOnPlane into 3D position, momentum and 6x6 covariance.
Determinstic Annealing Filter (DAF) implementation.
Event display designed to run with Genfit.
Exception class for error handling in GENFIT (provides storage for diagnostic information)
Measurement class implementing a measurement of all track parameters.
virtual SharedPlanePtr constructPlane(const StateOnPlane &state) const override
Construct (virtual) detector plane (use state's AbsTrackRep).
Collects information needed and produced by a AbsKalmanFitter implementations and is specific to one ...
Kalman filter implementation with linearization around a reference track.
Simple Kalman filter implementation.
#StateOnPlane with additional covariance matrix.
Measured coordinates on a plane.
Measurement class implementing a planar hit geometry (1 or 2D).
AbsTrackRep with 5D track parameterization in plane coordinates: (q/p, u', v', u, v)
Class for measurements implementing a space point hit geometry.
A state with arbitrary dimension defined in a DetPlane.
Object containing AbsMeasurement and AbsFitterInfo objects.
AbsFitterInfo * getFitterInfo(const AbsTrackRep *rep=nullptr) const
Get fitterInfo for rep. Per default, use cardinal rep.
Collection of TrackPoint objects, AbsTrackRep objects and FitStatus objects.
Class for measurements in wire detectors (Straw tubes and drift chambers) which can measure the coord...
double sqrt(double a)
sqrt for double
double eval(const std::vector< double > &spl, const std::vector< double > &vals, double x)
Evaluate spline (zero order or first order) in point x.
Defines for I/O streams used for error and debug printing.
std::shared_ptr< genfit::DetPlane > SharedPlanePtr
Shared Pointer to a DetPlane.
eMultipleMeasurementHandling
@ weightedAverage
weighted average between measurements; used by DAF