FEDRA emulsion software from the OPERA Collaboration
EdbEDAShowerTab Class Reference

#include <EdbEDAShowerTab.h>

Collaboration diagram for EdbEDAShowerTab:

Public Member Functions

void Button1 ()
 
void Button2 ()
 
void Button3 ()
 
double CalcdR_NoPropagation (EdbSegP *s1, EdbSegP *stest)
 
double CalcdR_WithPropagation (EdbSegP *s1, EdbSegP *stest)
 
double CalcdTheta (EdbSegP *s1, EdbSegP *s2)
 
void ChangeShowerRecParameters ()
 
void CheckAndSetPVRGeneric ()
 
void CheckBTDensity (EdbPVRec *pvr, TObjArray *tracks)
 
void CheckBTDensityCanv (EdbPVRec *pvr=NULL, TObjArray *tracks=NULL)
 
void CheckBTDensityUsingEdbPVRQuality ()
 
void CheckCosmicReco ()
 
void CheckInBT ()
 
Bool_t CheckPairDuplications (Int_t nmax, Int_t SegPID, Int_t SegID, Int_t Seg2PID, Int_t Seg2ID, TArrayI *SegmentPIDArray, TArrayI *SegmentIDArray, TArrayI *Segment2PIDArray, TArrayI *Segment2IDArray)
 
void CheckPVRec ()
 
void CheckQualitycut ()
 
TCanvasCreateCanvas (char *plot_name)
 
void DrawShower ()
 
void DrawShower (int iShowerNr, bool eClearShower=0)
 
void DrawShowerAll ()
 
void DrawShowerNext ()
 
void DrawShowerPrevious ()
 
void DrawShowers ()
 For the toggle Button! More...
 
void DrawSingleShower (EdbTrackP *shower)
 
 EdbEDAShowerTab ()
 
void FindPairings ()
 
void FindPairingsNewAlgo ()
 
void GO ()
 Make EVERY POSSIBLE OPTIONS!! ....still to be done... More...
 
void Help ()
 
void LoadShowerRecEdbPVRec ()
 
void MakeGUI ()
 
void MakeParameterWindowQualitySettings ()
 
void PlotShower ()
 
void PrintActualShower ()
 
void PrintQualityCutSetting ()
 
void PrintShower ()
 
void PrintShowerRecInBT ()
 
void PrintShowerRecParameters ()
 
void PrintShowerRecRecoShowerArray ()
 
void ReadComboBoxParameter ()
 
void ReadComboBoxParameterQualitySetting ()
 
void Reco ()
 
void Reset ()
 
void ResetShowerRecParameters ()
 
void SetInBTFromLinkedTracks ()
 
void SetInBTFromSelected ()
 
void SetNumberEntryShowerRecParameter (Double_t val)
 
void SetParaSetHighNBT ()
 
void SetParaSetHighPur ()
 
void SetPVRecFromFile ()
 

Private Attributes

TGCheckButton * eCheckButtonCheckQualitycut
 
TGCheckButton * eCheckButtonDrawShowers
 
TGCheckButton * eCheckButtonUseTSData
 
Bool_t eCheckQualitycut
 
TGComboBox * eComboBox_Parameter
 Buttons and so on: More...
 
TGComboBox * eComboBox_ParameterQualitySetting
 
Bool_t eDrawShowers
 
TObjArray * eInBTArray
 
Bool_t eIsAliLoaded
 
Bool_t eIspvrec_cleaned
 
Int_t eNShowers
 An own Shower array. Not the one of the ShowerRec object: More...
 
TGNumberEntry * eNumberEntry_ParaValue
 
TGMainFrame * eParamWindow
 
Int_t eQualityCutSetting
 
EdbTrackPeShower
 Actual Shower (current one, either for drawing or printing) More...
 
EdbShowerReceShowerRec
 The ShowerRec object: More...
 
TGNumberEntry * eTGNumberEntryDrawShowerNr
 
Bool_t eUseTSData
 Various Variables used for settings: More...
 
EdbPVRecpvrec_cleaned
 (==generic cleaned) More...
 
EdbPVRecpvrec_cleaned_cpfiles
 "fake" pvrec object with adapted BG quality cut in it. More...
 
EdbPVRecpvrec_cleaned_linkedtracks
 "fake" pvrec object with adapted BG quality cut in it. More...
 
EdbPVRecpvrec_cpfiles
 
EdbPVRecpvrec_generic
 
EdbPVRecpvrec_linkedtracks
 Two Pointers for the EdbPVRec object: More...
 
TObjArray * RecoShowerArray
 

Constructor & Destructor Documentation

◆ EdbEDAShowerTab()

EdbEDAShowerTab::EdbEDAShowerTab ( )
235 {
236
237 // Create ShowerRec Object
238 eShowerRec = new EdbShowerRec();
239
240 // Create Initiator BT array:
241 eInBTArray=new TObjArray();
242
243 eIsAliLoaded=kFALSE; // this means EdbPVrec object using ALL cp.roor Basetrack info!
244 eUseTSData=kTRUE;
245 eDrawShowers=kTRUE;
246 eCheckQualitycut=kTRUE;
248
249 RecoShowerArray= new TObjArray();
250
251 // Two Pointers for the EdbPVRec object:
255
256
257 // Reset actual shower:
258 eShower=0;
259
260 // Make the tab:
261 MakeGUI();
262
263}
TObjArray * RecoShowerArray
Definition: EdbEDAShowerTab.h:58
EdbTrackP * eShower
Actual Shower (current one, either for drawing or printing)
Definition: EdbEDAShowerTab.h:61
Bool_t eDrawShowers
Definition: EdbEDAShowerTab.h:50
TObjArray * eInBTArray
Definition: EdbEDAShowerTab.h:18
EdbPVRec * pvrec_cpfiles
Definition: EdbEDAShowerTab.h:22
EdbShowerRec * eShowerRec
The ShowerRec object:
Definition: EdbEDAShowerTab.h:14
EdbPVRec * pvrec_linkedtracks
Two Pointers for the EdbPVRec object:
Definition: EdbEDAShowerTab.h:21
EdbPVRec * pvrec_generic
Definition: EdbEDAShowerTab.h:31
Bool_t eUseTSData
Various Variables used for settings:
Definition: EdbEDAShowerTab.h:46
Bool_t eCheckQualitycut
Definition: EdbEDAShowerTab.h:51
Int_t eQualityCutSetting
Definition: EdbEDAShowerTab.h:53
Bool_t eIsAliLoaded
Definition: EdbEDAShowerTab.h:47
Definition: EdbShowerRec.h:25

Member Function Documentation

◆ Button1()

void EdbEDAShowerTab::Button1 ( )
inline
140 {
141 printf("**************** NOT IMPLEMETED *************************\n");
142 }

◆ Button2()

void EdbEDAShowerTab::Button2 ( )
inline
143 {
144 printf("Button2\n");
145 }

◆ Button3()

void EdbEDAShowerTab::Button3 ( )
inline
146 {
147 PrintShower();
148 printf("Print\n");
149 }
void PrintShower()
Definition: EdbEDAShowerTab.C:653

◆ CalcdR_NoPropagation()

double EdbEDAShowerTab::CalcdR_NoPropagation ( EdbSegP s1,
EdbSegP stest 
)
1222 {
1223 return TMath::Sqrt((s->X()-stest->X())*(s->X()-stest->X())+(s->Y()-stest->Y())*(s->Y()-stest->Y()));
1224}
Float_t X() const
Definition: EdbSegP.h:173
Float_t Y() const
Definition: EdbSegP.h:174
s
Definition: check_shower.C:55

◆ CalcdR_WithPropagation()

double EdbEDAShowerTab::CalcdR_WithPropagation ( EdbSegP s1,
EdbSegP stest 
)
1226 {
1227 if (s->Z()==stest->Z()) CalcdR_NoPropagation(s,stest);
1228 Double_t zorig;
1229 Double_t dR;
1230 zorig=s->Z();
1231 s->PropagateTo(stest->Z());
1232 dR=TMath::Sqrt( (s->X()-stest->X())*(s->X()-stest->X())+(s->Y()-stest->Y())*(s->Y()-stest->Y()) );
1233 s->PropagateTo(zorig);
1234 return dR;
1235}
double CalcdR_NoPropagation(EdbSegP *s1, EdbSegP *stest)
Definition: EdbEDAShowerTab.C:1222
Float_t Z() const
Definition: EdbSegP.h:153

◆ CalcdTheta()

double EdbEDAShowerTab::CalcdTheta ( EdbSegP s1,
EdbSegP s2 
)
1237 {
1238 Double_t tx1,tx2,ty1,ty2;
1239 tx1=s1->TX();
1240 tx2=s2->TX();
1241 ty1=s1->TY();
1242 ty2=s2->TY();
1243 Double_t dt= TMath::Sqrt( (tx1-tx2)*(tx1-tx2) + (ty1-ty2)*(ty1-ty2) );
1244 return dt;
1245}
Float_t TX() const
tangens = deltaX/deltaZ
Definition: EdbSegP.h:175
Float_t TY() const
tangens = deltaY/deltaZ
Definition: EdbSegP.h:176
EdbSegP * s1
Definition: tlg2couples.C:29
EdbSegP * s2
Definition: tlg2couples.C:30

◆ ChangeShowerRecParameters()

void EdbEDAShowerTab::ChangeShowerRecParameters ( )

Read Value from the value field and pass it over to eShowerRec

530 {
532 Double_t Val = eNumberEntry_ParaValue->GetNumber();
533 Double_t IdVal = eComboBox_Parameter->GetSelectedEntry()->EntryId();
534 eShowerRec->SetAlgoParameter(Val,IdVal);
536}
TGComboBox * eComboBox_Parameter
Buttons and so on:
Definition: EdbEDAShowerTab.h:35
TGNumberEntry * eNumberEntry_ParaValue
Definition: EdbEDAShowerTab.h:40
void SetAlgoParameter(Double_t paravalue, Int_t paranr)
Definition: EdbShowerRec.cxx:4818
void PrintParameters()
Definition: EdbShowerRec.cxx:4846

◆ CheckAndSetPVRGeneric()

