FEDRA emulsion software from the OPERA Collaboration
EdbScanTracking Class Reference

#include <EdbScanTracking.h>

Inheritance diagram for EdbScanTracking:
Collaboration diagram for EdbScanTracking:

Public Member Functions

 EdbScanTracking ()
 
void SaveHist (EdbID idset, EdbTrackAssembler &etra)
 
void TrackAli (EdbPVRec &ali, TEnv &env)
 
void TrackSetBT (EdbID id, TEnv &env)
 
virtual ~EdbScanTracking ()
 

Public Attributes

bool eDoRealign
 
int eNgapMax
 
int eNsegMin
 
EdbScanProceSproc
 

Constructor & Destructor Documentation

◆ EdbScanTracking()

EdbScanTracking::EdbScanTracking ( )
449{
450 eNsegMin=2;
451 eNgapMax=50;
452}
int eNsegMin
Definition: EdbScanTracking.h:84
int eNgapMax
Definition: EdbScanTracking.h:85

◆ ~EdbScanTracking()

virtual EdbScanTracking::~EdbScanTracking ( )
inlinevirtual
91{}

Member Function Documentation

◆ SaveHist()

void EdbScanTracking::SaveHist ( EdbID  idset,
EdbTrackAssembler etra 
)
594{
595 TFile f( Form("b%s.trk.root", idset.AsString()) ,"UPDATE");
596 if(etra.eHistNcnd) etra.eHistNcnd->Write();
597 if(etra.eHistProbAll) etra.eHistProbAll->Write();
598 if(etra.eHistProbBest) etra.eHistProbBest->Write();
599 TH1F *probrest = 0, *probPurity=0;
600 if(etra.eHistProbAll&&etra.eHistProbBest) {
601 probrest = (TH1F*)(etra.eHistProbAll->Clone("ProbRest"));
602 probrest->SetTitle("prob for the other candidates");
603 probrest->Add(etra.eHistProbBest,-1);
604 probrest->Write();
605 probPurity = (TH1F*)(probrest->Clone("ProbPurity"));
606 probPurity->SetTitle("Nother/Nall vs prob");
607 probPurity->Divide(etra.eHistProbAll);
608 probPurity->Write();
609 }
610 if(etra.eHistThetaAll) etra.eHistThetaAll->Write();
611 if(etra.eHistThetaBest) etra.eHistThetaBest->Write();
612 TH1F *thetarest = 0, *thetaPurity=0;
613 if(etra.eHistThetaAll&&etra.eHistThetaBest) {
614 thetarest = (TH1F*)(etra.eHistThetaAll->Clone("ThetaRest"));
615 thetarest->SetTitle("theta for other candidates");
616 thetarest->Add(etra.eHistThetaBest,-1);
617 thetarest->Write();
618 thetaPurity = (TH1F*)(thetarest->Clone("ThetaPurity"));
619 thetaPurity->SetTitle("Nother/Nall vs theta");
620 thetaPurity->Divide(etra.eHistThetaAll);
621 thetaPurity->Write();
622 }
623
624 gStyle->SetPalette(1);
625 gStyle->SetOptStat(1);
626 bool batch = gROOT->IsBatch();
627 gROOT->SetBatch();
628
629 TCanvas *c = new TCanvas("purity","tracking purity",900,800);
630 c->Divide(2,3);
631
632 c->cd(1)->SetLogy();
633 etra.eHistNcnd->SetAxisRange(0,10);
634 etra.eHistNcnd->Draw();
635
636 c->cd(2); probrest->SetLineColor(kBlue); probrest->Draw();
637
638 c->cd(3)->SetLogy();
639 etra.eHistProbAll->Draw();
640 etra.eHistProbBest->SetLineColor(kRed); etra.eHistProbBest->Draw("same");
641 probrest->SetLineColor(kBlue); probrest->Draw("same");
642
643 c->cd(5);
644 probPurity->Draw();
645
646 c->cd(4)->SetLogy();
647 etra.eHistThetaAll->Draw();
648 etra.eHistThetaBest->SetLineColor(kRed); etra.eHistThetaBest->Draw("same");
649 thetarest->SetLineColor(kBlue); thetarest->Draw("same");
650
651 c->cd(6);
652 thetaPurity->Draw();
653
654 c->Write();
655 SafeDelete(c);
656 gROOT->SetBatch(batch);
657 f.Close();
658}
FILE * f
Definition: RecDispMC.C:150
char * AsString() const
Definition: EdbID.cxx:26
TH1F * eHistNcnd
number of candidates after preliminary selection
Definition: EdbScanTracking.h:46
TH1F * eHistProbBest
prob of the best candidate
Definition: EdbScanTracking.h:42
TH1F * eHistThetaAll
theta of all candidate
Definition: EdbScanTracking.h:45
TH1F * eHistThetaBest
theta of the best candidate
Definition: EdbScanTracking.h:44
TH1F * eHistProbAll
prob of all candidate
Definition: EdbScanTracking.h:43
EdbID idset
Definition: emrec.cpp:35
new TCanvas()

