updates in the PWGPP/TRD
[u/mrichter/AliRoot.git] / PWGPP / TRD / info / AliTRDtriggerInfo.cxx
1 ////////////////////////////////////////////////////////////////////////////
2 //                                                                        //
3 //  Mergable Trigger List                                                 //
4 //                                                                        //
5 //  Authors:                                                              //
6 //    Alexandru Bercuci <A.Bercuci@gsi.de>                                //
7 //                                                                        //
8 ////////////////////////////////////////////////////////////////////////////
9
10
11
12 #include "TObjArray.h"
13 #include "TObjString.h"
14 #include "TH1.h"
15
16 #include "AliLog.h"
17
18 #include "AliTRDtriggerInfo.h"
19
20 ClassImp(AliTRDtriggerInfo)
21 //____________________________________________
22 AliTRDtriggerInfo::AliTRDtriggerInfo()
23   :TObject()
24   ,fTriggerList(NULL)
25   ,fHisto(NULL)
26 {
27 //  Constructor. Reset all fields.
28   memset(fTriggerSel, 0, kTriggerListSize *sizeof(Bool_t));
29   memset(fTriggerStat, 0, kTriggerListSize *sizeof(Int_t));
30 }
31
32 //____________________________________________
33 AliTRDtriggerInfo::~AliTRDtriggerInfo()
34 {
35 //  Destructor. 
36   if(fHisto) delete fHisto;
37   if(fTriggerList){ fTriggerList->Delete(); delete fTriggerList;}
38 }
39
40 //____________________________________________
41 Int_t AliTRDtriggerInfo::Add(const Char_t *trigger, Int_t nstat, Bool_t select)
42 {
43 // Add trigger named "trigger" to the list. If trigger is not on the list it is add.
44 // The statistics of this trigger is set to "nstat" [default=1]
45 // On return the position incremented
46
47   Int_t itrig(-1), nt(0);
48   if((nt=GetNTriggers()) && (itrig = GetTrigger(trigger))>=0){
49     fTriggerSel[itrig]  = select;
50     fTriggerStat[itrig]+= nstat;
51   } else {
52     if(!nt) fTriggerList = new TObjArray;
53     fTriggerList->AddAt(new TObjString(trigger), nt);
54     fTriggerStat[nt] = nstat;
55     fTriggerSel[nt]  = select;
56     nt++;
57   }
58   return nt;
59 }
60
61 //____________________________________________
62 void AliTRDtriggerInfo::Draw(Option_t *opt)
63 {
64 // Draw trigger statistics. Via parameter "opt" one can select a sublist of trigger names
65
66   Int_t nt(0);
67   if(!(nt = GetNTriggers())) return;
68   if(fHisto && fHisto->GetNbinsX() != nt){ delete fHisto; fHisto=NULL;}
69 //
70   if(!fHisto) fHisto = new TH1F("hTriggerStat", "Trigger Statistics;TRIGGER;Freq. [%]", nt, 0.5, nt+0.5);
71   TAxis *ax(fHisto->GetXaxis());  fHisto->Reset();
72   ax->SetTitleOffset(6.5); ax->CenterTitle();
73   fHisto->SetFillStyle(3001);fHisto->SetFillColor(kGreen);
74   fHisto->SetBarWidth(0.8);fHisto->SetBarOffset(0.1);
75 //
76   Int_t ibin(0);
77   if(strcmp(opt, "")!=0){
78     TString sopt(opt);
79     TObjArray *optTriggers = sopt.Tokenize(" ");
80     for(Int_t iopt(0); iopt<optTriggers->GetEntries(); iopt++){
81       const Char_t *trigger = ((TObjString*)(*optTriggers)[iopt])->GetName();
82       Int_t itrig(-1);
83       if((itrig = GetTrigger(trigger))<0) continue;
84       ibin++;
85       ax->SetBinLabel(ibin, trigger);
86       fHisto->SetBinContent(ibin, fTriggerStat[itrig]);
87     }
88     optTriggers->Delete(); delete optTriggers;
89   } else {
90     for(Int_t jtrig(nt); jtrig--;){
91       ibin++;
92       ax->SetBinLabel(ibin, ((TObjString*)(*fTriggerList)[jtrig])->GetName());
93       fHisto->SetBinContent(ibin, fTriggerStat[jtrig]);
94     }
95   }
96   ax->SetRange(1, ibin);
97   fHisto->Scale(1.e2/fHisto->Integral());
98   fHisto->Draw("hbar2");
99
100   return;
101 }
102
103 //____________________________________________
104 Int_t AliTRDtriggerInfo::GetNTriggers() const
105 {
106 // Return no of trigger classes
107   return fTriggerList?fTriggerList->GetEntries():0;
108 }
109
110 //____________________________________________
111 const char* AliTRDtriggerInfo::GetTrigger(Int_t it) const
112 {
113 // Return trigger by index. 0 in case it is not found
114   if(it<0||it>=GetNTriggers()) return 0;
115   return ((TObjString*)(*fTriggerList)[it])->GetName();
116 }
117
118 //____________________________________________
119 Int_t AliTRDtriggerInfo::GetTrigger(const char *trigger) const
120 {
121 // Return trigger index by name. -1 in case trigger not registred
122
123   for(Int_t jtrig(GetNTriggers()); jtrig--;){
124     if(((TObjString*)(*fTriggerList)[jtrig])->String().CompareTo(trigger)==0) return jtrig;
125   }
126   return -1;
127 }
128
129 //____________________________________________
130 Long64_t AliTRDtriggerInfo::Merge(TCollection* list)
131 {
132 // Merge list of trigger statistics. Add new triggers if neccessary
133
134   Int_t imerge(1);
135   AliTRDtriggerInfo *o(NULL);
136   TIterator *it(list->MakeIterator());
137   while((o  = (AliTRDtriggerInfo*)it->Next())){
138     for(Int_t jtrig(o->GetNTriggers()); jtrig--;) Add(((TObjString*)(*o->fTriggerList)[jtrig])->GetName(), o->fTriggerStat[jtrig], o->fTriggerSel[jtrig]);
139     imerge++;
140   }
141   return imerge;
142 }
143   
144 //____________________________________________
145 void AliTRDtriggerInfo::Print(Option_t *opt) const
146 {
147 // Print trigger statistics. If opt="a" detailed trigger/statistics is also provided
148
149   Int_t nt(GetNTriggers());
150   Int_t stat(0); for(Int_t jtrig(0); jtrig<nt; jtrig++) stat+=fTriggerStat[jtrig];
151   printf("TriggerClasses[%3d] Stat[%5d]\n", nt, stat);
152   if(!nt || strcmp(opt,"a")!=0) return;
153
154   for(Int_t jtrig(0); jtrig<nt; jtrig++) printf("  %3d \"%s\" [%5d]\n", jtrig, ((TObjString*)(*fTriggerList)[jtrig])->GetName(), fTriggerStat[jtrig]);
155 }
156
157