void EdbEDAShowerTab::CheckAndSetPVRGeneric ( )
1532 {
1533 cout << "EdbEDAShowerTab::CheckAndSetPVRGeneric()" << endl;
1534 cout << "EdbEDAShowerTab::CheckAndSetPVRGeneric() pvrec_cpfiles=" << pvrec_cpfiles << endl;
1535 cout << "EdbEDAShowerTab::CheckAndSetPVRGeneric() pvrec_linkedtracks=" << pvrec_linkedtracks << endl;
1536 cout << "EdbEDAShowerTab::CheckAndSetPVRGeneric() pvrec_cleaned_cpfiles=" << pvrec_cleaned_cpfiles << endl;
1537 cout << "EdbEDAShowerTab::CheckAndSetPVRGeneric() pvrec_cleaned_linkedtracks=" << pvrec_cleaned_linkedtracks << endl;
1538
1539 // Cases: LT/CP (eUseTSData);
1540 // Cases: Cleaned/Raw (eIspvrec_cleaned);
1541
1542
1543 if (eUseTSData) {
1546 }
1547 else {
1550 }
1551
1552 cout << "EdbEDAShowerTab::CheckAndSetPVRGeneric() Since eUseTSData="<< eUseTSData << " and eIspvrec_cleaned=" << eIspvrec_cleaned << endl;
1553 cout << "EdbEDAShowerTab::CheckAndSetPVRGeneric() We set as pvrec_generic = " << pvrec_generic << endl;
1554 cout << "EdbEDAShowerTab::CheckAndSetPVRGeneric() ...Done." << endl;
1555
1556 return;
1557}
Bool_t eIspvrec_cleaned
Definition: EdbEDAShowerTab.h:48
EdbPVRec * pvrec_cleaned_linkedtracks
"fake" pvrec object with adapted BG quality cut in it.
Definition: EdbEDAShowerTab.h:24
EdbPVRec * pvrec_cleaned_cpfiles
"fake" pvrec object with adapted BG quality cut in it.
Definition: EdbEDAShowerTab.h:25
#define NULL
Definition: nidaqmx.h:84

◆ CheckBTDensity()

void EdbEDAShowerTab::CheckBTDensity ( EdbPVRec pvr,
TObjArray *  tracks 
)

Check Basetrack density of totalscan for each plate

690 {
692 printf("EdbEDAShowerTab::CheckBTDensity...\n");
693 if (!eCheckQualitycut) return;
694
696 // Whole cp files data:
697 pvr= eShowerRec->GetEdbPVRec();
698 pvr=pvrec_generic;
699
700 //pvr->Print();
701 EdbPattern* pat = (EdbPattern*)pvr->GetPattern(0);
702 Int_t patN= pat->N();
703 Int_t pvrN= pvr->Npatterns()-1;
704 EdbSegP* seg=NULL;
705 Float_t BTpermm2;
706
707 TH2F* h_plate = new TH2F("All Plates","h_plate",125,0,125000,100,0,100000);
708 TH2F* BTDens = new TH2F("BTDens","BTDens",125,0,125000,100,0,100000);
709 TH2F* BTDensity[10];
710 TProfile* BTDensProfile = new TProfile("BTDensityProfile","BTDensityProfile",57,0,57,0,200);
711 TProfile* BTDensityProfile[10];
712 TH1F* h_BinContentDistribution = new TH1F("h_BinContentDistribution","h_BinContentDistribution",50,0,50);
713 TH1F* h_BinCont[58];
714
715 for (int l=0; l<10; ++l) {
716 BTDensity[l]=(TH2F*)BTDens->Clone();
717 BTDensityProfile[l]=(TProfile*)BTDensProfile->Clone();
718 BTDensityProfile[l]->SetLineColor(21+2*l);
719 }
720 for (int l=0; l<pvrN; ++l) {
721 h_BinCont[l]=(TH1F*)h_BinContentDistribution->Clone();
722 }
723
724 for (int j=0; j<pvrN; ++j) {
725 pat = (EdbPattern*)pvr->GetPattern(j);
726 patN= pat->N();
727 // Get all Segments of the pattern and Fill Area Histos:
728 for (int i=0; i<patN; ++i) {
729 seg=(EdbSegP*)pat->GetSegment(i);
730 h_plate->Fill(seg->X(),seg->Y());
731 // Fill for 4 different Chi2 Cut sets the for histograms
732 for (int l=0; l<10; ++l) {
733 float factor=0.15-(float)l*0.01;
734 if (seg->Chi2()<seg->W()*factor-1.0) BTDensity[l]->Fill(seg->X(),seg->Y());
735 }
736 }
737 // Get Bin Entries and fill profile Histos:
738 for (int ki=1; ki<h_plate->GetNbinsX(); ki++ ) {
739 for (int kj=1; kj<h_plate->GetNbinsY(); kj++ ) {
740 for (int l=0; l<10; ++l) {
741 BTpermm2=(Float_t)BTDensity[l]->GetBinContent(ki,kj); // ONLY WHEN BINAREA=1000*1000
742 if (BTpermm2==0) continue;
743 BTDensityProfile[l]->Fill(j,BTpermm2,1);
744 if (l==3) h_BinCont[j]->Fill(BTpermm2);
745 }
746 }
747 }
748 for (int l=0; l<10; ++l) {
749 BTDensity[l]->Reset();
750 }
751 }
752 // end of pattern loop
753
754 // Estimate Maximum by constant fit (just threshold determination).
755 Bool_t warning=kFALSE;
756 Int_t last_l=0;
757 float factor=0.15-(float)last_l*0.01;
758 TString Chi2Cut=Form("s.eChi2<s.eW*%f-1.0",factor);
759
760 for (int l=9; l>=0; --l) {
761 if (BTDensityProfile[l]->GetMaximum()<1) continue;
762 BTDensityProfile[l]->Fit("pol0","Q0");
763 TF1 *myfunc = BTDensityProfile[l]->GetFunction("pol0");
764 //cout << 0.15-(float)l*0.01 << " " << myfunc->GetParameter(0) << endl;
765 if (myfunc->GetParameter(0)>20 && !warning) {
766 warning=kTRUE;
767 last_l=l+1;
768 factor=0.15-(float)last_l*0.01;
769 Chi2Cut=Form("s.eChi2<s.eW*%f-1.0",factor);
770 cout << " WARNING BT DENSITY GREATER THAN 20 !! Appearing after maximum Chi2Cut of : " << Chi2Cut.Data() << endl;
771 }
772 }
773 factor=0.15-(float)last_l*0.01;
774 Chi2Cut=Form("s.eChi2<s.eW*%f-1.0",factor);
775 if (warning) cout << "Maximum Chi2Cut = " << Chi2Cut.Data() << endl;
776
777 if (eCheckQualitycut==kTRUE) {
778 cout << "----- Set Parameters for Chi2Cut --------"<<endl;
779 cout << "----- " << factor << " , 1.0 "<< endl;
780 eShowerRec->SetQualityCutValues(factor,1.0);
781 }
782 else {
783 cout << "----- No high BT density observed. --------"<<endl;
784 }
785 return;
786}
ParaSet Fill()
void LoadShowerRecEdbPVRec()
Definition: EdbEDAShowerTab.C:337
Definition: EdbPattern.h:273
Int_t Npatterns() const
Definition: EdbPattern.h:366
EdbPattern * GetPattern(int id) const
Definition: EdbPattern.cxx:1721
Definition: EdbSegP.h:21
Float_t Chi2() const
Definition: EdbSegP.h:157
Float_t W() const
Definition: EdbSegP.h:151
Int_t N() const
Definition: EdbPattern.h:86
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:66
void SetQualityCutValues(Double_t p0, Double_t p1)
Definition: EdbShowerRec.h:670
EdbPVRec * GetEdbPVRec() const
Definition: EdbShowerRec.h:681

◆ CheckBTDensityCanv()

void EdbEDAShowerTab::CheckBTDensityCanv ( EdbPVRec pvr = NULL,
TObjArray *  tracks = NULL 
)

Plot Basetrack density of totalscan for each plate
Using the tracks in the given PVRec object.

