SND@LHC Software
Loading...
Searching...
No Matches
EventDisplay.cc File Reference
#include "EventDisplay.h"
#include <assert.h>
#include <algorithm>
#include <cmath>
#include <exception>
#include <iostream>
#include <sys/time.h>
#include "AbsMeasurement.h"
#include "FullMeasurement.h"
#include "PlanarMeasurement.h"
#include "ProlateSpacepointMeasurement.h"
#include "SpacepointMeasurement.h"
#include "WireMeasurement.h"
#include "WirePointMeasurement.h"
#include "AbsTrackRep.h"
#include "ConstField.h"
#include "DetPlane.h"
#include "Exception.h"
#include "FieldManager.h"
#include "Tools.h"
#include "KalmanFitterInfo.h"
#include "KalmanFitter.h"
#include "DAF.h"
#include "KalmanFitterRefTrack.h"
#include "RKTrackRep.h"
#include <TApplication.h>
#include <TEveBrowser.h>
#include <TEveManager.h>
#include <TEveEventManager.h>
#include <TEveGeoNode.h>
#include <TEveGeoShape.h>
#include <TEveStraightLineSet.h>
#include <TEveTriangleSet.h>
#include <TDecompSVD.h>
#include <TGButton.h>
#include <TGLabel.h>
#include <TGNumberEntry.h>
#include <TGeoEltu.h>
#include <TGeoManager.h>
#include <TGeoMatrix.h>
#include <TGeoNode.h>
#include <TGeoSphere.h>
#include <TGeoTube.h>
#include <TMath.h>
#include <TMatrixT.h>
#include <TMatrixTSym.h>
#include <TMatrixDEigen.h>
#include <TROOT.h>
#include <TVector2.h>
#include <TVectorD.h>
#include <TSystem.h>
#include "boost/scoped_ptr.hpp"

Go to the source code of this file.

Functions

 ClassImp (genfit::EventDisplay) namespace genfit
 

Function Documentation

◆ ClassImp()

ClassImp ( genfit::EventDisplay  )

Definition at line 60 of file EventDisplay.cc.

