]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibTrigger.cxx
Adding helper functions for the trigger mask
[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       "pixel="<<hasPIXEL<<
187       "trd="<<hasTRD<<
188       "tof="<<hasTOF<<
189       "acorde="<<hasACORDE<<
190       "ntracks="<<ntracks<<
191       "\n";
192   }
193   //
194   xcont[8] = -1.;
195   for (Int_t itrack=0; itrack<ntracks; itrack++){
196     AliESDtrack *track=event->GetTrack(itrack);
197     if (!track) continue;
198     Float_t dca[2];
199     Double_t pxyz[3];
200     track->GetDZ(0.,0.,0.,event->GetMagneticField(),dca);
201     Bool_t status = track->GetPxPyPz(pxyz);
202     Double_t alpha = TMath::ATan2(pxyz[1],pxyz[0]);
203     xcont[1]=track->GetTPCNcls();
204     xcont[2]=dca[0];
205     xcont[3]=dca[1];
206     xcont[4]=alpha;
207     xcont[5]=track->GetParameter()[3];
208     xcont[6]=track->Pt();
209     xcont[7]=track->GetTPCsignal();    
210     histoAll->Fill(xcont);
211     histo->Fill(xcont);
212     //
213     //
214     if (cstream) {
215       Double_t mpt = track->GetParameter()[4];
216       (*cstream) << "Track" <<
217         "run="<<fRun<<
218         "time="<<fTime<<
219         "status="<<status<<
220         "ntracks="<<ntracks<<
221         "pixel="<<hasPIXEL<<
222         "trd="<<hasTRD<<
223         "tof="<<hasTOF<<
224         "acorde="<<hasACORDE<<
225         "ncl="<<xcont[1]<<
226         "dcaR="<<xcont[2]<<
227         "dcaZ="<<xcont[3]<<
228         "alpha="<<xcont[4]<<
229         "theta="<<xcont[5]<<
230         "pt="<<xcont[6]<<
231         "dEdx="<<xcont[7]<<     
232         "mpt="<<mpt<<
233         "\n";
234     }
235   }
236 }
237
238 THnSparse *AliTPCcalibTrigger::MakeHisto(const char* trigger){
239   //
240   // Make event/track histograms
241   // trigger histo 
242   //
243   //                 ntracks  nclMax dcaR dcaZ alpha   theta pt dEdx ev
244   Int_t    bins[9] = {50,     40,    20,  20,  18,     25,   25, 25, 2 };
245   //Int_t    bins[9] = {50*   20*    25*  25*  18*     25*   25* 25 };
246   Double_t xmin[9] = {0.,     0,      0,  -250,   -3.14, -1.5,  0, 0, -1.};
247   Double_t xmax[9] = {50,     160,  150,  250, 3.14,  1.5,   100, 100, 1.};
248   TString  axisName[9]={
249     "ntracks",
250     "ncl",
251     "dcaR",
252     "dcaZ",
253     "alpha",
254     "theta",
255     "pt",
256     "dEdx",
257     "ev"
258   };
259   TString  axisTitle[9]={
260     "Number of tracks",
261     "N_{cl}",
262     "dca_{R} (cm)",
263     "dca_{z} (cm)",
264     "alpha (mrad)",
265     "theta",
266     "p_{t} (GeV/c)",
267     "dEdx (a.u.)",
268     "ev"
269   };
270
271   
272   THnSparse *sparse = new THnSparseF(Form("his_%s",trigger), Form("his_%s",trigger), 9, bins, xmin, xmax);
273   for (Int_t iaxis=0; iaxis<9; iaxis++){
274     sparse->GetAxis(iaxis)->SetName(axisName[iaxis]);
275     sparse->GetAxis(iaxis)->SetTitle(axisTitle[iaxis]);
276   }
277   return sparse;
278 }
279
280 THnSparse * AliTPCcalibTrigger::GetHisto(const char *trigger) { 
281   //
282   // return histogram for given class
283   if (!fHisMap) fHisMap=new TMap;
284   return (THnSparse*) fHisMap->GetValue(trigger);
285 }
286
287 void   AliTPCcalibTrigger::AddHisto(const char *trigger, THnSparse *his) { 
288   if (!GetHisto(trigger)) {
289     TObjString *pstr = new TObjString(trigger);
290     fHisMap->Add(pstr,his);
291   }
292 }
293
294 TTree * AliTPCcalibTrigger::MakeTree(const char * fname){
295   //
296   //
297   //
298   TTreeSRedirector * sred = new TTreeSRedirector(fname);
299   TTreeStream &pcstream = (*sred)<<"Trigger";
300   //
301   //
302   TIterator* iterator = fHisMap->MakeIterator();
303   TObject * object=0;
304   //
305   while((object=iterator->Next())){
306     MakeTree(pcstream, object->GetName());
307   }
308   delete sred;
309   TFile *f = new TFile(fname);
310   TTree *tree = (TTree*)f->Get("Trigger");
311   return tree;
312 }
313
314
315 void AliTPCcalibTrigger::MakeTree(TTreeStream &pcstream, const char *tname){
316   //
317   //  TTreeSRedirector * sred = new TTreeSRedirector("trigger.root");
318   //  TTreeStream &pcstream = (*sred)<<"Trigger";
319   //
320   //AliTPCcalibTrigger *calibTrigger = this;
321   Double_t value;
322   THnSparse * his = GetHisto(tname);
323   if (!his) return;
324   //
325   Int_t *bins = new Int_t[100];
326   Int_t ndim = his->GetNdimensions();
327   Double_t position[10];
328   //
329   TObjString str(tname);
330   Bool_t isAll  = str.String().Contains("all");
331   Bool_t hasPIXEL=HasPIXEL(&str);
332   Bool_t hasTRD=HasTRD(&str);
333   Bool_t hasTOF=HasTOF(&str);
334   Bool_t hasACORDE=HasACORDE(&str);
335   for (Long64_t i = 0; i < his->GetNbins(); ++i) {
336     value = his->GetBinContent(i, bins);
337     pcstream<<"val="<<value;
338     pcstream<<"tname.="<<&str;
339     //
340     pcstream<<"all="<<isAll;
341     pcstream<<"pixel="<<hasPIXEL;
342     pcstream<<"trd="<<hasTRD;
343     pcstream<<"tof="<<hasTOF;
344     pcstream<<"acorde="<<hasACORDE;
345     //
346     for (Int_t idim = 0; idim < ndim; idim++) {
347       position[idim] = his->GetAxis(idim)->GetBinCenter(bins[idim]);
348       pcstream<<Form("%s=",his->GetAxis(idim)->GetName())<<position[idim];
349     }
350     pcstream<<"\n";
351   }
352 }
353
354
355 Bool_t AliTPCcalibTrigger::HasTOF(TObjString *tname){
356   //
357   Bool_t result = kFALSE;
358   result|=(tname->String().Contains("0OB")>0);
359   result|=(tname->String().Contains("0OC")>0);
360   return result;
361 }
362
363 Bool_t AliTPCcalibTrigger::HasACORDE(TObjString *tname){
364   Bool_t result = kFALSE;
365   result|=(tname->String().Contains("0ASL")>0);
366   result|=(tname->String().Contains("0AMU")>0);
367   result|=(tname->String().Contains("0ASC")>0);
368   return result;
369 }
370
371 Bool_t AliTPCcalibTrigger::HasPIXEL(TObjString *tname){
372   return (tname->String().Contains("0SCO")>0);
373 }
374
375 Bool_t AliTPCcalibTrigger::HasTRD(TObjString *tname){
376   return (tname->String().Contains("1H")>0);
377 }
378