790 {
793 printf("EdbEDAShowerTab::CheckBTDensityCanv...\n");
795
796
797 // Genergic edbPVrec data (means either cp or lt file, depends on how the user set it:)
798 // pvr= eShowerRec->GetEdbPVRec();
799 pvr=pvrec_generic;
800
801
802 // In Case this is set we directly use pvrec_cleaned of the EdbPVRQuality class.
804 cout << "EdbEDAShowerTab::CheckBTDensityCanv INFO : eQualityCutSetting==1 || eQualityCutSetting==2 We directly use pvrec_cleaned of the EdbPVRQuality class." << endl;
805 pvr=pvrec_cleaned;
806 }
807
808 // Check if exists:
809 if ( NULL == pvr ) {
810 cout << "EdbEDAShowerTab::CheckBTDensityCanv ERROR : NULL == pvr Return now." << endl;
811 return;
812 }
813
814 TCanvas *c1 = CreateCanvas("CheckBTDensity");
815 c1->Divide(2,2);
816 c1->cd(1);
817
818 EdbPattern* pat = (EdbPattern*)pvr->GetPattern(0);
819 Int_t patN= pat->N();
820 Int_t pvrN= pvr->Npatterns()-1;
821 Int_t pvrNMiddle= 0+(pvr->Npatterns()-1)/2;
822 EdbSegP* seg=NULL;
823 Float_t BTpermm2;
824
825 printf("EdbEDAShowerTab::Create Histos...\n");
826 TH2F* h_plate = new TH2F("All Plates","h_plate",125,0,125000,100,0,100000);
827 TH2F* BTDens = new TH2F("BTDens","BTDens",125,0,125000,100,0,100000);
828 TH2F* BTDensity[10];
829 TProfile* BTDensProfile = new TProfile("BTDensityProfile","BTDensityProfile",57,0,57,0,200);
830 TProfile* BTDensityProfile[10];
831 TH1F* h_BinContentDistribution = new TH1F("h_BinContentDistribution","h_BinContentDistribution",50,0,50);
832 TH1F* h_BinCont[58];
833
834 for (int l=0; l<10; ++l) {
835 BTDensity[l]=(TH2F*)BTDens->Clone();
836 BTDensityProfile[l]=(TProfile*)BTDensProfile->Clone();
837 BTDensityProfile[l]->SetLineColor(21+2*l);
838 }
839 for (int l=0; l<pvrN; ++l) {
840 h_BinCont[l]=(TH1F*)h_BinContentDistribution->Clone();
841 }
842
843 for (int j=0; j<pvrN; ++j) {
844 pat = (EdbPattern*)pvr->GetPattern(j);
845 patN= pat->N();
846 // Get all Segments of the pattern and Fill Area Histos:
847 for (int i=0; i<patN; ++i) {
848 seg=(EdbSegP*)pat->GetSegment(i);
849 h_plate->Fill(seg->X(),seg->Y());
850 // Fill for 4 different Chi2 Cut sets the for histograms
851 for (int l=0; l<10; ++l) {
852 float factor=0.15-(float)l*0.01;
853 if (seg->Chi2()<seg->W()*factor-1.0) BTDensity[l]->Fill(seg->X(),seg->Y());
854 }
855 }
856 // Get Bin Entries and fill profile Histos:
857 for (int ki=1; ki<h_plate->GetNbinsX(); ki++ ) {
858 for (int kj=1; kj<h_plate->GetNbinsY(); kj++ ) {
859 for (int l=0; l<10; ++l) {
860 BTpermm2=(Float_t)BTDensity[l]->GetBinContent(ki,kj); // ONLY WHEN BINAREA=1000*1000
861 if (BTpermm2==0) continue;
862 BTDensityProfile[l]->Fill(j,BTpermm2,1);
863 if (l==3) h_BinCont[j]->Fill(BTpermm2);
864 }
865 }
866 }
867 for (int l=0; l<10; ++l) {
868 BTDensity[l]->Reset();
869 }
870 }
871 // end of pattern loop.
872
873
874 printf("EdbEDAShowerTab::Draw Histos...\n");
875 // Draw Part:
876 c1->cd(1);
877 h_plate->Draw("colz");
878 c1->cd(2);
879 h_BinCont[0]->DrawNormalized();
880 h_BinCont[pvrNMiddle]->SetLineColor(2);
881 h_BinCont[pvrN-1]->SetLineColor(4);
882 h_BinCont[pvrNMiddle]->DrawNormalized("same");
883 h_BinCont[pvrN-1]->DrawNormalized("same");
884 TLegend* leg= new TLegend(0.77,0.6,0.98,0.8,"","brNDC");
885 leg->AddEntry(h_BinCont[0],"First plate","LP");
886 leg->AddEntry(h_BinCont[pvrNMiddle],"Middle plate","LP");
887 leg->AddEntry(h_BinCont[pvrN-1],"Last plate","LP");
888 leg->SetBorderSize(1);
889 leg->Draw();
890 c1->cd(3);
891 BTDensityProfile[0]->Draw("profileZ");
892 int xnpat=pvrN+5;
893 BTDensityProfile[0]->GetXaxis()->SetRangeUser(0,xnpat);
894 BTDensityProfile[3]->Draw("profileZsame");
895 BTDensityProfile[9]->Draw("profileZsame");
896 BTDensityProfile[3]->SetLineColor(2);
897 BTDensityProfile[3]->SetLineWidth(2);
898 BTDensityProfile[9]->SetLineColor(6);
899 BTDensityProfile[0]->GetYaxis()->SetRangeUser(0.0,1.2* BTDensityProfile[0]->GetMaximum());
900 TLegend* leg2= new TLegend(0.15,0.75,0.4,0.88,"","brNDC");
901 leg2->AddEntry(BTDensityProfile[0],"Chi2<W*0.15-1","LP");
902 leg2->AddEntry(BTDensityProfile[3],"Chi2<W*0.12-1","LP");
903 leg2->AddEntry(BTDensityProfile[9],"Chi2<W*0.06-1","LP");
904 leg2->SetBorderSize(1);
905 leg2->Draw();
906
907 TLegend* leg3= new TLegend(0.7,0.15,0.98,0.75,"","brNDC");
908 leg3->AddEntry(BTDensityProfile[0],"Chi2<W*0.15-1","LP");
909 leg3->AddEntry(BTDensityProfile[1],"Chi2<W*0.14-1","LP");
910 leg3->AddEntry(BTDensityProfile[2],"Chi2<W*0.13-1","LP");
911 leg3->AddEntry(BTDensityProfile[3],"Chi2<W*0.12-1","LP");
912 leg3->AddEntry(BTDensityProfile[4],"Chi2<W*0.11-1","LP");
913 leg3->AddEntry(BTDensityProfile[5],"Chi2<W*0.10-1","LP");
914 leg3->AddEntry(BTDensityProfile[6],"Chi2<W*0.09-1","LP");
915 leg3->AddEntry(BTDensityProfile[7],"Chi2<W*0.08-1","LP");
916 leg3->AddEntry(BTDensityProfile[8],"Chi2<W*0.07-1","LP");
917 leg3->AddEntry(BTDensityProfile[9],"Chi2<W*0.06-1","LP");
918 leg3->SetBorderSize(1);
919
920
921 c1->cd(4);
922 BTDensityProfile[0]->Draw("profileZ");
923 BTDensityProfile[0]->GetXaxis()->SetRangeUser(0,xnpat);
924 BTDensityProfile[0]->DrawCopy("profileZsame");
925 for (int l=0; l<10; ++l) {
926 BTDensityProfile[l]->DrawCopy("profileZsame");
927 }
928 leg3->Draw();
929 c1->Update();
930 // Draw Done....
931
932 printf("EdbEDAShowerTab::Estimate Maximum by constant fit...\n");
933 // Estimate Maximum by constant fit (just threshold determination).
934 Bool_t warning=kFALSE;
935 Int_t last_l=0;
936 float factor=0.15-(float)last_l*0.01;
937 TString Chi2Cut=Form("s.eChi2<s.eW*%f-1.0",factor);
938
939 for (int l=9; l>=0; --l) {
940 if (BTDensityProfile[l]->GetMaximum()<1) continue;
941 BTDensityProfile[l]->Fit("pol0","Q0");
942 TF1 *myfunc = BTDensityProfile[l]->GetFunction("pol0");
943 if (myfunc->GetParameter(0)>20 && !warning) {
944 warning=kTRUE;
945 last_l=l+1;
946 factor=0.15-(float)last_l*0.01;
947 Chi2Cut=Form("s.eChi2<s.eW*%f-1.0",factor);
948 cout << " WARNING BT DENSITY GREATER THAN 20 !! Appearing after maximum Chi2Cut of : " << Chi2Cut << endl;
949 }
950 }
951 factor=0.15-(float)last_l*0.01;
952 Chi2Cut=Form("s.eChi2<s.eW*%f-1.0",factor);
953 if (warning) {
954 cout << "Maximum Chi2Cut = " << Chi2Cut << endl;
955 }
956 else {
957 cout << "----- No high BT density observed. --------"<<endl;
958 }
959
960 return;
961}
TLegend * leg
Definition: Canv_SYSTEMATICS_ALLCOMBINED__RMSEnergy__vs__Energy__ELECTRON.C:122
EdbPVRec * pvrec_cleaned
(==generic cleaned)
Definition: EdbEDAShowerTab.h:26
TCanvas * CreateCanvas(char *plot_name)
Definition: EdbEDAShowerTab.h:164
TCanvas * c1
Definition: energy.C:13
new TCanvas()

◆ CheckBTDensityUsingEdbPVRQuality()

void EdbEDAShowerTab::CheckBTDensityUsingEdbPVRQuality ( )
1521 {
1523 pvrec_cleaned = qual->GetEdbPVRec(1);
1524 eIspvrec_cleaned=kTRUE;
1525 qual->Print();
1526}
Root Class Definition for EdbPVRQuality.
Definition: EdbPVRQuality.h:38
void Print()
Definition: EdbPVRQuality.cxx:1076
EdbPVRec * GetEdbPVRec()
Definition: EdbPVRQuality.h:232

◆ CheckCosmicReco()

void EdbEDAShowerTab::CheckCosmicReco ( )
1248 {
1249
1250 if (NULL == pvrec_generic ) {
1251 cout << "ERROR: NULL == pvrec_generic " << endl;
1252 return;
1253 }
1254
1255 cout << "" << endl;
1256
1257 // Get track set of TS, this is our Referecne Cosmic Track Set
1258 cout << "// Get track set of TS, this is our Referecne Cosmic Track Set" << endl;
1259
1260 cout << "gEDA->GetTrackSet(TS)->N() " << gEDA->GetTrackSet("TS")->N() << endl;
1261
1262 cout << "gEDA->GetTrackSet(TS)->GetPVRec()->Print(); :" << endl;
1263 gEDA->GetTrackSet("TS")->GetPVRec()->Print();
1264
1265 EdbEDATrackSet* reftrackset = gEDA->GetTrackSet("TS"); //snew EdbEDATrackSet("ReferencetrackSet");
1266 TObjArray* tracksArray=reftrackset->GetTracksBase();
1268 int maxTrackSegments=0;
1269
1270 // Loop to find longest track npl:
1271 cout << "// Loop to find longest track npl:" << endl;
1272 for (int i=0; i<tracksArray->GetEntries(); i++) {
1273 track=(EdbTrackP*) tracksArray->At(i);
1274 //track->Print();
1275 if (track->N()>maxTrackSegments) maxTrackSegments=track->N();
1276 }
1277 cout << "maxTrackSegments=" << maxTrackSegments << endl;
1278
1279
1280 if (NULL == pvrec_generic ) cout << "ERROR: NULL == pvrec_generic " << endl;
1281
1282 // Get number of patterns for edbpvrec objects
1283 cout << "// Get number of patterns for edbpvrec objects" << endl;
1284 int maxNpat=pvrec_generic->Npatterns();
1285 cout << "maxNpat=" << maxNpat << endl;
1286
1287 // Get number of patterns for edbpvrec objects
1288 TObjArray* TestRecoCosmicInBTArray = new TObjArray();
1289
1290 for (int i=0; i<tracksArray->GetEntries(); i++) {
1291 track=(EdbTrackP*) tracksArray->At(i);
1292 if (track->N()<maxNpat-6) continue;
1293 TestRecoCosmicInBTArray->Add(track);
1294 }
1295
1296 cout << "// Get number of TestRecoCosmicInBTArray" << TestRecoCosmicInBTArray->GetEntries() << endl;
1297
1298 eInBTArray=TestRecoCosmicInBTArray;
1299
1302
1303 cout << "" << endl;
1304
1305 cout << " void CheckCosmicReco() eShowerRec->PrintRecoShowerArray()..." << endl;
1307
1311
1312 if (eDrawShowers) {
1313 for (Int_t i=0; i<eNShowers; ++i) {
1314 DrawShower(i);
1315 }
1316 }
1317
1318
1319 cout << "// Get number of TestRecoCosmicInBTArray" << TestRecoCosmicInBTArray->GetEntries() << endl;
1320 cout << "// Get number of RecoShowerArray" << RecoShowerArray->GetEntries() << endl;
1321
1323 cout << "// Starting Tracks with Reconstructed Showers...:" << endl;
1324
1325 int nseg_track_minus_shower=0;
1326 TH1F* h_nseg_track_minus_shower=new TH1F("track_shower_difference","track_shower_difference",100,-50,50);
1327
1328 for (int i=0; i<TestRecoCosmicInBTArray->GetEntries(); i++) {
1329 track=(EdbTrackP*) TestRecoCosmicInBTArray->At(i);
1330// track->PrintNice();
1331 for (int j=i; j<RecoShowerArray->GetEntries(); j++) {
1333// showertrack->PrintNice();
1334 if (TMath::Abs(track->Z()-showertrack->Z())>10) continue;
1335 if (TMath::Abs(track->X()-showertrack->X())>10) continue;
1336 if (TMath::Abs(track->Y()-showertrack->Y())>10) continue;
1337
1338 cout << "Comparing Starting Track i=" << i << "with RecoShower j=" << j <<":";
1339 cout << "Nseg: " << track->N() << " <> " << showertrack->N() << endl;
1340
1341 nseg_track_minus_shower = track->N() - showertrack->N();
1342 h_nseg_track_minus_shower->Fill(nseg_track_minus_shower);
1343 }
1344 }
1345
1346 cout << "Drawing now histogram: positive Entries Mean More Tracks thank shower segemnts." << endl;
1347 cout << "Drawing now histogram: negative Entries Mean More shower segemnts than tracks." << endl;
1348 cout << "Drawing now histogram: In the ideal case, this histogram should peak at 0. If it has large negatie number than check, because the shower contains then toooo much tracks...." << endl;
1349
1350 TCanvas *c1 = CreateCanvas("CosmicCheck");
1351 c1->cd(1);
1352 h_nseg_track_minus_shower -> Draw();
1353 c1->cd();
1354
1355 if (h_nseg_track_minus_shower->GetMean()<0) cout << "WARNING !! MAYBE TOO MUCH BACKGROUND !!! " << endl;
1356}
EdbEDA * gEDA
Definition: EdbEDA.C:3
EdbTrackP * showertrack
Definition: Shower_E_FromShowerRoot.C:27
void DrawShower()
Definition: EdbEDAShowerTab.C:577
Int_t eNShowers
An own Shower array. Not the one of the ShowerRec object:
Definition: EdbEDAShowerTab.h:57
Definition: EdbEDATrackSet.h:178
int N()
Definition: EdbEDATrackSet.h:439
EdbPVRec * GetPVRec()
Definition: EdbEDATrackSet.h:320
TObjArray * GetTracksBase()
Definition: EdbEDATrackSet.h:275
EdbEDATrackSet * GetTrackSet(int i)
Definition: EdbEDA.h:617
void Print() const
Definition: EdbPattern.cxx:1693
void PrintRecoShowerArray()
Definition: EdbShowerRec.cxx:3477
void Execute()
Definition: EdbShowerRec.cxx:3937
void SetInBTArray(TObjArray *InBTArray)
Definition: EdbShowerRec.h:636
TObjArray * GetRecoShowerArray() const
Definition: EdbShowerRec.h:704
Int_t GetRecoShowerArrayN() const
Definition: EdbShowerRec.h:711
Definition: EdbPattern.h:113
Int_t N() const
Definition: EdbPattern.h:177
Definition: bitview.h:14
t Draw("sizeb>>NbBT", cut2)

