Belle II Software  release-08-01-10
EventDisplay.cc
1 
2 #include "EventDisplay.h"
3 
4 #include <assert.h>
5 #include <algorithm>
6 #include <cmath>
7 #include <exception>
8 #include <iostream>
9 #include <sys/time.h>
10 
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"
20 #include "DetPlane.h"
21 #include "Exception.h"
22 #include "FieldManager.h"
23 #include "Tools.h"
24 #include "KalmanFitterInfo.h"
25 #include "KalmanFitter.h"
26 #include "DAF.h"
27 #include "KalmanFitterRefTrack.h"
28 #include "RKTrackRep.h"
29 
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>
39 #include <TGButton.h>
40 #include <TGLabel.h>
41 #include <TGNumberEntry.h>
42 #include <TGeoEltu.h>
43 #include <TGeoManager.h>
44 #include <TGeoMatrix.h>
45 #include <TGeoNode.h>
46 #include <TGeoSphere.h>
47 #include <TGeoTube.h>
48 #include <TMath.h>
49 #include <TMatrixT.h>
50 #include <TMatrixTSym.h>
51 #include <TMatrixDSymEigen.h>
52 #include <TROOT.h>
53 #include <TVector2.h>
54 #include <TVectorD.h>
55 #include <TSystem.h>
56 
57 #include <memory>
58 
59 namespace genfit {
60 
61 
62 EventDisplay* EventDisplay::eventDisplay_ = nullptr;
63 
64 EventDisplay::EventDisplay() :
65  errorScale_(1.),
66  drawGeometry_(false),
67  drawDetectors_(true),
68  drawHits_(true),
69  drawErrors_(true),
70  drawPlanes_(true),
71  drawTrackMarkers_(true),
72  drawTrack_(true),
73  drawRefTrack_(true),
74  drawForward_(true),
75  drawBackward_(true),
76  drawAutoScale_(true),
77  drawScaleMan_(false),
78  drawSilent_(false),
79  drawCardinalRep_(true),
80  repId_(0),
81  drawAllTracks_(true),
82  trackId_(0),
83  refit_(false),
84  debugLvl_(0),
85  fitterId_(SimpleKalman),
86  mmHandling_(weightedAverage),
87  squareRootFormalism_(false),
88  dPVal_(1.E-3),
89  dRelChi2_(0.2),
90  dChi2Ref_(1.),
91  nMinIter_(2),
92  nMaxIter_(4),
93  nMaxFailed_(-1),
94  resort_(false)
95 {
96 
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;
101  }
102  if(!gEve) {
103  std::cout << "In EventDisplay ctor: gEve not found, creating..." << std::flush;
104  TEveManager::Create();
105  std::cout << "done!" << std::endl;
106  }
107 
108  eventId_ = 0;
109 
110 }
111 
112 void EventDisplay::setOptions(std::string opts) {
113 
114  if(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;
128  }
129  }
130 
131 }
132 
133 void EventDisplay::setErrScale(double errScale) { errorScale_ = errScale; }
134 
135 double EventDisplay::getErrScale() { return errorScale_; }
136 
137 EventDisplay* EventDisplay::getInstance() {
138 
139  if(eventDisplay_ == nullptr) {
140  eventDisplay_ = new EventDisplay();
141  }
142  return eventDisplay_;
143 
144 }
145 
146 EventDisplay::~EventDisplay() { reset(); }
147 
148 void EventDisplay::reset() {
149 
150  for(unsigned int i = 0; i < events_.size(); i++) {
151 
152  for(unsigned int j = 0; j < events_[i]->size(); j++) {
153 
154  delete events_[i]->at(j);
155 
156  }
157  delete events_[i];
158  }
159 
160  events_.clear();
161 }
162 
163 
164 void EventDisplay::addEvent(std::vector<Track*>& tracks) {
165 
166  std::vector<Track*>* vec = new std::vector<Track*>;
167 
168  for(unsigned int i = 0; i < tracks.size(); i++) {
169  vec->push_back(new Track(*(tracks[i])));
170  }
171 
172  events_.push_back(vec);
173 }
174 
175 
176 void EventDisplay::addEvent(std::vector<const Track*>& tracks) {
177 
178  std::vector<Track*>* vec = new std::vector<Track*>;
179 
180  for(unsigned int i = 0; i < tracks.size(); i++) {
181  vec->push_back(new Track(*(tracks[i])));
182  }
183 
184  events_.push_back(vec);
185 }
186 
187 
188 void EventDisplay::addEvent(const Track* tr) {
189 
190  std::vector<Track*>* vec = new std::vector<Track*>;
191  vec->push_back(new Track(*tr));
192  events_.push_back(vec);
193 }
194 
195 
196 void EventDisplay::next(unsigned int stp) {
197 
198  gotoEvent(eventId_ + stp);
199 
200 }
201 
202 void EventDisplay::prev(unsigned int stp) {
203 
204  if(events_.size() == 0) return;
205  if(eventId_ < stp) {
206  gotoEvent(0);
207  } else {
208  gotoEvent(eventId_ - stp);
209  }
210 
211 }
212 
213 int EventDisplay::getNEvents() { return events_.size(); }
214 
215 
216 void EventDisplay::gotoEvent(unsigned int id) {
217 
218  if (events_.size() == 0)
219  return;
220  else if(id >= events_.size())
221  id = events_.size() - 1;
222 
223  bool resetCam = true;
224 
225  if (id == eventId_)
226  resetCam = false;
227 
228  eventId_ = id;
229 
230  std::cout << "At event " << id << std::endl;
231  if (gEve->GetCurrentEvent()) {
232  gEve->GetCurrentEvent()->DestroyElements();
233  }
234  double old_error_scale = errorScale_;
235  drawEvent(eventId_, resetCam);
236  if(old_error_scale != errorScale_) {
237  if (gEve->GetCurrentEvent()) {
238  gEve->GetCurrentEvent()->DestroyElements();
239  }
240  drawEvent(eventId_, resetCam); // if autoscaling changed the error, draw again.
241  }
242  errorScale_ = old_error_scale;
243 
244 }
245 
246 void EventDisplay::open() {
247 
248  std::cout << "EventDisplay::open(); " << getNEvents() << " events loaded" << std::endl;
249 
250  if(getNEvents() > 0) {
251  double old_error_scale = errorScale_;
252  drawEvent(0);
253  if(old_error_scale != errorScale_) {
254  std::cout << "autoscaling changed the error, draw again." << std::endl;
255  gotoEvent(0); // if autoscaling changed the error, draw again.
256  }
257  errorScale_ = old_error_scale;
258  }
259 
260 
261  if(!drawSilent_) {
262  makeGui();
263  gApplication->Run(kTRUE);
264  }
265 
266  std::cout << "opened" << std::endl;
267 
268 }
269 
270 
271 void EventDisplay::drawEvent(unsigned int id, bool resetCam) {
272 
273  std::cout << "EventDisplay::drawEvent(" << id << ")" << std::endl;
274 
275 
276  // draw the geometry, does not really work yet. If it's fixed, the docu in the header file should be changed.
277  if(drawGeometry_) {
278  TGeoNode* top_node = gGeoManager->GetTopNode();
279  assert(top_node != nullptr);
280 
281  //Set transparency & color of geometry
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);
288  }
289 
290  TEveGeoTopNode* eve_top_node = new TEveGeoTopNode(gGeoManager, top_node);
291  eve_top_node->IncDenyDestroy();
292  gEve->AddElement(eve_top_node);
293  }
294 
295 
296  for(unsigned int i = 0; i < events_.at(id)->size(); i++) { // loop over all tracks in an event
297 
298  if (!drawAllTracks_ && trackId_ != i)
299  continue;
300 
301  Track* track = events_[id]->at(i);
302  try {
303  track->checkConsistency();
304  } catch (genfit::Exception& e) {
305  std::cerr<< e.getExcString() <<std::endl;
306  continue;
307  }
308 
309  std::unique_ptr<Track> refittedTrack(nullptr);
310  if (refit_) {
311 
312  std::cout << "Refit track:" << std::endl;
313 
314  std::unique_ptr<AbsKalmanFitter> fitter;
315  switch (fitterId_) {
316  case SimpleKalman:
317  fitter.reset(new KalmanFitter(nMaxIter_, dPVal_));
318  fitter->setMultipleMeasurementHandling(mmHandling_);
319  (static_cast<KalmanFitter*>(fitter.get()))->useSquareRootFormalism(squareRootFormalism_);
320  break;
321 
322  case RefKalman:
323  fitter.reset(new KalmanFitterRefTrack(nMaxIter_, dPVal_));
324  fitter->setMultipleMeasurementHandling(mmHandling_);
325  static_cast<KalmanFitterRefTrack*>(fitter.get())->setDeltaChi2Ref(dChi2Ref_);
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  ( static_cast<KalmanFitterRefTrack*>( (static_cast<DAF*>(fitter.get()))->getKalman() ) )->setDeltaChi2Ref(dChi2Ref_);
335  break;
336 
337  }
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_);
343 
344 
345  refittedTrack.reset(new Track(*track));
346  refittedTrack->deleteFitterInfo();
347 
348  if (debugLvl_>0)
349  refittedTrack->Print("C");
350 
351  timeval startcputime, endcputime;
352 
353  try{
354  gettimeofday(&startcputime, nullptr);
355  fitter->processTrack(refittedTrack.get(), resort_);
356  gettimeofday(&endcputime, nullptr);
357  }
358  catch(genfit::Exception& e){
359  std::cerr << e.what();
360  std::cerr << "Exception, could not refit track" << std::endl;
361  continue;
362  }
363 
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";
366 
367  try {
368  refittedTrack->checkConsistency();
369  } catch (genfit::Exception& e) {
370  std::cerr<< e.getExcString() <<std::endl;
371  continue;
372  }
373 
374  track = refittedTrack.get();
375  }
376 
377 
378 
379 
380  AbsTrackRep* rep;
381 
382  if (drawCardinalRep_) {
383  rep = track->getCardinalRep();
384  std::cout << "Draw cardinal rep" << std::endl;
385  }
386  else {
387  if (repId_ >= track->getNumReps())
388  repId_ = track->getNumReps() - 1;
389  rep = track->getTrackRep(repId_);
390  std::cout << "Draw rep" << repId_ << std::endl;
391  }
392 
393  if (debugLvl_>0) {
394  std::cout << "track " << i << std::endl;
395  //track->Print();
396  track->Print("C");
397  track->getFitStatus(rep)->Print();
398 
399  if (track->getFitStatus(rep)->isFitted()) {
400  try {
401  std::cout << "fitted state: \n";
402  track->getFittedState().Print();
403  }
404  catch (Exception& e) {
405  std::cerr << e.what();
406  }
407  }
408  }
409 
410 
411 
412  rep->setPropDir(0);
413 
414  unsigned int numhits = track->getNumPointsWithMeasurement();
415 
416  KalmanFitterInfo* fi;
417  KalmanFitterInfo* prevFi = 0;
418  const MeasuredStateOnPlane* fittedState(nullptr);
419  const MeasuredStateOnPlane* prevFittedState(nullptr);
420 
421  for(unsigned int j = 0; j < numhits; j++) { // loop over all hits in the track
422 
423  fittedState = nullptr;
424 
425  TrackPoint* tp = track->getPointWithMeasurement(j);
426  if (! tp->hasRawMeasurements()) {
427  std::cerr<<"trackPoint has no raw measurements"<<std::endl;
428  continue;
429  }
430 
431  const AbsMeasurement* m = tp->getRawMeasurement();
432  int hit_coords_dim = m->getDim();
433 
434  // check if multiple AbsMeasurements are of same type
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))
440  sameTypes = false;
441  }
442  if (!sameTypes) {
443  std::cerr<<"cannot draw trackpoint containing multiple Measurements of differend types"<<std::endl;
444  continue;
445  }
446  }
447 
448 
449 
450  // get the fitter infos ------------------------------------------------------------------
451  if (! tp->hasFitterInfo(rep)) {
452  std::cerr<<"trackPoint has no fitterInfo for rep"<<std::endl;
453  continue;
454  }
455 
456  AbsFitterInfo* fitterInfo = tp->getFitterInfo(rep);
457 
458  fi = dynamic_cast<KalmanFitterInfo*>(fitterInfo);
459  if(fi == nullptr) {
460  std::cerr<<"can only display KalmanFitterInfo"<<std::endl;
461  continue;
462  }
463  if (! fi->hasPredictionsAndUpdates()) {
464  std::cerr<<"KalmanFitterInfo does not have all predictions and updates"<<std::endl;
465  //continue;
466  }
467  else {
468  try {
469  fittedState = &(fi->getFittedState(true));
470  }
471  catch (Exception& e) {
472  std::cerr << e.what();
473  std::cerr<<"can not get fitted state"<<std::endl;
474  fittedState = nullptr;
475  prevFi = fi;
476  prevFittedState = fittedState;
477  continue;
478  }
479  }
480 
481  if (fittedState == nullptr) {
482  if (fi->hasForwardUpdate()) {
483  fittedState = fi->getForwardUpdate();
484  }
485  else if (fi->hasBackwardUpdate()) {
486  fittedState = fi->getBackwardUpdate();
487  }
488  else if (fi->hasForwardPrediction()) {
489  fittedState = fi->getForwardPrediction();
490  }
491  else if (fi->hasBackwardPrediction()) {
492  fittedState = fi->getBackwardPrediction();
493  }
494  }
495 
496  if (fittedState == nullptr) {
497  std::cout << "cannot get any state from fitterInfo, continue.\n";
498  prevFi = fi;
499  prevFittedState = fittedState;
500  continue;
501  }
502 
503  TVector3 track_pos = fittedState->getPos();
504  double charge = fittedState->getCharge();
505 
506  //std::cout << "trackPos: "; track_pos.Print();
507 
508 
509  // determine measurement type
510  bool full_hit = (dynamic_cast<const FullMeasurement*>(m) != nullptr);
511  bool planar_hit = (dynamic_cast<const PlanarMeasurement*>(m) != nullptr);
512  bool planar_pixel_hit = planar_hit && hit_coords_dim == 2;
513  bool space_hit = (dynamic_cast<const SpacepointMeasurement*>(m) != nullptr);
514  bool wire_hit = m && m->isLeftRightMeasurement();
515  bool wirepoint_hit = wire_hit && (dynamic_cast<const WirePointMeasurement*>(m) != nullptr);
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;
518  continue;
519  }
520 
521 
522  // loop over MeasurementOnPlanes
523  unsigned int nMeas = fi->getNumMeasurements();
524  for (unsigned int iMeas = 0; iMeas < nMeas; ++iMeas) {
525 
526  if (iMeas > 0 && wire_hit)
527  break;
528 
529  const MeasurementOnPlane* mop = fi->getMeasurementOnPlane(iMeas);
530  const TVectorT<double>& hit_coords = mop->getState();
531  const TMatrixTSym<double>& hit_cov = mop->getCov();
532 
533  // finished getting the hit infos -----------------------------------------------------
534 
535  // sort hit infos into variables ------------------------------------------------------
536  TVector3 o = fittedState->getPlane()->getO();
537  TVector3 u = fittedState->getPlane()->getU();
538  TVector3 v = fittedState->getPlane()->getV();
539 
540  double_t hit_u = 0;
541  double_t hit_v = 0;
542  double_t plane_size = 4;
543  TVector2 stripDir(1,0);
544 
545  if(planar_hit) {
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));
550  }
551  hit_u = hit_coords(0);
552  } else {
553  hit_u = hit_coords(0);
554  hit_v = hit_coords(1);
555  }
556  } else if (wire_hit) {
557  hit_u = fabs(hit_coords(0));
558  hit_v = v*(track_pos-o); // move the covariance tube so that the track goes through it
559  if (wirepoint_hit) {
560  hit_v = hit_coords(1);
561  }
562  }
563 
564  if(plane_size < 4) plane_size = 4;
565  // finished setting variables ---------------------------------------------------------
566 
567  // draw planes if corresponding option is set -----------------------------------------
568  if(iMeas == 0 &&
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)); // move the plane along the wire until the track goes through it
573  TEveBox* box = boxCreator(o + move, u, v, plane_size, plane_size, 0.01);
574  if (drawDetectors_ && planar_hit) {
575  box->SetMainColor(kCyan);
576  } else {
577  box->SetMainColor(kGray);
578  }
579  box->SetMainTransparency(50);
580  gEve->AddElement(box);
581  }
582  // finished drawing planes ------------------------------------------------------------
583 
584  // draw track if corresponding option is set ------------------------------------------
585  try {
586  if (j == 0) {
587  if (drawBackward_) {
588  MeasuredStateOnPlane update ( *fi->getBackwardUpdate() );
589  update.extrapolateBy(-3.);
590  makeLines(&update, fi->getBackwardUpdate(), rep, kMagenta, 1, drawTrackMarkers_, drawErrors_, 1);
591  }
592  }
593  if (j > 0 && prevFi != nullptr) {
594  if(drawTrack_) {
595  makeLines(prevFittedState, fittedState, rep, charge > 0 ? kRed : kBlue, 1, drawTrackMarkers_, drawErrors_, 3);
596  if (drawErrors_) { // make sure to draw errors in both directions
597  makeLines(prevFittedState, fittedState, rep, charge > 0 ? kRed : kBlue, 1, false, drawErrors_, 0, 0);
598  }
599  }
600  if (drawForward_) {
601  makeLines(prevFi->getForwardUpdate(), fi->getForwardPrediction(), rep, kCyan, 1, drawTrackMarkers_, drawErrors_, 1, 0);
602  if (j == numhits-1) {
603  MeasuredStateOnPlane update ( *fi->getForwardUpdate() );
604  update.extrapolateBy(3.);
605  makeLines(fi->getForwardUpdate(), &update, rep, kCyan, 1, drawTrackMarkers_, drawErrors_, 1, 0);
606  }
607  }
608  if (drawBackward_) {
609  makeLines(prevFi->getBackwardPrediction(), fi->getBackwardUpdate(), rep, kMagenta, 1, drawTrackMarkers_, drawErrors_, 1);
610  }
611  // draw reference track if corresponding option is set ------------------------------------------
612  if(drawRefTrack_ && fi->hasReferenceState() && prevFi->hasReferenceState())
613  makeLines(prevFi->getReferenceState(), fi->getReferenceState(), rep, charge > 0 ? kRed + 2 : kBlue + 2, 2, drawTrackMarkers_, false, 3);
614  }
615  else if (j > 0 && prevFi == nullptr) {
616  std::cout << "previous FitterInfo == nullptr \n";
617  }
618  }
619  catch (Exception& e) {
620  std::cerr << "extrapolation failed, cannot draw track" << std::endl;
621  std::cerr << e.what();
622  }
623 
624  // draw detectors if option is set, only important for wire hits ----------------------
625  if(drawDetectors_) {
626 
627  if(wire_hit) {
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));
631 
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()); // move the tube to the right place and rotate it correctly
636  TVector3 move = v*(v*(track_pos-o)); // move the tube along the wire until the track goes through it
637  TGeoCombiTrans det_trans(o(0) + move.X(),
638  o(1) + move.Y(),
639  o(2) + move.Z(),
640  &det_rot);
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);
646  }
647  }
648 
649  }
650  // finished drawing detectors ---------------------------------------------------------
651 
652  if(drawHits_) {
653 
654  // draw planar hits, with distinction between strip and pixel hits ----------------
655  if (full_hit) {
656 
657  StateOnPlane dummy(rep);
658  StateOnPlane dummy2(TVectorD(rep->getDim()), static_cast<const FullMeasurement*>(m)->constructPlane(dummy), rep);
659  MeasuredStateOnPlane sop = *(static_cast<const FullMeasurement*>(m)->constructMeasurementsOnPlane(dummy2)[0]);
660  sop.getCov()*=errorScale_;
661 
662  MeasuredStateOnPlane prevSop(sop);
663  prevSop.extrapolateBy(-3);
664  makeLines(&sop, &prevSop, rep, kYellow, 1, false, true, 0, 0);
665 
666  prevSop = sop;
667  prevSop.extrapolateBy(3);
668  makeLines(&sop, &prevSop, rep, kYellow, 1, false, true, 0, 0);
669  }
670 
671  if(planar_hit) {
672  if(!planar_pixel_hit) {
673  TEveBox* hit_box;
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);
681  } else {
682  // calculate eigenvalues to draw error-ellipse ----------------------------
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));
690  // finished calcluating, got the values -----------------------------------
691 
692  // do autoscaling if necessary --------------------------------------------
693  if(drawAutoScale_) {
694  double min_cov = std::min(pseudo_res_0,pseudo_res_1);
695  if(min_cov < 1e-5) {
696  std::cout << "Track " << i << ", Hit " << j << ": Invalid covariance matrix (Eigenvalue < 1e-5), autoscaling not possible!" << std::endl;
697  } else {
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;
701  errorScale_ *= cor;
702  pseudo_res_0 *= cor;
703  pseudo_res_1 *= cor;
704  std::cout << " to " << errorScale_ << std::endl;
705  }
706  }
707  }
708  // finished autoscaling ---------------------------------------------------
709 
710  // calculate the semiaxis of the error ellipse ----------------------------
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);
716  // finished calculating ---------------------------------------------------
717 
718  // rotate and translate everything correctly ------------------------------
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);
724  // finished rotating and translating --------------------------------------
725 
726  cov_shape->SetMainColor(kYellow);
727  cov_shape->SetMainTransparency(0);
728  gEve->AddElement(cov_shape);
729  }
730  }
731  // finished drawing planar hits ---------------------------------------------------
732 
733  // draw spacepoint hits -----------------------------------------------------------
734  if(space_hit) {
735  {
736  // get eigenvalues of covariance to know how to draw the ellipsoid ------------
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);
747  // got everything we need -----------------------------------------------------
748 
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);
753 
754  if (! det_rot.IsValid()){
755  // hackish fix if eigenvectors are not orthonogonal
756  if (fabs(eVec2*eVec3) > 1.e-10)
757  eVec3 = eVec1.Cross(eVec2);
758 
759  det_rot.SetAngles(eVec1.Theta()*radDeg, eVec1.Phi()*radDeg,
760  eVec2.Theta()*radDeg, eVec2.Phi()*radDeg,
761  eVec3.Theta()*radDeg, eVec3.Phi()*radDeg);
762  }
763 
764  // set the scaled eigenvalues -------------------------------------------------
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));
768  if(drawScaleMan_) { // override again if necessary
769  pseudo_res_0 = errorScale_*0.5;
770  pseudo_res_1 = errorScale_*0.5;
771  pseudo_res_2 = errorScale_*0.5;
772  }
773  // finished scaling -----------------------------------------------------------
774 
775  // autoscale if necessary -----------------------------------------------------
776  if(drawAutoScale_) {
777  double min_cov = std::min(pseudo_res_0,std::min(pseudo_res_1,pseudo_res_2));
778  if(min_cov < 1e-5) {
779  std::cout << "Track " << i << ", Hit " << j << ": Invalid covariance matrix (Eigenvalue < 1e-5), autoscaling not possible!" << std::endl;
780  } else {
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;
784  errorScale_ *= cor;
785  pseudo_res_0 *= cor;
786  pseudo_res_1 *= cor;
787  pseudo_res_2 *= cor;
788  std::cout << " to " << errorScale_ << std::endl;
789 
790  }
791  }
792  }
793  // finished autoscaling -------------------------------------------------------
794 
795  // rotate and translate -------------------------------------------------------
796  TGeoGenTrans det_trans(o(0),o(1),o(2),
797  //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
798  //1/(pseudo_res_0),1/(pseudo_res_1),1/(pseudo_res_2),
799  pseudo_res_0, pseudo_res_1, pseudo_res_2,
800  &det_rot);
801  cov_shape->SetTransMatrix(det_trans);
802  // finished rotating and translating ------------------------------------------
803 
804  cov_shape->SetMainColor(kYellow);
805  cov_shape->SetMainTransparency(10);
806  gEve->AddElement(cov_shape);
807  }
808 
809 
810  {
811  // calculate eigenvalues to draw error-ellipse ----------------------------
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));
819  // finished calcluating, got the values -----------------------------------
820 
821  // do autoscaling if necessary --------------------------------------------
822  if(drawAutoScale_) {
823  double min_cov = std::min(pseudo_res_0,pseudo_res_1);
824  if(min_cov < 1e-5) {
825  std::cout << "Track " << i << ", Hit " << j << ": Invalid covariance matrix (Eigenvalue < 1e-5), autoscaling not possible!" << std::endl;
826  } else {
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;
830  errorScale_ *= cor;
831  pseudo_res_0 *= cor;
832  pseudo_res_1 *= cor;
833  std::cout << " to " << errorScale_ << std::endl;
834  }
835  }
836  }
837  // finished autoscaling ---------------------------------------------------
838 
839  // calculate the semiaxis of the error ellipse ----------------------------
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);
845  // finished calculating ---------------------------------------------------
846 
847  // rotate and translate everything correctly ------------------------------
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);
852  /*if (! det_rot.IsValid()){
853  u_semiaxis.Print();
854  v_semiaxis.Print();
855  norm.Print();
856  }*/
857  TGeoCombiTrans det_trans(pix_pos(0),pix_pos(1),pix_pos(2), &det_rot);
858  cov_shape->SetTransMatrix(det_trans);
859  // finished rotating and translating --------------------------------------
860 
861  cov_shape->SetMainColor(kYellow);
862  cov_shape->SetMainTransparency(0);
863  gEve->AddElement(cov_shape);
864  }
865  }
866  // finished drawing spacepoint hits -----------------------------------------------
867 
868  // draw wire hits -----------------------------------------------------------------
869  if(wire_hit) {
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));
875 
876  // autoscale if necessary -----------------------------------------------------
877  if(drawAutoScale_) {
878  if(pseudo_res_0 < 1e-5) {
879  std::cout << "Track " << i << ", Hit " << j << ": Invalid wire resolution (< 1e-5), autoscaling not possible!" << std::endl;
880  } else {
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;
884  errorScale_ *= cor;
885  pseudo_res_0 *= cor;
886  std::cout << " to " << errorScale_ << std::endl;
887  }
888  }
889 
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;
892  } else {
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;
896  errorScale_ *= cor;
897  pseudo_res_1 *= cor;
898  std::cout << " to " << errorScale_ << std::endl;
899  }
900  }
901  }
902  // finished autoscaling -------------------------------------------------------
903 
904  TEveBox* hit_box;
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);
910 
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);
915  }
916  // finished drawing wire hits -----------------------------------------------------
917 
918  } // finished drawing hits
919 
920  } // finished looping over MeasurmentOnPlanes
921 
922 
923  prevFi = fi;
924  prevFittedState = fittedState;
925 
926  }
927 
928  }
929 
930  gEve->Redraw3D(resetCam);
931 
932 }
933 
934 
935 
936 
937 TEveBox* EventDisplay::boxCreator(TVector3 o, TVector3 u, TVector3 v, float ud, float vd, float depth) {
938 
939  TEveBox* box = new TEveBox("detPlane_shape");
940  float vertices[24];
941 
942  TVector3 norm = u.Cross(v);
943  u *= (0.5*ud);
944  v *= (0.5*vd);
945  norm *= (0.5*depth);
946 
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);
950 
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);
954 
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);
958 
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);
962 
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);
966 
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);
970 
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);
974 
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);
978 
979 
980  for(int k = 0; k < 24; k += 3) box->SetVertex((k/3), vertices[k], vertices[k+1], vertices[k+2]);
981 
982  return box;
983 
984 }
985 
986 
987 void EventDisplay::makeLines(const StateOnPlane* prevState, const StateOnPlane* state, const AbsTrackRep* rep,
988  const Color_t& color, const Style_t& style, bool drawMarkers, bool drawErrors, double lineWidth, int markerPos)
989 {
990  if (prevState == nullptr || state == nullptr) {
991  std::cerr << "prevState == nullptr || state == nullptr\n";
992  return;
993  }
994 
995  TVector3 pos, dir, oldPos, oldDir;
996  rep->getPosDir(*state, pos, dir);
997  rep->getPosDir(*prevState, oldPos, oldDir);
998 
999  double distA = (pos-oldPos).Mag();
1000  double distB = distA;
1001  if ((pos-oldPos)*oldDir < 0)
1002  distA *= -1.;
1003  if ((pos-oldPos)*dir < 0)
1004  distB *= -1.;
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);
1014  if (drawMarkers) {
1015  if (markerPos == 0)
1016  lineSet->AddMarker(oldPos(0), oldPos(1), oldPos(2));
1017  else
1018  lineSet->AddMarker(pos(0), pos(1), pos(2));
1019  }
1020 
1021  if (lineWidth > 0)
1022  gEve->AddElement(lineSet);
1023 
1024 
1025  if (drawErrors) {
1026  const MeasuredStateOnPlane* measuredState;
1027  if (markerPos == 0)
1028  measuredState = dynamic_cast<const MeasuredStateOnPlane*>(prevState);
1029  else
1030  measuredState = dynamic_cast<const MeasuredStateOnPlane*>(state);
1031 
1032  if (measuredState != nullptr) {
1033 
1034  // step for evaluate at a distance from the original plane
1035  TVector3 eval;
1036  if (markerPos == 0) {
1037  if (fabs(distA) < 1.) {
1038  distA < 0 ? distA = -1 : distA = 1;
1039  }
1040  eval = 0.2 * distA * oldDir;
1041  }
1042  else {
1043  if (fabs(distB) < 1.) {
1044  distB < 0 ? distB = -1 : distB = 1;
1045  }
1046  eval = -0.2 * distB * dir;
1047  }
1048 
1049 
1050  // get cov at first plane
1051  TMatrixDSym cov;
1052  TVector3 position, direction;
1053  rep->getPosMomCov(*measuredState, position, direction, cov);
1054 
1055  // get eigenvalues & -vectors
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;
1060  // limit
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);
1065 
1066  // get two largest eigenvalues/-vectors
1067  if (ev0 < ev1 && ev0 < ev2) {
1068  eVec1.SetXYZ(eVec(0,1),eVec(1,1),eVec(2,1));
1069  eVec1 *= sqrt(ev1);
1070  eVec2.SetXYZ(eVec(0,2),eVec(1,2),eVec(2,2));
1071  eVec2 *= sqrt(ev2);
1072  }
1073  else if (ev1 < ev0 && ev1 < ev2) {
1074  eVec1.SetXYZ(eVec(0,0),eVec(1,0),eVec(2,0));
1075  eVec1 *= sqrt(ev0);
1076  eVec2.SetXYZ(eVec(0,2),eVec(1,2),eVec(2,2));
1077  eVec2 *= sqrt(ev2);
1078  }
1079  else {
1080  eVec1.SetXYZ(eVec(0,0),eVec(1,0),eVec(2,0));
1081  eVec1 *= sqrt(ev0);
1082  eVec2.SetXYZ(eVec(0,1),eVec(1,1),eVec(2,1));
1083  eVec2 *= sqrt(ev1);
1084  }
1085 
1086  if (eVec1.Cross(eVec2)*eval < 0)
1087  eVec2 *= -1;
1088  //assert(eVec1.Cross(eVec2)*eval > 0);
1089 
1090  const TVector3 oldEVec1(eVec1);
1091  const TVector3 oldEVec2(eVec2);
1092 
1093  const int nEdges = 24;
1094  std::vector<TVector3> vertices;
1095 
1096  vertices.push_back(position);
1097 
1098  // vertices at plane
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);
1102  }
1103 
1104 
1105 
1106  DetPlane* newPlane = new DetPlane(*(measuredState->getPlane()));
1107  newPlane->setO(position + eval);
1108 
1109  MeasuredStateOnPlane stateCopy(*measuredState);
1110  try{
1111  rep->extrapolateToPlane(stateCopy, SharedPlanePtr(newPlane));
1112  }
1113  catch(Exception& e){
1114  std::cerr<<e.what();
1115  return;
1116  }
1117 
1118  // get cov at 2nd plane
1119  rep->getPosMomCov(stateCopy, position, direction, cov);
1120 
1121  // get eigenvalues & -vectors
1122  TMatrixDSymEigen eigen_values2(cov.GetSub(0,2, 0,2));
1123  ev = eigen_values2.GetEigenValues();
1124  eVec = eigen_values2.GetEigenVectors();
1125  // limit
1126  ev0 = std::min(ev(0), maxErr);
1127  ev1 = std::min(ev(1), maxErr);
1128  ev2 = std::min(ev(2), maxErr);
1129 
1130  // get two largest eigenvalues/-vectors
1131  if (ev0 < ev1 && ev0 < ev2) {
1132  eVec1.SetXYZ(eVec(0,1),eVec(1,1),eVec(2,1));
1133  eVec1 *= sqrt(ev1);
1134  eVec2.SetXYZ(eVec(0,2),eVec(1,2),eVec(2,2));
1135  eVec2 *= sqrt(ev2);
1136  }
1137  else if (ev1 < ev0 && ev1 < ev2) {
1138  eVec1.SetXYZ(eVec(0,0),eVec(1,0),eVec(2,0));
1139  eVec1 *= sqrt(ev0);
1140  eVec2.SetXYZ(eVec(0,2),eVec(1,2),eVec(2,2));
1141  eVec2 *= sqrt(ev2);
1142  }
1143  else {
1144  eVec1.SetXYZ(eVec(0,0),eVec(1,0),eVec(2,0));
1145  eVec1 *= sqrt(ev0);
1146  eVec2.SetXYZ(eVec(0,1),eVec(1,1),eVec(2,1));
1147  eVec2 *= sqrt(ev1);
1148  }
1149 
1150  if (eVec1.Cross(eVec2)*eval < 0)
1151  eVec2 *= -1;
1152  //assert(eVec1.Cross(eVec2)*eval > 0);
1153 
1154  if (oldEVec1*eVec1 < 0) {
1155  eVec1 *= -1;
1156  eVec2 *= -1;
1157  }
1158 
1159  // vertices at 2nd plane
1160  double angle0 = eVec1.Angle(oldEVec1);
1161  if (eVec1*(eval.Cross(oldEVec1)) < 0)
1162  angle0 *= -1;
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);
1166  }
1167 
1168  vertices.push_back(position);
1169 
1170 
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());
1174  }
1175 
1176  assert(vertices.size() == 2*nEdges+2);
1177 
1178  int iTri(0);
1179  for (int i=0; i<nEdges; ++i) {
1180  //error_shape->SetTriangle(iTri++, 0, i+1, (i+1)%nEdges+1);
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);
1183  //error_shape->SetTriangle(iTri++, 2*nEdges+1, i+1+nEdges, (i+1)%nEdges+1+nEdges);
1184  }
1185 
1186  //assert(iTri == nEdges*4);
1187 
1188  error_shape->SetMainColor(color);
1189  error_shape->SetMainTransparency(25);
1190  gEve->AddElement(error_shape);
1191  }
1192  }
1193 }
1194 
1195 
1196 void EventDisplay::makeGui() {
1197 
1198  TEveBrowser* browser = gEve->GetBrowser();
1199  browser->StartEmbedding(TRootBrowser::kLeft);
1200 
1201  TGMainFrame* frmMain = new TGMainFrame(gClient->GetRoot(), 1000, 600);
1202  frmMain->SetWindowName("XX GUI");
1203  frmMain->SetCleanup(kDeepCleanup);
1204 
1205  TGLabel* lbl = 0;
1206  TGTextButton* tb = 0;
1207  EventDisplay* fh = EventDisplay::getInstance();
1208 
1209  TGHorizontalFrame* hf = new TGHorizontalFrame(frmMain); {
1210  // evt number entry
1211  lbl = new TGLabel(hf, "Go to event: ");
1212  hf->AddFrame(lbl);
1213  guiEvent = new TGNumberEntry(hf, 0, 9,999, TGNumberFormat::kNESInteger,
1214  TGNumberFormat::kNEANonNegative,
1215  TGNumberFormat::kNELLimitMinMax,
1216  0, 99999);
1217  hf->AddFrame(guiEvent);
1218  guiEvent->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiGoto()");
1219 
1220  // redraw button
1221  tb = new TGTextButton(hf, "Redraw Event");
1222  hf->AddFrame(tb);
1223  tb->Connect("Clicked()", "genfit::EventDisplay", fh, "guiGoto()");
1224  }
1225  frmMain->AddFrame(hf);
1226 
1227  // draw options
1228  hf = new TGHorizontalFrame(frmMain); {
1229  lbl = new TGLabel(hf, "\n Draw Options");
1230  hf->AddFrame(lbl);
1231  }
1232  frmMain->AddFrame(hf);
1233 
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()");
1239  }
1240  frmMain->AddFrame(hf);
1241 
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()");
1247  }
1248  frmMain->AddFrame(hf);
1249 
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()");
1255  }
1256  frmMain->AddFrame(hf);
1257 
1258 
1259 
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()");
1265  }
1266  frmMain->AddFrame(hf);
1267 
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()");
1273  }
1274  frmMain->AddFrame(hf);
1275 
1276 
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()");
1282  }
1283  frmMain->AddFrame(hf);
1284 
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()");
1290  }
1291  frmMain->AddFrame(hf);
1292 
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()");
1298  }
1299  frmMain->AddFrame(hf);
1300 
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()");
1306  }
1307  frmMain->AddFrame(hf);
1308 
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()");
1314  }
1315  frmMain->AddFrame(hf);
1316 
1317 
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()");
1323  }
1324  frmMain->AddFrame(hf);
1325 
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()");
1331  }
1332  frmMain->AddFrame(hf);
1333 
1334  hf = new TGHorizontalFrame(frmMain); {
1335  guiErrorScale_ = new TGNumberEntry(hf, errorScale_, 6,999, TGNumberFormat::kNESReal,
1336  TGNumberFormat::kNEANonNegative,
1337  TGNumberFormat::kNELLimitMinMax,
1338  1.E-4, 1.E5);
1339  hf->AddFrame(guiErrorScale_);
1340  guiErrorScale_->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1341  lbl = new TGLabel(hf, "Error scale");
1342  hf->AddFrame(lbl);
1343  }
1344  frmMain->AddFrame(hf);
1345 
1346 
1347 
1348  hf = new TGHorizontalFrame(frmMain); {
1349  lbl = new TGLabel(hf, "\n TrackRep options");
1350  hf->AddFrame(lbl);
1351  }
1352  frmMain->AddFrame(hf);
1353 
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()");
1359  }
1360  frmMain->AddFrame(hf);
1361 
1362  hf = new TGHorizontalFrame(frmMain); {
1363  guiRepId_ = new TGNumberEntry(hf, repId_, 6,999, TGNumberFormat::kNESInteger,
1364  TGNumberFormat::kNEANonNegative,
1365  TGNumberFormat::kNELLimitMinMax,
1366  0, 99);
1367  hf->AddFrame(guiRepId_);
1368  guiRepId_->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1369  lbl = new TGLabel(hf, "Else draw rep with id");
1370  hf->AddFrame(lbl);
1371  }
1372  frmMain->AddFrame(hf);
1373 
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()");
1379  }
1380  frmMain->AddFrame(hf);
1381 
1382  hf = new TGHorizontalFrame(frmMain); {
1383  guiTrackId_ = new TGNumberEntry(hf, trackId_, 6,999, TGNumberFormat::kNESInteger,
1384  TGNumberFormat::kNEANonNegative,
1385  TGNumberFormat::kNELLimitMinMax,
1386  0, 99);
1387  hf->AddFrame(guiTrackId_);
1388  guiTrackId_->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1389  lbl = new TGLabel(hf, "Else draw track nr. ");
1390  hf->AddFrame(lbl);
1391  }
1392  frmMain->AddFrame(hf);
1393 
1394 
1395 
1396  frmMain->MapSubwindows();
1397  frmMain->Resize();
1398  frmMain->MapWindow();
1399 
1400  browser->StopEmbedding();
1401  browser->SetTabTitle("Draw Control", 0);
1402 
1403 
1404  browser->StartEmbedding(TRootBrowser::kLeft);
1405  TGMainFrame* frmMain2 = new TGMainFrame(gClient->GetRoot(), 1000, 600);
1406  frmMain2->SetWindowName("XX GUI");
1407  frmMain2->SetCleanup(kDeepCleanup);
1408 
1409  hf = new TGHorizontalFrame(frmMain2); {
1410  // evt number entry
1411  lbl = new TGLabel(hf, "Go to event: ");
1412  hf->AddFrame(lbl);
1413  guiEvent2 = new TGNumberEntry(hf, 0, 9,999, TGNumberFormat::kNESInteger,
1414  TGNumberFormat::kNEANonNegative,
1415  TGNumberFormat::kNELLimitMinMax,
1416  0, 99999);
1417  hf->AddFrame(guiEvent2);
1418  guiEvent2->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiGoto2()");
1419 
1420  // redraw button
1421  tb = new TGTextButton(hf, "Redraw Event");
1422  hf->AddFrame(tb);
1423  tb->Connect("Clicked()", "genfit::EventDisplay", fh, "guiGoto()");
1424  }
1425  frmMain2->AddFrame(hf);
1426 
1427  hf = new TGHorizontalFrame(frmMain2); {
1428  lbl = new TGLabel(hf, "\n Fitting options");
1429  hf->AddFrame(lbl);
1430  }
1431  frmMain2->AddFrame(hf);
1432 
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()");
1438  }
1439  frmMain2->AddFrame(hf);
1440 
1441  hf = new TGHorizontalFrame(frmMain2); {
1442  guiDebugLvl_ = new TGNumberEntry(hf, debugLvl_, 6,999, TGNumberFormat::kNESInteger,
1443  TGNumberFormat::kNEANonNegative,
1444  TGNumberFormat::kNELLimitMinMax,
1445  0, 999);
1446  hf->AddFrame(guiDebugLvl_);
1447  guiDebugLvl_->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1448  lbl = new TGLabel(hf, "debug level");
1449  hf->AddFrame(lbl);
1450  }
1451  frmMain2->AddFrame(hf);
1452 
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();
1463  }
1464  frmMain2->AddFrame(hf);
1465 
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();
1482  }
1483  frmMain2->AddFrame(hf);
1484 
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()");
1490  }
1491  frmMain2->AddFrame(hf);
1492 
1493  hf = new TGHorizontalFrame(frmMain2); {
1494  guiDPVal_ = new TGNumberEntry(hf, dPVal_, 6,9999, TGNumberFormat::kNESReal,
1495  TGNumberFormat::kNEANonNegative,
1496  TGNumberFormat::kNELLimitMinMax,
1497  0, 999);
1498  hf->AddFrame(guiDPVal_);
1499  guiDPVal_->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1500  lbl = new TGLabel(hf, "delta pVal (convergence criterium)");
1501  hf->AddFrame(lbl);
1502  }
1503  frmMain2->AddFrame(hf);
1504 
1505  hf = new TGHorizontalFrame(frmMain2); {
1506  guiRelChi2_ = new TGNumberEntry(hf, dRelChi2_, 6,9999, TGNumberFormat::kNESReal,
1507  TGNumberFormat::kNEANonNegative,
1508  TGNumberFormat::kNELLimitMinMax,
1509  0, 999);
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)");
1513  hf->AddFrame(lbl);
1514  }
1515  frmMain2->AddFrame(hf);
1516 
1517  hf = new TGHorizontalFrame(frmMain2); {
1518  guiDChi2Ref_ = new TGNumberEntry(hf, dChi2Ref_, 6,9999, TGNumberFormat::kNESReal,
1519  TGNumberFormat::kNEANonNegative,
1520  TGNumberFormat::kNELLimitMinMax,
1521  0, 999);
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)");
1525  hf->AddFrame(lbl);
1526  }
1527  frmMain2->AddFrame(hf);
1528 
1529  hf = new TGHorizontalFrame(frmMain2); {
1530  guiNMinIter_ = new TGNumberEntry(hf, nMinIter_, 6,999, TGNumberFormat::kNESInteger,
1531  TGNumberFormat::kNEANonNegative,
1532  TGNumberFormat::kNELLimitMinMax,
1533  1, 100);
1534  hf->AddFrame(guiNMinIter_);
1535  guiNMinIter_->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1536  lbl = new TGLabel(hf, "Minimum nr of iterations");
1537  hf->AddFrame(lbl);
1538  }
1539  frmMain2->AddFrame(hf);
1540 
1541  hf = new TGHorizontalFrame(frmMain2); {
1542  guiNMaxIter_ = new TGNumberEntry(hf, nMaxIter_, 6,999, TGNumberFormat::kNESInteger,
1543  TGNumberFormat::kNEANonNegative,
1544  TGNumberFormat::kNELLimitMinMax,
1545  1, 100);
1546  hf->AddFrame(guiNMaxIter_);
1547  guiNMaxIter_->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1548  lbl = new TGLabel(hf, "Maximum nr of iterations");
1549  hf->AddFrame(lbl);
1550  }
1551  frmMain2->AddFrame(hf);
1552 
1553  hf = new TGHorizontalFrame(frmMain2); {
1554  guiNMaxFailed_ = new TGNumberEntry(hf, nMaxFailed_, 6,999, TGNumberFormat::kNESInteger,
1555  TGNumberFormat::kNEAAnyNumber,
1556  TGNumberFormat::kNELLimitMinMax,
1557  -1, 1000);
1558  hf->AddFrame(guiNMaxFailed_);
1559  guiNMaxFailed_->Connect("ValueSet(Long_t)", "genfit::EventDisplay", fh, "guiSetDrawParams()");
1560  lbl = new TGLabel(hf, "Maximum nr of failed hits");
1561  hf->AddFrame(lbl);
1562  }
1563  frmMain2->AddFrame(hf);
1564 
1565 
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()");
1571  }
1572  frmMain2->AddFrame(hf);
1573 
1574 
1575 
1576 
1577  frmMain2->MapSubwindows();
1578  frmMain2->Resize();
1579  frmMain2->MapWindow();
1580 
1581  browser->StopEmbedding();
1582  browser->SetTabTitle("Refit Control", 0);
1583 }
1584 
1585 
1586 void EventDisplay::guiGoto(){
1587  Long_t n = guiEvent->GetNumberEntry()->GetIntNumber();
1588  guiEvent2->SetIntNumber(n);
1589  gotoEvent(n);
1590 }
1591 
1592 void EventDisplay::guiGoto2(){
1593  Long_t n = guiEvent2->GetNumberEntry()->GetIntNumber();
1594  guiEvent->SetIntNumber(n);
1595  gotoEvent(n);
1596 }
1597 
1598 
1599 void EventDisplay::guiSetDrawParams(){
1600 
1601  drawGeometry_ = guiDrawGeometry_->IsOn();
1602  drawDetectors_ = guiDrawDetectors_->IsOn();
1603  drawHits_ = guiDrawHits_->IsOn();
1604  drawErrors_ = guiDrawErrors_->IsOn();
1605 
1606  drawPlanes_ = guiDrawPlanes_->IsOn();
1607  drawTrackMarkers_ = guiDrawTrackMarkers_->IsOn();
1608  drawTrack_ = guiDrawTrack_->IsOn();
1609  drawRefTrack_ = guiDrawRefTrack_->IsOn();
1610  drawForward_ = guiDrawForward_->IsOn();
1611  drawBackward_ = guiDrawBackward_->IsOn();
1612 
1613  drawAutoScale_ = guiDrawAutoScale_->IsOn();
1614  drawScaleMan_ = guiDrawScaleMan_->IsOn();
1615 
1616  errorScale_ = guiErrorScale_->GetNumberEntry()->GetNumber();
1617 
1618  drawCardinalRep_ = guiDrawCardinalRep_->IsOn();
1619  repId_ = guiRepId_->GetNumberEntry()->GetNumber();
1620 
1621  drawAllTracks_ = guiDrawAllTracks_->IsOn();
1622  trackId_ = guiTrackId_->GetNumberEntry()->GetNumber();
1623 
1624 
1625  refit_ = guiRefit_->IsOn();
1626  debugLvl_ = guiDebugLvl_->GetNumberEntry()->GetNumber();
1627 
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();
1636 
1637  gotoEvent(eventId_);
1638 }
1639 
1640 
1641 void EventDisplay::guiSelectFitterId(int val){
1642  fitterId_ = eFitterType(val-1);
1643  gotoEvent(eventId_);
1644 }
1645 
1646 void EventDisplay::guiSelectMmHandling(int val){
1647  mmHandling_ = eMultipleMeasurementHandling(val-1);
1648  gotoEvent(eventId_);
1649 }
1650 
1651 
1652 } // end of namespace genfit
This class collects all information needed and produced by a specific AbsFitter and is specific to on...
Definition: AbsFitterInfo.h:42
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.
Definition: AbsTrackRep.h:66
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.
Definition: AbsTrackRep.h:252
virtual double extrapolateToPlane(StateOnPlane &state, const genfit::SharedPlanePtr &plane, bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0
Extrapolates the state to plane, and returns the extrapolation length and, via reference,...
void setPropDir(int dir)
Set propagation direction. (-1, 0, 1) -> (backward, auto, forward).
Definition: AbsTrackRep.h:336
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.
Definition: DAF.h:49
Event display designed to run with Genfit.
Definition: EventDisplay.h:59
Exception class for error handling in GENFIT (provides storage for diagnostic information)
Definition: Exception.h:48
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.
Definition: KalmanFitter.h:40
#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)
Definition: RKTrackRep.h:72
Class for measurements implementing a space point hit geometry.
A state with arbitrary dimension defined in a DetPlane.
Definition: StateOnPlane.h:47
Object containing AbsMeasurement and AbsFitterInfo objects.
Definition: TrackPoint.h:46
AbsFitterInfo * getFitterInfo(const AbsTrackRep *rep=nullptr) const
Get fitterInfo for rep. Per default, use cardinal rep.
Definition: TrackPoint.cc:170
Collection of TrackPoint objects, AbsTrackRep objects and FitStatus objects.
Definition: Track.h:71
Class for measurements in wire detectors (Straw tubes and drift chambers) which can measure the coord...
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
double eval(const std::vector< double > &spl, const std::vector< double > &vals, double x)
Evaluate spline (zero order or first order) in point x.
Definition: tools.h:115
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