]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibTrigger.cxx
Fix of coverity warning (sent by email, but not seen in the coverity-webpage and...
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibTrigger.cxx
1
2 /**************************************************************************
3  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  *                                                                        *
5  * Author: The ALICE Off-line Project.                                    *
6  * Contributors are mentioned in the code where appropriate.              *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 /*
18   // Load libraries
19   gSystem->Load("libANALYSIS");
20   gSystem->Load("libTPCcalib");
21   
22   
23   .x ~/NimStyle.C
24   gSystem->Load("libANALYSIS");
25   gSystem->Load("libTPCcalib");
26
27   TFile f("CalibObjects.root");
28   AliTPCcalibTrigger *calibTrigger = (AliTPCcalibTrigger *)f->Get("TPCCalib")->FindObject("calibTrigger");
29
30
31   gSystem->AddIncludePath("-I$ALICE_ROOT/TPC/macros");
32   gROOT->LoadMacro("$ALICE_ROOT/TPC/macros/AliXRDPROOFtoolkit.cxx+")
33   AliXRDPROOFtoolkit tool;
34   TChain * chainTrack = tool.MakeChain("trigger.txt","Track",0,10200);
35   chainTrack->Lookup();
36   TChain * chainEvent = tool.MakeChain("trigger.txt","Event",0,10200);
37   chainEvent->Lookup();
38
39
40 */
41
42 #include "Riostream.h"
43 #include "TChain.h"
44 #include "TTree.h"
45 #include "TH1F.h"
46 #include "TH2F.h"
47 #include "TH3F.h"
48 #include "THnSparse.h"
49 #include "TList.h"
50 #include "TMath.h"
51 #include "TCanvas.h"
52 #include "TFile.h"
53 #include "TF1.h"
54 #include "TVectorD.h"
55 #include "TProfile.h"
56 #include "TGraphErrors.h"
57 #include "TCanvas.h"
58
59 #include "AliTPCclusterMI.h"
60 #include "AliTPCseed.h"
61 #include "AliESDVertex.h"
62 #include "AliESDEvent.h"
63 #include "AliESDfriend.h"
64 #include "AliESDInputHandler.h"
65 #include "AliAnalysisManager.h"
66
67 #include "AliTracker.h"
68 #include "AliMagF.h"
69 #include "AliTPCCalROC.h"
70
71 #include "AliLog.h"
72
73 #include "AliTPCcalibTrigger.h"
74
75 #include "TTreeStream.h"
76 #include "AliTPCTracklet.h"
77 #include "TTimeStamp.h"
78 #include "AliTPCcalibDB.h"
79 #include "AliTPCcalibLaser.h"
80 #include "AliDCSSensorArray.h"
81 #include "AliDCSSensor.h"
82
83 ClassImp(AliTPCcalibTrigger)
84
85 AliTPCcalibTrigger::AliTPCcalibTrigger():
86   AliTPCcalibBase("calibTrigger","calibTrigger"),
87   fHisMap(0)
88 {
89
90 }
91
92 AliTPCcalibTrigger::AliTPCcalibTrigger(const char * name, const char * title):
93   AliTPCcalibBase(name,title),
94   fHisMap(0)
95 {
96   //
97   //
98   //
99   fHisMap = new TMap;
100 }
101
102 AliTPCcalibTrigger::~AliTPCcalibTrigger(){
103   //
104   // delete histograms
105   // class is owner of all histograms
106   //
107   if (!fHisMap) return;
108   fHisMap->SetOwner(kTRUE);
109   fHisMap->DeleteAll();
110   delete fHisMap;
111 }
112
113
114 Long64_t AliTPCcalibTrigger::Merge(TCollection *li) {
115   //
116   // Merge histograms
117   //
118   TIterator* iter = li->MakeIterator();
119   AliTPCcalibTrigger* cal = 0;
120
121   while ((cal = (AliTPCcalibTrigger*)iter->Next())) {
122     if (!cal->InheritsFrom(AliTPCcalibTrigger::Class())) {
123       Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
124       return -1;
125     }
126     TMap * addMap=(cal->fHisMap);
127     if(!addMap) return 0;
128     TIterator* iterator = addMap->MakeIterator();
129     iterator->Reset();
130     //    TPair* addPair=0;
131     TObject *object=0;
132     //
133     while((object=iterator->Next())){
134       THnSparse* his1 = dynamic_cast<THnSparseF*>(cal->fHisMap->GetValue(object->GetName()));
135       if (!his1) continue;      
136       his1->Print();
137       THnSparse* his0 = dynamic_cast<THnSparseF*>(fHisMap->GetValue(object->GetName()));
138
139       if(!his0){
140         his0=MakeHisto(object->GetName());
141         fHisMap->Add(new TObjString(object->GetName()),his0);
142       }
143       his0->Add(his1);
144     }
145   }
146   return 0;
147 }
148
149
150
151 void AliTPCcalibTrigger::Process(AliESDEvent *event){
152   //
153   //
154   //
155   if (!event) return;
156   const TString &trigger = event->GetFiredTriggerClasses();
157   TTreeSRedirector * cstream =  GetDebugStreamer();
158   //
159   TObjString str(event->GetFiredTriggerClasses());
160   Bool_t hasPIXEL=HasPIXEL(&str);
161   Int_t  hasTRD=HasTRD(&str);
162   Bool_t hasTOF=HasTOF(&str);
163   Bool_t hasACORDE=HasACORDE(&str);
164   //
165   if (!GetHisto(trigger.Data())){
166     AddHisto(trigger.Data(),MakeHisto(trigger.Data()));
167   }
168   if (!GetHisto("all")){
169     AddHisto("all",MakeHisto("all"));
170   }
171
172   THnSparse *histoAll = GetHisto("all");
173   THnSparse *histo = GetHisto(trigger.Data());
174   Double_t xcont[9]={0,0,0,0,0,0,0,0,0};
175   
176   Int_t ntracks = event->GetNumberOfTracks();
177   xcont[0] = ntracks;
178   xcont[8] = 1;
179   //
180   // GetLongest track
181   //  
182   AliESDtrack * longest=0;
183   Int_t nclmax=0;
184   for (Int_t itrack=0; itrack<ntracks; itrack++){
185     AliESDtrack *track=event->GetTrack(itrack);
186     if (!track) continue;
187     if (track->GetTPCNcls()<=nclmax) continue;
188     nclmax = track->GetTPCNcls();
189     longest= track;
190   }
191   xcont[1]  =nclmax;
192   histoAll->Fill(xcont);
193   histo->Fill(xcont);
194   if (cstream) {
195     (*cstream) << "Event" <<
196       "run="<<fRun<<
197       "time="<<fTime<<
198       "tname.="<<&str<<
199       "pixel="<<hasPIXEL<<
200       "trd="<<hasTRD<<
201       "tof="<<hasTOF<<
202       "acorde="<<hasACORDE<<
203       "ntracks="<<ntracks<<
204       "\n";
205   }
206   //
207   xcont[8] = -1.;
208   for (Int_t itrack=0; itrack<ntracks; itrack++){
209     AliESDtrack *track=event->GetTrack(itrack);
210     if (!track) continue;
211     Float_t dca[2];
212     Double_t pxyz[3];
213     track->GetDZ(0.,0.,0.,event->GetMagneticField(),dca);
214     Bool_t status = track->GetPxPyPz(pxyz);
215     Double_t alpha = TMath::ATan2(pxyz[1],pxyz[0]);
216     xcont[1]=track->GetTPCNcls();
217     xcont[2]=dca[0];
218     xcont[3]=dca[1];
219     xcont[4]=alpha;
220     xcont[5]=track->GetParameter()[3];
221     xcont[6]=track->Pt();
222     xcont[7]=track->GetTPCsignal();    
223     histoAll->Fill(xcont);
224     histo->Fill(xcont);
225     //
226     //
227     if (cstream) {
228       Double_t mpt = track->GetParameter()[4];
229       Int_t kokot[1000];
230       Int_t nclITS=track->GetITSclusters(kokot);
231       Int_t nclTPC=track->GetTPCNcls();
232       Int_t nclTRD=track->GetTRDclusters(kokot);
233       Int_t ntlTRD=track->GetTRDntracklets();
234       ULong_t tstatus = track->GetStatus();
235       (*cstream) << "Track" <<
236         "run="<<fRun<<
237         "time="<<fTime<<
238         "tname.="<<&str<<
239         "status="<<status<<     
240         "tstatus="<<tstatus<<   
241         //
242         "ntracks="<<ntracks<<
243         "tstatus="<<status<<
244         "nclITS="<<nclITS<<
245         "nclTPC="<<nclTPC<<
246         "nclTRD="<<nclTRD<<
247   "ntlTRD="<<ntlTRD<<
248         //
249         "pixel="<<hasPIXEL<<
250         "trd="<<hasTRD<<
251         "tof="<<hasTOF<<
252         "acorde="<<hasACORDE<<
253         "ncl="<<xcont[1]<<
254         "dcaR="<<xcont[2]<<
255         "dcaZ="<<xcont[3]<<
256         "alpha="<<xcont[4]<<
257         "theta="<<xcont[5]<<
258         "pt="<<xcont[6]<<
259         "dEdx="<<xcont[7]<<     
260         "mpt="<<mpt<<
261         "\n";
262     }
263   }
264 }
265
266 THnSparse *AliTPCcalibTrigger::MakeHisto(const char* trigger){
267   //
268   // Make event/track histograms
269   // trigger histo 
270   //
271   //                 ntracks  nclMax dcaR dcaZ alpha   theta pt dEdx ev
272   Int_t    bins[9] = {50,     40,    20,  20,  18,     25,   25, 25, 2 };
273   //Int_t    bins[9] = {50*   20*    25*  25*  18*     25*   25* 25 };
274   Double_t xmin[9] = {0.,     0,      0,  -250,   -3.14, -1.5,  0, 0, -1.};
275   Double_t xmax[9] = {50,     160,  150,  250, 3.14,  1.5,   100, 100, 1.};
276   TString  axisName[9]={
277     "ntracks",
278     "ncl",
279     "dcaR",
280     "dcaZ",
281     "alpha",
282     "theta",
283     "pt",
284     "dEdx",
285     "ev"
286   };
287   TString  axisTitle[9]={
288     "Number of tracks",
289     "N_{cl}",
290     "dca_{R} (cm)",
291     "dca_{z} (cm)",
292     "alpha (mrad)",
293     "theta",
294     "p_{t} (GeV/c)",
295     "dEdx (a.u.)",
296     "ev"
297   };
298
299   
300   THnSparse *sparse = new THnSparseF(Form("his_%s",trigger), Form("his_%s",trigger), 9, bins, xmin, xmax);
301   for (Int_t iaxis=0; iaxis<9; iaxis++){
302     sparse->GetAxis(iaxis)->SetName(axisName[iaxis]);
303     sparse->GetAxis(iaxis)->SetTitle(axisTitle[iaxis]);
304   }
305   return sparse;
306 }
307
308 THnSparse * AliTPCcalibTrigger::GetHisto(const char *trigger) { 
309   //
310   // return histogram for given class
311   if (!fHisMap) fHisMap=new TMap;
312   return (THnSparse*) fHisMap->GetValue(trigger);
313 }
314
315 void   AliTPCcalibTrigger::AddHisto(const char *trigger, THnSparse *his) { 
316   if (!GetHisto(trigger)) {
317     TObjString *pstr = new TObjString(trigger);
318     fHisMap->Add(pstr,his);
319   }
320 }
321
322 TTree * AliTPCcalibTrigger::MakeTree(const char * fname){
323   //
324   //
325   //
326   TTreeSRedirector * sred = new TTreeSRedirector(fname);
327   TTreeStream &pcstream = (*sred)<<"Trigger";
328   //
329   //
330   TIterator* iterator = fHisMap->MakeIterator();
331   TObject * object=0;
332   //
333   while((object=iterator->Next())){
334     MakeTree(pcstream, object->GetName());
335   }
336   delete sred;
337   TFile *f = new TFile(fname);
338   TTree *tree = (TTree*)f->Get("Trigger");
339   return tree;
340 }
341
342
343 void AliTPCcalibTrigger::MakeTree(TTreeStream &pcstream, const char *tname){
344   //
345   //  TTreeSRedirector * sred = new TTreeSRedirector("trigger.root");
346   //  TTreeStream &pcstream = (*sred)<<"Trigger";
347   //
348   //AliTPCcalibTrigger *calibTrigger = this;
349   Double_t value;
350   THnSparse * his = GetHisto(tname);
351   if (!his) return;
352   //
353   Int_t bins[1000];
354   Int_t ndim = his->GetNdimensions();
355   Double_t position[10];
356   //
357   TObjString str(tname);
358   Bool_t isAll  = str.String().Contains("all");
359   Bool_t hasPIXEL=HasPIXEL(&str);
360   Int_t  hasTRD=HasTRD(&str);
361   Bool_t hasTOF=HasTOF(&str);
362   Bool_t hasACORDE=HasACORDE(&str);
363   for (Long64_t i = 0; i < his->GetNbins(); ++i) {
364     value = his->GetBinContent(i, bins);
365     pcstream<<"val="<<value;
366     pcstream<<"tname.="<<&str;
367     //
368     pcstream<<"all="<<isAll;
369     pcstream<<"pixel="<<hasPIXEL;
370     pcstream<<"trd="<<hasTRD;
371     pcstream<<"tof="<<hasTOF;
372     pcstream<<"acorde="<<hasACORDE;
373     //
374     for (Int_t idim = 0; idim < ndim; idim++) {
375       position[idim] = his->GetAxis(idim)->GetBinCenter(bins[idim]);
376       pcstream<<Form("%s=",his->GetAxis(idim)->GetName())<<position[idim];
377     }
378     pcstream<<"\n";
379   }
380 }
381
382
383 Bool_t AliTPCcalibTrigger::HasTOF(TObjString *tname){
384   //
385   Bool_t result = kFALSE;
386   result|=(tname->String().Contains("0OB")>0);
387   result|=(tname->String().Contains("0OC")>0);
388   return result;
389 }
390
391 Bool_t AliTPCcalibTrigger::HasACORDE(TObjString *tname){
392   Bool_t result = kFALSE;
393   result|=(tname->String().Contains("0ASL")>0);
394   result|=(tname->String().Contains("0AMU")>0);
395   result|=(tname->String().Contains("0ASC")>0);
396   return result;
397 }
398
399 Bool_t AliTPCcalibTrigger::HasPIXEL(TObjString *tname){
400   return (tname->String().Contains("0SCO")>0);
401 }
402
403 Int_t AliTPCcalibTrigger::HasTRD(TObjString *tname){
404   //
405   // Returns a mask containing TRD trigger information
406   // 0: No TRD trigger fired
407   // 1: TRD L1 fired
408   // 2: TRD L0 (krypton trigger) fired
409   //
410   Int_t result = 0;
411   if(tname->String().Contains("TRD")) result = 1;     // Normal TRD L1 name
412   if(tname->String().Contains("0HPT1")) result = 1;   // Old TRD L1 name
413   if(tname->String().Contains("0HWU") && !tname->String().Contains("TRD")) result = 2;  // pretrigger always input for L1
414   return result;
415 }
416