◆ CheckInBT()

void EdbEDAShowerTab::CheckInBT ( )
447 {
448 if (eShowerRec->GetInBTArrayN()==0) {
449 new TGMsgBox(gClient->GetRoot(), gEve->GetMainWindow(),"EDA", "No Initiator Basetrack selected!");//, 4, 4);
450 return;
451 }
452}
Int_t GetInBTArrayN() const
Definition: EdbShowerRec.h:708

◆ CheckPairDuplications()

Bool_t EdbEDAShowerTab::CheckPairDuplications ( Int_t  nmax,
Int_t  SegPID,
Int_t  SegID,
Int_t  Seg2PID,
Int_t  Seg2ID,
TArrayI *  SegmentPIDArray,
TArrayI *  SegmentIDArray,
TArrayI *  Segment2PIDArray,
TArrayI *  Segment2IDArray 
)
1209{
1210
1211 for (Int_t i=0; i<nmax; i++) {
1212 // PID and ID of Seg and Seg2 to be exchanged for duplications
1213 if ( SegPID==Segment2PIDArray->At(i) && Seg2PID==SegmentPIDArray->At(i) && SegID==Segment2IDArray->At(i) && Seg2ID==SegmentIDArray->At(i)) {
1214 //cout << "Found duplocation for... return true"<<endl;
1215 return kTRUE;
1216 }
1217 }
1218 return kFALSE;
1219}

◆ CheckPVRec()

void EdbEDAShowerTab::CheckPVRec ( )

Log(2, "EdbEDAShowerTab", "--- CheckPVRec() eShowerRec->SetUsePattern_ID(kTRUE)"); if (PIDallSame) eShowerRec->SetUsePattern_ID(kTRUE);

292 {
293
294 Log(2, "EdbEDAShowerTab", "--- CheckPVRec()");
295
296 cout << " ----- Why we need this? During the data storage, sometimes the PID of the " << endl;
297 cout << " ----- EdbPVRec object is always zero. For example, when using o2root downloaded database data. " << endl;
298 cout << " ----- (Not when using simulation data with ORFEO)" << endl;
299 cout << " ----- In this case we have to tell the shower library to use EdbPVRec::ID() instead of EdbPVRec::PID()" << endl;
300
302
303 cout << "EdbEDAShowerTab::CheckPVRec() pvrec_generic = " << pvrec_generic << endl;
304
306 cout << pvrec_generic->Npatterns() << endl;
307 Bool_t PIDallSame=kFALSE;
308
309 // Quick check:
310 // Check if all are of the same value....
311 for (int i=1; i<npat; i++) {
313 PIDallSame=kTRUE;
314 break;
315 }
316 }
317 if (PIDallSame) cout << "PIDallSame ... !!! .... ONe should take rather ID() instead of PID()" << endl;
318 if (!PIDallSame) cout << "PID are different. One can take PIDs for Reconstruction/eAliSub creation..." << endl;
319
322
323
324 // Check PVRec with EdbPVRQualityCut:
325// EdbPVRQuality* qual = new EdbPVRQuality(pvrec_generic);
326// cout << "EdbEDAShowerTab::CheckPVRec() pvrec_generic = " pvrec_generic << endl;
327
328
329
330 Log(2, "EdbEDAShowerTab", "--- CheckPVRec()...done.");
331 return;
332}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
Int_t npat
Definition: Xi2HatStartScript.C:33
int PID() const
Definition: EdbPattern.h:320

◆ CheckQualitycut()

void EdbEDAShowerTab::CheckQualitycut ( )
454 {
456 cout << "EdbEDAShowerTab::CheckQualitycut() Set To: "<< eCheckQualitycut << endl;
457 return;
458}
TGCheckButton * eCheckButtonCheckQualitycut
Definition: EdbEDAShowerTab.h:39

◆ CreateCanvas()

TCanvas * EdbEDAShowerTab::CreateCanvas ( char *  plot_name)
inline

----------------------------------------------------------------------------------------—
Usually already decleared in EdbEdaPlotTab, but since this is for now only of interest
of shower tab we redeclarate it here....

164 {
165 // --- Create an embedded canvas
166 gEve->GetBrowser()->StartEmbedding(1);
167 gROOT->ProcessLineFast("new TCanvas");
168 TCanvas *c1 = (TCanvas*) gPad;
169 gEve->GetBrowser()->StopEmbedding(plot_name);
170 return c1;
171 }

◆ DrawShower() [1/2]

void EdbEDAShowerTab::DrawShower ( )
577 {
578 int nr = eTGNumberEntryDrawShowerNr->GetNumber();
579 DrawShower(nr,1);
580}
TGNumberEntry * eTGNumberEntryDrawShowerNr
Definition: EdbEDAShowerTab.h:41

◆ DrawShower() [2/2]

void EdbEDAShowerTab::DrawShower ( int  iShowerNr,
bool  eClearShower = 0 
)
600 {
601 cout << "void EdbEDAShowerTab::DrawShower( " << iShowerNr << endl;
602
603 // Set Actual Shower Array:
604 TObjArray *ShowerArray = RecoShowerArray;
605
606 // Set Actual Shower:
607 EdbTrackP* Shower = (EdbTrackP*)ShowerArray->At(iShowerNr);
608 eShower= Shower;
609
610 // Print It:
612
613 // Create Segment array for the TrackSet:
614 TObjArray *segments = new TObjArray();
615 for (int j=0; j<Shower->N(); j++) {
616 //cout << "j= " << j << endl;
617 EdbSegP *s = Shower->GetSegment(j);
618 //s->PrintNice();
619 segments->Add( s);
620 }
621
622 cout << "EdbEDATrackSet *set = gEDA->GetTrackSet(BT);"<<endl;
624 set->SetExtendMode(kExtendUpDown);
625 cout << "set->N(); " << set->N() << endl;
626 //set->Print();
627
628 // If Checkbox Overlay is Clicked, the old shower will not be erased,
629 // the new one will be overlayn!
630 if (set->N()>0 && eClearShower) set->Clear();
631 cout << "set->Clear();"<<endl;
632 set->AddSegments(segments);
633
634 // Draw tracks
635 gEDA->Redraw();
636
637 // Draw Shower starting points as blue circles:
638 DrawSingleShower( Shower );
639
640 return;
641}
void PrintActualShower()
Definition: EdbEDAShowerTab.C:661
void DrawSingleShower(EdbTrackP *shower)
Definition: EdbEDAShowerTab.C:670
void Redraw()
Definition: EdbEDA.h:680
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:195
EdbScanSet * set
Definition: emtraceback.cpp:14
@ kExtendUpDown
Definition: EdbEDAUtil.h:91

◆ DrawShowerAll()

void EdbEDAShowerTab::DrawShowerAll ( )
582 {
583 for (int i=0; i<eNShowers; ++i) DrawShower(i,0);
584}

◆ DrawShowerNext()

void EdbEDAShowerTab::DrawShowerNext ( )
593 {
594 int nr = eTGNumberEntryDrawShowerNr->GetNumber()+1;
595 if (nr>=eNShowers) return;
596 DrawShower(nr,1);
597 eTGNumberEntryDrawShowerNr->SetNumber(nr);
598}

◆ DrawShowerPrevious()

void EdbEDAShowerTab::DrawShowerPrevious ( )
586 {
587 int nr = eTGNumberEntryDrawShowerNr->GetNumber()-1;
588 if (nr<0) return;
589 DrawShower(nr,1);
590 eTGNumberEntryDrawShowerNr->SetNumber(nr);
591}

◆ DrawShowers()

void EdbEDAShowerTab::DrawShowers ( )

For the toggle Button!

572 {
573 cout << "void EdbEDAShowerTab::DrawShowers() eCheckButtonDrawShowers->GetState() " << eCheckButtonDrawShowers->GetState()<<endl;
575}
TGCheckButton * eCheckButtonDrawShowers
Definition: EdbEDAShowerTab.h:38

◆ DrawSingleShower()

void EdbEDAShowerTab::DrawSingleShower ( EdbTrackP shower)
670 {
671
672 Log(3, "EdbEDAShowerTab", "--- void DrawSingleShower()");
673
674 EdbSegP *s = shower->GetSegmentFirst();
675 // Draw single vertex.
676 TEvePointSet *ps = new TEvePointSet();
677 ps->SetElementName(Form("Shower (Id, nseg, X,Y,Z) %3d %d %8.1lf %8.1lf %8.1lf", shower->ID(), shower->N(), s->X(), s->Y(), s->Z()));
678 ps->SetMarkerStyle(25);
679 ps->SetNextPoint(s->X(), s->Y(), s->Z()*gEDA->GetScaleZ());
680 ps->SetMarkerColor(kBlue);
681 ps->SetMarkerSize(ps->GetMarkerSize()*1.8);
682 ps->SetUserData(shower);
683 gEve->AddElement(ps);
684 gEve->Redraw3D();
685 return;
686}
double GetScaleZ()
Definition: EdbEDA.h:732
Int_t ID() const
Definition: EdbSegP.h:147
EdbSegP * GetSegmentFirst() const
Definition: EdbPattern.h:189

◆ FindPairings()

void EdbEDAShowerTab::FindPairings ( )
1037 {
1038
1039 cout << "void EdbEDAShowerTab::FindPairings() " << endl;
1040
1041 if (NULL == pvrec_generic ) {
1042 cout << "ERROR: NULL == pvrec_generic " << endl;
1043 return;
1044 }
1045
1046 EdbSegP* InBT=NULL;
1047 EdbSegP* Segment=NULL;
1048 EdbSegP* Segment2=NULL;
1049 EdbSegP* Segment_Sum=new EdbSegP(0,0,0,0,0,0);
1050 Float_t IP_Pair_To_InBT=0;
1051 Float_t x_av,y_av,z_av,tx_av,ty_av,distZ;
1052 Int_t RecoShowerArrayN=0;
1053 TArrayI* SegmentPIDArray = new TArrayI(99999999);
1054 TArrayI* SegmentIDArray = new TArrayI(99999999);
1055 TArrayI* Segment2PIDArray = new TArrayI(99999999);
1056 TArrayI* Segment2IDArray = new TArrayI(99999999);
1057
1058 TObjArray* InBTArray=gEDA->GetSelected();
1059 Int_t eInBTArrayN=InBTArray->GetEntries();
1060 cout << eInBTArrayN << endl;
1061
1062 if (eInBTArrayN==0) {
1063 cout << " NO Vertex or BT Selected! Do nothing and return. "<< endl;
1064 return;
1065 }
1066
1067// For Pairings Search we want to have the full cp files data:
1068// LoadShowerRecEdbPVRec();
1069// eShowerRec->GetEdbPVRec()->Print();
1070
1071// gEDA->GetPVRec() doesnt exist anymore, now its splitted via the different
1072// track sets:
1073// but ususally this has only entries of linked tracks (even if labelled as "TS")
1074// EdbPVRec* eAli_Sub=gEDA->GetTrackSet("TS")->GetPVRec();
1075// Whole cp files data:
1076// EdbPVRec* eAli_Sub=eShowerRec->GetEdbPVRec();
1077// eAli_Sub->Print();
1078// return;
1079
1080
1081 EdbPVRec* eAli_Sub=pvrec_generic;
1082
1083 for (Int_t i=eInBTArrayN-1; i>=0; --i) {
1084
1085 // Get InitiatorBT from eInBTArray InBT
1086 InBT=(EdbSegP*)InBTArray->At(i);
1087 if (gEDBDEBUGLEVEL>1) InBT->PrintNice();
1088
1089 Int_t npat=eAli_Sub->Npatterns()-1;
1090 Int_t pat_one_bt_cnt_max,pat_two_bt_cnt_max=0;
1091 EdbPattern* pat_one=0;
1092 EdbPattern* pat_two=0;
1093
1094 for (Int_t pat_one_cnt=0; pat_one_cnt<npat; ++pat_one_cnt) {
1095 cout << "Looping over pattern one, counter=" << pat_one_cnt << endl;
1096 for (Int_t pat_two_cnt=0; pat_two_cnt<npat; ++pat_two_cnt) {
1097 cout << " Looping over pattern two, counter=" << pat_two_cnt << endl;
1098
1099 // Now apply cut conditions: GS GAMMA SEARCH Alg:
1100 if (TMath::Abs(pat_one_cnt-pat_two_cnt)>1) continue;
1101
1102 pat_one=(EdbPattern*)eAli_Sub->GetPattern(pat_one_cnt);
1103 pat_two=(EdbPattern*)eAli_Sub->GetPattern(pat_two_cnt);
1104 pat_one_bt_cnt_max=eAli_Sub->GetPattern(pat_one_cnt)->GetN();
1105 pat_two_bt_cnt_max=eAli_Sub->GetPattern(pat_two_cnt)->GetN();
1106
1107 if (pat_one->Z()-InBT->Z()>25000) continue;
1108 if (pat_one->Z()-InBT->Z()<0) continue;
1109 if (pat_two->Z()-InBT->Z()>25000) continue;
1110 if (pat_two->Z()-InBT->Z()<0) continue;
1111
1112 for (Int_t pat_one_bt_cnt=0; pat_one_bt_cnt<pat_one_bt_cnt_max; ++pat_one_bt_cnt) {
1113 for (Int_t pat_two_bt_cnt=0; pat_two_bt_cnt<pat_two_bt_cnt_max; ++pat_two_bt_cnt) {
1114 // cout << pat_one_cnt <<" " << pat_two_cnt <<" " << pat_one_bt_cnt <<" " << pat_two_bt_cnt <<" " << endl;
1115
1116 Segment = (EdbSegP*)pat_one->GetSegment(pat_one_bt_cnt);
1117 Segment2 = (EdbSegP*)pat_two->GetSegment(pat_two_bt_cnt);
1118 if (Segment2==Segment) continue;
1119
1120 // Now apply cut conditions: GS GAMMA SEARCH Alg --------------------
1121 // At last: Remove // better: check for duplicated pairings:
1122 if (CheckPairDuplications(RecoShowerArrayN,Segment->PID(),Segment->ID(),Segment2->PID(),Segment2->ID(), SegmentPIDArray,SegmentIDArray,Segment2PIDArray,Segment2IDArray)) continue;
1123 // d) Check if dist Z to vtx (BT) is ok:
1124 // TestSegment has to have greater Z than Vertex/InBT
1125 distZ=Segment->Z()-InBT->Z();
1126 if (distZ>25000) continue;
1127 // a) Check dR between tracks:
1128 if (CalcdR_NoPropagation(Segment,Segment2)>100) continue;
1129 // b) Check dT between tracks:
1130 if (CalcdTheta(Segment,Segment2)>0.1) continue;
1131 // c) Check dMinDist between tracks:
1132 if (CalcDmin(Segment,Segment2)>30.0) continue;
1133
1134 x_av=Segment2->X()+(Segment->X()-Segment2->X())/2.0;
1135 y_av=Segment2->Y()+(Segment->Y()-Segment2->Y())/2.0;
1136 z_av=Segment2->Z()+(Segment->Z()-Segment2->Z())/2.0;
1137 tx_av=Segment2->TX()+(Segment->TX()-Segment2->TX())/2.0;
1138 ty_av=Segment2->TY()+(Segment->TY()-Segment2->TY())/2.0;
1139 Segment_Sum->SetX(x_av);
1140 Segment_Sum->SetY(y_av);
1141 Segment_Sum->SetTX(tx_av);
1142 Segment_Sum->SetTY(ty_av);
1143 Segment_Sum->SetZ(z_av);
1144
1145 // d) Check if IP to vtx (BT) is ok:
1146 IP_Pair_To_InBT=CalcIP(Segment_Sum, InBT->X(),InBT->Y(),InBT->Z());
1147 if (IP_Pair_To_InBT>200) continue;
1148
1149 SegmentPIDArray->AddAt(Segment->PID(),RecoShowerArrayN);
1150 SegmentIDArray->AddAt(Segment->ID(),RecoShowerArrayN);
1151 Segment2PIDArray->AddAt(Segment2->PID(),RecoShowerArrayN);
1152 Segment2IDArray->AddAt(Segment2->ID(),RecoShowerArrayN);
1153 // end of cut conditions: GS GAMMA SEARCH Alg --------------------
1154
1155
1156 // Create new EdbShowerP Object for storage;
1157// // See EdbShowerP why I have to call the Constructor as "unique" ideficable value
1158 EdbTrackP* RecoShower = new EdbTrackP();
1159// // Using this algorithm we search only pairs, so shower has by definition only two segments:
1160 RecoShower -> AddSegment(Segment);
1161 RecoShower -> AddSegment(Segment2);
1162// // Obligatory when Shower Reconstruction is finished:
1163// RecoShower ->Update();
1164// // Print:
1165 RecoShower ->PrintNice();
1166// RecoShower ->PrintSegments();
1167// Segment_Sum->PrintNice();
1168// cout <<"------------"<< endl;
1169//
1170
1171 // Add Shower to Array:
1172 RecoShowerArray->Add(RecoShower);
1173 RecoShowerArrayN++;
1174
1175 }
1176 }
1177
1178 }
1179 }
1180
1181
1182 } // Loop over InBT
1183
1184 // Delete unnecessary objects:
1185 delete SegmentPIDArray;
1186 delete Segment2PIDArray;
1187 delete SegmentIDArray;
1188 delete Segment2IDArray;
1189
1190
1191 InBT->PrintNice();
1192 cout << "EdbShowAlg_GS::Execute()...In total we found (raw) Gamma Candidates: " << RecoShowerArrayN << endl;
1193 cout << "EdbShowAlg_GS::Execute()...done." << endl;
1194
1195 cout << "...TO DO... Return the eRecoShowerArray and Draw it...." << endl;
1196
1197 if (eDrawShowers) {
1198 cout << " void Reco() Draw All Reconstructed Showers" << endl;
1199 for (Int_t i=0; i<RecoShowerArrayN; ++i) {
1200 cout << " void Reco() Draw Reconstructed Showers i= " << i << endl;
1201 DrawShower(i);
1202 }
1203 }
1204 return;
1205}
Double_t CalcIP(EdbSegP *s, EdbVertex *v)
Definition: ShowRec.cpp:8872
TObjArray * GetSelected(void)
Definition: EdbEDA.h:405
Bool_t CheckPairDuplications(Int_t nmax, Int_t SegPID, Int_t SegID, Int_t Seg2PID, Int_t Seg2ID, TArrayI *SegmentPIDArray, TArrayI *SegmentIDArray, TArrayI *Segment2PIDArray, TArrayI *Segment2IDArray)
Definition: EdbEDAShowerTab.C:1208
double CalcdTheta(EdbSegP *s1, EdbSegP *s2)
Definition: EdbEDAShowerTab.C:1237
Definition: EdbPVRec.h:148
void SetY(Float_t y)
Definition: EdbSegP.h:178
void SetTX(Float_t tx)
Definition: EdbSegP.h:179
void SetX(Float_t x)
Definition: EdbSegP.h:177
void SetZ(float z)
Definition: EdbSegP.h:125
void SetTY(Float_t ty)
other functions
Definition: EdbSegP.h:180
void PrintNice() const
Definition: EdbSegP.cxx:392
Int_t PID() const
Definition: EdbSegP.h:148
Float_t Z() const
Definition: EdbPattern.h:84
Int_t GetN() const
Definition: EdbPattern.h:65
void PrintNice()
Definition: EdbPattern.cxx:1169
gEDBDEBUGLEVEL
Definition: energy.C:7
double CalcDmin(EdbSegP *seg1, EdbSegP *seg2, double *dminz=NULL)
Definition: EdbEDAUtil.C:239

◆ FindPairingsNewAlgo()

void EdbEDAShowerTab::FindPairingsNewAlgo ( )
968 {
969 cout << "void EdbEDAShowerTab::FindPairings() " << endl;
970
971 if (NULL == pvrec_generic ) {
972 cout << "ERROR: NULL == pvrec_generic " << endl;
973 return;
974 }
975
976 /*
977 // Get possible Initiator BT Array.
978 TObjArray* InBTArray=gEDA->GetSelected();
979 Int_t eInBTArrayN=InBTArray->GetEntries();
980 cout << eInBTArrayN << endl;
981
982 if (eInBTArrayN==0) {
983 cout << " NO Vertex or BT Selected! Do nothing and return. "<< endl;
984 return;
985 }
986
987 // Get possible selected Vertex BT Array.
988 EdbVertex* eVertex=gEDA->GetSelectedVertex();
989 if (eVertex==NULL) {
990 cout << "WARNING:::gEDA->GetSelectedVertex() is NULL POINTER .... DONT USE A VERTEX AS INPUT. Do nothing and return." << endl;
991 return;
992 }
993
994 // return;
995
996
997 // For Pairings Search we want to have the full cp files data:
998 // LoadShowerRecEdbPVRec();
999 // eShowerRec->GetEdbPVRec()->Print();
1000 // EdbPVRec* eAli_Sub=eShowerRec->GetEdbPVRec();
1001 // eAli_Sub->Print();
1002
1003
1004 EdbPVRec* eAli_Sub=pvrec_generic;
1005 eNShowers=0;
1006 RecoShowerArray=NULL;
1007 // eVertex->Print();
1008 // return;
1009
1010
1011 EdbShowAlg_GS* Alg = new EdbShowAlg_GS();
1012 Alg->SetEdbPVRec(eAli_Sub);
1013 if (eVertex!=NULL) Alg->SetInVtx(eVertex);
1014 else Alg->SetInVtxArray(InBTArray);
1015 cout << Alg->GetInVtxArrayN() << endl;
1016 // gEDBDEBUGLEVEL=3;
1017 Alg->Execute();
1018
1019 eNShowers=Alg->GetRecoShowerArrayN();
1020 RecoShowerArray=Alg->GetRecoShowerArray();
1021
1022 if (eDrawShowers) {
1023 //cout << " void Reco() Draw All Reconstructed Showers" << endl;
1024 for (Int_t i=0; i<eNShowers; ++i) {
1025 //cout << " void Reco() Draw Reconstructed Showers i= " << i << endl;
1026 DrawShower(i);
1027 }
1028 }
1029
1030 */
1031}

◆ GO()

void EdbEDAShowerTab::GO ( )

Make EVERY POSSIBLE OPTIONS!! ....still to be done...

1452 {
1453 cout << "-----------------------------------------" << endl;
1454 cout << endl << endl;
1455 cout << "EdbEDAShowerTab::GO()... Not implemented yet." << endl;
1456 cout << "-----------------------------------------" << endl;
1457}

◆ Help()

void EdbEDAShowerTab::Help ( )
1460 {
1461 cout << "-----------------------------------------" << endl;
1462 cout << endl << endl;
1463 cout << "EdbEDAShowerTab::Help() Last date updated: 12th July, 2010; frank.meisel@lhep.unibe.ch " << endl;
1464 cout << endl << endl;
1465 cout << "HowTo Use this tab and generally:" << endl;
1466 cout << "* Click a BaseTrack and then click RecoButton." << endl;
1467 cout << "" << endl;
1468 cout << "* UseTS: Normally we have the linked_tracks dataset and the full cp files." << endl;
1469 cout << "* UseTS: Whereas standard track search is done on the volume of the linked_tracks, other..." << endl;
1470 cout << "* UseTS: operations like microtrack search or shower reconstruction can (and shall) operate..." << endl;
1471 cout << "* UseTS: on full cp dataset. To load the EdbPVrec object, toggle this button. Loading full..." << endl;
1472 cout << "* UseTS: EdbPVrec object may take a while, but it needs to be done only one time during session." << endl;
1473 cout << "" << endl;
1474 cout << "* UseQuality: Normally the standard shower reco algorithm is backgound stable up to a level..." << endl;
1475 cout << "* UseQuality: of 20 BT/mm2. In case the BG level is higher then here the readaption of ..." << endl;
1476 cout << "* UseQuality: the Chi2 vs Grain quality cut is adapted. For now it is working by constant BT density fit ..." << endl;
1477 cout << "* UseQuality: over the whole EdbPVRecObject (It may be extended Plate/by/Plate later)." << endl;
1478 cout << "" << endl;
1479 cout << "* Draw: If toggled, then reconstructed showers will be stored in EDA TrackSet _BT_ ." << endl;
1480 cout << "* Plot: Plot interesting quantities for the showers stored in the internal shower array." << endl;
1481 cout << "* Plot: For the time being, it plots the quantities for all reconstructed showers at once." << endl;
1482 cout << "* Nr: Draw Shower Number ... ." << endl;
1483 cout << "* Prev/All/Next: Draw Shower: Previous, Next, or All again.." << endl;
1484 cout << "" << endl;
1485 cout << "* CheckPVR: ... Check the EdbPVrec Object." << endl;
1486 cout << "" << endl;
1487 cout << "* CheckBTDens: Check BaseTrack Density of the generic PVRec object. It is strongly recommended that you" << endl;
1488 cout << "* CheckBTDens: click this button and watch the density not going to strong over 20-30 BT/mm2." << endl;
1489 cout << "" << endl;
1490 cout << "* CheckCosmicReco: This is used for a analysis wheter the background is acceptable. We start here from" << endl;
1491 cout << "* CheckCosmicReco: passing through cosmics and then run the reco on them. In the best case, " << endl;
1492 cout << "* CheckCosmicReco: the shower should contain only out the tracks (since they are cosmics) and therefore N_shower <= N_track." << endl;
1493 cout << "* CheckCosmicReco: In the histogram which appears, the difference in number of segments is plotted." << endl;
1494 cout << "* CheckCosmicReco: If its largely in the negative region this is an indication that there is too much BG." << endl;
1495 cout << "" << endl;
1496 cout << "* SetPVR: Set the EdbPVrec Object from a (root)file. ..." << endl;
1497 cout << "" << endl;
1498 cout << "* RecoParams: This box is expert status, here you can change parameters of the reco algorithm itself." << endl;
1499 cout << "* RecoParams: Select Parameter, change to your desired value (units are mrad,microns) and click change button." << endl;
1500 cout << "* RecoParams: Three predefined sets are there: Standard; " << endl;
1501 cout << "* RecoParams: Three predefined sets are there: HighNBT (large acceptance, only good in case of loooow background.)"<< endl;
1502 cout << "* RecoParams: Three predefined sets are there: HighPur (tight cuts, good in case of large background.)"<< endl;
1503 cout << "" << endl;
1504
1505 new TGMsgBox(gClient->GetRoot(), gEve->GetMainWindow(),"EDA", "HELP BOX! Please be so kind and read the shell output. Thank you.");
1506 cout << "-----------------------------------------" << endl;
1507 return;
1508}

◆ LoadShowerRecEdbPVRec()

void EdbEDAShowerTab::LoadShowerRecEdbPVRec ( )

Either the TS data is already loaded for EDA Object, or we have to reload it for the ShowerRec object:

