File indexing completed on 2025-08-05 08:18:23
0001
0002 #include "EventDisplay.h"
0003
0004 #include <assert.h>
0005 #include <algorithm>
0006 #include <cmath>
0007 #include <exception>
0008 #include <iostream>
0009 #include <sys/time.h>
0010
0011 #include "AbsMeasurement.h"
0012 #include "FullMeasurement.h"
0013 #include "PlanarMeasurement.h"
0014 #include "ProlateSpacepointMeasurement.h"
0015 #include "SpacepointMeasurement.h"
0016 #include "WireMeasurement.h"
0017 #include "WirePointMeasurement.h"
0018 #include "AbsTrackRep.h"
0019 #include "ConstField.h"
0020 #include "DetPlane.h"
0021 #include "Exception.h"
0022 #include "FieldManager.h"
0023 #include "Tools.h"
0024 #include "KalmanFitterInfo.h"
0025 #include "KalmanFitter.h"
0026 #include "DAF.h"
0027 #include "KalmanFitterRefTrack.h"
0028 #include "RKTrackRep.h"
0029
0030 #include <TApplication.h>
0031 #include <TEveBrowser.h>
0032 #include <TEveManager.h>
0033 #include <TEveEventManager.h>
0034 #include <TEveGeoNode.h>
0035 #include <TEveGeoShape.h>
0036 #include <TEveStraightLineSet.h>
0037 #include <TEveTriangleSet.h>
0038 #include <TDecompSVD.h>
0039 #include <TGButton.h>
0040 #include <TGLabel.h>
0041 #include <TGNumberEntry.h>
0042 #include <TGeoEltu.h>
0043 #include <TGeoManager.h>
0044 #include <TGeoMatrix.h>
0045 #include <TGeoNode.h>
0046 #include <TGeoSphere.h>
0047 #include <TGeoTube.h>
0048 #include <TMath.h>
0049 #include <TMatrixT.h>
0050 #include <TMatrixTSym.h>
0051 #include <TMatrixDSymEigen.h>
0052 #include <TROOT.h>
0053 #include <TVector2.h>
0054 #include <TVectorD.h>
0055 #include <TSystem.h>
0056
0057 #include <memory>
0058
0059 namespace genfit {
0060
0061
0062 EventDisplay* EventDisplay::eventDisplay_ = nullptr;
0063
0064 EventDisplay::EventDisplay() :
0065 errorScale_(1.),
0066 drawGeometry_(false),
0067 drawDetectors_(true),
0068 drawHits_(true),
0069 drawErrors_(true),
0070 drawPlanes_(true),
0071 drawTrackMarkers_(true),
0072 drawTrack_(true),
0073 drawRefTrack_(true),
0074 drawForward_(true),
0075 drawBackward_(true),
0076 drawAutoScale_(true),
0077 drawScaleMan_(false),
0078 drawSilent_(false),
0079 drawCardinalRep_(true),
0080 repId_(0),
0081 drawAllTracks_(true),
0082 trackId_(0),
0083 refit_(false),
0084 debugLvl_(0),
0085 fitterId_(SimpleKalman),
0086 mmHandling_(weightedAverage),
0087 squareRootFormalism_(false),
0088 dPVal_(1.E-3),
0089 dRelChi2_(0.2),
0090 dChi2Ref_(1.),
0091 nMinIter_(2),
0092 nMaxIter_(4),
0093 nMaxFailed_(-1),
0094 resort_(false)
0095 {
0096
0097 if((!gApplication) || (gApplication && gApplication->TestBit(TApplication::kDefaultApplication))) {
0098 std::cout << "In EventDisplay ctor: gApplication not found, creating..." << std::flush;
0099 new TApplication("ROOT_application", 0, 0);
0100 std::cout << "done!" << std::endl;
0101 }
0102 if(!gEve) {
0103 std::cout << "In EventDisplay ctor: gEve not found, creating..." << std::flush;
0104 TEveManager::Create();
0105 std::cout << "done!" << std::endl;
0106 }
0107
0108 eventId_ = 0;
0109
0110 }
0111
0112 void EventDisplay::setOptions(std::string opts) {
0113
0114 if(opts != "") {
0115 for(size_t i = 0; i < opts.length(); ++i) {
0116 if(opts[i] == 'A') drawAutoScale_ = true;
0117 if(opts[i] == 'B') drawBackward_ = true;
0118 if(opts[i] == 'D') drawDetectors_ = true;
0119 if(opts[i] == 'E') drawErrors_ = true;
0120 if(opts[i] == 'F') drawForward_ = true;
0121 if(opts[i] == 'H') drawHits_ = true;
0122 if(opts[i] == 'M') drawTrackMarkers_ = true;
0123 if(opts[i] == 'P') drawPlanes_ = true;
0124 if(opts[i] == 'S') drawScaleMan_ = true;
0125 if(opts[i] == 'T') drawTrack_ = true;
0126 if(opts[i] == 'X') drawSilent_ = true;
0127 if(opts[i] == 'G') drawGeometry_ = true;
0128 }
0129 }
0130
0131 }
0132
0133 void EventDisplay::setErrScale(double errScale) { errorScale_ = errScale; }
0134
0135 double EventDisplay::getErrScale() { return errorScale_; }
0136
0137 EventDisplay* EventDisplay::getInstance() {
0138
0139 if(eventDisplay_ == nullptr) {
0140 eventDisplay_ = new EventDisplay();
0141 }
0142 return eventDisplay_;
0143
0144 }
0145
0146 EventDisplay::~EventDisplay() { reset(); }
0147
0148 void EventDisplay::reset() {
0149
0150 for(unsigned int i = 0; i < events_.size(); i++) {
0151
0152 for(unsigned int j = 0; j < events_[i]->size(); j++) {
0153
0154 delete events_[i]->at(j);
0155
0156 }
0157 delete events_[i];
0158 }
0159
0160 events_.clear();
0161 }
0162
0163
0164 void EventDisplay::addEvent(std::vector<Track*>& tracks) {
0165
0166 std::vector<Track*>* vec = new std::vector<Track*>;
0167
0168 for(unsigned int i = 0; i < tracks.size(); i++) {
0169 vec->push_back(new Track(*(tracks[i])));
0170 }
0171
0172 events_.push_back(vec);
0173 }
0174
0175
0176 void EventDisplay::addEvent(std::vector<const Track*>& tracks) {
0177
0178 std::vector<Track*>* vec = new std::vector<Track*>;
0179
0180 for(unsigned int i = 0; i < tracks.size(); i++) {
0181 vec->push_back(new Track(*(tracks[i])));
0182 }
0183
0184 events_.push_back(vec);
0185 }
0186
0187
0188 void EventDisplay::addEvent(const Track* tr) {
0189
0190 std::vector<Track*>* vec = new std::vector<Track*>;
0191 vec->push_back(new Track(*tr));
0192 events_.push_back(vec);
0193 }
0194
0195
0196 void EventDisplay::next(unsigned int stp) {
0197
0198 gotoEvent(eventId_ + stp);
0199
0200 }
0201
0202 void EventDisplay::prev(unsigned int stp) {
0203
0204 if(events_.size() == 0) return;
0205 if(eventId_ < stp) {
0206 gotoEvent(0);
0207 } else {
0208 gotoEvent(eventId_ - stp);
0209 }
0210
0211 }
0212
0213 int EventDisplay::getNEvents() { return events_.size(); }
0214
0215
0216 void EventDisplay::gotoEvent(unsigned int id) {
0217
0218 if (events_.size() == 0)
0219 return;
0220 else if(id >= events_.size())
0221 id = events_.size() - 1;
0222
0223 bool resetCam = true;
0224
0225 if (id == eventId_)
0226 resetCam = false;
0227
0228 eventId_ = id;
0229
0230 std::cout << "At event " << id << std::endl;
0231 if (gEve->GetCurrentEvent()) {
0232 gEve->GetCurrentEvent()->DestroyElements();
0233 }
0234 double old_error_scale = errorScale_;
0235 drawEvent(eventId_, resetCam);
0236 if(old_error_scale != errorScale_) {
0237 if (gEve->GetCurrentEvent()) {
0238 gEve->GetCurrentEvent()->DestroyElements();
0239 }
0240 drawEvent(eventId_, resetCam);
0241 }
0242 errorScale_ = old_error_scale;
0243
0244 }
0245
0246 void EventDisplay::open() {
0247
0248 std::cout << "EventDisplay::open(); " << getNEvents() << " events loaded" << std::endl;
0249
0250 if(getNEvents() > 0) {
0251 double old_error_scale = errorScale_;
0252 drawEvent(0);
0253 if(old_error_scale != errorScale_) {
0254 std::cout << "autoscaling changed the error, draw again." << std::endl;
0255 gotoEvent(0);
0256 }
0257 errorScale_ = old_error_scale;
0258 }
0259
0260
0261 if(!drawSilent_) {
0262 makeGui();
0263 gApplication->Run(kTRUE);
0264 }
0265
0266 std::cout << "opened" << std::endl;
0267
0268 }
0269
0270
0271 void EventDisplay::drawEvent(unsigned int id, bool resetCam) {
0272
0273 std::cout << "EventDisplay::drawEvent(" << id << ")" << std::endl;
0274
0275
0276
0277 if(drawGeometry_) {
0278 TGeoNode* top_node = gGeoManager->GetTopNode();
0279 assert(top_node != nullptr);
0280
0281
0282 TObjArray* volumes = gGeoManager->GetListOfVolumes();
0283 for(int i = 0; i < volumes->GetEntriesFast(); i++) {
0284 TGeoVolume* volume = dynamic_cast<TGeoVolume*>(volumes->At(i));
0285 assert(volume != nullptr);
0286 volume->SetLineColor(12);
0287 volume->SetTransparency(50);
0288 }
0289
0290 TEveGeoTopNode* eve_top_node = new TEveGeoTopNode(gGeoManager, top_node);
0291 eve_top_node->IncDenyDestroy();
0292 gEve->AddElement(eve_top_node);
0293 }
0294
0295
0296 for(unsigned int i = 0; i < events_.at(id)->size(); i++) {
0297
0298 if (!drawAllTracks_ && trackId_ != i)
0299 continue;
0300
0301 Track* track = events_[id]->at(i);
0302 try {
0303 track->checkConsistency();
0304 } catch (genfit::Exception& e) {
0305 std::cerr<< e.getExcString() <<std::endl;
0306 continue;
0307 }
0308
0309 std::unique_ptr<Track> refittedTrack(nullptr);
0310 if (refit_) {
0311
0312 std::cout << "Refit track:" << std::endl;
0313
0314 std::unique_ptr<AbsKalmanFitter> fitter;
0315 switch (fitterId_) {
0316 case SimpleKalman:
0317 fitter.reset(new KalmanFitter(nMaxIter_, dPVal_));
0318 fitter->setMultipleMeasurementHandling(mmHandling_);
0319 (static_cast<KalmanFitter*>(fitter.get()))->useSquareRootFormalism(squareRootFormalism_);
0320 break;
0321
0322 case RefKalman:
0323 fitter.reset(new KalmanFitterRefTrack(nMaxIter_, dPVal_));
0324 fitter->setMultipleMeasurementHandling(mmHandling_);
0325 static_cast<KalmanFitterRefTrack*>(fitter.get())->setDeltaChi2Ref(dChi2Ref_);
0326 break;
0327
0328 case DafSimple:
0329 fitter.reset(new DAF(false));
0330 ( static_cast<KalmanFitter*>( (static_cast<DAF*>(fitter.get()))->getKalman() ) )->useSquareRootFormalism(squareRootFormalism_);
0331 break;
0332 case DafRef:
0333 fitter.reset(new DAF());
0334 ( static_cast<KalmanFitterRefTrack*>( (static_cast<DAF*>(fitter.get()))->getKalman() ) )->setDeltaChi2Ref(dChi2Ref_);
0335 break;
0336
0337 }
0338 fitter->setDebugLvl(std::max(0, (int)debugLvl_-1));
0339 fitter->setMinIterations(nMinIter_);
0340 fitter->setMaxIterations(nMaxIter_);
0341 fitter->setRelChi2Change(dRelChi2_);
0342 fitter->setMaxFailedHits(nMaxFailed_);
0343
0344
0345 refittedTrack.reset(new Track(*track));
0346 refittedTrack->deleteFitterInfo();
0347
0348 if (debugLvl_>0)
0349 refittedTrack->Print("C");
0350
0351 timeval startcputime, endcputime;
0352
0353 try{
0354 gettimeofday(&startcputime, nullptr);
0355 fitter->processTrack(refittedTrack.get(), resort_);
0356 gettimeofday(&endcputime, nullptr);
0357 }
0358 catch(genfit::Exception& e){
0359 std::cerr << e.what();
0360 std::cerr << "Exception, could not refit track" << std::endl;
0361 continue;
0362 }
0363
0364 int microseconds = 1000000*(endcputime.tv_sec - startcputime.tv_sec) + (endcputime.tv_usec - startcputime.tv_usec);
0365 std::cout << "it took " << double(microseconds) / 1000 << " ms of CPU to fit the track\n";
0366
0367 try {
0368 refittedTrack->checkConsistency();
0369 } catch (genfit::Exception& e) {
0370 std::cerr<< e.getExcString() <<std::endl;
0371 continue;
0372 }
0373
0374 track = refittedTrack.get();
0375 }
0376
0377
0378
0379
0380 AbsTrackRep* rep;
0381
0382 if (drawCardinalRep_) {
0383 rep = track->getCardinalRep();
0384 std::cout << "Draw cardinal rep" << std::endl;
0385 }
0386 else {
0387 if (repId_ >= track->getNumReps())
0388 repId_ = track->getNumReps() - 1;
0389 rep = track->getTrackRep(repId_);
0390 std::cout << "Draw rep" << repId_ << std::endl;
0391 }
0392
0393 if (debugLvl_>0) {
0394 std::cout << "track " << i << std::endl;
0395
0396 track->Print("C");
0397 track->getFitStatus(rep)->Print();
0398
0399 if (track->getFitStatus(rep)->isFitted()) {
0400 try {
0401 std::cout << "fitted state: \n";
0402 track->getFittedState().Print();
0403 }
0404 catch (Exception& e) {
0405 std::cerr << e.what();
0406 }
0407 }
0408 }
0409
0410
0411
0412 rep->setPropDir(0);
0413
0414 unsigned int numhits = track->getNumPointsWithMeasurement();
0415
0416 KalmanFitterInfo* fi;
0417 KalmanFitterInfo* prevFi = 0;
0418 const MeasuredStateOnPlane* fittedState(nullptr);
0419 const MeasuredStateOnPlane* prevFittedState(nullptr);
0420
0421 for(unsigned int j = 0; j < numhits; j++) {
0422
0423 fittedState = nullptr;
0424
0425 TrackPoint* tp = track->getPointWithMeasurement(j);
0426 if (! tp->hasRawMeasurements()) {
0427 std::cerr<<"trackPoint has no raw measurements"<<std::endl;
0428 continue;
0429 }
0430
0431 const AbsMeasurement* m = tp->getRawMeasurement();
0432 int hit_coords_dim = m->getDim();
0433
0434
0435 if (tp->getNumRawMeasurements() > 1) {
0436 bool sameTypes(true);
0437 for (unsigned int iM=1; iM<tp->getNumRawMeasurements(); ++iM) {
0438 auto& rawMeasurement = *(tp->getRawMeasurement(iM));
0439 if (typeid(rawMeasurement) != typeid(*m))
0440 sameTypes = false;
0441 }
0442 if (!sameTypes) {
0443 std::cerr<<"cannot draw trackpoint containing multiple Measurements of differend types"<<std::endl;
0444 continue;
0445 }
0446 }
0447
0448
0449
0450
0451 if (! tp->hasFitterInfo(rep)) {
0452 std::cerr<<"trackPoint has no fitterInfo for rep"<<std::endl;
0453 continue;
0454 }
0455
0456 AbsFitterInfo* fitterInfo = tp->getFitterInfo(rep);
0457
0458 fi = dynamic_cast<KalmanFitterInfo*>(fitterInfo);
0459 if(fi == nullptr) {
0460 std::cerr<<"can only display KalmanFitterInfo"<<std::endl;
0461 continue;
0462 }
0463 if (! fi->hasPredictionsAndUpdates()) {
0464 std::cerr<<"KalmanFitterInfo does not have all predictions and updates"<<std::endl;
0465
0466 }
0467 else {
0468 try {
0469 fittedState = &(fi->getFittedState(true));
0470 }
0471 catch (Exception& e) {
0472 std::cerr << e.what();
0473 std::cerr<<"can not get fitted state"<<std::endl;
0474 fittedState = nullptr;
0475 prevFi = fi;
0476 prevFittedState = fittedState;
0477 continue;
0478 }
0479 }
0480
0481 if (fittedState == nullptr) {
0482 if (fi->hasForwardUpdate()) {
0483 fittedState = fi->getForwardUpdate();
0484 }
0485 else if (fi->hasBackwardUpdate()) {
0486 fittedState = fi->getBackwardUpdate();
0487 }
0488 else if (fi->hasForwardPrediction()) {
0489 fittedState = fi->getForwardPrediction();
0490 }
0491 else if (fi->hasBackwardPrediction()) {
0492 fittedState = fi->getBackwardPrediction();
0493 }
0494 }
0495
0496 if (fittedState == nullptr) {
0497 std::cout << "cannot get any state from fitterInfo, continue.\n";
0498 prevFi = fi;
0499 prevFittedState = fittedState;
0500 continue;
0501 }
0502
0503 TVector3 track_pos = fittedState->getPos();
0504 double charge = fittedState->getCharge();
0505
0506
0507
0508
0509
0510 bool full_hit = (dynamic_cast<const FullMeasurement*>(m) != nullptr);
0511 bool planar_hit = (dynamic_cast<const PlanarMeasurement*>(m) != nullptr);
0512 bool planar_pixel_hit = planar_hit && hit_coords_dim == 2;
0513 bool space_hit = (dynamic_cast<const SpacepointMeasurement*>(m) != nullptr);
0514 bool wire_hit = m && m->isLeftRightMeasurement();
0515 bool wirepoint_hit = wire_hit && (dynamic_cast<const WirePointMeasurement*>(m) != nullptr);
0516 if (!full_hit && !planar_hit && !planar_pixel_hit && !space_hit && !wire_hit && !wirepoint_hit) {
0517 std::cout << "Track " << i << ", Hit " << j << ": Unknown measurement type: skipping hit!" << std::endl;
0518 continue;
0519 }
0520
0521
0522
0523 unsigned int nMeas = fi->getNumMeasurements();
0524 for (unsigned int iMeas = 0; iMeas < nMeas; ++iMeas) {
0525
0526 if (iMeas > 0 && wire_hit)
0527 break;
0528
0529 const MeasurementOnPlane* mop = fi->getMeasurementOnPlane(iMeas);
0530 const TVectorT<double>& hit_coords = mop->getState();
0531 const TMatrixTSym<double>& hit_cov = mop->getCov();
0532
0533
0534
0535
0536 TVector3 o = fittedState->getPlane()->getO();
0537 TVector3 u = fittedState->getPlane()->getU();
0538 TVector3 v = fittedState->getPlane()->getV();
0539
0540 double_t hit_u = 0;
0541 double_t hit_v = 0;
0542 double_t plane_size = 4;
0543 TVector2 stripDir(1,0);
0544
0545 if(planar_hit) {
0546 if(!planar_pixel_hit) {
0547 if (dynamic_cast<RKTrackRep*>(rep) != nullptr) {
0548 const TMatrixD& H = mop->getHMatrix()->getMatrix();
0549 stripDir.Set(H(0,3), H(0,4));
0550 }
0551 hit_u = hit_coords(0);
0552 } else {
0553 hit_u = hit_coords(0);
0554 hit_v = hit_coords(1);
0555 }
0556 } else if (wire_hit) {
0557 hit_u = fabs(hit_coords(0));
0558 hit_v = v*(track_pos-o);
0559 if (wirepoint_hit) {
0560 hit_v = hit_coords(1);
0561 }
0562 }
0563
0564 if(plane_size < 4) plane_size = 4;
0565
0566
0567
0568 if(iMeas == 0 &&
0569 (drawPlanes_ || (drawDetectors_ && planar_hit))) {
0570 TVector3 move(0,0,0);
0571 if (planar_hit) move = track_pos-o;
0572 if (wire_hit) move = v*(v*(track_pos-o));
0573 TEveBox* box = boxCreator(o + move, u, v, plane_size, plane_size, 0.01);
0574 if (drawDetectors_ && planar_hit) {
0575 box->SetMainColor(kCyan);
0576 } else {
0577 box->SetMainColor(kGray);
0578 }
0579 box->SetMainTransparency(50);
0580 gEve->AddElement(box);
0581 }
0582
0583
0584
0585 try {
0586 if (j == 0) {
0587 if (drawBackward_) {
0588 MeasuredStateOnPlane update ( *fi->getBackwardUpdate() );
0589 update.extrapolateBy(-3.);
0590 makeLines(&update, fi->getBackwardUpdate(), rep, kMagenta, 1, drawTrackMarkers_, drawErrors_, 1);
0591 }
0592 }
0593 if (j > 0 && prevFi != nullptr) {
0594 if(drawTrack_) {
0595 makeLines(prevFittedState, fittedState, rep, charge > 0 ? kRed : kBlue, 1, drawTrackMarkers_, drawErrors_, 3);
0596 if (drawErrors_) {
0597 makeLines(prevFittedState, fittedState, rep, charge > 0 ? kRed : kBlue, 1, false, drawErrors_, 0, 0);
0598 }
0599 }
0600 if (drawForward_) {
0601 makeLines(prevFi->getForwardUpdate(), fi->getForwardPrediction(), rep, kCyan, 1, drawTrackMarkers_, drawErrors_, 1, 0);
0602 if (j == numhits-1) {
0603 MeasuredStateOnPlane update ( *fi->getForwardUpdate() );
0604 update.extrapolateBy(3.);
0605 makeLines(fi->getForwardUpdate(), &update, rep, kCyan, 1, drawTrackMarkers_, drawErrors_, 1, 0);
0606 }
0607 }
0608 if (drawBackward_) {
0609 makeLines(prevFi->getBackwardPrediction(), fi->getBackwardUpdate(), rep, kMagenta, 1, drawTrackMarkers_, drawErrors_, 1);
0610 }
0611
0612 if(drawRefTrack_ && fi->hasReferenceState() && prevFi->hasReferenceState())
0613 makeLines(prevFi->getReferenceState(), fi->getReferenceState(), rep, charge > 0 ? kRed + 2 : kBlue + 2, 2, drawTrackMarkers_, false, 3);
0614 }
0615 else if (j > 0 && prevFi == nullptr) {
0616 std::cout << "previous FitterInfo == nullptr \n";
0617 }
0618 }
0619 catch (Exception& e) {
0620 std::cerr << "extrapolation failed, cannot draw track" << std::endl;
0621 std::cerr << e.what();
0622 }
0623
0624
0625 if(drawDetectors_) {
0626
0627 if(wire_hit) {
0628 TEveGeoShape* det_shape = new TEveGeoShape("det_shape");
0629 det_shape->IncDenyDestroy();
0630 det_shape->SetShape(new TGeoTube(std::max(0., (double)(hit_u-0.0105/2.)), hit_u+0.0105/2., plane_size));
0631
0632 TVector3 norm = u.Cross(v);
0633 TGeoRotation det_rot("det_rot", (u.Theta()*180)/TMath::Pi(), (u.Phi()*180)/TMath::Pi(),
0634 (norm.Theta()*180)/TMath::Pi(), (norm.Phi()*180)/TMath::Pi(),
0635 (v.Theta()*180)/TMath::Pi(), (v.Phi()*180)/TMath::Pi());
0636 TVector3 move = v*(v*(track_pos-o));
0637 TGeoCombiTrans det_trans(o(0) + move.X(),
0638 o(1) + move.Y(),
0639 o(2) + move.Z(),
0640 &det_rot);
0641 det_shape->SetTransMatrix(det_trans);
0642 det_shape->SetMainColor(kCyan);
0643 det_shape->SetMainTransparency(25);
0644 if((drawHits_ && (hit_u+0.0105/2 > 0)) || !drawHits_) {
0645 gEve->AddElement(det_shape);
0646 }
0647 }
0648
0649 }
0650
0651
0652 if(drawHits_) {
0653
0654
0655 if (full_hit) {
0656
0657 StateOnPlane dummy(rep);
0658 StateOnPlane dummy2(TVectorD(rep->getDim()), static_cast<const FullMeasurement*>(m)->constructPlane(dummy), rep);
0659 MeasuredStateOnPlane sop = *(static_cast<const FullMeasurement*>(m)->constructMeasurementsOnPlane(dummy2)[0]);
0660 sop.getCov()*=errorScale_;
0661
0662 MeasuredStateOnPlane prevSop(sop);
0663 prevSop.extrapolateBy(-3);
0664 makeLines(&sop, &prevSop, rep, kYellow, 1, false, true, 0, 0);
0665
0666 prevSop = sop;
0667 prevSop.extrapolateBy(3);
0668 makeLines(&sop, &prevSop, rep, kYellow, 1, false, true, 0, 0);
0669 }
0670
0671 if(planar_hit) {
0672 if(!planar_pixel_hit) {
0673 TEveBox* hit_box;
0674 TVector3 stripDir3 = stripDir.X()*u + stripDir.Y()*v;
0675 TVector3 stripDir3perp = stripDir.Y()*u - stripDir.X()*v;
0676 TVector3 move = stripDir3perp*(stripDir3perp*(track_pos-o));
0677 hit_box = boxCreator((o + move + hit_u*stripDir3), stripDir3, stripDir3perp, errorScale_*std::sqrt(hit_cov(0,0)), plane_size, 0.0105);
0678 hit_box->SetMainColor(kYellow);
0679 hit_box->SetMainTransparency(0);
0680 gEve->AddElement(hit_box);
0681 } else {
0682
0683 TMatrixDSymEigen eigen_values(hit_cov);
0684 TEveGeoShape* cov_shape = new TEveGeoShape("cov_shape");
0685 cov_shape->IncDenyDestroy();
0686 TVectorT<double> ev = eigen_values.GetEigenValues();
0687 TMatrixT<double> eVec = eigen_values.GetEigenVectors();
0688 double pseudo_res_0 = errorScale_*std::sqrt(ev(0));
0689 double pseudo_res_1 = errorScale_*std::sqrt(ev(1));
0690
0691
0692
0693 if(drawAutoScale_) {
0694 double min_cov = std::min(pseudo_res_0,pseudo_res_1);
0695 if(min_cov < 1e-5) {
0696 std::cout << "Track " << i << ", Hit " << j << ": Invalid covariance matrix (Eigenvalue < 1e-5), autoscaling not possible!" << std::endl;
0697 } else {
0698 if(min_cov < 0.049) {
0699 double cor = 0.05 / min_cov;
0700 std::cout << "Track " << i << ", Hit " << j << ": Pixel covariance too small, rescaling by " << cor;
0701 errorScale_ *= cor;
0702 pseudo_res_0 *= cor;
0703 pseudo_res_1 *= cor;
0704 std::cout << " to " << errorScale_ << std::endl;
0705 }
0706 }
0707 }
0708
0709
0710
0711 cov_shape->SetShape(new TGeoEltu(pseudo_res_0, pseudo_res_1, 0.0105));
0712 TVector3 pix_pos = o + hit_u*u + hit_v*v;
0713 TVector3 u_semiaxis = (pix_pos + eVec(0,0)*u + eVec(1,0)*v)-pix_pos;
0714 TVector3 v_semiaxis = (pix_pos + eVec(0,1)*u + eVec(1,1)*v)-pix_pos;
0715 TVector3 norm = u.Cross(v);
0716
0717
0718
0719 TGeoRotation det_rot("det_rot", (u_semiaxis.Theta()*180)/TMath::Pi(), (u_semiaxis.Phi()*180)/TMath::Pi(),
0720 (v_semiaxis.Theta()*180)/TMath::Pi(), (v_semiaxis.Phi()*180)/TMath::Pi(),
0721 (norm.Theta()*180)/TMath::Pi(), (norm.Phi()*180)/TMath::Pi());
0722 TGeoCombiTrans det_trans(pix_pos(0),pix_pos(1),pix_pos(2), &det_rot);
0723 cov_shape->SetTransMatrix(det_trans);
0724
0725
0726 cov_shape->SetMainColor(kYellow);
0727 cov_shape->SetMainTransparency(0);
0728 gEve->AddElement(cov_shape);
0729 }
0730 }
0731
0732
0733
0734 if(space_hit) {
0735 {
0736
0737 TMatrixDSymEigen eigen_values(m->getRawHitCov());
0738 TEveGeoShape* cov_shape = new TEveGeoShape("cov_shape");
0739 cov_shape->IncDenyDestroy();
0740 cov_shape->SetShape(new TGeoSphere(0.,1.));
0741 TVectorT<double> ev = eigen_values.GetEigenValues();
0742 TMatrixT<double> eVec = eigen_values.GetEigenVectors();
0743 TVector3 eVec1(eVec(0,0),eVec(1,0),eVec(2,0));
0744 TVector3 eVec2(eVec(0,1),eVec(1,1),eVec(2,1));
0745 TVector3 eVec3(eVec(0,2),eVec(1,2),eVec(2,2));
0746 const TVector3 norm = u.Cross(v);
0747
0748
0749 static const double radDeg(180./TMath::Pi());
0750 TGeoRotation det_rot("det_rot", eVec1.Theta()*radDeg, eVec1.Phi()*radDeg,
0751 eVec2.Theta()*radDeg, eVec2.Phi()*radDeg,
0752 eVec3.Theta()*radDeg, eVec3.Phi()*radDeg);
0753
0754 if (! det_rot.IsValid()){
0755
0756 if (fabs(eVec2*eVec3) > 1.e-10)
0757 eVec3 = eVec1.Cross(eVec2);
0758
0759 det_rot.SetAngles(eVec1.Theta()*radDeg, eVec1.Phi()*radDeg,
0760 eVec2.Theta()*radDeg, eVec2.Phi()*radDeg,
0761 eVec3.Theta()*radDeg, eVec3.Phi()*radDeg);
0762 }
0763
0764
0765 double pseudo_res_0 = errorScale_*std::sqrt(ev(0));
0766 double pseudo_res_1 = errorScale_*std::sqrt(ev(1));
0767 double pseudo_res_2 = errorScale_*std::sqrt(ev(2));
0768 if(drawScaleMan_) {
0769 pseudo_res_0 = errorScale_*0.5;
0770 pseudo_res_1 = errorScale_*0.5;
0771 pseudo_res_2 = errorScale_*0.5;
0772 }
0773
0774
0775
0776 if(drawAutoScale_) {
0777 double min_cov = std::min(pseudo_res_0,std::min(pseudo_res_1,pseudo_res_2));
0778 if(min_cov < 1e-5) {
0779 std::cout << "Track " << i << ", Hit " << j << ": Invalid covariance matrix (Eigenvalue < 1e-5), autoscaling not possible!" << std::endl;
0780 } else {
0781 if(min_cov <= 0.149) {
0782 double cor = 0.15 / min_cov;
0783 std::cout << "Track " << i << ", Hit " << j << ": Space hit covariance too small, rescaling by " << cor;
0784 errorScale_ *= cor;
0785 pseudo_res_0 *= cor;
0786 pseudo_res_1 *= cor;
0787 pseudo_res_2 *= cor;
0788 std::cout << " to " << errorScale_ << std::endl;
0789
0790 }
0791 }
0792 }
0793
0794
0795
0796 TGeoGenTrans det_trans(o(0),o(1),o(2),
0797
0798
0799 pseudo_res_0, pseudo_res_1, pseudo_res_2,
0800 &det_rot);
0801 cov_shape->SetTransMatrix(det_trans);
0802
0803
0804 cov_shape->SetMainColor(kYellow);
0805 cov_shape->SetMainTransparency(10);
0806 gEve->AddElement(cov_shape);
0807 }
0808
0809
0810 {
0811
0812 TMatrixDSymEigen eigen_values(hit_cov);
0813 TEveGeoShape* cov_shape = new TEveGeoShape("cov_shape");
0814 cov_shape->IncDenyDestroy();
0815 TVectorT<double> ev = eigen_values.GetEigenValues();
0816 TMatrixT<double> eVec = eigen_values.GetEigenVectors();
0817 double pseudo_res_0 = errorScale_*std::sqrt(ev(0));
0818 double pseudo_res_1 = errorScale_*std::sqrt(ev(1));
0819
0820
0821
0822 if(drawAutoScale_) {
0823 double min_cov = std::min(pseudo_res_0,pseudo_res_1);
0824 if(min_cov < 1e-5) {
0825 std::cout << "Track " << i << ", Hit " << j << ": Invalid covariance matrix (Eigenvalue < 1e-5), autoscaling not possible!" << std::endl;
0826 } else {
0827 if(min_cov < 0.049) {
0828 double cor = 0.05 / min_cov;
0829 std::cout << "Track " << i << ", Hit " << j << ": Pixel covariance too small, rescaling by " << cor;
0830 errorScale_ *= cor;
0831 pseudo_res_0 *= cor;
0832 pseudo_res_1 *= cor;
0833 std::cout << " to " << errorScale_ << std::endl;
0834 }
0835 }
0836 }
0837
0838
0839
0840 cov_shape->SetShape(new TGeoEltu(pseudo_res_0, pseudo_res_1, 0.0105));
0841 TVector3 pix_pos = o + hit_u*u + hit_v*v;
0842 TVector3 u_semiaxis = (pix_pos + eVec(0,0)*u + eVec(1,0)*v)-pix_pos;
0843 TVector3 v_semiaxis = (pix_pos + eVec(0,1)*u + eVec(1,1)*v)-pix_pos;
0844 TVector3 norm = u.Cross(v);
0845
0846
0847
0848 static const double radDeg(180./TMath::Pi());
0849 TGeoRotation det_rot("det_rot", u_semiaxis.Theta()*radDeg, u_semiaxis.Phi()*radDeg,
0850 v_semiaxis.Theta()*radDeg, v_semiaxis.Phi()*radDeg,
0851 norm.Theta()*radDeg, norm.Phi()*radDeg);
0852
0853
0854
0855
0856
0857 TGeoCombiTrans det_trans(pix_pos(0),pix_pos(1),pix_pos(2), &det_rot);
0858 cov_shape->SetTransMatrix(det_trans);
0859
0860
0861 cov_shape->SetMainColor(kYellow);
0862 cov_shape->SetMainTransparency(0);
0863 gEve->AddElement(cov_shape);
0864 }
0865 }
0866
0867
0868
0869 if(wire_hit) {
0870 TEveGeoShape* cov_shape = new TEveGeoShape("cov_shape");
0871 cov_shape->IncDenyDestroy();
0872 double pseudo_res_0 = errorScale_*std::sqrt(hit_cov(0,0));
0873 double pseudo_res_1 = plane_size;
0874 if (wirepoint_hit) pseudo_res_1 = errorScale_*std::sqrt(hit_cov(1,1));
0875
0876
0877 if(drawAutoScale_) {
0878 if(pseudo_res_0 < 1e-5) {
0879 std::cout << "Track " << i << ", Hit " << j << ": Invalid wire resolution (< 1e-5), autoscaling not possible!" << std::endl;
0880 } else {
0881 if(pseudo_res_0 < 0.0049) {
0882 double cor = 0.005 / pseudo_res_0;
0883 std::cout << "Track " << i << ", Hit " << j << ": Wire covariance too small, rescaling by " << cor;
0884 errorScale_ *= cor;
0885 pseudo_res_0 *= cor;
0886 std::cout << " to " << errorScale_ << std::endl;
0887 }
0888 }
0889
0890 if(wirepoint_hit && pseudo_res_1 < 1e-5) {
0891 std::cout << "Track " << i << ", Hit " << j << ": Invalid wire resolution (< 1e-5), autoscaling not possible!" << std::endl;
0892 } else {
0893 if(pseudo_res_1 < 0.0049) {
0894 double cor = 0.005 / pseudo_res_1;
0895 std::cout << "Track " << i << ", Hit " << j << ": Wire covariance too small, rescaling by " << cor;
0896 errorScale_ *= cor;
0897 pseudo_res_1 *= cor;
0898 std::cout << " to " << errorScale_ << std::endl;
0899 }
0900 }
0901 }
0902
0903
0904 TEveBox* hit_box;
0905 TVector3 move = v*(v*(track_pos-o));
0906 hit_box = boxCreator((o + move + hit_u*u), u, v, errorScale_*std::sqrt(hit_cov(0,0)), pseudo_res_1, 0.0105);
0907 hit_box->SetMainColor(kYellow);
0908 hit_box->SetMainTransparency(0);
0909 gEve->AddElement(hit_box);
0910
0911 hit_box = boxCreator((o + move - hit_u*u), u, v, errorScale_*std::sqrt(hit_cov(0,0)), pseudo_res_1, 0.0105);
0912 hit_box->SetMainColor(kYellow);
0913 hit_box->SetMainTransparency(0);
0914 gEve->AddElement(hit_box);
0915 }
0916
0917
0918 }
0919
0920 }
0921
0922
0923 prevFi = fi;
0924 prevFittedState = fittedState;
0925
0926 }
0927
0928 }
0929
0930 gEve->Redraw3D(resetCam);
0931
0932 }
0933
0934
0935
0936
0937 TEveBox* EventDisplay::boxCreator(TVector3 o, TVector3 u, TVector3 v, float ud, float vd, float depth) {
0938
0939 TEveBox* box = new TEveBox("detPlane_shape");
0940 float vertices[24];
0941
0942 TVector3 norm = u.Cross(v);
0943 u *= (0.5*ud);
0944 v *= (0.5*vd);
0945 norm *= (0.5*depth);
0946
0947 vertices[0] = o(0) - u(0) - v(0) - norm(0);
0948 vertices[1] = o(1) - u(1) - v(1) - norm(1);
0949 vertices[2] = o(2) - u(2) - v(2) - norm(2);
0950
0951 vertices[3] = o(0) + u(0) - v(0) - norm(0);
0952 vertices[4] = o(1) + u(1) - v(1) - norm(1);
0953 vertices[5] = o(2) + u(2) - v(2) - norm(2);
0954
0955 vertices[6] = o(0) + u(0) - v(0) + norm(0);
0956 vertices[7] = o(1) + u(1) - v(1) + norm(1);
0957 vertices[8] = o(2) + u(2) - v(2) + norm(2);
0958
0959 vertices[9] = o(0) - u(0) - v(0) + norm(0);
0960 vertices[10] = o(1) - u(1) - v(1) + norm(1);
0961 vertices[11] = o(2) - u(2) - v(2) + norm(2);
0962
0963 vertices[12] = o(0) - u(0) + v(0) - norm(0);
0964 vertices[13] = o(1) - u(1) + v(1) - norm(1);
0965 vertices[14] = o(2) - u(2) + v(2) - norm(2);
0966
0967 vertices[15] = o(0) + u(0) + v(0) - norm(0);
0968 vertices[16] = o(1) + u(1) + v(1) - norm(1);
0969 vertices[17] = o(2) + u(2) + v(2) - norm(2);
0970
0971 vertices[18] = o(0) + u(0) + v(0) + norm(0);
0972 vertices[19] = o(1) + u(1) + v(1) + norm(1);
0973 vertices[20] = o(2) + u(2) + v(2) + norm(2);
0974
0975 vertices[21] = o(0) - u(0) + v(0) + norm(0);
0976 vertices[22] = o(1) - u(1) + v(1) + norm(1);
0977 vertices[23] = o(2) - u(2) + v(2) + norm(2);
0978
0979
0980 for(int k = 0; k < 24; k += 3) box->SetVertex((k/3), vertices[k], vertices[k+1], vertices[k+2]);
0981
0982 return box;
0983
0984 }
0985
0986
0987 void EventDisplay::makeLines(const StateOnPlane* prevState, const StateOnPlane* state, const AbsTrackRep* rep,
0988 const Color_t& color, const Style_t& style, bool drawMarkers, bool drawErrors, double lineWidth, int markerPos)
0989 {
0990 if (prevState == nullptr || state == nullptr) {
0991 std::cerr << "prevState == nullptr || state == nullptr\n";
0992 return;
0993 }
0994
0995 TVector3 pos, dir, oldPos, oldDir;
0996 rep->getPosDir(*state, pos, dir);
0997 rep->getPosDir(*prevState, oldPos, oldDir);
0998
0999 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
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
1051 TMatrixDSym cov;
1052 TVector3 position, direction;
1053 rep->getPosMomCov(*measuredState, position, direction, cov);
1054
1055
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
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
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
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
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
1119 rep->getPosMomCov(stateCopy, position, direction, cov);
1120
1121
1122 TMatrixDSymEigen eigen_values2(cov.GetSub(0,2, 0,2));
1123 ev = eigen_values2.GetEigenValues();
1124 eVec = eigen_values2.GetEigenVectors();
1125
1126 ev0 = std::min(ev(0), maxErr);
1127 ev1 = std::min(ev(1), maxErr);
1128 ev2 = std::min(ev(2), maxErr);
1129
1130
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
1153
1154 if (oldEVec1*eVec1 < 0) {
1155 eVec1 *= -1;
1156 eVec2 *= -1;
1157 }
1158
1159
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
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
1184 }
1185
1186
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
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
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
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
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
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 }