]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/trigger/AliTriggerTask.cxx
adding plot to histogram spd online trigger bits
[u/mrichter/AliRoot.git] / PWG0 / trigger / AliTriggerTask.cxx
1 /* $Id: AliTriggerTask.cxx 35782 2009-10-22 11:54:31Z jgrosseo $ */
2
3 #include "AliTriggerTask.h"
4
5 #include <TCanvas.h>
6 #include <TFile.h>
7 #include <TChain.h>
8 #include <TH1F.h>
9 #include <TH2F.h>
10 #include <TH3F.h>
11
12 #include <AliLog.h>
13 #include <AliESDEvent.h>
14 #include <AliHeader.h>
15 #include <AliAnalysisManager.h>
16 #include <AliESDInputHandler.h>
17 #include <AliESDHeader.h>
18 #include <AliTriggerAnalysis.h>
19
20 ClassImp(AliTriggerTask)
21
22 AliTriggerTask::AliTriggerTask(const char* opt) :
23   AliAnalysisTask("AliTriggerTask", ""),
24   fESD(0),
25   fOutput(0),
26   fOption(opt),
27   fStartTime(0),
28   fEndTime(0),
29   fUseOrbits(kFALSE),
30   fFirstOrbit(0),
31   fLastOrbit(0),
32   fNTriggers(0),
33   fTriggerList(0),
34   fStats(0),
35   fTrigger(0)
36 {
37   //
38   // Constructor. Initialization of pointers
39   //
40
41   // Define input and output slots here
42   DefineInput(0, TChain::Class());
43   DefineOutput(0, TList::Class());
44   
45   fNTriggers = 14;
46   
47   static AliTriggerAnalysis::Trigger triggerList[] = { AliTriggerAnalysis::kAcceptAll, AliTriggerAnalysis::kFPANY, AliTriggerAnalysis::kMB1, AliTriggerAnalysis::kMB2, AliTriggerAnalysis::kMB3, AliTriggerAnalysis::kSPDGFO, AliTriggerAnalysis::kSPDGFOBits, AliTriggerAnalysis::kV0A, AliTriggerAnalysis::kV0C, AliTriggerAnalysis::kZDC, AliTriggerAnalysis::kZDCA, AliTriggerAnalysis::kZDCC, AliTriggerAnalysis::kFMDA, AliTriggerAnalysis::kFMDC };
48   fTriggerList = triggerList;
49   
50   fStats = new TH1*[fNTriggers];
51   
52   fTrigger = new AliTriggerAnalysis;
53   fTrigger->EnableHistograms();
54   
55   AliLog::SetClassDebugLevel("AliTriggerTask", AliLog::kWarning);
56 }
57
58 AliTriggerTask::~AliTriggerTask()
59 {
60   //
61   // Destructor
62   //
63
64   // histograms are in the output list and deleted when the output
65   // list is deleted by the TSelector dtor
66
67   if (fOutput) {
68     delete fOutput;
69     fOutput = 0;
70   }
71 }
72
73 //________________________________________________________________________
74 void AliTriggerTask::ConnectInputData(Option_t *)
75 {
76   // Connect ESD
77   // Called once
78
79   Printf("AliTriggerTask::ConnectInputData called");
80
81   AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
82
83   if (!esdH) {
84     Printf("ERROR: Could not get ESDInputHandler");
85   } else {
86     fESD = esdH->GetEvent();
87     
88     TString branches("AliESDHeader Vertex AliMultiplicity ALIESDVZERO ALIESDZDC FMD");
89
90     // Enable only the needed branches
91     esdH->SetActiveBranches(branches);
92   }
93 }
94
95 void AliTriggerTask::CreateOutputObjects()
96 {
97   // create result objects and add to output list
98
99   Printf("AliTriggerTask::CreateOutputObjects");
100
101   fOutput = new TList;
102   fOutput->SetOwner();
103   
104   if (fStartTime == fEndTime)
105     AliWarning("Start and endtime not set. Automatic binning will be used. This does not work in parallel systems");
106
107   Int_t nBins = 1000;
108   if (fEndTime - fStartTime > 0)
109     nBins = fEndTime - fStartTime + 1;
110   if (nBins > 10000)
111     nBins = 10000;
112   
113   Int_t start = 0;
114   Int_t end = fEndTime - fStartTime;
115   
116   if (fUseOrbits)
117   {
118     start = fStartTime;
119     end = fEndTime;
120   }
121   
122   for (Int_t i=0; i<fNTriggers; i++)
123   {
124     fStats[i] = new TH1F(Form("fStats_%d", i), Form("%s;%s;counts", AliTriggerAnalysis::GetTriggerName(fTriggerList[i]), (fUseOrbits) ? "orbit number" : "time"), nBins, start - 0.5, end + 0.5);
125     fOutput->Add(fStats[i]);
126   }
127   
128   fFirstOrbit = new TParameter<Long_t> ("fFirstOrbit", 0);
129   fLastOrbit = new TParameter<Long_t> ("fLastOrbit", 0);
130   fOutput->Add(fFirstOrbit);
131   fOutput->Add(fLastOrbit);
132   
133   fOutput->Add(fTrigger);
134 }
135
136 void AliTriggerTask::Exec(Option_t*)
137 {
138   // process the event
139
140   // post the data already here
141   PostData(0, fOutput);
142
143   if (!fESD)
144   {
145     AliError("ESD branch not available");
146     return;
147   }
148   
149   // check event type (should be PHYSICS = 7)
150   AliESDHeader* esdHeader = fESD->GetHeader();
151   if (!esdHeader)
152   {
153     Printf("ERROR: esdHeader could not be retrieved");
154     return;
155   }
156   
157   /*
158   UInt_t eventType = esdHeader->GetEventType();
159   if (eventType != 7)
160   {
161     Printf("Skipping event because it is of type %d", eventType);
162     return;
163   }
164   */
165   //Printf("Trigger classes: %s:", fESD->GetFiredTriggerClasses().Data());
166   
167   fTrigger->FillHistograms(fESD);
168   
169   UInt_t timeStamp = 0;
170   if (fUseOrbits)
171     timeStamp = fESD->GetOrbitNumber();
172   else
173     timeStamp = fESD->GetTimeStamp() - fStartTime;
174   //Printf("%d", timeStamp);
175   
176   for (Int_t i = 0; i < fNTriggers; i++)
177   {
178     Bool_t triggered = fTrigger->IsOfflineTriggerFired(fESD, fTriggerList[i]);
179     if (triggered)
180       fStats[i]->Fill(timeStamp);
181     //Printf("%s: %d", AliTriggerAnalysis::GetTriggerName(fTriggerList[i]), triggered);
182   }
183   
184   if (fFirstOrbit->GetVal() == 0)
185     fFirstOrbit->SetVal(esdHeader->GetOrbitNumber());
186   else
187     fFirstOrbit->SetVal(TMath::Min(fFirstOrbit->GetVal(), (Long_t) esdHeader->GetOrbitNumber()));
188   
189   fLastOrbit->SetVal(TMath::Max(fLastOrbit->GetVal(), (Long_t) esdHeader->GetOrbitNumber()));
190 }
191
192 void AliTriggerTask::Terminate(Option_t *)
193 {
194   // The Terminate() function is the last function to be called during
195   // a query. It always runs on the client, it can be used to present
196   // the results graphically or save the results to file.
197
198   fOutput = dynamic_cast<TList*> (GetOutputData(0));
199   if (!fOutput)
200     Printf("ERROR: fOutput not available");
201     
202   fOutput->Print();
203
204   if (fOutput)
205   {
206     for (Int_t i=0; i<fNTriggers; i++)
207       fStats[i] = dynamic_cast<TH1*> (fOutput->FindObject(Form("fStats_%d", i)));
208     fTrigger = dynamic_cast<AliTriggerAnalysis*> (fOutput->FindObject("AliTriggerAnalysis"));
209   }
210
211   TFile* fout = new TFile("trigger.root", "RECREATE");
212
213   for (Int_t i=0; i<fNTriggers; i++)
214     if (fStats[i])
215       fStats[i]->Write();
216   if (fTrigger)
217     fTrigger->WriteHistograms();
218     
219   if (fFirstOrbit)
220     fFirstOrbit->Dump();
221   if (fLastOrbit)
222     fLastOrbit->Dump();
223     
224   fout->Write();
225   fout->Close();
226   
227   Int_t nX = (Int_t) TMath::Sqrt(fNTriggers);
228   Int_t nY = nX;
229   
230   while (nX * nY < fNTriggers)
231   {
232     if (nX == nY)
233       nX++;
234     else
235       nY++;
236   }
237   
238   TCanvas* c = new TCanvas("c", "c", 800, 800);
239   c->Divide(nX, nY);
240   
241   Printf("+++++++++ TRIGGER STATS:");
242
243   Int_t base = 1;
244   if (fStats[0])
245     base = (Int_t) fStats[0]->Integral();
246
247   Int_t length = fEndTime - fStartTime;
248   
249   for (Int_t i=0; i<fNTriggers; i++)
250     if (fStats[i])
251     {
252       c->cd(i+1);
253       fStats[i]->Draw();
254       Printf("%s: %d triggers | %f %% of all triggered | Rate: %f Hz", AliTriggerAnalysis::GetTriggerName(fTriggerList[i]), (UInt_t) fStats[i]->Integral(), 100.0 * fStats[i]->Integral() / base, (length > 0) ? (fStats[i]->Integral() / length) : -1);
255     }
256
257   Printf("Writting result to trigger.root");
258 }