337 {
340 Log(2, "EdbEDAShowerTab", "--- LoadShowerRecEdbPVRec()");
341
342 Bool_t UseTSData=eCheckButtonUseTSData->GetState();
343 Bool_t found_pvrec_linkedtracks=kFALSE;
344 Bool_t found_pvrec_cpfiles=kFALSE;
345
346 if (gEDBDEBUGLEVEL>2) {
347 cout << "UseTSData = " << UseTSData ;
348 cout << ", gEDA->GetDataType() = " << gEDA->GetDataType();
349 cout << ", eShowerRec->GetAliLoaded(); " << eShowerRec->GetAliLoaded();
350 cout << ", eIsAliLoaded " << eIsAliLoaded << endl;
351
352 cout << "--- Adresses for the PVrec objects: pvrec_cpfiles: " << pvrec_cpfiles << " pvrec_linkedtracks: " << pvrec_linkedtracks << endl;
353 cout << "gEDA->GetDataType() " << gEDA->GetDataType() << endl;
354 cout << "gEDA->GetTrackSet(TS)->N() " << gEDA->GetTrackSet("TS")->N() << endl;
355
356 cout << "gEDA->GetTrackSet(TS)->GetPVRec()->Print(); :" << endl;
357 gEDA->GetTrackSet("TS")->GetPVRec()->Print();
358 }
359
360
361 if (gEDA->GetDataType()==100 && UseTSData==kTRUE && eIsAliLoaded==kFALSE) {
362 Log(2, "EdbEDAShowerTab::LoadShowerRecEdbPVRec", "--- Read in EdbPVRec from cp.root files. This may take a while. PLEASE BE PATIENT...");
364 eIsAliLoaded=kTRUE;
365 // This is the EdbPVRec object from cp.root files:
367 Log(2, "EdbEDAShowerTab::LoadShowerRecEdbPVRec", "--- ...ok. Set: eShowerRec->LoadEdbPVRec() to use cp files for reconstruction");
368 }
369 else {
370
371 // Check if the PVR for "TS" exists:
372 if (gEDA->GetTrackSet("TS")->N()==0) {
373 cout << "if (gEDA->GetTrackSet(TS)->N()==0) {"<<endl;
374 cout << "now try to get pvrec_linkedtracks from the trackset of the selected track" << endl;
375 cout << "gEDA->GetSelectedTrack() = " << gEDA->GetSelectedTrack() << endl;
376 if (NULL == gEDA->GetSelectedTrack()) {
377 cout << "EdbEDAShowerTab::LoadShowerRecEdbPVRec() ERROR No Selected Track!!! " << endl;
378 cout << "EdbEDAShowerTab::LoadShowerRecEdbPVRec() ERROR Return Now: Check Your TrackSet Data Please!" << endl;
382 return;
383 }
385 cout << "done: pvrec_linkedtracks=gEDA->GetTrackSet(gEDA->GetSelectedTrack())->GetPVRec();" << endl;
386 //gEDA->GetSelectedTrack()->Print();
388 cout << gEDA->GetTrackSet(gEDA->GetSelectedTrack())->N() << endl;
389 }
390 else {
391 cout << "Trying to get PVRecObject now: "<<endl;
393 cout << "pvrec_linkedtracks " << pvrec_linkedtracks << endl;
395 cout << "pvrec_linkedtracks is NULL. Obviously, the line " << endl;
396 cout << "pvrec_linkedtracks=gEDA->GetTrackSet(TS)->GetPVRec();"<< endl;
397 cout << "did NOT work, try somethign different now:" << endl;
398 cout << "pvrec_linkedtracks=gEDA->GetPVR();"<< endl;
400 cout << "pvrec_linkedtracks " << pvrec_linkedtracks << endl;
401 }
402 }
403
404 // Check now which have found.
405 if (pvrec_linkedtracks!=NULL) found_pvrec_linkedtracks=kTRUE;
406 if (pvrec_cpfiles!=NULL) found_pvrec_cpfiles=kTRUE;
407
408
409 if (UseTSData==kTRUE && found_pvrec_cpfiles) {
410 cout << "UseTSData==kTRUE && found_pvrec_cpfiles" << endl;
411 Log(2, "EdbEDAShowerTab::LoadShowerRecEdbPVRec", "--- Set: eShowerRec->pvrec_cpfiles() to use cp files for reconstruction:");
414 Log(2, "EdbEDAShowerTab::LoadShowerRecEdbPVRec", "--- Set: eShowerRec->pvrec_cpfiles() to use cp files for reconstruction...done.");
415 }
416 else {
417 cout << "neither UseTSData==kTRUE nor found_pvrec_cpfiles. take pvrec_linkedtracks !" << endl;
418 Log(2, "EdbEDAShowerTab::LoadShowerRecEdbPVRec", "--- Set: eShowerRec->pvrec_linkedtracks() to use cp files for reconstruction");
421 Log(2, "EdbEDAShowerTab::LoadShowerRecEdbPVRec", "--- Set: eShowerRec->pvrec_linkedtracks() to use cp files for reconstruction...don.");
422 }
423
424
425 } // of else from // (gEDA->GetDataType()==100 && UseTSData==kTRUE && eIsAliLoaded==kFALSE)
426
427
428 if (gEDBDEBUGLEVEL>2) {
429 cout << "Print Once Again the eShowerRec->GetEdbPVRec()->Print();"<<endl;
431 cout << "--- Adresses for the PVrec objects: pvrec_cpfiles: " << pvrec_cpfiles << " pvrec_linkedtracks: " << pvrec_linkedtracks << endl;
432 cout << "--- We take " << eShowerRec->GetEdbPVRec() << " for pvrec_generic " << endl;
433 }
434
436
437
438 Log(2, "EdbEDAShowerTab", "--- LoadShowerRecEdbPVRec()...Done.");
439 return;
440}
int GetDataType()
Definition: EdbEDA.h:234
EdbPVRec * GetPVR()
Definition: EdbEDA.h:225
EdbTrackP * GetSelectedTrack(int i=-1)
Definition: EdbEDA.h:421
TGCheckButton * eCheckButtonUseTSData
Definition: EdbEDAShowerTab.h:37
void Print()
Definition: EdbEDATrackSet.h:258
void SetEdbPVRec(EdbPVRec *Ali)
Definition: EdbShowerRec.h:629
void LoadEdbPVRec()
Definition: EdbShowerRec.cxx:3689
Bool_t GetAliLoaded() const
Definition: EdbShowerRec.h:690

◆ MakeGUI()

void EdbEDAShowerTab::MakeGUI ( )

◆ MakeParameterWindowQualitySettings()

void EdbEDAShowerTab::MakeParameterWindowQualitySettings ( )
1387 {
1388
1389 // main frame
1390 TGMainFrame *fMainFrame = new TGMainFrame(gClient->GetRoot(),10,300,kMainFrame | kVerticalFrame);
1391 fMainFrame->SetLayoutBroken(kTRUE);
1392
1393 eParamWindow = fMainFrame;
1394
1395 int posy=10;
1396 int posx=10;
1397 int dx;
1398
1399 TGLabel* fLabel = new TGLabel(fMainFrame,"QualitySettings Parameters");
1400 fMainFrame->AddFrame(fLabel);
1401 fLabel->MoveResize(posx,posy,dx=140,20);
1402 posx+=dx+10;
1403
1405 posy+=20;
1406 posx=5;
1407
1408 eComboBox_ParameterQualitySetting = new TGComboBox(fMainFrame,"eQualitySetting");
1409 eComboBox_ParameterQualitySetting->AddEntry("No Qual Cut",-1);
1410 eComboBox_ParameterQualitySetting->AddEntry("Global Qual Cut",0);
1411 eComboBox_ParameterQualitySetting->AddEntry("PlateByPlate, Constant BT Density",1);
1412 eComboBox_ParameterQualitySetting->AddEntry("PlateByPlate, Constant Quality",2);
1414 eComboBox_ParameterQualitySetting->MoveResize(posx,posy,dx=140,20);
1415 fMainFrame->AddFrame(eComboBox_ParameterQualitySetting);
1416 eComboBox_ParameterQualitySetting->Connect("Selected(const char *)", "EdbEDAShowerTab", this, "ReadComboBoxParameterQualitySetting()");
1417
1418
1419 posx+=dx+5;
1420 TGTextButton *eButton = new TGTextButton(fMainFrame,"Check BT density directly");
1421 eButton->MoveResize(posx,posy,dx=120,20);
1422 eButton->Connect("Clicked()","EdbEDAShowerTab", this,"CheckBTDensityCanv()");
1423
1424
1425
1426 fMainFrame->SetMWMHints(kMWMDecorAll,kMWMFuncAll,kMWMInputModeless);
1427 fMainFrame->MapSubwindows();
1428 fMainFrame->Resize(fMainFrame->GetDefaultSize());
1429 fMainFrame->MapWindow();
1430 fMainFrame->Resize(800,280);
1431
1432 gClient->WaitFor(fMainFrame);
1433
1434}
TGMainFrame * eParamWindow
Definition: EdbEDAShowerTab.h:43
TGComboBox * eComboBox_ParameterQualitySetting
Definition: EdbEDAShowerTab.h:36

◆ PlotShower()

void EdbEDAShowerTab::PlotShower ( )

Decide wheter to plot the histograms for one shower, or for whole Shower Array: ....

1362 {
1363
1367
1368 TH1F* h_nseg= new TH1F("h_nseg","h_nseg",100,0,100);
1369 TH1F* h_npl= new TH1F("h_npl","h_npl",58,0,58);
1370 for (int i=0; i<RecoShowerArray->GetEntries(); i++) {
1371 EdbTrackP* sh = (EdbTrackP*)RecoShowerArray->At(i);
1372 h_nseg->Fill(sh->N());
1373 h_npl->Fill(sh->Npl());
1374 }
1375 TCanvas *c1 = CreateCanvas("ShowerPlots");
1376 c1->Divide(1,2);
1377 c1->cd(1);
1378 h_nseg->Draw();
1379 c1->cd(2);
1380 h_npl->Draw();
1381 c1->cd();
1382}
Int_t Npl() const
Definition: EdbPattern.h:171

◆ PrintActualShower()

void EdbEDAShowerTab::PrintActualShower ( )
661 {
662 if (NULL == eShower) {
663 cout << "EdbEDAShowerTab::PrintActualShower() WARNING Actual Shower: eShower == NULLPOINTER. Do nothing. " << endl;
664 return;
665 }
667}

◆ PrintQualityCutSetting()

void EdbEDAShowerTab::PrintQualityCutSetting ( )

Print important information about the QualityCut Setting:

1439 {
1441 cout << "EdbEDAShowerTab::PrintQualityCutSetting()" << endl;
1442 cout << "EdbEDAShowerTab::PrintQualityCutSetting() eQualityCutSetting= " << eQualityCutSetting << endl;
1443 cout << "EdbEDAShowerTab::PrintQualityCutSetting() ----- For information: -----" << endl;
1444 cout << "EdbEDAShowerTab::PrintQualityCutSetting() eQualityCutSetting==-1 : No Quality Cut" << endl;
1445 cout << "EdbEDAShowerTab::PrintQualityCutSetting() eQualityCutSetting==0 : Global Quality Cut by constant BT density." << endl;
1446 cout << "EdbEDAShowerTab::PrintQualityCutSetting() eQualityCutSetting==1 : Plate by Plate Cut by constant BT density." << endl;
1447 cout << "EdbEDAShowerTab::PrintQualityCutSetting() eQualityCutSetting==2 : Plate by Plate Cut by constant Quality." << endl;
1448}

◆ PrintShower()

void EdbEDAShowerTab::PrintShower ( )

! QUICK and DIRTY solution !!!

653 {
655 TFile* file=TFile::Open("Shower.root");
656 TTree* ShowerTree = (TTree *)file->Get("treebranch");
657 ShowerTree->Scan("xb[0]:yb[0]:zb[0]:txb[0]:tyb[0]:sizeb:eProb1:eProb90:Energy:EnergySigma:trackdensb");
658 return;
659}
TFile * file
Definition: write_pvr.C:3

◆ PrintShowerRecInBT()

void EdbEDAShowerTab::PrintShowerRecInBT ( )
inline
100 {
101 cout << "DBG1"<<endl;
103 cout << "DBG2"<<endl;
104 return;
105 }
void PrintInitiatorBTs()
Definition: EdbShowerRec.cxx:2424

◆ PrintShowerRecParameters()

void EdbEDAShowerTab::PrintShowerRecParameters ( )
inline
112 {
114 return;
115 }

◆ PrintShowerRecRecoShowerArray()

void EdbEDAShowerTab::PrintShowerRecRecoShowerArray ( )
inline
106 {
108 return;
109 }

◆ ReadComboBoxParameter()

void EdbEDAShowerTab::ReadComboBoxParameter ( )

Read Name from the Combobox and put the corresponging value in value field:

538 {
540 cout << "eComboBox_Parameter->GetSelectedEntry->EntryId()" << endl;
541 Double_t IdVal = eComboBox_Parameter->GetSelectedEntry()->EntryId();
542 cout << "IdVal = " << IdVal << endl;
543 Double_t Val = eShowerRec->GetAlgoParameter(IdVal);
544 cout << "Val = " << Val << endl;
546}
void SetNumberEntryShowerRecParameter(Double_t val)
Definition: EdbEDAShowerTab.C:550
Double_t GetAlgoParameter(Int_t paranr)
Definition: EdbShowerRec.cxx:4832

◆ ReadComboBoxParameterQualitySetting()

void EdbEDAShowerTab::ReadComboBoxParameterQualitySetting ( )

Read Name from the ReadComboBoxParameterQualitySetting and put the corresponging value in value field:

1512 {
1514 cout << "eComboBox_Parameter->GetSelectedEntry->EntryId()" << endl;
1515 Double_t IdVal = eComboBox_ParameterQualitySetting->GetSelectedEntry()->EntryId();
1516 eQualityCutSetting=IdVal;
1517 cout << "eQualityCutSetting = " << eQualityCutSetting << endl;
1518}

◆ Reco()

void EdbEDAShowerTab::Reco ( )
463 {
464 cout << " void Reco()" << endl;
465
466// Reset eShowerRec Arrays: InBTArray and RecoShowerArray....
469
470
471 cout << " void Reco() Clear gEDA->GetTrackSet(BT); Array..." << endl;
473 set->Clear();
474
475 cout << " void Reco() Set InBTArray..." << endl;
476 cout << gEDA->GetSelected() << " gEDA->GetSelected() " << endl;
477
478 eInBTArray=(TObjArray *) gEDA->GetSelected()->Clone();
479 TObjArray* arr=(TObjArray*)eInBTArray->Clone();
480//arr->Print();
482
483 cout << " void Reco() Check For CheckInBT..." << endl;
484 CheckInBT();
485
486
487// Recheck if anything is selected:
488 if (eShowerRec->GetInBTArrayN()==0) return;
489
490 if (gEDBDEBUGLEVEL>2) cout << " void Reco() PrintShowerRecInBT..." << endl;
492
493 cout << " void Reco() Check For gAli..." << endl;
494 if (eIsAliLoaded==kFALSE) LoadShowerRecEdbPVRec();
495
496 cout << " void Reco() CheckPVRec..." << endl;
497 CheckPVRec();
498
499 cout << " void Reco() Check For suited BT density cut if needed..." << endl;
501
502 cout << " void Reco() eShowerRec->SetUseAliSub(0)..." << endl;
504
505 cout << " void Reco() eShowerRec->Execute()..." << endl;
507
508
509 if (gEDBDEBUGLEVEL>2) cout << " void Reco() eShowerRec->PrintRecoShowerArray()..." << endl;
511
514
515 if (eDrawShowers) {
516 for (Int_t i=0; i<eNShowers; ++i) {
517 DrawShower(i);
518 }
519 }
520
521// Set Factor values back:
522 cout << "// Set Factor values back: eShowerRec->SetQualityCutValues(0.15,1.0); " << endl;
524
525 return;
526}
void PrintShowerRecInBT()
Definition: EdbEDAShowerTab.h:100
void CheckBTDensity(EdbPVRec *pvr, TObjArray *tracks)
Definition: EdbEDAShowerTab.C:690
void CheckPVRec()
Definition: EdbEDAShowerTab.C:292
void CheckInBT()
Definition: EdbEDAShowerTab.C:447
void SetUseAliSub(Bool_t UseAliSub)
Definition: EdbShowerRec.h:625
void ResetInBTArray()
Definition: EdbShowerRec.h:719
void ResetRecoShowerArray()
Definition: EdbShowerRec.h:723

◆ Reset()

void EdbEDAShowerTab::Reset ( )
643 {
644
645 cout << "----- void EdbEDAShowerTab::Reset() -----"<<endl;
648 set->Clear();
649 gEDA->Redraw();
650 return;
651}
void ResetShowerRecParameters()
Definition: EdbEDAShowerTab.h:116

◆ ResetShowerRecParameters()

void EdbEDAShowerTab::ResetShowerRecParameters ( )
inline
116 {
118 return;
119 }
void ResetAlgoParameters()
Definition: EdbShowerRec.cxx:4806

◆ SetInBTFromLinkedTracks()

void EdbEDAShowerTab::SetInBTFromLinkedTracks ( )
inline
87 {
88 cout << "EdbEDAShowerTab::WARNING: SetInBTFromLinkedTracks() not yet implemented." << endl;
89 //eShowerRec->SetInBTArray(NULL); /// TODO
90 }

◆ SetInBTFromSelected()

void EdbEDAShowerTab::SetInBTFromSelected ( )
442 {
443 new TGMsgBox(gClient->GetRoot(), gEve->GetMainWindow(),"EDA", "Not yet supported! Here will come possibility to set Inititator BaseTracks either from reconstructed track, scanback tracks, BTs attaced to vertex, or similar...");//, 4, 4);
444 return;
445}

◆ SetNumberEntryShowerRecParameter()

void EdbEDAShowerTab::SetNumberEntryShowerRecParameter ( Double_t  val)
550 {
551 eNumberEntry_ParaValue->SetNumber(val);
552}

◆ SetParaSetHighNBT()

void EdbEDAShowerTab::SetParaSetHighNBT ( )

◆ SetParaSetHighPur()

void EdbEDAShowerTab::SetParaSetHighPur ( )

◆ SetPVRecFromFile()

void EdbEDAShowerTab::SetPVRecFromFile ( )
268 {
269
270 Log(2, "EdbEDAShowerTab", "--- SetPVRecFromFile()");
271 TFile* file = new TFile("File_EdbPVRec.root");
272 pvrec_cpfiles = (EdbPVRec*)file->Get("EdbPVRec");
273 if (NULL==pvrec_cpfiles) return;
274
275 if (pvrec_cpfiles->eTracks == NULL) cout << "EdbEDAShowerTab::SetPVRecFromFile() WARNING This file has no eTracks! "<<endl;
276 if (pvrec_cpfiles->eVTX == NULL) cout << "EdbEDAShowerTab::SetPVRecFromFile() WARNING This file has no eVTX! "<<endl;
277
280
281 eIsAliLoaded=kTRUE;
282
283 // Directly set it for the ShowerRec:
285
286 Log(2, "EdbEDAShowerTab", "--- SetPVRecFromFile()...done.");
287}
TObjArray * eVTX
array of vertex
Definition: EdbPVRec.h:162
TObjArray * eTracks
Definition: EdbPVRec.h:161

Member Data Documentation

◆ eCheckButtonCheckQualitycut

TGCheckButton* EdbEDAShowerTab::eCheckButtonCheckQualitycut
private

◆ eCheckButtonDrawShowers

TGCheckButton* EdbEDAShowerTab::eCheckButtonDrawShowers
private

◆ eCheckButtonUseTSData

TGCheckButton* EdbEDAShowerTab::eCheckButtonUseTSData
private

◆ eCheckQualitycut

Bool_t EdbEDAShowerTab::eCheckQualitycut
private

◆ eComboBox_Parameter

TGComboBox* EdbEDAShowerTab::eComboBox_Parameter
private

Buttons and so on:

◆ eComboBox_ParameterQualitySetting

TGComboBox* EdbEDAShowerTab::eComboBox_ParameterQualitySetting
private

◆ eDrawShowers

Bool_t EdbEDAShowerTab::eDrawShowers
private

◆ eInBTArray

TObjArray* EdbEDAShowerTab::eInBTArray
private

The InitiatorBT array. Has to be variable otherwise it will be out of scope:

◆ eIsAliLoaded

Bool_t EdbEDAShowerTab::eIsAliLoaded
private

◆ eIspvrec_cleaned

Bool_t EdbEDAShowerTab::eIspvrec_cleaned
private

◆ eNShowers

Int_t EdbEDAShowerTab::eNShowers
private

An own Shower array. Not the one of the ShowerRec object:

◆ eNumberEntry_ParaValue

TGNumberEntry* EdbEDAShowerTab::eNumberEntry_ParaValue
private

◆ eParamWindow

TGMainFrame* EdbEDAShowerTab::eParamWindow
private

◆ eQualityCutSetting

Int_t EdbEDAShowerTab::eQualityCutSetting
private

◆ eShower

EdbTrackP* EdbEDAShowerTab::eShower
private

Actual Shower (current one, either for drawing or printing)

◆ eShowerRec

EdbShowerRec* EdbEDAShowerTab::eShowerRec
private

The ShowerRec object:

◆ eTGNumberEntryDrawShowerNr

TGNumberEntry* EdbEDAShowerTab::eTGNumberEntryDrawShowerNr
private

◆ eUseTSData

Bool_t EdbEDAShowerTab::eUseTSData
private

Various Variables used for settings:

◆ pvrec_cleaned

EdbPVRec* EdbEDAShowerTab::pvrec_cleaned
private

(==generic cleaned)

◆ pvrec_cleaned_cpfiles

EdbPVRec* EdbEDAShowerTab::pvrec_cleaned_cpfiles
private

"fake" pvrec object with adapted BG quality cut in it.

◆ pvrec_cleaned_linkedtracks

EdbPVRec* EdbEDAShowerTab::pvrec_cleaned_linkedtracks
private

"fake" pvrec object with adapted BG quality cut in it.

◆ pvrec_cpfiles

EdbPVRec* EdbEDAShowerTab::pvrec_cpfiles
private

◆ pvrec_generic

EdbPVRec* EdbEDAShowerTab::pvrec_generic
private

One Pointer for the generic EdbPVRec object:
This is the object where ALWASY the shower reco is performed!

◆ pvrec_linkedtracks

EdbPVRec* EdbEDAShowerTab::pvrec_linkedtracks
private

Two Pointers for the EdbPVRec object:

◆ RecoShowerArray

TObjArray* EdbEDAShowerTab::RecoShowerArray
private

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