]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/JetTasks/AliPhiCorrelationsQATask.cxx
update to trigger correlation study
[u/mrichter/AliRoot.git] / PWG4 / JetTasks / AliPhiCorrelationsQATask.cxx
1 /* $Id: AliPhiCorrelationsQATask.cxx 47416 2011-02-15 13:55:09Z jgrosseo $ */
2
3 #include "AliPhiCorrelationsQATask.h"
4
5 #include <TH2F.h>
6 #include <TParticle.h>
7 #include <TFile.h>
8
9 #include <AliLog.h>
10 #include <AliESDVertex.h>
11 #include <AliESDEvent.h>
12 #include <AliMCEvent.h>
13 #include <AliStack.h>
14 #include <AliHeader.h>
15 #include <AliAnalysisManager.h>
16 #include <AliCentrality.h>
17
18 #include "AliESDtrackCuts.h"
19
20 ClassImp(AliPhiCorrelationsQATask)
21
22 AliPhiCorrelationsQATask::AliPhiCorrelationsQATask(const char* opt) :
23   AliAnalysisTaskSE("AliPhiCorrelationsQATask"),
24   fOutput(0),
25   fOption(opt),
26   fEsdTrackCuts(0),
27   fEsdTrackCuts2(0),
28   fCheckITS(0),
29   fGlobalTracks(0),
30   fCentralityCorrelation(0),
31   fDCAPrimaries(0),
32   fDCASecondaries(0),
33   fUseUncheckedCentrality(kFALSE)
34 {
35   //
36   // Constructor. Initialization of pointers
37   //
38
39   // Define input and output slots here
40   DefineOutput(1, TList::Class());
41   
42   AliLog::SetClassDebugLevel("AliPhiCorrelationsQATask", AliLog::kWarning);
43 }
44
45 AliPhiCorrelationsQATask::~AliPhiCorrelationsQATask()
46 {
47   //
48   // Destructor
49   //
50
51   // histograms are in the output list and deleted when the output
52   // list is deleted by the TSelector dtor
53
54   if (fOutput && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) {
55     delete fOutput;
56     fOutput = 0;
57   }
58 }
59
60 void AliPhiCorrelationsQATask::UserCreateOutputObjects()
61 {
62   // create result objects and add to output list
63
64   Printf("AliPhiCorrelationsQATask::CreateOutputObjects");
65   //AliLog::SetClassDebugLevel("AliPhysicsSelection", AliLog::kDebug);
66
67   fOutput = new TList;
68   fOutput->SetOwner();
69
70   fCentralityCorrelation = new TH2F("fCentralityCorrelation", ";v0 centr;spd centr", 100, 0, 100.001, 100, 0, 100.001);
71   fOutput->Add(fCentralityCorrelation);
72
73   fEsdTrackCuts = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
74   fEsdTrackCuts->SetMinNClustersTPC(70);
75   // disable DCA cut
76   fEsdTrackCuts->SetMaxDCAToVertexZ();
77   fEsdTrackCuts->SetMaxDCAToVertexXY();
78   fEsdTrackCuts->SetDCAToVertex2D(kFALSE);
79   fEsdTrackCuts->SetName("cuts_quality_only");
80   fEsdTrackCuts->DefineHistograms();
81   fOutput->Add(fEsdTrackCuts);
82   
83   fEsdTrackCuts2 = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
84   fEsdTrackCuts2->SetMinNClustersTPC(70);
85   fEsdTrackCuts2->SetName("cuts_quality_dca");
86   fEsdTrackCuts2->DefineHistograms();
87   fOutput->Add(fEsdTrackCuts2);
88   
89   fGlobalTracks = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010();
90   fGlobalTracks->SetName("global_cuts");
91   fGlobalTracks->DefineHistograms();
92   fOutput->Add(fGlobalTracks);
93         
94   fCheckITS = new AliESDtrackCuts;
95   fCheckITS->SetRequireITSRefit(kTRUE);
96   fCheckITS->SetClusterRequirementITS(AliESDtrackCuts::kSPD,
97                                       AliESDtrackCuts::kAny);
98   fCheckITS->SetDCAToVertex2D(kFALSE);
99   fCheckITS->SetRequireSigmaToVertex(kFALSE);
100   fCheckITS->SetName("check_its");
101   fCheckITS->DefineHistograms();
102   fOutput->Add(fCheckITS);
103   
104   fDCAPrimaries = new TH2F("fDCAPrimaries", ";dca_xy;dca_z", 1000, -5, 5, 1000, -5, 5);
105   fDCASecondaries = (TH2F*) fDCAPrimaries->Clone("fDCASecondaries");
106   
107   fOutput->Add(fDCAPrimaries);
108   fOutput->Add(fDCASecondaries);
109 }
110
111 void AliPhiCorrelationsQATask::UserExec(Option_t*)
112 {
113   // process the event
114   
115   //Printf("AliPhiCorrelationsQATask::UserExec");
116
117   AliESDEvent* esd = (AliESDEvent*) fInputEvent;
118   if (!esd)
119     AliFatal("No input event");
120     
121   // post the data already here
122   PostData(1, fOutput);
123
124   // vertex cut
125   const AliESDVertex* vtxESD = esd->GetPrimaryVertex();
126   const AliESDVertex *vtxSPD = esd->GetPrimaryVertexSPD();
127   if (!vtxESD)
128     return;
129   if (!vtxSPD)
130     return;
131     
132   Double_t vtx[3];
133   vtxESD->GetXYZ(vtx);
134   if (TMath::Abs(vtx[2]) > 10)
135     return;
136
137   AliCentrality *centralityObj = esd->GetCentrality();
138   if (!centralityObj)
139     AliFatal("Centrality object is 0");
140     
141   Float_t v0Centrality = -1;
142   Float_t spdCentrality = -1; 
143     
144   if (fUseUncheckedCentrality)
145   {
146     v0Centrality = centralityObj->GetCentralityPercentileUnchecked("V0M");
147     spdCentrality = centralityObj->GetCentralityPercentileUnchecked("CL1");
148   }
149   else
150   {
151     v0Centrality = centralityObj->GetCentralityPercentile("V0M");
152     spdCentrality = centralityObj->GetCentralityPercentile("CL1");
153   }
154   
155   fCentralityCorrelation->Fill(v0Centrality, spdCentrality);
156       
157   AliStack* stack = 0;
158   if (fMCEvent)
159   {
160     stack = fMCEvent->Stack();
161     if (!stack)
162       AliFatal("Stack is 0");
163   }
164   
165   Int_t count = 0;
166     
167   for (Int_t i=0; i<esd->GetNumberOfTracks(); i++)
168   {
169     AliESDtrack* esdTrack = esd->GetTrack(i); 
170     
171     if (0)
172     {
173       if (!fGlobalTracks->AcceptTrack(esdTrack))
174         continue;
175     }
176     
177     if (!fEsdTrackCuts->AcceptTrack(esdTrack))
178       continue;
179       
180     count++;
181     
182     //if (v0Centrality < 0)
183     //  continue;
184       
185     Float_t b[2];
186     Float_t bCov[3];
187     esdTrack->GetImpactParameters(b, bCov);
188     
189     Bool_t primary = kTRUE;
190     if (stack)
191       primary = stack->IsPhysicalPrimary(TMath::Abs(esdTrack->GetLabel()));
192     
193     if (primary)
194       fDCAPrimaries->Fill(b[0], b[1]);
195     else
196       fDCASecondaries->Fill(b[0], b[1]);
197       
198     // fill histograms of second object
199     if (!fEsdTrackCuts2->AcceptTrack(esdTrack))
200       continue;
201       
202     if (fCheckITS->AcceptTrack(esdTrack))
203     {
204       if (!fGlobalTracks->AcceptTrack(esdTrack))
205         continue;
206     }
207     
208     continue;
209     
210     // create a tpc only track
211     AliESDtrack *track = AliESDtrackCuts::GetTPCOnlyTrack(esd,esdTrack->GetID());
212     if (!track) 
213       continue;
214
215     if (track->Pt()>0.)
216     {
217       // only constrain tracks above threshold
218       AliExternalTrackParam exParam;
219       // take the B-feild from the ESD, no 3D fieldMap available at this point
220       Bool_t relate = false;
221       relate = track->RelateToVertex(vtxSPD,esd->GetMagneticField(),kVeryBig,&exParam);
222       if(!relate){
223         delete track;
224         continue;
225       }
226       track->Set(exParam.GetX(),exParam.GetAlpha(),exParam.GetParameter(),exParam.GetCovariance());
227     }
228     
229     Float_t b2[2];
230     track->GetImpactParameters(b2,bCov);
231     
232     //Printf("%.2f %.2f %.2f %.2f", b[0], b2[0], b[1], b2[1]);
233     
234     delete track;
235   }      
236   
237   //Printf("%.2f: %d %d tracks (out of %d)", v0Centrality, centralityObj->GetQuality(), count, esd->GetNumberOfTracks());
238 }
239
240 void AliPhiCorrelationsQATask::Terminate(Option_t *)
241 {
242   // terminate
243   
244   Printf("In AliPhiCorrelationsQATask::Terminate...");
245   
246   fOutput = dynamic_cast<TList*> (GetOutputData(1));
247   if (!fOutput)
248   {
249     Printf("ERROR: fOutput not available");
250     return;
251   }
252
253   fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (fOutput->FindObject("cuts_quality_only"));
254   if (!fEsdTrackCuts)
255   {
256     Printf("ERROR: fEsdTrackCuts not available");
257     return;
258   }
259   
260   fEsdTrackCuts2 = dynamic_cast<AliESDtrackCuts*> (fOutput->FindObject("cuts_quality_dca"));
261   if (!fEsdTrackCuts2)
262   {
263     Printf("ERROR: fEsdTrackCuts2 not available");
264     return;
265   }
266   
267   TFile* file = TFile::Open("track_cuts.root", "RECREATE");
268   fEsdTrackCuts->SaveHistograms();
269   fEsdTrackCuts2->SaveHistograms();
270   file->Close();
271 }