◆ TrackAli()

void EdbScanTracking::TrackAli ( EdbPVRec ali,
TEnv &  env 
)
662{
664
665 etra.eCond.SetSigma0( env.GetValue("fedra.track.Sigma0" , "3 3 0.005 0.005") );
666 etra.eCond.SetPulsRamp0( env.GetValue("fedra.track.PulsRamp0" , "15 20") );
667 etra.eCond.SetPulsRamp04( env.GetValue("fedra.track.PulsRamp04" , "15 20") );
668 etra.eCond.SetDegrad( env.GetValue("fedra.track.Degrad" , 4) );
669 etra.eCond.SetRadX0( env.GetValue("fedra.track.RadX0" , 5810.) );
670
671 etra.eDTmax = env.GetValue("fedra.track.DTmax" , 0.07 );
672 etra.eDRmax = env.GetValue("fedra.track.DRmax" , 45. );
673 etra.eDZGapMax = env.GetValue("fedra.track.DZGapMax" , 5000. );
674 etra.eProbMin = env.GetValue("fedra.track.probmin" , 0.001 );
675
676 bool do_misalign = env.GetValue("fedra.track.do_misalign" , 0 );
677 int npass = env.GetValue("fedra.track.npass" , 1 );
678 float misalign_offset = env.GetValue("fedra.track.misalign_offset", 500. );
679 //bool do_local_corr = env.GetValue("fedra.track.do_local_corr" , 1 );
680 bool eDoRealign = env.GetValue("fedra.track.do_realign" , 0 );
681 bool do_comb = env.GetValue("fedra.track.do_comb" , 0 );
682 eNsegMin = env.GetValue("fedra.track.NsegMin" , 2 );
683 float momentum = env.GetValue("fedra.track.momentum", 2. );
684 etra.InitTrZMap( env.GetValue("fedra.track.TrZmap", "2400 0 120000 2000 0 100000 30" ) );
685
686 EdbAffine2D misalign[60];
687 if(do_misalign) {
688 // 1 2 3 4 5 6 7 8 9
689 int dx[9] = {0,0,1,1,1,0,-1,-1,-1};
690 int dy[9] = {0,1,1,0,-1,-1,-1,0,1};
691 for(int i=0; i<60; i++) {
692 misalign[i].ShiftX( dx[i%9] * misalign_offset );
693 misalign[i].ShiftY( dy[i%9] * misalign_offset );
694 if(gEDBDEBUGLEVEL>1) printf("%d | %d %d\n",i, dx[i%9], dy[i%9]);
695 }
696 }
697
698 int npl = ali.Npatterns();
699
700 // read segments and use them for tracking
701 for(int ipass=0; ipass<npass; ipass++) {
702 if(gEDBDEBUGLEVEL>1) printf("\n\n*************** ipass=%d ************\n",ipass);
703 etra.eCollisionsRate=0;
704 for(int i=0; i<npl; i++) {
705
706 EdbPattern &p = *ali.GetPattern(i);
707
708 //p.SetZ(plate->Z());
709 p.SetSegmentsZ();
710 p.SetID(i);
711 p.SetPID(i);
712 p.SetSegmentsPID();
713// p.SetSegmentsPlate(id->ePlate);
714 if(gEDBDEBUGLEVEL>1) printf("pattern with z: %f\n", p.Z());
715
716 if(do_misalign) {
717 p.Transform(&misalign[i]);
718 Log(2,"EdbScanTracking::TrackSetBT","apply misalignment of %f",misalign_offset);
719 }
720
721 if(i>0) etra.ExtrapolateTracksToZ(p.Z());
722 //if( eDoRealign && i==1 ) etra.CheckPatternAlignment(p,1);
723 //if( eDoRealign && i>1 ) etra.CheckPatternAlignment(p,2);
724 etra.FillTrZMap();
725 etra.AddPattern(p);
726 }
727 }
728
729 int ntr = etra.Tracks().GetEntries();
730
731 for(int i=0; i<ntr; i++) {
732 EdbTrackP *t = (EdbTrackP *)(etra.Tracks().At(i));
733 if(t->N()<eNsegMin) t->SetFlag(-10);
734 }
735
736 etra.SetSegmentsErrors();
737 etra.FitTracks();
738
739 TObjArray selectedTracks(ntr);
740 if(do_comb) {
741 etra.CombTracks(selectedTracks);
742 } else {
743 int cnt=0;
744 for(int i=0; i<ntr; i++) {
745 EdbTrackP *t = (EdbTrackP *)(etra.Tracks().At(i));
746 if(t->Flag()!=-10) {
747 t->SetID(cnt++);
748 t->SetCounters();
749 t->SetSegmentsTrack();
750 t->SetP(momentum);
751 selectedTracks.Add(t);
752 }
753 }
754 }
755
756 EdbID idset;
757 EdbDataProc::MakeTracksTree( selectedTracks, 0., 0., Form("b%s.trk.root", idset.AsString()) );
758 TFile f( Form("b%s.trk.root", idset.AsString()) ,"UPDATE");
759 env.Write();
760 f.Close();
761
762 SaveHist(idset,etra);
763}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
Definition: EdbAffine.h:17
void ShiftX(float d)
Definition: EdbAffine.h:64
void ShiftY(float d)
Definition: EdbAffine.h:65
static int MakeTracksTree(EdbPVRec *ali=0, const char *file="linked_tracks.root")
Definition: EdbDataSet.cxx:2506
Definition: EdbID.h:7
Definition: EdbPattern.h:273
Int_t Npatterns() const
Definition: EdbPattern.h:366
EdbPattern * GetPattern(int id) const
Definition: EdbPattern.cxx:1721
void SetPulsRamp0(float p1, float p2)
Definition: EdbScanCond.h:74
void SetDegrad(float d)
Definition: EdbScanCond.h:71
void SetSigma0(float x, float y, float tx, float ty)
Definition: EdbScanCond.h:62
void SetPulsRamp04(float p1, float p2)
Definition: EdbScanCond.h:75
void SetRadX0(float x0)
Definition: EdbScanCond.h:57
bool eDoRealign
Definition: EdbScanTracking.h:87
void SaveHist(EdbID idset, EdbTrackAssembler &etra)
Definition: EdbScanTracking.cxx:593
generic class for the tracks assembling from segments
Definition: EdbScanTracking.h:17
TObjArray & Tracks()
Definition: EdbScanTracking.h:75
void FillTrZMap()
Definition: EdbScanTracking.cxx:290
EdbScanCond eCond
Definition: EdbScanTracking.h:40
void SetSegmentsErrors()
Definition: EdbScanTracking.cxx:336
float eDZGapMax
maxgap acceptance for the fast preselection
Definition: EdbScanTracking.h:35
void ExtrapolateTracksToZ(float z, int nsegmin=0)
Definition: EdbScanTracking.cxx:259
void InitTrZMap(const char *str)
Definition: EdbScanTracking.cxx:301
void AddPattern(EdbPattern &p)
Definition: EdbScanTracking.cxx:105
float eProbMin
min acceptable probability for segments preselection
Definition: EdbScanTracking.h:36
void FitTracks()
Definition: EdbScanTracking.cxx:354
float eDRmax
position acceptance for the fast preselection
Definition: EdbScanTracking.h:34
float eDTmax
angular acceptance for the fast preselection
Definition: EdbScanTracking.h:33
void CombTracks(TObjArray &selected)
Definition: EdbScanTracking.cxx:370
int eCollisionsRate
Definition: EdbScanTracking.h:39
Definition: EdbPattern.h:113
TTree * t
Definition: check_shower.C:4
gEDBDEBUGLEVEL
Definition: energy.C:7
EdbPVRec * ali
Definition: test_oracle.C:9
float momentum
Definition: check_vertex.C:20
p
Definition: testBGReduction_AllMethods.C:8