62 {
63
64
65EventDisplay* EventDisplay::eventDisplay_ = NULL;
66
67EventDisplay::EventDisplay() :
68 errorScale_(1.),
69 drawGeometry_(false),
70 drawDetectors_(true),
71 drawHits_(true),
72 drawErrors_(true),
73 drawPlanes_(true),
74 drawTrackMarkers_(true),
75 drawTrack_(true),
76 drawRefTrack_(true),
77 drawForward_(true),
78 drawBackward_(true),
79 drawAutoScale_(true),
80 drawScaleMan_(false),
81 drawSilent_(false),
82 drawCardinalRep_(true),
83 repId_(0),
84 drawAllTracks_(true),
85 trackId_(0),
86 refit_(false),
87 debugLvl_(0),
88 fitterId_(SimpleKalman),
89 mmHandling_(weightedAverage),
90 squareRootFormalism_(false),
91 dPVal_(1.E-3),
92 dRelChi2_(0.2),
93 nMinIter_(2),
94 nMaxIter_(4),
95 nMaxFailed_(-1),
96 resort_(false)
97{
98
99 if((!gApplication) || (gApplication && gApplication->TestBit(TApplication::kDefaultApplication))) {
100 std::cout << "In EventDisplay ctor: gApplication not found, creating..." << std::flush;
101 new TApplication("ROOT_application", 0, 0);
102 std::cout << "done!" << std::endl;
103 }
104 if(!gEve) {
105 std::cout << "In EventDisplay ctor: gEve not found, creating..." << std::flush;
106 TEveManager::Create();
107 std::cout << "done!" << std::endl;
108 }
109
110 eventId_ = 0;
111
112}
113
114void EventDisplay::setOptions(std::string opts) {
115
116 if(opts != "") {
117 for(size_t i = 0; i < opts.length(); ++i) {
118 if(opts[i] == 'A') drawAutoScale_ = true;
119 if(opts[i] == 'B') drawBackward_ = true;
120 if(opts[i] == 'D') drawDetectors_ = true;
121 if(opts[i] == 'E') drawErrors_ = true;
122 if(opts[i] == 'F') drawForward_ = true;
123 if(opts[i] == 'H') drawHits_ = true;
124 if(opts[i] == 'M') drawTrackMarkers_ = true;
125 if(opts[i] == 'P') drawPlanes_ = true;
126 if(opts[i] == 'S') drawScaleMan_ = true;
127 if(opts[i] == 'T') drawTrack_ = true;
128 if(opts[i] == 'X') drawSilent_ = true;
129 if(opts[i] == 'G') drawGeometry_ = true;
130 }
131 }
132
133}
134
135void EventDisplay::setErrScale(double errScale) { errorScale_ = errScale; }
136
137double EventDisplay::getErrScale() { return errorScale_; }
138
139EventDisplay* EventDisplay::getInstance() {
140
141 if(eventDisplay_ == NULL) {
142 eventDisplay_ = new EventDisplay();
143 }
144 return eventDisplay_;
145
146}
147
148EventDisplay::~EventDisplay() { reset(); }
149
150void EventDisplay::reset() {
151
152 for(unsigned int i = 0; i < events_.size(); i++) {
153
154 for(unsigned int j = 0; j < events_[i]->size(); j++) {
155
156 delete events_[i]->at(j);
157
158 }
159 delete events_[i];
160 }
161
162 events_.clear();
163}
164
165
166void EventDisplay::addEvent(std::vector<Track*>& tracks) {
167
168 std::vector<Track*>* vec = new std::vector<Track*>;
169
170 for(unsigned int i = 0; i < tracks.size(); i++) {
171 vec->push_back(new Track(*(tracks[i])));
172 }
173
174 events_.push_back(vec);
175}
176
177
178void EventDisplay::addEvent(std::vector<const Track*>& tracks) {
179
180 std::vector<Track*>* vec = new std::vector<Track*>;
181
182 for(unsigned int i = 0; i < tracks.size(); i++) {
183 vec->push_back(new Track(*(tracks[i])));
184 }
185
186 events_.push_back(vec);
187}
188
189
190void EventDisplay::addEvent(const Track* tr) {
191
192 std::vector<Track*>* vec = new std::vector<Track*>;
193 vec->push_back(new Track(*tr));
194 events_.push_back(vec);
195}
196
197
198void EventDisplay::next(unsigned int stp) {
199
200 gotoEvent(eventId_ + stp);
201
202}
203
204void EventDisplay::prev(unsigned int stp) {
205
206 if(events_.size() == 0) return;
207 if(eventId_ < stp) {
208 gotoEvent(0);
209 } else {
210 gotoEvent(eventId_ - stp);
211 }
212
213}
214
215int EventDisplay::getNEvents() { return events_.size(); }
216
217
218void EventDisplay::gotoEvent(unsigned int id) {
219
220 if (events_.size() == 0)
221 return;
222 else if(id >= events_.size())
223 id = events_.size() - 1;
224
225 bool resetCam = true;
226
227 if (id == eventId_)
228 resetCam = false;
229
230 eventId_ = id;
231
232 std::cout << "At event " << id << std::endl;
233 if (gEve->GetCurrentEvent()) {
234 gEve->GetCurrentEvent()->DestroyElements();
235 }
236 double old_error_scale = errorScale_;
237 drawEvent(eventId_, resetCam);
238 if(old_error_scale != errorScale_) {
239 if (gEve->GetCurrentEvent()) {
240 gEve->GetCurrentEvent()->DestroyElements();
241 }
242 drawEvent(eventId_, resetCam); // if autoscaling changed the error, draw again.
243 }
244 errorScale_ = old_error_scale;
245
246}
247
248void EventDisplay::open() {
249
250 std::cout << "EventDisplay::open(); " << getNEvents() << " events loaded" << std::endl;
251
252 if(getNEvents() > 0) {
253 double old_error_scale = errorScale_;
254 drawEvent(0);
255 if(old_error_scale != errorScale_) {
256 std::cout << "autoscaling changed the error, draw again." << std::endl;
257 gotoEvent(0); // if autoscaling changed the error, draw again.
258 }
259 errorScale_ = old_error_scale;
260 }
261
262
263 if(!drawSilent_) {
264 makeGui();
265 gApplication->Run(kTRUE);
266 }
267
268 std::cout << "opened" << std::endl;
269
270}
271
272
273void EventDisplay::drawEvent(unsigned int id, bool resetCam) {
274
275 std::cout << "EventDisplay::drawEvent(" << id << ")" << std::endl;
276
277
278 // draw the geometry, does not really work yet. If it's fixed, the docu in the header file should be changed.
279 if(drawGeometry_) {
280 TGeoNode* top_node = gGeoManager->GetTopNode();
281 assert(top_node != NULL);
282
283 //Set transparency & color of geometry
284 TObjArray* volumes = gGeoManager->GetListOfVolumes();
285 for(int i = 0; i < volumes->GetEntriesFast(); i++) {
286 TGeoVolume* volume = dynamic_cast<TGeoVolume*>(volumes->At(i));
287 assert(volume != NULL);
288 volume->SetLineColor(12);
289 volume->SetTransparency(50);
290 }
291
292 TEveGeoTopNode* eve_top_node = new TEveGeoTopNode(gGeoManager, top_node);
293 eve_top_node->IncDenyDestroy();
294 gEve->AddElement(eve_top_node);
295 }
296
297
298 for(unsigned int i = 0; i < events_.at(id)->size(); i++) { // loop over all tracks in an event
299
300 if (!drawAllTracks_ && trackId_ != i)
301 continue;
302
303 Track* track = events_[id]->at(i);
304 if (! track->checkConsistency()){
305 std::cerr<<"track is not consistent"<<std::endl;
306 continue;
307 }
308
309
310 boost::scoped_ptr<Track> refittedTrack(NULL);
311 if (refit_) {
312
313 std::cout << "Refit track:" << std::endl;
314
315 boost::scoped_ptr<AbsKalmanFitter> fitter;
316 switch (fitterId_) {
317 case SimpleKalman:
318 fitter.reset(new KalmanFitter(nMaxIter_, dPVal_));
319 fitter->setMultipleMeasurementHandling(mmHandling_);
320 (static_cast<KalmanFitter*>(fitter.get()))->useSquareRootFormalism(squareRootFormalism_);
321 break;
322
323 case RefKalman:
324 fitter.reset(new KalmanFitterRefTrack(nMaxIter_, dPVal_));
325 fitter->setMultipleMeasurementHandling(mmHandling_);
326 break;
327
328 case DafSimple:
329 fitter.reset(new DAF(false));
330 ( static_cast<KalmanFitter*>( (static_cast<DAF*>(fitter.get()))->getKalman() ) )->useSquareRootFormalism(squareRootFormalism_);
331 break;
332 case DafRef:
333 fitter.reset(new DAF());
334 break;
335
336 }
337 fitter->setDebugLvl(std::max(0, (int)debugLvl_-1));
338 fitter->setMinIterations(nMinIter_);
339 fitter->setMaxIterations(nMaxIter_);
340 fitter->setRelChi2Change(dRelChi2_);
341 fitter->setMaxFailedHits(nMaxFailed_);
342
343
344 refittedTrack.reset(new Track(*track));
345 refittedTrack->deleteFitterInfo();
346
347 if (debugLvl_>0)
348 refittedTrack->Print("C");
349
350 timeval startcputime, endcputime;
351
352 try{
353 gettimeofday(&startcputime, NULL);
354 fitter->processTrack(refittedTrack.get(), resort_);
355 gettimeofday(&endcputime, NULL);
356 }
357 catch(genfit::Exception& e){
358 std::cerr << e.what();
359 std::cerr << "Exception, could not refit track" << std::endl;
360 continue;
361 }
362
363 int microseconds = 1000000*(endcputime.tv_sec - startcputime.tv_sec) + (endcputime.tv_usec - startcputime.tv_usec);
364 std::cout << "it took " << double(microseconds) / 1000 << " ms of CPU to fit the track\n";
365
366 if (! refittedTrack->checkConsistency()){
367 std::cerr<<"refittedTrack is not consistent"<<std::endl;
368 continue;
369 }
370
371 track = refittedTrack.get();
372 }
373
374
375
376
377 AbsTrackRep* rep;
378
379 if (drawCardinalRep_) {
380 rep = track->getCardinalRep();
381 std::cout << "Draw cardinal rep" << std::endl;
382 }
383 else {
384 if (repId_ >= track->getNumReps())
385 repId_ = track->getNumReps() - 1;
386 rep = track->getTrackRep(repId_);
387 std::cout << "Draw rep" << repId_ << std::endl;
388 }
389
390 if (debugLvl_>0) {
391 //track->Print();
392 track->Print("C");
393 track->getFitStatus(rep)->Print();
394
395 if (track->getFitStatus(rep)->isFitted()) {
396 try {
397 std::cout << "fitted state: \n";
398 track->getFittedState().Print();
399 }
400 catch (Exception& e) {
401 std::cerr << e.what();
402 }
403 }
404 }
405
406
407
408 rep->setPropDir(0);
409
410 unsigned int numhits = track->getNumPointsWithMeasurement();
411
412 KalmanFitterInfo* fi;
413 KalmanFitterInfo* prevFi = 0;
414 const MeasuredStateOnPlane* fittedState(NULL);
415 const MeasuredStateOnPlane* prevFittedState(NULL);
416
417 for(unsigned int j = 0; j < numhits; j++) { // loop over all hits in the track
418
419 fittedState = NULL;
420
421 TrackPoint* tp = track->getPointWithMeasurement(j);
422 if (! tp->hasRawMeasurements()) {
423 std::cerr<<"trackPoint has no raw measurements"<<std::endl;
424 continue;
425 }
426
427 const AbsMeasurement* m = tp->getRawMeasurement();
428 int hit_coords_dim = m->getDim();
429
430 // check if multiple AbsMeasurements are of same type
431 if (tp->getNumRawMeasurements() > 1) {
432 bool sameTypes(true);
433 for (unsigned int iM=1; iM<tp->getNumRawMeasurements(); ++iM) {
434 AbsMeasurement* tmp = tp->getRawMeasurement(iM);
435 if (typeid(*tmp) != typeid(*m))
436 sameTypes = false;
437 }
438 if (!sameTypes) {
439 std::cerr<<"cannot draw trackpoint containing multiple Measurements of differend types"<<std::endl;
440 continue;
441 }
442 }
443
444
445
446 // get the fitter infos ------------------------------------------------------------------
447 if (! tp->hasFitterInfo(rep)) {
448 std::cerr<<"trackPoint has no fitterInfo for rep"<<std::endl;
449 continue;
450 }
451
452 AbsFitterInfo* fitterInfo = tp->getFitterInfo(rep);
453
454 fi = dynamic_cast<KalmanFitterInfo*>(fitterInfo);
455 if(fi == NULL) {
456 std::cerr<<"can only display KalmanFitterInfo"<<std::endl;
457 continue;
458 }
459 if (! fi->hasPredictionsAndUpdates()) {
460 std::cerr<<"KalmanFitterInfo does not have all predictions and updates"<<std::endl;
461 //continue;
462 }
463 else {
464 try {
465 fittedState = &(fi->getFittedState(true));
466 }
467 catch (Exception& e) {
468 std::cerr << e.what();
469 std::cerr<<"can not get fitted state"<<std::endl;
470 fittedState = NULL;
471 prevFi = fi;
472 prevFittedState = fittedState;
473 continue;
474 }
475 }
476
477 if (fittedState == NULL) {
478 if (fi->hasForwardUpdate()) {
479 fittedState = fi->getForwardUpdate();
480 }
481 else if (fi->hasBackwardUpdate()) {
482 fittedState = fi->getBackwardUpdate();
483 }
484 else if (fi->hasForwardPrediction()) {
485 fittedState = fi->getForwardPrediction();
486 }
487 else if (fi->hasBackwardPrediction()) {
488 fittedState = fi->getBackwardPrediction();
489 }
490 }
491
492 if (fittedState == NULL) {
493 std::cout << "canot get any state from fitterInfo, continue.\n";
494 prevFi = fi;
495 prevFittedState = fittedState;
496 continue;
497 }
498
499 TVector3 track_pos = fittedState->getPos();
500 double charge = fittedState->getCharge();
501
502 //std::cout << "trackPos: "; track_pos.Print();
503
504
505 // determine measurement type
506 bool full_hit = false;
507 bool planar_hit = false;
508 bool planar_pixel_hit = false;
509 bool space_hit = false;
510 bool wire_hit = false;
511 bool wirepoint_hit = false;
512 if (dynamic_cast<const FullMeasurement*>(m) != NULL) {
513 full_hit = true;
514 }
515 else if(dynamic_cast<const PlanarMeasurement*>(m) != NULL) {
516 planar_hit = true;
517 if(hit_coords_dim == 2) {
518 planar_pixel_hit = true;
519 }
520 } else if (dynamic_cast<const SpacepointMeasurement*>(m) != NULL) {
521 space_hit = true;
522 } else if (dynamic_cast<const WireMeasurement*>(m) != NULL) {
523 wire_hit = true;
524 if (dynamic_cast<const WirePointMeasurement*>(m) != NULL) {
525 wirepoint_hit = true;
526 }
527 } else {
528 std::cout << "Track " << i << ", Hit " << j << ": Unknown measurement type: skipping hit!" << std::endl;
529 continue;
530 }
531
532
533 // loop over MeasurementOnPlanes
534 for (unsigned int iMeas = 0; iMeas < fi->getNumMeasurements(); ++iMeas) {
535
536 if (iMeas > 0 && wire_hit)
537 break;
538
539 const MeasurementOnPlane* mop = fi->getMeasurementOnPlane(iMeas);
540 const TVectorT<double>& hit_coords = mop->getState();
541 const TMatrixTSym<double>& hit_cov = mop->getCov();
542
543 // finished getting the hit infos -----------------------------------------------------
544
545 // sort hit infos into variables ------------------------------------------------------
546 TVector3 o = fittedState->getPlane()->getO();
547 TVector3 u = fittedState->getPlane()->getU();
548 TVector3 v = fittedState->getPlane()->getV();
549
550 double_t hit_u = 0;
551 double_t hit_v = 0;
552 double_t plane_size = 4;
553 TVector2 stripDir(1,0);
554
555 if(planar_hit) {
556 if(!planar_pixel_hit) {
557 if (dynamic_cast<RKTrackRep*>(rep) != NULL) {
558 const TMatrixD& H = mop->getHMatrix()->getMatrix();
559 stripDir.Set(H(0,3), H(0,4));
560 }
561 hit_u = hit_coords(0);
562 } else {
563 hit_u = hit_coords(0);
564 hit_v = hit_coords(1);
565 }
566 } else if (wire_hit) {
567 hit_u = fabs(hit_coords(0));
568 hit_v = v*(track_pos-o); // move the covariance tube so that the track goes through it
569 if (wirepoint_hit) {
570 hit_v = hit_coords(1);
571 }
572 }
573
574 if(plane_size < 4) plane_size = 4;
575 // finished setting variables ---------------------------------------------------------
576
577 // draw planes if corresponding option is set -----------------------------------------
578 if(iMeas == 0 &&
579 (drawPlanes_ || (drawDetectors_ && planar_hit))) {
580 TVector3 move(0,0,0);
581 if (planar_hit) move = track_pos-o;
582 if (wire_hit) move = v*(v*(track_pos-o)); // move the plane along the wire until the track goes through it
583 TEveBox* box = boxCreator(o + move, u, v, plane_size, plane_size, 0.01);
584 if (drawDetectors_ && planar_hit) {
585 box->SetMainColor(kCyan);
586 } else {
587 box->SetMainColor(kGray);
588 }
589 box->SetMainTransparency(50);
590 gEve->AddElement(box);
591 }
592 // finished drawing planes ------------------------------------------------------------
593
594 // draw track if corresponding option is set ------------------------------------------
595 if (j > 0 && prevFi != NULL) {
596 if(drawTrack_) {
597 makeLines(prevFittedState, fittedState, rep, charge > 0 ? kRed : kBlue, 1, drawTrackMarkers_, drawErrors_, 3);
598 if (drawErrors_) { // make sure to draw errors in both directions
599 makeLines(prevFittedState, fittedState, rep, charge > 0 ? kRed : kBlue, 1, false, drawErrors_, 0, 0);
600 }
601 }
602 if (drawForward_)
603 makeLines(prevFi->getForwardUpdate(), fi->getForwardPrediction(), rep, kCyan, 1, drawTrackMarkers_, drawErrors_, 1, 0);
604 if (drawBackward_)
605 makeLines(prevFi->getBackwardPrediction(), fi->getBackwardUpdate(), rep, kMagenta, 1, drawTrackMarkers_, drawErrors_, 1);
606 // draw reference track if corresponding option is set ------------------------------------------
607 if(drawRefTrack_ && fi->hasReferenceState() && prevFi->hasReferenceState())
608 makeLines(prevFi->getReferenceState(), fi->getReferenceState(), rep, charge > 0 ? kRed + 2 : kBlue + 2, 2, drawTrackMarkers_, false, 3);
609 }
610 else if (j > 0 && prevFi == NULL) {
611 std::cout << "previous FitterInfo == NULL \n";
612 }
613
614 // draw detectors if option is set, only important for wire hits ----------------------
615 if(drawDetectors_) {
616
617 if(wire_hit) {
618 TEveGeoShape* det_shape = new TEveGeoShape("det_shape");
619 det_shape->IncDenyDestroy();
620 det_shape->SetShape(new TGeoTube(std::max(0., (double)(hit_u-0.0105/2.)), hit_u+0.0105/2., plane_size));
621
622 TVector3 norm = u.Cross(v);
623 TGeoRotation det_rot("det_rot", (u.Theta()*180)/TMath::Pi(), (u.Phi()*180)/TMath::Pi(),
624 (norm.Theta()*180)/TMath::Pi(), (norm.Phi()*180)/TMath::Pi(),
625 (v.Theta()*180)/TMath::Pi(), (v.Phi()*180)/TMath::Pi()); // move the tube to the right place and rotate it correctly
626 TVector3 move = v*(v*(track_pos-o)); // move the tube along the wire until the track goes through it
627 TGeoCombiTrans det_trans(o(0) + move.X(),
628 o(1) + move.Y(),
629 o(2) + move.Z(),
630 &det_rot);
631 det_shape->SetTransMatrix(det_trans);
632 det_shape->SetMainColor(kCyan);
633 det_shape->SetMainTransparency(25);
634 if((drawHits_ && (hit_u+0.0105/2 > 0)) || !drawHits_) {
635 gEve->AddElement(det_shape);
636 }
637 }
638
639 }
640 // finished drawing detectors ---------------------------------------------------------
641
642 if(drawHits_) {
643
644 // draw planar hits, with distinction between strip and pixel hits ----------------
645 if (full_hit) {
646
647 StateOnPlane dummy(rep);
648 StateOnPlane dummy2(TVectorD(rep->getDim()), static_cast<const FullMeasurement*>(m)->constructPlane(dummy), rep);
649 MeasuredStateOnPlane sop = *(static_cast<const FullMeasurement*>(m)->constructMeasurementsOnPlane(dummy2)[0]);
650 sop.getCov()*=errorScale_;
651
652 MeasuredStateOnPlane prevSop(sop);
653 prevSop.extrapolateBy(-3);
654 makeLines(&sop, &prevSop, rep, kYellow, 1, false, true, 0, 0);
655
656 prevSop = sop;
657 prevSop.extrapolateBy(3);
658 makeLines(&sop, &prevSop, rep, kYellow, 1, false, true, 0, 0);
659 }
660
661 if(planar_hit) {
662 if(!planar_pixel_hit) {
663 TEveBox* hit_box;
664 TVector3 stripDir3 = stripDir.X()*u + stripDir.Y()*v;
665 TVector3 stripDir3perp = stripDir.Y()*u - stripDir.X()*v;
666 TVector3 move = stripDir3perp*(stripDir3perp*(track_pos-o));
667 hit_box = boxCreator((o + move + hit_u*stripDir3), stripDir3, stripDir3perp, errorScale_*std::sqrt(hit_cov(0,0)), plane_size, 0.0105);
668 hit_box->SetMainColor(kYellow);
669 hit_box->SetMainTransparency(0);
670 gEve->AddElement(hit_box);
671 } else {
672 // calculate eigenvalues to draw error-ellipse ----------------------------
673 TMatrixDEigen eigen_values(hit_cov);
674 TEveGeoShape* cov_shape = new TEveGeoShape("cov_shape");
675 cov_shape->IncDenyDestroy();
676 TMatrixT<double> ev = eigen_values.GetEigenValues();
677 TMatrixT<double> eVec = eigen_values.GetEigenVectors();
678 double pseudo_res_0 = errorScale_*std::sqrt(ev(0,0));
679 double pseudo_res_1 = errorScale_*std::sqrt(ev(1,1));
680 // finished calcluating, got the values -----------------------------------
681
682 // do autoscaling if necessary --------------------------------------------
683 if(drawAutoScale_) {
684 double min_cov = std::min(pseudo_res_0,pseudo_res_1);
685 if(min_cov < 1e-5) {
686 std::cout << "Track " << i << ", Hit " << j << ": Invalid covariance matrix (Eigenvalue < 1e-5), autoscaling not possible!" << std::endl;
687 } else {
688 if(min_cov < 0.049) {
689 double cor = 0.05 / min_cov;
690 std::cout << "Track " << i << ", Hit " << j << ": Pixel covariance too small, rescaling by " << cor;
691 errorScale_ *= cor;
692 pseudo_res_0 *= cor;
693 pseudo_res_1 *= cor;
694 std::cout << " to " << errorScale_ << std::endl;
695 }
696 }
697 }
698 // finished autoscaling ---------------------------------------------------
699
700 // calculate the semiaxis of the error ellipse ----------------------------
701 cov_shape->SetShape(new TGeoEltu(pseudo_res_0, pseudo_res_1, 0.0105));
702 TVector3 pix_pos = o + hit_u*u + hit_v*v;
703 TVector3 u_semiaxis = (pix_pos + eVec(0,0)*u + eVec(1,0)*v)-pix_pos;
704 TVector3 v_semiaxis = (pix_pos + eVec(0,1)*u + eVec(1,1)*v)-pix_pos;
705 TVector3 norm = u.Cross(v);
706 // finished calculating ---------------------------------------------------
707
708 // rotate and translate everything correctly ------------------------------
709 TGeoRotation det_rot("det_rot", (u_semiaxis.Theta()*180)/TMath::Pi(), (u_semiaxis.Phi()*180)/TMath::Pi(),
710 (v_semiaxis.Theta()*180)/TMath::Pi(), (v_semiaxis.Phi()*180)/TMath::Pi(),
711 (norm.Theta()*180)/TMath::Pi(), (norm.Phi()*180)/TMath::Pi());
712 TGeoCombiTrans det_trans(pix_pos(0),pix_pos(1),pix_pos(2), &det_rot);
713 cov_shape->SetTransMatrix(det_trans);
714 // finished rotating and translating --------------------------------------
715
716 cov_shape->SetMainColor(kYellow);
717 cov_shape->SetMainTransparency(0);
718 gEve->AddElement(cov_shape);
719 }
720 }
721 // finished drawing planar hits ---------------------------------------------------
722
723 // draw spacepoint hits -----------------------------------------------------------
724 if(space_hit) {
725 {
726 // get eigenvalues of covariance to know how to draw the ellipsoid ------------
727 TMatrixDEigen eigen_values(m->getRawHitCov());
728 TEveGeoShape* cov_shape = new TEveGeoShape("cov_shape");
729 cov_shape->IncDenyDestroy();
730 cov_shape->SetShape(new TGeoSphere(0.,1.));
731 TMatrixT<double> ev = eigen_values.GetEigenValues();
732 TMatrixT<double> eVec = eigen_values.GetEigenVectors();
733 TVector3 eVec1(eVec(0,0),eVec(1,0),eVec(2,0));
734 TVector3 eVec2(eVec(0,1),eVec(1,1),eVec(2,1));
735 TVector3 eVec3(eVec(0,2),eVec(1,2),eVec(2,2));
736 const TVector3 norm = u.Cross(v);
737 // got everything we need -----------------------------------------------------
738
739 static const double radDeg(180./TMath::Pi());
740 TGeoRotation det_rot("det_rot", eVec1.Theta()*radDeg, eVec1.Phi()*radDeg,
741 eVec2.Theta()*radDeg, eVec2.Phi()*radDeg,
742 eVec3.Theta()*radDeg, eVec3.Phi()*radDeg);
743
744 if (! det_rot.IsValid()){
745 // hackish fix if eigenvectors are not orthonogonal
746 if (fabs(eVec2*eVec3) > 1.e-10)
747 eVec3 = eVec1.Cross(eVec2);
748
749 det_rot.SetAngles(eVec1.Theta()*radDeg, eVec1.Phi()*radDeg,
750 eVec2.Theta()*radDeg, eVec2.Phi()*radDeg,
751 eVec3.Theta()*radDeg, eVec3.Phi()*radDeg);
752 }
753
754 // set the scaled eigenvalues -------------------------------------------------
755 double pseudo_res_0 = errorScale_*std::sqrt(ev(0,0));
756 double pseudo_res_1 = errorScale_*std::sqrt(ev(1,1));
757 double pseudo_res_2 = errorScale_*std::sqrt(ev(2,2));
758 if(drawScaleMan_) { // override again if necessary
759 pseudo_res_0 = errorScale_*0.5;
760 pseudo_res_1 = errorScale_*0.5;
761 pseudo_res_2 = errorScale_*0.5;
762 }
763 // finished scaling -----------------------------------------------------------
764
765 // autoscale if necessary -----------------------------------------------------
766 if(drawAutoScale_) {
767 double min_cov = std::min(pseudo_res_0,std::min(pseudo_res_1,pseudo_res_2));
768 if(min_cov < 1e-5) {
769 std::cout << "Track " << i << ", Hit " << j << ": Invalid covariance matrix (Eigenvalue < 1e-5), autoscaling not possible!" << std::endl;
770 } else {
771 if(min_cov <= 0.149) {
772 double cor = 0.15 / min_cov;
773 std::cout << "Track " << i << ", Hit " << j << ": Space hit covariance too small, rescaling by " << cor;
774 errorScale_ *= cor;
775 pseudo_res_0 *= cor;
776 pseudo_res_1 *= cor;
777 pseudo_res_2 *= cor;
778 std::cout << " to " << errorScale_ << std::endl;
779
780 }
781 }
782 }
783 // finished autoscaling -------------------------------------------------------
784
785 // rotate and translate -------------------------------------------------------
786 TGeoGenTrans det_trans(o(0),o(1),o(2),
787 //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
788 //1/(pseudo_res_0),1/(pseudo_res_1),1/(pseudo_res_2),
789 pseudo_res_0, pseudo_res_1, pseudo_res_2,
790 &det_rot);
791 cov_shape->SetTransMatrix(det_trans);
792 // finished rotating and translating ------------------------------------------
793
794 cov_shape->SetMainColor(kYellow);
795 cov_shape->SetMainTransparency(10);
796 gEve->AddElement(cov_shape);
797 }
798
799
800 {
801 // calculate eigenvalues to draw error-ellipse ----------------------------
802 TMatrixDEigen eigen_values(hit_cov);
803 TEveGeoShape* cov_shape = new TEveGeoShape("cov_shape");
804 cov_shape->IncDenyDestroy();
805 TMatrixT<double> ev = eigen_values.GetEigenValues();
806 TMatrixT<double> eVec = eigen_values.GetEigenVectors();
807 double pseudo_res_0 = errorScale_*std::sqrt(ev(0,0));
808 double pseudo_res_1 = errorScale_*std::sqrt(ev(1,1));
809 // finished calcluating, got the values -----------------------------------
810
811 // do autoscaling if necessary --------------------------------------------
812 if(drawAutoScale_) {
813 double min_cov = std::min(pseudo_res_0,pseudo_res_1);
814 if(min_cov < 1e-5) {
815 std::cout << "Track " << i << ", Hit " << j << ": Invalid covariance matrix (Eigenvalue < 1e-5), autoscaling not possible!" << std::endl;
816 } else {
817 if(min_cov < 0.049) {
818 double cor = 0.05 / min_cov;
819 std::cout << "Track " << i << ", Hit " << j << ": Pixel covariance too small, rescaling by " << cor;
820 errorScale_ *= cor;
821 pseudo_res_0 *= cor;
822 pseudo_res_1 *= cor;
823 std::cout << " to " << errorScale_ << std::endl;
824 }
825 }
826 }
827 // finished autoscaling ---------------------------------------------------
828
829 // calculate the semiaxis of the error ellipse ----------------------------
830 cov_shape->SetShape(new TGeoEltu(pseudo_res_0, pseudo_res_1, 0.0105));
831 TVector3 pix_pos = o + hit_u*u + hit_v*v;
832 TVector3 u_semiaxis = (pix_pos + eVec(0,0)*u + eVec(1,0)*v)-pix_pos;
833 TVector3 v_semiaxis = (pix_pos + eVec(0,1)*u + eVec(1,1)*v)-pix_pos;
834 TVector3 norm = u.Cross(v);
835 // finished calculating ---------------------------------------------------
836
837 // rotate and translate everything correctly ------------------------------
838 static const double radDeg(180./TMath::Pi());
839 TGeoRotation det_rot("det_rot", u_semiaxis.Theta()*radDeg, u_semiaxis.Phi()*radDeg,
840 v_semiaxis.Theta()*radDeg, v_semiaxis.Phi()*radDeg,
841 norm.Theta()*radDeg, norm.Phi()*radDeg);
842 /*if (! det_rot.IsValid()){
843 u_semiaxis.Print();
844 v_semiaxis.Print();
845 norm.Print();
846 }*/
847 TGeoCombiTrans det_trans(pix_pos(0),pix_pos(1),pix_pos(2), &det_rot);
848 cov_shape->SetTransMatrix(det_trans);
849 // finished rotating and translating --------------------------------------
850
851 cov_shape->SetMainColor(kYellow);
852 cov_shape->SetMainTransparency(0);
853 gEve->AddElement(cov_shape);
854 }
855 }
856 // finished drawing spacepoint hits -----------------------------------------------
857
858 // draw wire hits -----------------------------------------------------------------
859 if(wire_hit) {
860 TEveGeoShape* cov_shape = new TEveGeoShape("cov_shape");
861 cov_shape->IncDenyDestroy();
862 double pseudo_res_0 = errorScale_*std::sqrt(hit_cov(0,0));
863 double pseudo_res_1 = plane_size;
864 if (wirepoint_hit) pseudo_res_1 = errorScale_*std::sqrt(hit_cov(1,1));
865
866 // autoscale if necessary -----------------------------------------------------
867 if(drawAutoScale_) {
868 if(pseudo_res_0 < 1e-5) {
869 std::cout << "Track " << i << ", Hit " << j << ": Invalid wire resolution (< 1e-5), autoscaling not possible!" << std::endl;
870 } else {
871 if(pseudo_res_0 < 0.0049) {
872 double cor = 0.005 / pseudo_res_0;
873 std::cout << "Track " << i << ", Hit " << j << ": Wire covariance too small, rescaling by " << cor;
874 errorScale_ *= cor;
875 pseudo_res_0 *= cor;
876 std::cout << " to " << errorScale_ << std::endl;
877 }
878 }
879
880 if(wirepoint_hit && pseudo_res_1 < 1e-5) {
881 std::cout << "Track " << i << ", Hit " << j << ": Invalid wire resolution (< 1e-5), autoscaling not possible!" << std::endl;
882 } else {
883 if(pseudo_res_1 < 0.0049) {
884 double cor = 0.005 / pseudo_res_1;
885 std::cout << "Track " << i << ", Hit " << j << ": Wire covariance too small, rescaling by " << cor;
886 errorScale_ *= cor;
887 pseudo_res_1 *= cor;
888 std::cout << " to " << errorScale_ << std::endl;
889 }
890 }
891 }
892 // finished autoscaling -------------------------------------------------------
893
894 TEveBox* hit_box;
895 TVector3 move = v*(v*(track_pos-o));
896 hit_box = boxCreator((o + move + hit_u*u), u, v, errorScale_*std::sqrt(hit_cov(0,0)), pseudo_res_1, 0.0105);
897 hit_box->SetMainColor(kYellow);
898 hit_box->SetMainTransparency(0);
899 gEve->AddElement(hit_box);
900
901 hit_box = boxCreator((o + move - hit_u*u), u, v, errorScale_*std::sqrt(hit_cov(0,0)), pseudo_res_1, 0.0105);
902 hit_box->SetMainColor(kYellow);
903 hit_box->SetMainTransparency(0);
904 gEve->AddElement(hit_box);
905 }
906 // finished drawing wire hits -----------------------------------------------------
907
908 } // finished drawing hits
909
910 } // finished looping over MeasurmentOnPlanes
911
912
913 prevFi = fi;
914 prevFittedState = fittedState;
915
916 }
917
918 }
919
920 gEve->Redraw3D(resetCam);
921
922}
923
924
925
926
927TEveBox* EventDisplay::boxCreator(TVector3 o, TVector3 u, TVector3 v, float ud, float vd, float depth) {
928
929 TEveBox* box = new TEveBox("detPlane_shape");
930 float vertices[24];
931
932 TVector3 norm = u.Cross(v);
933 u *= (0.5*ud);
934 v *= (0.5*vd);
935 norm *= (0.5*depth);
936
937 vertices[0] = o(0) - u(0) - v(0) - norm(0);
938 vertices[1] = o(1) - u(1) - v(1) - norm(1);
939 vertices[2] = o(2) - u(2) - v(2) - norm(2);
940
941 vertices[3] = o(0) + u(0) - v(0) - norm(0);
942 vertices[4] = o(1) + u(1) - v(1) - norm(1);
943 vertices[5] = o(2) + u(2) - v(2) - norm(2);
944
945 vertices[6] = o(0) + u(0) - v(0) + norm(0);
946 vertices[7] = o(1) + u(1) - v(1) + norm(1);
947 vertices[8] = o(2) + u(2) - v(2) + norm(2);
948
949 vertices[9] = o(0) - u(0) - v(0) + norm(0);
950 vertices[10] = o(1) - u(1) - v(1) + norm(1);
951 vertices[11] = o(2) - u(2) - v(2) + norm(2);
952
953 vertices[12] = o(0) - u(0) + v(0) - norm(0);
954 vertices[13] = o(1) - u(1) + v(1) - norm(1);
955 vertices[14] = o(2) - u(2) + v(2) - norm(2);
956
957 vertices[15] = o(0) + u(0) + v(0) - norm(0);
958 vertices[16] = o(1) + u(1) + v(1) - norm(1);
959 vertices[17] = o(2) + u(2) + v(2) - norm(2);
960
961 vertices[18] = o(0) + u(0) + v(0) + norm(0);
962 vertices[19] = o(1) + u(1) + v(1) + norm(1);
963 vertices[20] = o(2) + u(2) + v(2) + norm(2);
964
965 vertices[21] = o(0) - u(0) + v(0) + norm(0);
966 vertices[22] = o(1) - u(1) + v(1) + norm(1);
967 vertices[23] = o(2) - u(2) + v(2) + norm(2);
968
969
970 for(int k = 0; k < 24; k += 3) box->SetVertex((k/3), vertices[k], vertices[k+1], vertices[k+2]);
971
972 return box;
973
974}
975
976
977void EventDisplay::makeLines(const StateOnPlane* prevState, const StateOnPlane* state, const AbsTrackRep* rep,
978 const Color_t& color, const Style_t& style, bool drawMarkers, bool drawErrors, double lineWidth, int markerPos)
979{
980 if (prevState == NULL || state == NULL) {
981 std::cerr << "prevState == NULL || state == NULL\n";
982 return;
983 }
984
985 TVector3 pos, dir, oldPos, oldDir;
986 rep->getPosDir(*state, pos, dir);
987 rep->getPosDir(*prevState, oldPos, oldDir);
988
989 double distA = (pos-oldPos).Mag();
990 double distB = distA;
991 if ((pos-oldPos)*oldDir < 0)
992 distA *= -1.;
993 if ((pos-oldPos)*dir < 0)
994 distB *= -1.;
995 TVector3 intermediate1 = oldPos + 0.3 * distA * oldDir;
996 TVector3 intermediate2 = pos - 0.3 * distB * dir;
997 TEveStraightLineSet* lineSet = new TEveStraightLineSet;
998 lineSet->AddLine(oldPos(0), oldPos(1), oldPos(2), intermediate1(0), intermediate1(1), intermediate1(2));
999 lineSet->AddLine(intermediate1(0), intermediate1(1), intermediate1(2), intermediate2(0), intermediate2(1), intermediate2(2));
1000 lineSet->AddLine(intermediate2(0), intermediate2(1), intermediate2(2), pos(0), pos(1), pos(2));
1001 lineSet->SetLineColor(color);
1002 lineSet->SetLineStyle(style);
1003 lineSet->SetLineWidth(lineWidth);
1004 if (drawMarkers) {
1005 if (markerPos == 0)
1006 lineSet->AddMarker(oldPos(0), oldPos(1), oldPos(2));
1007 else
1008 lineSet->AddMarker(pos(0), pos(1), pos(2));
1009 }
1010
1011 if (lineWidth > 0)
1012 gEve->AddElement(lineSet);
1013
1014
1015 if (drawErrors) {
1016 const MeasuredStateOnPlane* measuredState;
1017 if (markerPos == 0)
1018 measuredState = dynamic_cast<const MeasuredStateOnPlane*>(prevState);
1019 else
1020 measuredState = dynamic_cast<const MeasuredStateOnPlane*>(state);
1021
1022 if (measuredState != NULL) {
1023
1024 // step for evaluate at a distance from the original plane
1025 TVector3 eval;
1026 if (markerPos == 0) {
1027 if (fabs(distA) < 1.) {
1028 distA < 0 ? distA = -1 : distA = 1;
1029 }
1030 eval = 0.2 * distA * oldDir;
1031 }
1032 else {
1033 if (fabs(distB) < 1.) {
1034 distB < 0 ? distB = -1 : distB = 1;
1035 }
1036 eval = -0.2 * distB * dir;
1037 }
1038
1039
1040 // get cov at first plane
1041 TMatrixDSym cov;
1042 TVector3 position, direction;
1043 rep->getPosMomCov(*measuredState, position, direction, cov);
1044
1045 // get eigenvalues & -vectors
1046 TMatrixDEigen eigen_values(cov.GetSub(0,2, 0,2));
1047 TMatrixT<double> ev = eigen_values.GetEigenValues();
1048 TMatrixT<double> eVec = eigen_values.GetEigenVectors();
1049 TVector3 eVec1, eVec2;
1050 // limit
1051 static const double maxErr = 1000.;
1052 double ev0 = std::min(ev(0,0), maxErr);
1053 double ev1 = std::min(ev(1,1), maxErr);
1054 double ev2 = std::min(ev(2,2), maxErr);
1055
1056 // get two largest eigenvalues/-vectors
1057 if (ev0 < ev1 && ev0 < ev2) {
1058 eVec1.SetXYZ(eVec(0,1),eVec(1,1),eVec(2,1));
1059 eVec1 *= sqrt(ev1);
1060 eVec2.SetXYZ(eVec(0,2),eVec(1,2),eVec(2,2));
1061 eVec2 *= sqrt(ev2);
1062 }
1063 else if (ev1 < ev0 && ev1 < ev2) {
1064 eVec1.SetXYZ(eVec(0,0),eVec(1,0),eVec(2,0));
1065 eVec1 *= sqrt(ev0);
1066 eVec2.SetXYZ(eVec(0,2),eVec(1,2),eVec(2,2));
1067 eVec2 *= sqrt(ev2);
1068 }
1069 else {
1070 eVec1.SetXYZ(eVec(0,0),eVec(1,0),eVec(2,0));
1071 eVec1 *= sqrt(ev0);
1072 eVec2.SetXYZ(eVec(0,1),eVec(1,1),eVec(2,1));
1073 eVec2 *= sqrt(ev1);
1074 }
1075
1076 if (eVec1.Cross(eVec2)*eval < 0)
1077 eVec2 *= -1;
1078 //assert(eVec1.Cross(eVec2)*eval > 0);
1079
1080 const TVector3 oldEVec1(eVec1);
1081 const TVector3 oldEVec2(eVec2);
1082
1083 const int nEdges = 24;
1084 std::vector<TVector3> vertices;
1085
1086 vertices.push_back(position);
1087
1088 // vertices at plane
1089 for (int i=0; i<nEdges; ++i) {
1090 const double angle = 2*TMath::Pi()/nEdges * i;
1091 vertices.push_back(position + cos(angle)*eVec1 + sin(angle)*eVec2);
1092 }
1093
1094
1095
1096 DetPlane* newPlane = new DetPlane(*(measuredState->getPlane()));
1097 newPlane->setO(position + eval);
1098
1099 MeasuredStateOnPlane stateCopy(*measuredState);
1100 try{
1101 rep->extrapolateToPlane(stateCopy, SharedPlanePtr(newPlane));
1102 }
1103 catch(Exception& e){
1104 std::cerr<<e.what();
1105 return;
1106 }
1107
1108 // get cov at 2nd plane
1109 rep->getPosMomCov(stateCopy, position, direction, cov);
1110
1111 // get eigenvalues & -vectors
1112 TMatrixDEigen eigen_values2(cov.GetSub(0,2, 0,2));
1113 ev = eigen_values2.GetEigenValues();
1114 eVec = eigen_values2.GetEigenVectors();
1115 // limit
1116 ev0 = std::min(ev(0,0), maxErr);
1117 ev1 = std::min(ev(1,1), maxErr);
1118 ev2 = std::min(ev(2,2), maxErr);
1119
1120 // get two largest eigenvalues/-vectors
1121 if (ev0 < ev1 && ev0 < ev2) {
1122 eVec1.SetXYZ(eVec(0,1),eVec(1,1),eVec(2,1));
1123 eVec1 *= sqrt(ev1);
1124 eVec2.SetXYZ(eVec(0,2),eVec(1,2),eVec(2,2));
1125 eVec2 *= sqrt(ev2);
1126 }
1127 else if (ev1 < ev0 && ev1 < ev2) {
1128 eVec1.SetXYZ(eVec(0,0),eVec(1,0),eVec(2,0));
1129 eVec1 *= sqrt(ev0);
1130 eVec2.SetXYZ(eVec(0,2),eVec(1,2),eVec(2,2));
1131 eVec2 *= sqrt(ev2);
1132 }
1133 else {
1134 eVec1.SetXYZ(eVec(0,0),eVec(1,0),eVec(2,0));
1135 eVec1 *= sqrt(ev0);
1136 eVec2.SetXYZ(eVec(0,1),eVec(1,1),eVec(2,1));
1137 eVec2 *= sqrt(ev1);
1138 }
1139
1140 if (eVec1.Cross(eVec2)*eval < 0)
1141 eVec2 *= -1;
1142 //assert(eVec1.Cross(eVec2)*eval > 0);
1143
1144 if (oldEVec1*eVec1 < 0) {
1145 eVec1 *= -1;
1146 eVec2 *= -1;
1147 }
1148
1149 // vertices at 2nd plane
1150 double angle0 = eVec1.Angle(oldEVec1);
1151 if (eVec1*(eval.Cross(oldEVec1)) < 0)
1152 angle0 *= -1;
1153 for (int i=0; i<nEdges; ++i) {
1154 const double angle = 2*TMath::Pi()/nEdges * i - angle0;
1155 vertices.push_back(position + cos(angle)*eVec1 + sin(angle)*eVec2);
1156 }
1157
1158 vertices.push_back(position);
1159
1160
1161 TEveTriangleSet* error_shape = new TEveTriangleSet(vertices.size(), nEdges*2);
1162 for(unsigned int k = 0; k < vertices.size(); ++k) {
1163 error_shape->SetVertex(k, vertices[k].X(), vertices[k].Y(), vertices[k].Z());
1164 }
1165
1166 assert(vertices.size() == 2*nEdges+2);
1167
1168 int iTri(0);
1169 for (int i=0; i<nEdges; ++i) {
1170 //error_shape->SetTriangle(iTri++, 0, i+1, (i+1)%nEdges+1);
1171 error_shape->SetTriangle(iTri++, i+1, i+1+nEdges, (i+1)%nEdges+1);
1172 error_shape->SetTriangle(iTri++, (i+1)%nEdges+1, i+1+nEdges, (i+1)%nEdges+1+nEdges);
1173 //error_shape->SetTriangle(iTri++, 2*nEdges+1, i+1+nEdges, (i+1)%nEdges+1+nEdges);
1174 }
1175
1176 //assert(iTri == nEdges*4);
1177
1178 error_shape->SetMainColor(color);
1179 error_shape->SetMainTransparency(25);
1180 gEve->AddElement(error_shape);
1181 }
1182 }
1183}
1184
1185
1186void EventDisplay::makeGui() {
1187
1188 TEveBrowser* browser = gEve->GetBrowser();
1189 browser->StartEmbedding(TRootBrowser::kLeft);
1190
1191 TGMainFrame* frmMain = new TGMainFrame(gClient->GetRoot(), 1000, 600);
1192 frmMain->SetWindowName("XX GUI");
1193 frmMain->SetCleanup(kDeepCleanup);
1194
1195 TGLabel* lbl = 0;
1196 TGTextButton* tb = 0;
1197 EventDisplay* fh = EventDisplay::getInstance();
1198
1199 TGHorizontalFrame* hf = new TGHorizontalFrame(frmMain); {
1200 // evt number entry
1201 lbl = new TGLabel(hf, "Go to event: ");
1202 hf->AddFrame(lbl);
1203 guiEvent = new TGNumberEntry(hf, 0, 9,999, TGNumberFormat::kNESInteger,
1204 TGNumberFormat::kNEANonNegative,
1205 TGNumberFormat::kNELLimitMinMax,
1206 0, 99999);
1207 hf->AddFrame(guiEvent);
1208 guiEvent->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiGoto()");
1209
1210 // redraw button
1211 tb = new TGTextButton(hf, "Redraw Event");
1212 hf->AddFrame(tb);
1213 tb->Connect("Clicked()", "genfit::EventDisplay", fh, "guiGoto()");
1214 }
1215 frmMain->AddFrame(hf);
1216
1217 // draw options
1218 hf = new TGHorizontalFrame(frmMain); {
1219 lbl = new TGLabel(hf, "\n Draw Options");
1220 hf->AddFrame(lbl);
1221 }
1222 frmMain->AddFrame(hf);
1223
1224 hf = new TGHorizontalFrame(frmMain); {
1225 guiDrawGeometry_ = new TGCheckButton(hf, "Draw geometry");
1226 if(drawGeometry_) guiDrawGeometry_->Toggle();
1227 hf->AddFrame(guiDrawGeometry_);
1228 guiDrawGeometry_->Connect("Toggled(Bool_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1229 }
1230 frmMain->AddFrame(hf);
1231
1232 hf = new TGHorizontalFrame(frmMain); {
1233 guiDrawDetectors_ = new TGCheckButton(hf, "Draw detectors");
1234 if(drawDetectors_) guiDrawDetectors_->Toggle();
1235 hf->AddFrame(guiDrawDetectors_);
1236 guiDrawDetectors_->Connect("Toggled(Bool_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1237 }
1238 frmMain->AddFrame(hf);
1239
1240 hf = new TGHorizontalFrame(frmMain); {
1241 guiDrawHits_ = new TGCheckButton(hf, "Draw hits");
1242 if(drawHits_) guiDrawHits_->Toggle();
1243 hf->AddFrame(guiDrawHits_);
1244 guiDrawHits_->Connect("Toggled(Bool_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1245 }
1246 frmMain->AddFrame(hf);
1247
1248
1249
1250 hf = new TGHorizontalFrame(frmMain); {
1251 guiDrawPlanes_ = new TGCheckButton(hf, "Draw planes");
1252 if(drawPlanes_) guiDrawPlanes_->Toggle();
1253 hf->AddFrame(guiDrawPlanes_);
1254 guiDrawPlanes_->Connect("Toggled(Bool_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1255 }
1256 frmMain->AddFrame(hf);
1257
1258 hf = new TGHorizontalFrame(frmMain); {
1259 guiDrawTrackMarkers_ = new TGCheckButton(hf, "Draw track markers");
1260 if(drawTrackMarkers_) guiDrawTrackMarkers_->Toggle();
1261 hf->AddFrame(guiDrawTrackMarkers_);
1262 guiDrawTrackMarkers_->Connect("Toggled(Bool_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1263 }
1264 frmMain->AddFrame(hf);
1265
1266
1267 hf = new TGHorizontalFrame(frmMain); {
1268 guiDrawTrack_ = new TGCheckButton(hf, "Draw track");
1269 if(drawTrack_) guiDrawTrack_->Toggle();
1270 hf->AddFrame(guiDrawTrack_);
1271 guiDrawTrack_->Connect("Toggled(Bool_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1272 }
1273 frmMain->AddFrame(hf);
1274
1275 hf = new TGHorizontalFrame(frmMain); {
1276 guiDrawRefTrack_ = new TGCheckButton(hf, "Draw reference track");
1277 if(drawRefTrack_) guiDrawRefTrack_->Toggle();
1278 hf->AddFrame(guiDrawRefTrack_);
1279 guiDrawRefTrack_->Connect("Toggled(Bool_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1280 }
1281 frmMain->AddFrame(hf);
1282
1283 hf = new TGHorizontalFrame(frmMain); {
1284 guiDrawErrors_ = new TGCheckButton(hf, "Draw track errors");
1285 if(drawErrors_) guiDrawErrors_->Toggle();
1286 hf->AddFrame(guiDrawErrors_);
1287 guiDrawErrors_->Connect("Toggled(Bool_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1288 }
1289 frmMain->AddFrame(hf);
1290
1291 hf = new TGHorizontalFrame(frmMain); {
1292 guiDrawForward_ = new TGCheckButton(hf, "Draw forward fit");
1293 if(drawForward_) guiDrawForward_->Toggle();
1294 hf->AddFrame(guiDrawForward_);
1295 guiDrawForward_->Connect("Toggled(Bool_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1296 }
1297 frmMain->AddFrame(hf);
1298
1299 hf = new TGHorizontalFrame(frmMain); {
1300 guiDrawBackward_ = new TGCheckButton(hf, "Draw backward fit");
1301 if(drawBackward_) guiDrawBackward_->Toggle();
1302 hf->AddFrame(guiDrawBackward_);
1303 guiDrawBackward_->Connect("Toggled(Bool_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1304 }
1305 frmMain->AddFrame(hf);
1306
1307
1308 hf = new TGHorizontalFrame(frmMain); {
1309 guiDrawAutoScale_ = new TGCheckButton(hf, "Auto-scale errors");
1310 if(drawAutoScale_) guiDrawAutoScale_->Toggle();
1311 hf->AddFrame(guiDrawAutoScale_);
1312 guiDrawAutoScale_->Connect("Toggled(Bool_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1313 }
1314 frmMain->AddFrame(hf);
1315
1316 hf = new TGHorizontalFrame(frmMain); {
1317 guiDrawScaleMan_ = new TGCheckButton(hf, "Manually scale errors");
1318 if(drawScaleMan_) guiDrawScaleMan_->Toggle();
1319 hf->AddFrame(guiDrawScaleMan_);
1320 guiDrawScaleMan_->Connect("Toggled(Bool_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1321 }
1322 frmMain->AddFrame(hf);
1323
1324 hf = new TGHorizontalFrame(frmMain); {
1325 guiErrorScale_ = new TGNumberEntry(hf, errorScale_, 6,999, TGNumberFormat::kNESReal,
1326 TGNumberFormat::kNEANonNegative,
1327 TGNumberFormat::kNELLimitMinMax,
1328 1.E-4, 1.E5);
1329 hf->AddFrame(guiErrorScale_);
1330 guiErrorScale_->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1331 lbl = new TGLabel(hf, "Error scale");
1332 hf->AddFrame(lbl);
1333 }
1334 frmMain->AddFrame(hf);
1335
1336
1337
1338 hf = new TGHorizontalFrame(frmMain); {
1339 lbl = new TGLabel(hf, "\n TrackRep options");
1340 hf->AddFrame(lbl);
1341 }
1342 frmMain->AddFrame(hf);
1343
1344 hf = new TGHorizontalFrame(frmMain); {
1345 guiDrawCardinalRep_ = new TGCheckButton(hf, "Draw cardinal rep");
1346 if(drawCardinalRep_) guiDrawCardinalRep_->Toggle();
1347 hf->AddFrame(guiDrawCardinalRep_);
1348 guiDrawCardinalRep_->Connect("Toggled(Bool_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1349 }
1350 frmMain->AddFrame(hf);
1351
1352 hf = new TGHorizontalFrame(frmMain); {
1353 guiRepId_ = new TGNumberEntry(hf, repId_, 6,999, TGNumberFormat::kNESInteger,
1354 TGNumberFormat::kNEANonNegative,
1355 TGNumberFormat::kNELLimitMinMax,
1356 0, 99);
1357 hf->AddFrame(guiRepId_);
1358 guiRepId_->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1359 lbl = new TGLabel(hf, "Else draw rep with id");
1360 hf->AddFrame(lbl);
1361 }
1362 frmMain->AddFrame(hf);
1363
1364 hf = new TGHorizontalFrame(frmMain); {
1365 guiDrawAllTracks_ = new TGCheckButton(hf, "Draw all tracks");
1366 if(drawAllTracks_) guiDrawAllTracks_->Toggle();
1367 hf->AddFrame(guiDrawAllTracks_);
1368 guiDrawAllTracks_->Connect("Toggled(Bool_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1369 }
1370 frmMain->AddFrame(hf);
1371
1372 hf = new TGHorizontalFrame(frmMain); {
1373 guiTrackId_ = new TGNumberEntry(hf, trackId_, 6,999, TGNumberFormat::kNESInteger,
1374 TGNumberFormat::kNEANonNegative,
1375 TGNumberFormat::kNELLimitMinMax,
1376 0, 99);
1377 hf->AddFrame(guiTrackId_);
1378 guiTrackId_->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1379 lbl = new TGLabel(hf, "Else draw track nr. ");
1380 hf->AddFrame(lbl);
1381 }
1382 frmMain->AddFrame(hf);
1383
1384
1385
1386 frmMain->MapSubwindows();
1387 frmMain->Resize();
1388 frmMain->MapWindow();
1389
1390 browser->StopEmbedding();
1391 browser->SetTabTitle("Draw Control", 0);
1392
1393
1394 browser->StartEmbedding(TRootBrowser::kLeft);
1395 TGMainFrame* frmMain2 = new TGMainFrame(gClient->GetRoot(), 1000, 600);
1396 frmMain2->SetWindowName("XX GUI");
1397 frmMain2->SetCleanup(kDeepCleanup);
1398
1399 hf = new TGHorizontalFrame(frmMain2); {
1400 // evt number entry
1401 lbl = new TGLabel(hf, "Go to event: ");
1402 hf->AddFrame(lbl);
1403 guiEvent2 = new TGNumberEntry(hf, 0, 9,999, TGNumberFormat::kNESInteger,
1404 TGNumberFormat::kNEANonNegative,
1405 TGNumberFormat::kNELLimitMinMax,
1406 0, 99999);
1407 hf->AddFrame(guiEvent2);
1408 guiEvent2->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiGoto2()");
1409
1410 // redraw button
1411 tb = new TGTextButton(hf, "Redraw Event");
1412 hf->AddFrame(tb);
1413 tb->Connect("Clicked()", "genfit::EventDisplay", fh, "guiGoto()");
1414 }
1415 frmMain2->AddFrame(hf);
1416
1417 hf = new TGHorizontalFrame(frmMain2); {
1418 lbl = new TGLabel(hf, "\n Fitting options");
1419 hf->AddFrame(lbl);
1420 }
1421 frmMain2->AddFrame(hf);
1422
1423 hf = new TGHorizontalFrame(frmMain2); {
1424 guiRefit_ = new TGCheckButton(hf, "Refit");
1425 if(refit_) guiRefit_->Toggle();
1426 hf->AddFrame(guiRefit_);
1427 guiRefit_->Connect("Toggled(Bool_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1428 }
1429 frmMain2->AddFrame(hf);
1430
1431 hf = new TGHorizontalFrame(frmMain2); {
1432 guiDebugLvl_ = new TGNumberEntry(hf, debugLvl_, 6,999, TGNumberFormat::kNESInteger,
1433 TGNumberFormat::kNEANonNegative,
1434 TGNumberFormat::kNELLimitMinMax,
1435 0, 999);
1436 hf->AddFrame(guiDebugLvl_);
1437 guiDebugLvl_->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1438 lbl = new TGLabel(hf, "debug level");
1439 hf->AddFrame(lbl);
1440 }
1441 frmMain2->AddFrame(hf);
1442
1443 hf = new TGHorizontalFrame(frmMain2); {
1444 guiFitterId_ = new TGButtonGroup(hf,"Fitter type:");
1445 guiFitterId_->Connect("Clicked(Int_t)","genfit::EventDisplay", fh, "guiSelectFitterId(int)");
1446 hf->AddFrame(guiFitterId_, new TGLayoutHints(kLHintsTop));
1447 TGRadioButton* fitterId_button = new TGRadioButton(guiFitterId_, "Simple Kalman");
1448 new TGRadioButton(guiFitterId_, "Reference Kalman");
1449 new TGRadioButton(guiFitterId_, "DAF w/ simple Kalman");
1450 new TGRadioButton(guiFitterId_, "DAF w/ reference Kalman");
1451 fitterId_button->SetDown(true, false);
1452 guiFitterId_->Show();
1453 }
1454 frmMain2->AddFrame(hf);
1455
1456 hf = new TGHorizontalFrame(frmMain2); {
1457 guiMmHandling_ = new TGButtonGroup(hf,"Multiple measurement handling in Kalman:");
1458 guiMmHandling_->Connect("Clicked(Int_t)","genfit::EventDisplay", fh, "guiSelectMmHandling(int)");
1459 hf->AddFrame(guiMmHandling_, new TGLayoutHints(kLHintsTop));
1460 TGRadioButton* mmHandling_button = new TGRadioButton(guiMmHandling_, "weighted average");
1461 new TGRadioButton(guiMmHandling_, "unweighted average");
1462 new TGRadioButton(guiMmHandling_, "weighted, closest to reference");
1463 new TGRadioButton(guiMmHandling_, "unweighted, closest to reference");
1464 new TGRadioButton(guiMmHandling_, "weighted, closest to prediction");
1465 new TGRadioButton(guiMmHandling_, "unweighted, closest to prediction");
1466 new TGRadioButton(guiMmHandling_, "weighted, closest to reference for WireMeasurements, weighted average else");
1467 new TGRadioButton(guiMmHandling_, "unweighted, closest to reference for WireMeasurements, unweighted average else");
1468 new TGRadioButton(guiMmHandling_, "weighted, closest to prediction for WireMeasurements, weighted average else");
1469 new TGRadioButton(guiMmHandling_, "unweighted, closest to prediction for WireMeasurements, unweighted average else");
1470 mmHandling_button->SetDown(true, false);
1471 guiMmHandling_->Show();
1472 }
1473 frmMain2->AddFrame(hf);
1474
1475 hf = new TGHorizontalFrame(frmMain2); {
1476 guiSquareRootFormalism_ = new TGCheckButton(hf, "Use square root formalism (simple Kalman/simple DAF)");
1477 if(squareRootFormalism_) guiSquareRootFormalism_->Toggle();
1478 hf->AddFrame(guiSquareRootFormalism_);
1479 guiSquareRootFormalism_->Connect("Toggled(Bool_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1480 }
1481 frmMain2->AddFrame(hf);
1482
1483 hf = new TGHorizontalFrame(frmMain2); {
1484 guiDPVal_ = new TGNumberEntry(hf, dPVal_, 6,9999, TGNumberFormat::kNESReal,
1485 TGNumberFormat::kNEANonNegative,
1486 TGNumberFormat::kNELLimitMinMax,
1487 0, 999);
1488 hf->AddFrame(guiDPVal_);
1489 guiDPVal_->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1490 lbl = new TGLabel(hf, "delta pVal (convergence criterium)");
1491 hf->AddFrame(lbl);
1492 }
1493 frmMain2->AddFrame(hf);
1494
1495 hf = new TGHorizontalFrame(frmMain2); {
1496 guiRelChi2_ = new TGNumberEntry(hf, dRelChi2_, 6,9999, TGNumberFormat::kNESReal,
1497 TGNumberFormat::kNEANonNegative,
1498 TGNumberFormat::kNELLimitMinMax,
1499 0, 999);
1500 hf->AddFrame(guiRelChi2_);
1501 guiRelChi2_->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1502 lbl = new TGLabel(hf, "rel chi^2 change (non-convergence criterium)");
1503 hf->AddFrame(lbl);
1504 }
1505 frmMain2->AddFrame(hf);
1506
1507 hf = new TGHorizontalFrame(frmMain2); {
1508 guiNMinIter_ = new TGNumberEntry(hf, nMinIter_, 6,999, TGNumberFormat::kNESInteger,
1509 TGNumberFormat::kNEANonNegative,
1510 TGNumberFormat::kNELLimitMinMax,
1511 1, 100);
1512 hf->AddFrame(guiNMinIter_);
1513 guiNMinIter_->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1514 lbl = new TGLabel(hf, "Minimum nr of iterations");
1515 hf->AddFrame(lbl);
1516 }
1517 frmMain2->AddFrame(hf);
1518
1519 hf = new TGHorizontalFrame(frmMain2); {
1520 guiNMaxIter_ = new TGNumberEntry(hf, nMaxIter_, 6,999, TGNumberFormat::kNESInteger,
1521 TGNumberFormat::kNEANonNegative,
1522 TGNumberFormat::kNELLimitMinMax,
1523 1, 100);
1524 hf->AddFrame(guiNMaxIter_);
1525 guiNMaxIter_->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1526 lbl = new TGLabel(hf, "Maximum nr of iterations");
1527 hf->AddFrame(lbl);
1528 }
1529 frmMain2->AddFrame(hf);
1530
1531 hf = new TGHorizontalFrame(frmMain2); {
1532 guiNMaxFailed_ = new TGNumberEntry(hf, nMaxFailed_, 6,999, TGNumberFormat::kNESInteger,
1533 TGNumberFormat::kNEAAnyNumber,
1534 TGNumberFormat::kNELLimitMinMax,
1535 -1, 1000);
1536 hf->AddFrame(guiNMaxFailed_);
1537 guiNMaxFailed_->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1538 lbl = new TGLabel(hf, "Maximum nr of failed hits");
1539 hf->AddFrame(lbl);
1540 }
1541 frmMain2->AddFrame(hf);
1542
1543
1544 hf = new TGHorizontalFrame(frmMain2); {
1545 guiResort_ = new TGCheckButton(hf, "Resort track");
1546 if(resort_) guiResort_->Toggle();
1547 hf->AddFrame(guiResort_);
1548 guiResort_->Connect("Toggled(Bool_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1549 }
1550 frmMain2->AddFrame(hf);
1551
1552
1553
1554
1555 frmMain2->MapSubwindows();
1556 frmMain2->Resize();
1557 frmMain2->MapWindow();
1558
1559 browser->StopEmbedding();
1560 browser->SetTabTitle("Refit Control", 0);
1561}
1562
1563
1564void EventDisplay::guiGoto(){
1565 Long_t n = guiEvent->GetNumberEntry()->GetIntNumber();
1566 guiEvent2->SetIntNumber(n);
1567 gotoEvent(n);
1568}
1569
1570void EventDisplay::guiGoto2(){
1571 Long_t n = guiEvent2->GetNumberEntry()->GetIntNumber();
1572 guiEvent->SetIntNumber(n);
1573 gotoEvent(n);
1574}
1575
1576
1577void EventDisplay::guiSetDrawParams(){
1578
1579 drawGeometry_ = guiDrawGeometry_->IsOn();
1580 drawDetectors_ = guiDrawDetectors_->IsOn();
1581 drawHits_ = guiDrawHits_->IsOn();
1582 drawErrors_ = guiDrawErrors_->IsOn();
1583
1584 drawPlanes_ = guiDrawPlanes_->IsOn();
1585 drawTrackMarkers_ = guiDrawTrackMarkers_->IsOn();
1586 drawTrack_ = guiDrawTrack_->IsOn();
1587 drawRefTrack_ = guiDrawRefTrack_->IsOn();
1588 drawForward_ = guiDrawForward_->IsOn();
1589 drawBackward_ = guiDrawBackward_->IsOn();
1590
1591 drawAutoScale_ = guiDrawAutoScale_->IsOn();
1592 drawScaleMan_ = guiDrawScaleMan_->IsOn();
1593
1594 errorScale_ = guiErrorScale_->GetNumberEntry()->GetNumber();
1595
1596 drawCardinalRep_ = guiDrawCardinalRep_->IsOn();
1597 repId_ = guiRepId_->GetNumberEntry()->GetNumber();
1598
1599 drawAllTracks_ = guiDrawAllTracks_->IsOn();
1600 trackId_ = guiTrackId_->GetNumberEntry()->GetNumber();
1601
1602
1603 refit_ = guiRefit_->IsOn();
1604 debugLvl_ = guiDebugLvl_->GetNumberEntry()->GetNumber();
1605
1606 squareRootFormalism_ = guiSquareRootFormalism_->IsOn();
1607 dPVal_ = guiDPVal_->GetNumberEntry()->GetNumber();
1608 dRelChi2_ = guiRelChi2_->GetNumberEntry()->GetNumber();
1609 nMinIter_ = guiNMinIter_->GetNumberEntry()->GetNumber();
1610 nMaxIter_ = guiNMaxIter_->GetNumberEntry()->GetNumber();
1611 nMaxFailed_ = guiNMaxFailed_->GetNumberEntry()->GetNumber();
1612 resort_ = guiResort_->IsOn();
1613
1614 gotoEvent(eventId_);
1615}
1616
1617
1618void EventDisplay::guiSelectFitterId(int val){
1619 fitterId_ = eFitterType(val-1);
1620 gotoEvent(eventId_);
1621}
1622
1623void EventDisplay::guiSelectMmHandling(int val){
1624 mmHandling_ = eMultipleMeasurementHandling(val-1);
1625 gotoEvent(eventId_);
1626}
1627
1628
1629} // end of namespace genfit
Double_t m
Exception class for error handling in GENFIT (provides storage for diagnostic information)
Definition Exception.h:48
virtual const char * what() const
Standard error message handling for exceptions. use like "std::cerr << e.what();".
Definition Exception.cc:51
int i
Definition ShipAna.py:86
eMultipleMeasurementHandling
@ SimpleKalman
boost::shared_ptr< genfit::DetPlane > SharedPlanePtr
Shared Pointer to a DetPlane.
integer, parameter double
Definition mpdef.f90:38