]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGPP/global/AliGlobalFBFqa.cxx
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGPP / global / AliGlobalFBFqa.cxx
CommitLineData
81b9e9f2 1/**************************************************************************
2 * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16//////////////////////////////////////////////////////
17// Basic QA task to monitor the observables related to
18// flow and balance function analysis. It is intended
19// for both pp and PbPb data.
20//
21// Mainteiners:
22// Carlos Perez (cperez@cern.ch)
23// Alis Rodriguez (alisrm@nikhef.nl)
24//////////////////////////////////////////////////////
25
26#include "TChain.h"
27#include "TList.h"
28#include "TH2D.h"
29
30#include "AliAnalysisTaskSE.h"
31#include "AliAnalysisManager.h"
32
33#include "AliCentrality.h"
34
35#include "AliESDEvent.h"
36#include "AliESDInputHandler.h"
37#include "AliESDVZERO.h"
38#include "AliESDVertex.h"
39#include "AliESDtrack.h"
40#include "AliESDtrackCuts.h"
41
42#include "AliGlobalFBFqa.h"
43
44ClassImp(AliGlobalFBFqa)
45
46// C O N S T R U C T O R =======================================================
47AliGlobalFBFqa::AliGlobalFBFqa() :
48 AliAnalysisTaskSE(), fDebugger(kFALSE), fOutputList(NULL), fEvents(NULL),
49 fPhiRinVZERO(NULL), fPhiEtaTPC50(NULL), fPhiEtaITSSA(NULL),
50 fPsiCenVZERO(NULL), fPsiCenTPC50(NULL), fPsiCenITSSA(NULL) {
51 // default constructor
52 if(fDebugger) printf("AliGlobalFBFqa:: Default Constructor\n");
53}
54// C O N S T R U C T O R =======================================================
55AliGlobalFBFqa::AliGlobalFBFqa(const char *name) :
56 AliAnalysisTaskSE(name), fDebugger(kFALSE), fOutputList(NULL), fEvents(NULL),
57 fPhiRinVZERO(NULL), fPhiEtaTPC50(NULL), fPhiEtaITSSA(NULL),
58 fPsiCenVZERO(NULL), fPsiCenTPC50(NULL), fPsiCenITSSA(NULL) {
59 // named constructor
60 if(fDebugger) printf("AliGlobalFBFqa:: Named Constructor\n");
61 DefineInput( 0,TChain::Class());
62 DefineOutput(1,TList::Class());
63}
64// D E S T R U C T O R =========================================================
65AliGlobalFBFqa::~AliGlobalFBFqa() {
66 // destructor
67 if(fDebugger) printf("AliGlobalFBFqa:: Destructor\n");
68 if(fOutputList) delete fOutputList;
69}
70// U S E R C R E A T E O U T P U T O B J E C T S =========================
71void AliGlobalFBFqa::UserCreateOutputObjects() {
72 // user create output object
73 if(fDebugger) printf("AliGlobalFBFqa:: UserCreateOutputObjects\n");
74 fOutputList = new TList();
75 fOutputList->SetOwner();
76 // number of events
77 TList *tQAEvents = new TList();
78 tQAEvents->SetName("Events");
79 tQAEvents->SetOwner();
80 fEvents = new TH2D("Events","Events;TRK;V0M", 20,0,100, 20,0,100);
81 tQAEvents->Add(fEvents);
82 fOutputList->Add(tQAEvents);
83 TList *tQAPhaseSpace = new TList();
84 tQAPhaseSpace->SetOwner();
85 tQAPhaseSpace->SetName("PhaseSpace");
86 fPhiRinVZERO = new TH2D("PhiRinVZERO","PhiRinVZERO;#phi (rad);Ring Number",60,0,TMath::Pi(),8,0,8);
87 tQAPhaseSpace->Add(fPhiRinVZERO);
88 fPhiEtaTPC50 = new TH2D("PhiEtaTPC50","PhiEtaTPC50;#phi (rad);#eta",60,0,TMath::Pi(),16,-1.6,1.6);
89 tQAPhaseSpace->Add(fPhiEtaTPC50);
90 fPhiEtaITSSA = new TH2D("PhiEtaITSSA","PhiEtaITSSA;#phi (rad);#eta",60,0,TMath::Pi(),20,-2,2);
91 tQAPhaseSpace->Add(fPhiEtaITSSA);
92 fOutputList->Add(tQAPhaseSpace);
93 TList *tQAEventPlane = new TList();
94 tQAEventPlane->SetOwner();
95 tQAEventPlane->SetName("EventPlane");
96 fPsiCenVZERO = new TH2D("PsiCenVZERO","PsiCenVZERO;#phi (rad);Centrality TRK",60,0,TMath::Pi(),20,0,100);
97 tQAEventPlane->Add(fPsiCenVZERO);
98 fPsiCenTPC50 = new TH2D("PsiTPCCen50","PsiCenTPC50;#phi (rad);Centrality V0M",60,0,TMath::Pi(),20,0,100);
99 tQAEventPlane->Add(fPsiCenTPC50);
100 fPsiCenITSSA = new TH2D("PsiCenITSSA","PsiCenITSSA;#phi (rad);Centrality V0M",60,0,TMath::Pi(),20,0,100);
101 tQAEventPlane->Add(fPsiCenITSSA);
102 fOutputList->Add(tQAEventPlane);
103 //
104 PostData(1,fOutputList);
105}
106// U S E R E X E C ===========================================================
107void AliGlobalFBFqa::UserExec(Option_t *) {
108 // user exec
109 if(fDebugger) printf("AliGlobalFBFqa:: UserExec\n");
110 AliESDEvent *myESD = dynamic_cast<AliESDEvent*>(InputEvent());
111 if(!myESD) {
112 if(fDebugger) printf("AliGlobalFBFqa:: UserExec ESDEvent not found\n");
113 return;
114 }
115 AliESDVZERO *myVZero = myESD->GetVZEROData();
116 if(!myVZero) {
117 if(fDebugger) printf("AliGlobalFBFqa:: UserExec VZERO not found\n");
118 return;
119 }
120 AliESDVertex *myVertex = (AliESDVertex*) myESD->GetPrimaryVertex();
121 if(!myVertex) {
122 if(fDebugger) printf("AliGlobalFBFqa:: UserExec Vertex not found\n");
123 return;
124 }
125 if(myVertex->GetNContributors() < 2) {
126 if(fDebugger) printf("AliGlobalFBFqa:: UserExec poor vertex\n");
127 return;
128 }
129 Double_t ccTRK=-1, ccV0M=-1;
130 AliCentrality *myCentrality = myESD->GetCentrality();
131 if(myCentrality) {
132 ccTRK = myCentrality->GetCentralityPercentileUnchecked("TRK");
133 ccV0M = myCentrality->GetCentralityPercentileUnchecked("V0M");
134 } else
135 if(fDebugger) printf("AliGlobalFBFqa:: UserExec no centrality object\n");
136 //============================================================================
137 if(fDebugger)
138 printf("AliGlobalFBFqa:: UserExec cc says %.1f and %.1f\n", ccTRK,ccV0M);
139 fEvents->Fill( ccTRK, ccV0M );
140 Double_t dPhi, dEta;
141 Double_t cosTPC=0, sinTPC=0;
142 Double_t cosITS=0, sinITS=0;
143 Double_t cosVZE=0, sinVZE=0;
144 // Global tracks QA
145 if( (ccV0M>5)&&(ccV0M<90) ) { // cutting out edges of multiplicity
146 Int_t nTracks = myESD->GetNumberOfTracks();
147 AliESDtrack *track;
148 AliESDtrackCuts *myTPC=CutsTPC50Generic();
149 AliESDtrackCuts *myITS=CutsITSSAGeneric();
150 for(Int_t i=0; i!=nTracks; ++i) {
151 track = (AliESDtrack*) myESD->GetTrack( i );
152 dPhi = track->Phi();
153 dEta = track->Eta();
154 if( myTPC->IsSelected(track) ) { // TPC50
155 cosTPC += TMath::Cos(2*dPhi);
156 sinTPC += TMath::Sin(2*dPhi);
157 fPhiEtaTPC50->Fill( dPhi, dEta );
158 }
159 if( myITS->IsSelected(track) ) { // ITSSA
160 cosITS += TMath::Cos(2*dPhi);
161 sinITS += TMath::Sin(2*dPhi);
162 fPhiEtaITSSA->Fill( dPhi, dEta );
163 }
164 }
165 double dPsiTPC = 0.5*TMath::ATan2( sinTPC, cosTPC ) + TMath::PiOver2();
166 double dPsiITS = 0.5*TMath::ATan2( sinITS, cosITS ) + TMath::PiOver2();
167 fPsiCenTPC50->Fill( dPsiTPC, ccV0M );
168 fPsiCenITSSA->Fill( dPsiITS, ccV0M );
169 }
170 // VZERO QA
171 if( (ccTRK>5)&&(ccTRK<90) ) { // cutting out edges of multiplicity
172 for(int i=0;i!=64;++i) {
173 dPhi = TMath::Pi()*(i%8)/4;
174 dEta = i/8;
175 cosVZE += myVZero->GetMultiplicity(i)*TMath::Cos(2*dPhi);
176 sinVZE += myVZero->GetMultiplicity(i)*TMath::Sin(2*dPhi);
177 fPhiRinVZERO->Fill( dPhi, dEta );
178 }
179 double dPsiVZE = 0.5*TMath::ATan2( sinVZE, cosVZE ) + TMath::PiOver2();
180 fPsiCenVZERO->Fill( dPsiVZE, ccTRK );
181 }
182 PostData(1,fOutputList);
183 return;
184}
185
186AliESDtrackCuts* AliGlobalFBFqa::CutsITSSAGeneric() {
187 AliESDtrackCuts* esdTrackCuts = new AliESDtrackCuts;
188 esdTrackCuts->SetRequireITSStandAlone(kTRUE);
189 esdTrackCuts->SetRequireITSPureStandAlone(kFALSE);
190 esdTrackCuts->SetRequireITSRefit(kTRUE);
191 esdTrackCuts->SetMinNClustersITS(4);
192 esdTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD, AliESDtrackCuts::kAny);
193 esdTrackCuts->SetMaxChi2PerClusterITS(2.5);
194 // 7*(0.0033+0.0045/pt^1.3)
195 esdTrackCuts->SetMaxDCAToVertexXYPtDep("0.0231+0.0315/pt^1.3");
196 esdTrackCuts->SetRequireITSPid(kTRUE);
197 esdTrackCuts->SetMaxNOfMissingITSPoints(1);
198 return esdTrackCuts;
199}
200
201AliESDtrackCuts* AliGlobalFBFqa::CutsTPC50Generic() {
202 AliESDtrackCuts* esdTrackCuts = new AliESDtrackCuts;
203 esdTrackCuts->SetMinNClustersTPC(50);
204 esdTrackCuts->SetMaxChi2PerClusterTPC(4);
205 esdTrackCuts->SetAcceptKinkDaughters(kFALSE);
206 esdTrackCuts->SetMaxDCAToVertexZ(3.2);
207 esdTrackCuts->SetMaxDCAToVertexXY(2.4);
208 esdTrackCuts->SetDCAToVertex2D(kTRUE);
209 return esdTrackCuts;
210}
211