◆ TrackSetBT()

void EdbScanTracking::TrackSetBT ( EdbID  id,
TEnv &  env 
)
456{
457
458 // read scanset object
460 if(!ss) { Log(1,"EdbScanTracking::TrackSetBT",
461 "Error! set for %s do not found",idset.AsString()); return; }
462
463 int npl = ss->eIDS.GetEntries();
464 if(npl<2) { Log(1,"EdbScanTracking::TrackSetBT", "Warning! npl<2 : %d stop tracking!",npl); return; }
465
466 // create and init tracking object
468
469 etra.eCond.SetSigma0( env.GetValue("fedra.track.Sigma0" , "3 3 0.005 0.005") );
470 etra.eCond.SetPulsRamp0( env.GetValue("fedra.track.PulsRamp0" , "15 20") );
471 etra.eCond.SetPulsRamp04( env.GetValue("fedra.track.PulsRamp04" , "15 20") );
472 etra.eCond.SetDegrad( env.GetValue("fedra.track.Degrad" , 4) );
473 etra.SetRadLength( env.GetValue("fedra.track.RadX0" , 5810.) );
474 etra.eDoUseMCS = env.GetValue("fedra.track.do_use_mcs" , 0 );
475
476 etra.eDTmax = env.GetValue("fedra.track.DTmax" , 0.07 );
477 etra.eDRmax = env.GetValue("fedra.track.DRmax" , 45. );
478 etra.eDZGapMax = env.GetValue("fedra.track.DZGapMax" , 5000. );
479 etra.eProbMin = env.GetValue("fedra.track.probmin" , 0.001 );
480 bool do_erase = env.GetValue("fedra.track.erase" , false );
481 const char *cut = env.GetValue("fedra.readCPcut" , "1" );
482 bool do_misalign = env.GetValue("fedra.track.do_misalign" , 0 );
483 int npass = env.GetValue("fedra.track.npass" , 1 );
484 float misalign_offset = env.GetValue("fedra.track.misalign_offset", 500. );
485 bool do_local_corr = env.GetValue("fedra.track.do_local_corr" , 1 );
486 bool eDoRealign = env.GetValue("fedra.track.do_realign" , 0 );
487 bool do_comb = env.GetValue("fedra.track.do_comb" , 0 );
488 eNsegMin = env.GetValue("fedra.track.NsegMin" , 2 );
489 float momentum = env.GetValue("fedra.track.momentum" , 2. );
490 etra.SetMomentum (momentum);
491
492 etra.InitTrZMap( env.GetValue("fedra.track.TrZmap", "2400 0 120000 2000 0 100000 30" ) );
493
494 //EdbPattern p;
495
496 EdbAffine2D misalign[60];
497 if(do_misalign) {
498 // 1 2 3 4 5 6 7 8 9
499 int dx[9] = {0,0,1,1,1,0,-1,-1,-1};
500 int dy[9] = {0,1,1,0,-1,-1,-1,0,1};
501 for(int i=0; i<60; i++) {
502 misalign[i].ShiftX( dx[i%9] * misalign_offset );
503 misalign[i].ShiftY( dy[i%9] * misalign_offset );
504 if(gEDBDEBUGLEVEL>1) printf("%d | %d %d\n",i, dx[i%9], dy[i%9]);
505 }
506 }
507
508 // read segments and use them for tracking
509 for(int ipass=0; ipass<npass; ipass++) {
510 if(gEDBDEBUGLEVEL>1) printf("\n\n*************** ipass=%d ************\n",ipass);
511 etra.eCollisionsRate=0;
512 for(int i=0; i<npl; i++) {
513 EdbID *id = ss->GetID(i);
514
515 EdbPlateP *plate = ss->GetPlate(id->ePlate);
516
518 eSproc->ReadPatCPnopar(p,*id, cut, do_erase);
519 p.SetZ(plate->Z());
520 p.SetSegmentsZ();
521 p.SetID(i);
522 p.SetPID(i);
523 p.SetSegmentsPID();
524 //plate->Print();
525 p.Transform( plate->GetAffineXY() );
526 p.TransformShr( plate->Shr() );
527 p.TransformA( plate->GetAffineTXTY() );
528 p.SetSegmentsPlate(id->ePlate);
529
530 if(do_local_corr) {
531 int nseg = p.N();
532 for(int j=0; j<nseg; j++) {
533 EdbSegP *s = p.GetSegment(j);
534 plate->CorrectSegLocal(*s);
535 }
536 }
537
538
539 if(do_misalign) {
540 p.Transform(&misalign[i]);
541 Log(2,"EdbScanTracking::TrackSetBT","apply misalignment of %f",misalign_offset);
542 //misalign[i].Print();
543 }
544
545 if(i>0) etra.ExtrapolateTracksToZ(p.Z());
546 if(eDoRealign)
547 {
548 if( i==1 ) etra.CheckPatternAlignment(p,*plate,1);
549 if( i>1 ) etra.CheckPatternAlignment(p,*plate,1);
550 }
551 etra.FillTrZMap();
552 etra.AddPattern(p);
553 }
554 }
555
557
558 int ntr = etra.Tracks().GetEntries();
559
560 for(int i=0; i<ntr; i++) {
561 EdbTrackP *t = (EdbTrackP *)(etra.Tracks().At(i));
562 if(t->N()<eNsegMin) t->SetFlag(-10);
563 }
564
565 etra.SetSegmentsErrors();
566 etra.FitTracks();
567
568 TObjArray selectedTracks(ntr);
569 if(do_comb) {
570 etra.CombTracks(selectedTracks);
571 } else {
572 int cnt=0;
573 for(int i=0; i<ntr; i++) {
574 EdbTrackP *t = (EdbTrackP *)(etra.Tracks().At(i));
575 if(t->Flag()!=-10) {
576 t->SetID(cnt++);
577 t->SetCounters();
578 t->SetSegmentsTrack();
579 t->SetP(momentum);
580 selectedTracks.Add(t);
581 }
582 }
583 }
584 EdbDataProc::MakeTracksTree( selectedTracks, 0., 0., Form("b%s.trk.root", idset.AsString()) );
585 TFile f( Form("b%s.trk.root", idset.AsString()) ,"UPDATE");
586 env.Write();
587 f.Close();
588
589 SaveHist(idset,etra);
590}
Definition: EdbBrick.h:14
int ReadPatCPnopar(EdbPattern &pat, EdbID id, TCut cut="1", bool do_erase=false, bool read_mt=false)
Definition: EdbScanProc.cxx:681
int WriteScanSet(EdbID id, EdbScanSet &ss)
Definition: EdbScanProc.cxx:1428
EdbScanSet * ReadScanSet(EdbID id)
Definition: EdbScanProc.cxx:1482
Definition: EdbScanSet.h:11
EdbScanProc * eSproc
Definition: EdbScanTracking.h:86
Definition: EdbSegP.h:21
void CheckPatternAlignment(EdbPattern &p, EdbPlateP &plate, int nsegmin)
Definition: EdbScanTracking.cxx:74
void SetMomentum(float p)
Definition: EdbScanTracking.h:52
bool eDoUseMCS
flag to use MultipleScattering addition for chi2
Definition: EdbScanTracking.h:37
void SetRadLength(float x0)
Definition: EdbScanTracking.h:53
TCut cut
Definition: check_shower.C:6
s
Definition: check_shower.C:55
ss
Definition: energy.C:62
Int_t plate
Definition: merge_Energy_SytematicSources_Electron.C:1
UInt_t id
Definition: tlg2couples.C:117

Member Data Documentation

◆ eDoRealign

bool EdbScanTracking::eDoRealign

◆ eNgapMax

int EdbScanTracking::eNgapMax

◆ eNsegMin

int EdbScanTracking::eNsegMin

◆ eSproc

EdbScanProc* EdbScanTracking::eSproc

The documentation for this class was generated from the following files: