]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/ITS/AliAnalysisTaskdEdxSSDQA.cxx
Additional histograms
[u/mrichter/AliRoot.git] / PWG1 / ITS / AliAnalysisTaskdEdxSSDQA.cxx
1 #include "TChain.h"
2 #include "TTree.h"
3 #include "TH1F.h"
4 #include "TH2F.h"
5 #include "TFile.h"
6
7 #include "AliAnalysisTaskdEdxSSDQA.h"
8 #include "AliAnalysisManager.h"
9
10 #include "AliESDEvent.h"
11 #include "AliESDInputHandler.h"
12
13
14 #include "AliESDtrack.h"
15 #include "AliExternalTrackParam.h"
16 #include "Riostream.h"
17
18 #include "AliESDfriend.h"
19 #include "AliESDfriendTrack.h"
20 #include "AliTrackPointArray.h"
21 #include "AliCDBManager.h"
22 #include "AliGeomManager.h"
23 #include "TMath.h"
24
25
26 using namespace std;
27 ClassImp(AliAnalysisTaskdEdxSSDQA)
28
29 //________________________________________________________________________
30 AliAnalysisTaskdEdxSSDQA::AliAnalysisTaskdEdxSSDQA(const char *name) 
31 : AliAnalysisTaskSE(name), fHist1(0),fHist2(0),fHist3(0),fListOfHistos(0),fPcut(0.0)
32 {
33   // Constructor
34         // Define input and output slots here
35   // Input slot #1 works with a TChain
36   DefineOutput(1, TList::Class());
37 }
38
39 //________________________________________________________________________
40 void AliAnalysisTaskdEdxSSDQA::UserCreateOutputObjects()
41 {
42   // Create histograms
43   // Called once
44         
45         fListOfHistos = new TList();
46         fHist1 =new TH2F("QAChargeRatio","QAChargeRatio;Module;CR",1698,-0.5,1697.5,80,-1.0,1.0);
47         fListOfHistos->Add(fHist1);
48         fHist2=new TH2F("QACharge","QACharge;Module;Q",1698,-0.5,1697.5,150,0,300);
49         fListOfHistos->Add(fHist2);
50         fHist3=new TH2F("ChargeRatiovCharge","CRvQ;Q;CR",200,0,2000,100,-1.0,1.0);
51         fListOfHistos->Add(fHist3);
52 }
53 //________________________________________________________________________
54 void AliAnalysisTaskdEdxSSDQA::LocalInit() 
55 {
56         
57         Printf("end of LocalInit");
58 }
59
60 //________________________________________________________________________
61
62 void AliAnalysisTaskdEdxSSDQA::UserExec(Option_t *) 
63 {
64     // Main loop
65     // Called for each event
66     AliESDEvent* esd = dynamic_cast<AliESDEvent*> (InputEvent());
67     if(!esd) 
68     {
69         Printf("ERROR: Input ESD Event not available");
70         return;
71     }
72     
73     Int_t nTracks=esd->GetNumberOfTracks();
74     
75     
76
77         if(!ESDfriend())
78         {
79                 Printf("problem with friend");
80                 return;
81         }
82         
83         AliTrackPointArray*  trackar=0x0;
84         Bool_t l5;
85         Bool_t l6;
86         Int_t npoints;
87         AliTrackPoint point;
88         Int_t nlayer=0;
89         Int_t id=0;
90         Float_t chargeratio=0.0;
91         Float_t charge=0.0;
92         for(int itr=0;itr<nTracks;itr++)
93         {
94                 AliESDtrack* track= esd->GetTrack(itr);
95
96                 if(TMath::Abs(track->Eta())>0.9)
97                         continue;
98                 if(track->GetP()>10.0)
99                         continue;       
100                 if (track->IsOn(AliESDtrack::kITSrefit)&&track->IsOn(AliESDtrack::kTPCrefit))
101                 {
102                          l5=track->HasPointOnITSLayer(4);
103                          l6=track->HasPointOnITSLayer(5);
104                         if (!(l5||l6))//only tracks with SSD point
105                                 continue;
106                         Double_t tmpQESD[4]={-1.0,-1.0,-1.0,-1.0};
107                         track->GetITSdEdxSamples(tmpQESD);
108                         trackar=(AliTrackPointArray*)track->GetTrackPointArray(); 
109                         if(!trackar)
110                                 continue;
111                         npoints=trackar->GetNPoints();                  
112                         for(int itnp=0;itnp<npoints;itnp++)             
113                         {                                                                       
114                                 if(trackar->GetPoint(point,itnp))
115                                         nlayer=AliGeomManager::VolUIDToLayer(point.GetVolumeID(), id);//layer number
116                                 else
117                                         continue;               
118                                 if(nlayer==5||nlayer==6)
119                                 {
120                                         if(point.GetCharge()>0.0&&point.IsExtra()==kFALSE)//do not use additional clusters
121                                         {               
122                                                  chargeratio=point.GetChargeRatio();    
123                                                  charge=point.GetCharge();                                       
124                                                 if(nlayer==5&&tmpQESD[2]>0.0)
125                                                 {
126                                                         fHist1->Fill(id,chargeratio);
127                                                         if(track->GetP()>fPcut)
128                                                                 fHist2->Fill(id,tmpQESD[2]);
129                                                 }       
130                                                 if(nlayer==6&&tmpQESD[3]>0.0)
131                                                 {
132                                                         fHist1->Fill(id+748,chargeratio);
133                                                         if(track->GetP()>fPcut)
134                                                                 fHist2->Fill(id+748,tmpQESD[3]);
135                                                 }
136                                                 fHist3->Fill(charge,chargeratio);       
137                                         }       
138                                 }
139                         }       
140                 }
141         }
142         
143         // Post output data.
144         PostData(1,  fListOfHistos);
145 }      
146
147 //________________________________________________________________________
148 void AliAnalysisTaskdEdxSSDQA::Terminate(Option_t *) 
149 {
150
151         Printf("end of Terminate");
152 }
153 //_____________________________________________________________________________
154
155
156