]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibTrigger.cxx
New task to check TPC-ITS track prolongation eff with cosmics
[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 Long64_t AliTPCcalibTrigger::Merge(TCollection *li) {
103   //
104   // Merge histograms
105   //
106   TIterator* iter = li->MakeIterator();
107   AliTPCcalibTrigger* cal = 0;
108
109   while ((cal = (AliTPCcalibTrigger*)iter->Next())) {
110     if (!cal->InheritsFrom(AliTPCcalibTrigger::Class())) {
111       Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
112       return -1;
113     }
114     TMap * addMap=(cal->fHisMap);
115     if(!addMap) return 0;
116     TIterator* iterator = addMap->MakeIterator();
117     iterator->Reset();
118     //    TPair* addPair=0;
119     TObject *object=0;
120     //
121     while((object=iterator->Next())){
122       THnSparse* his1 = dynamic_cast<THnSparseF*>(cal->fHisMap->GetValue(object->GetName()));
123       if (!his1) continue;      
124       his1->Print();
125       THnSparse* his0 = dynamic_cast<THnSparseF*>(fHisMap->GetValue(object->GetName()));
126
127       if(!his0){
128         his0=MakeHisto(object->GetName());
129         fHisMap->Add(new TObjString(object->GetName()),his0);
130       }
131       his0->Add(his1);
132     }
133   }
134   return 0;
135 }
136
137
138
139 void AliTPCcalibTrigger::Process(AliESDEvent *event){
140   //
141   //
142   //
143   if (!event) return;
144   const TString &trigger = event->GetFiredTriggerClasses();
145   TTreeSRedirector * cstream =  GetDebugStreamer();
146   //
147   TObjString str(event->GetFiredTriggerClasses());
148   Bool_t hasPIXEL=HasPIXEL(&str);
149   Bool_t hasTRD=HasTRD(&str);
150   Bool_t hasTOF=HasTOF(&str);
151   Bool_t hasACORDE=HasACORDE(&str);
152   //
153   if (!GetHisto(trigger.Data())){
154     AddHisto(trigger.Data(),MakeHisto(trigger.Data()));
155   }
156   if (!GetHisto("all")){
157     AddHisto("all",MakeHisto("all"));
158   }
159
160   THnSparse *histoAll = GetHisto("all");
161   THnSparse *histo = GetHisto(trigger.Data());
162   Double_t xcont[9]={0,0,0,0,0,0,0,0,0};
163   
164   Int_t ntracks = event->GetNumberOfTracks();
165   xcont[0] = ntracks;
166   xcont[8] = 1;
167   //
168   // GetLongest track
169   //  
170   AliESDtrack * longest=0;
171   Int_t nclmax=0;
172   for (Int_t itrack=0; itrack<ntracks; itrack++){
173     AliESDtrack *track=event->GetTrack(itrack);
174     if (!track) continue;
175     if (track->GetTPCNcls()<=nclmax) continue;
176     nclmax = track->GetTPCNcls();
177     longest= track;
178   }
179   xcont[1]  =nclmax;
180   histoAll->Fill(xcont);
181   histo->Fill(xcont);
182   if (cstream) {
183     (*cstream) << "Event" <<
184       "run="<<fRun<<
185       "time="<<fTime<<
186       "tname.="<<&str<<
187       "pixel="<<hasPIXEL<<
188       "trd="<<hasTRD<<
189       "tof="<<hasTOF<<
190       "acorde="<<hasACORDE<<
191       "ntracks="<<ntracks<<
192       "\n";
193   }
194   //
195   xcont[8] = -1.;
196   for (Int_t itrack=0; itrack<ntracks; itrack++){
197     AliESDtrack *track=event->GetTrack(itrack);
198     if (!track) continue;
199     Float_t dca[2];
200     Double_t pxyz[3];
201     track->GetDZ(0.,0.,0.,event->GetMagneticField(),dca);
202     Bool_t status = track->GetPxPyPz(pxyz);
203     Double_t alpha = TMath::ATan2(pxyz[1],pxyz[0]);
204     xcont[1]=track->GetTPCNcls();
205     xcont[2]=dca[0];
206     xcont[3]=dca[1];
207     xcont[4]=alpha;
208     xcont[5]=track->GetParameter()[3];
209     xcont[6]=track->Pt();
210     xcont[7]=track->GetTPCsignal();    
211     histoAll->Fill(xcont);
212     histo->Fill(xcont);
213     //
214     //
215     if (cstream) {
216       Double_t mpt = track->GetParameter()[4];
217       Int_t kokot[1000];
218       Int_t nclITS=track->GetITSclusters(kokot);
219       Int_t nclTPC=track->GetTPCNcls();
220       Int_t nclTRD=track->GetTRDclusters(kokot);
221       ULong_t tstatus = track->GetStatus();
222       (*cstream) << "Track" <<
223         "run="<<fRun<<
224         "time="<<fTime<<
225         "tname.="<<&str<<
226         "status="<<status<<     
227         //
228         "ntracks="<<ntracks<<
229         "tstatus="<<status<<
230         "nclITS="<<nclITS<<
231         "nclTPC="<<nclTPC<<
232         "nclTRD="<<nclTRD<<
233         //
234         "pixel="<<hasPIXEL<<
235         "trd="<<hasTRD<<
236         "tof="<<hasTOF<<
237         "acorde="<<hasACORDE<<
238         "ncl="<<xcont[1]<<
239         "dcaR="<<xcont[2]<<
240         "dcaZ="<<xcont[3]<<
241         "alpha="<<xcont[4]<<
242         "theta="<<xcont[5]<<
243         "pt="<<xcont[6]<<
244         "dEdx="<<xcont[7]<<     
245         "mpt="<<mpt<<
246         "\n";
247     }
248   }
249 }
250
251 THnSparse *AliTPCcalibTrigger::MakeHisto(const char* trigger){
252   //
253   // Make event/track histograms
254   // trigger histo 
255   //
256   //                 ntracks  nclMax dcaR dcaZ alpha   theta pt dEdx ev
257   Int_t    bins[9] = {50,     40,    20,  20,  18,     25,   25, 25, 2 };
258   //Int_t    bins[9] = {50*   20*    25*  25*  18*     25*   25* 25 };
259   Double_t xmin[9] = {0.,     0,      0,  -250,   -3.14, -1.5,  0, 0, -1.};
260   Double_t xmax[9] = {50,     160,  150,  250, 3.14,  1.5,   100, 100, 1.};
261   TString  axisName[9]={
262     "ntracks",
263     "ncl",
264     "dcaR",
265     "dcaZ",
266     "alpha",
267     "theta",
268     "pt",
269     "dEdx",
270     "ev"
271   };
272   TString  axisTitle[9]={
273     "Number of tracks",
274     "N_{cl}",
275     "dca_{R} (cm)",
276     "dca_{z} (cm)",
277     "alpha (mrad)",
278     "theta",
279     "p_{t} (GeV/c)",
280     "dEdx (a.u.)",
281     "ev"
282   };
283
284   
285   THnSparse *sparse = new THnSparseF(Form("his_%s",trigger), Form("his_%s",trigger), 9, bins, xmin, xmax);
286   for (Int_t iaxis=0; iaxis<9; iaxis++){
287     sparse->GetAxis(iaxis)->SetName(axisName[iaxis]);
288     sparse->GetAxis(iaxis)->SetTitle(axisTitle[iaxis]);
289   }
290   return sparse;
291 }
292
293 THnSparse * AliTPCcalibTrigger::GetHisto(const char *trigger) { 
294   //
295   // return histogram for given class
296   if (!fHisMap) fHisMap=new TMap;
297   return (THnSparse*) fHisMap->GetValue(trigger);
298 }
299
300 void   AliTPCcalibTrigger::AddHisto(const char *trigger, THnSparse *his) { 
301   if (!GetHisto(trigger)) {
302     TObjString *pstr = new TObjString(trigger);
303     fHisMap->Add(pstr,his);
304   }
305 }
306
307 TTree * AliTPCcalibTrigger::MakeTree(const char * fname){
308   //
309   //
310   //
311   TTreeSRedirector * sred = new TTreeSRedirector(fname);
312   TTreeStream &pcstream = (*sred)<<"Trigger";
313   //
314   //
315   TIterator* iterator = fHisMap->MakeIterator();
316   TObject * object=0;
317   //
318   while((object=iterator->Next())){
319     MakeTree(pcstream, object->GetName());
320   }
321   delete sred;
322   TFile *f = new TFile(fname);
323   TTree *tree = (TTree*)f->Get("Trigger");
324   return tree;
325 }
326
327
328 void AliTPCcalibTrigger::MakeTree(TTreeStream &pcstream, const char *tname){
329   //
330   //  TTreeSRedirector * sred = new TTreeSRedirector("trigger.root");
331   //  TTreeStream &pcstream = (*sred)<<"Trigger";
332   //
333   //AliTPCcalibTrigger *calibTrigger = this;
334   Double_t value;
335   THnSparse * his = GetHisto(tname);
336   if (!his) return;
337   //
338   Int_t *bins = new Int_t[100];
339   Int_t ndim = his->GetNdimensions();
340   Double_t position[10];
341   //
342   TObjString str(tname);
343   Bool_t isAll  = str.String().Contains("all");
344   Bool_t hasPIXEL=HasPIXEL(&str);
345   Bool_t hasTRD=HasTRD(&str);
346   Bool_t hasTOF=HasTOF(&str);
347   Bool_t hasACORDE=HasACORDE(&str);
348   for (Long64_t i = 0; i < his->GetNbins(); ++i) {
349     value = his->GetBinContent(i, bins);
350     pcstream<<"val="<<value;
351     pcstream<<"tname.="<<&str;
352     //
353     pcstream<<"all="<<isAll;
354     pcstream<<"pixel="<<hasPIXEL;
355     pcstream<<"trd="<<hasTRD;
356     pcstream<<"tof="<<hasTOF;
357     pcstream<<"acorde="<<hasACORDE;
358     //
359     for (Int_t idim = 0; idim < ndim; idim++) {
360       position[idim] = his->GetAxis(idim)->GetBinCenter(bins[idim]);
361       pcstream<<Form("%s=",his->GetAxis(idim)->GetName())<<position[idim];
362     }
363     pcstream<<"\n";
364   }
365 }
366
367
368 Bool_t AliTPCcalibTrigger::HasTOF(TObjString *tname){
369   //
370   Bool_t result = kFALSE;
371   result|=(tname->String().Contains("0OB")>0);
372   result|=(tname->String().Contains("0OC")>0);
373   return result;
374 }
375
376 Bool_t AliTPCcalibTrigger::HasACORDE(TObjString *tname){
377   Bool_t result = kFALSE;
378   result|=(tname->String().Contains("0ASL")>0);
379   result|=(tname->String().Contains("0AMU")>0);
380   result|=(tname->String().Contains("0ASC")>0);
381   return result;
382 }
383
384 Bool_t AliTPCcalibTrigger::HasPIXEL(TObjString *tname){
385   return (tname->String().Contains("0SCO")>0);
386 }
387
388 Bool_t AliTPCcalibTrigger::HasTRD(TObjString *tname){
389   Bool_t result = kFALSE;
390   result|=(tname->String().Contains("TRD")>0);
391   result|=(tname->String().Contains("1H")>0);
392   return result;
393 }
394