1 /**************************************************************************
\r
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
\r
4 * Author: The ALICE Off-line Project. *
\r
5 * Contributors are mentioned in the code where appropriate. *
\r
7 * Permission to use, copy, modify and distribute this software and its *
\r
8 * documentation strictly for non-commercial purposes is hereby granted *
\r
9 * without fee, provided that the above copyright notice appears in all *
\r
10 * copies and that both the copyright notice and this permission notice *
\r
11 * appear in the supporting documentation. The authors make no claims *
\r
12 * about the suitability of this software for any purpose. It is *
\r
13 * provided "as is" without express or implied warranty. *
\r
14 **************************************************************************/
\r
16 //------------------------------------------------------------------------------
\r
17 // Implementation of the AliPerformanceTask class. It checks reconstruction performance
\r
18 // for the reconstructed vs MC particle tracks under several conditions. For real data
\r
19 // the control QA histograms are filled.
\r
21 // The comparison output objects deriving from AliPerformanceObject
\r
22 // (e.g. AliPerformanceRes, AliPerformanceEff, AliPerformanceDEdx, AliPerformanceDCA ...)
\r
23 // are stored in the output file (details in description of these classes).
\r
25 // Author: J.Otwinowski 01/04/2009
\r
26 // Changes by M.Knichel 15/10/2010
\r
27 //------------------------------------------------------------------------------
\r
34 #include "TCanvas.h"
\r
37 #include "TSystem.h"
\r
39 #include "AliAnalysisTask.h"
\r
40 #include "AliAnalysisManager.h"
\r
41 #include "AliESDEvent.h"
\r
42 #include "AliESDfriend.h"
\r
43 #include "AliMCEvent.h"
\r
44 #include "AliESDInputHandler.h"
\r
45 #include "AliMCEventHandler.h"
\r
46 #include "AliESDVertex.h"
\r
47 #include "AliMagF.h"
\r
48 #include "AliTracker.h"
\r
49 #include "AliGeomManager.h"
\r
50 #include "AliCDBManager.h"
\r
52 #include "AliCentrality.h"
\r
53 #include "AliESDVZERO.h"
\r
54 #include "AliMultiplicity.h"
\r
56 #include "AliMCInfo.h"
\r
57 #include "AliESDRecInfo.h"
\r
58 #include "AliMCInfoCuts.h"
\r
59 #include "AliRecInfoCuts.h"
\r
60 #include "AliPerformanceObject.h"
\r
61 #include "AliTPCPerformanceSummary.h"
\r
62 #include "AliPerformanceTPC.h"
\r
63 #include "AliPerformanceDEdx.h"
\r
64 #include "AliPerformanceMatch.h"
\r
65 #include "AliPerformanceTask.h"
\r
68 using namespace std;
\r
70 ClassImp(AliPerformanceTask)
\r
72 //_____________________________________________________________________________
\r
73 AliPerformanceTask::AliPerformanceTask()
\r
74 : AliAnalysisTaskSE("Performance")
\r
82 , fUseMCInfo(kFALSE)
\r
83 , fUseESDfriend(kFALSE)
\r
85 , fUseTerminate(kTRUE)
\r
88 , fUseCentralityBin(0)
\r
90 // Dummy Constructor
\r
91 // should not be used
\r
94 //_____________________________________________________________________________
\r
95 AliPerformanceTask::AliPerformanceTask(const char *name, const char */*title*/)
\r
96 : AliAnalysisTaskSE(name)
\r
101 , fOutputSummary(0)
\r
104 , fUseMCInfo(kFALSE)
\r
105 , fUseESDfriend(kFALSE)
\r
107 , fUseTerminate(kTRUE)
\r
108 , fUseCentrality(0)
\r
110 , fUseCentralityBin(0)
\r
114 // Define input and output slots here
\r
115 DefineOutput(0, TTree::Class());
\r
116 DefineOutput(1, TList::Class());
\r
118 // create the list for comparison objects
\r
119 fCompList = new TList;
\r
122 //_____________________________________________________________________________
\r
123 AliPerformanceTask::~AliPerformanceTask()
\r
125 if (!AliAnalysisManager::GetAnalysisManager()->IsProofMode()) {
\r
126 if (fOutput) delete fOutput; fOutput = 0;
\r
127 if (fOutputSummary) delete fOutputSummary; fOutputSummary = 0;
\r
128 if (fCompList) delete fCompList; fCompList = 0;
\r
132 //_____________________________________________________________________________
\r
133 Bool_t AliPerformanceTask::AddPerformanceObject(AliPerformanceObject *pObj)
\r
135 // add comparison object to the list
\r
137 Printf("ERROR: Could not add comparison object");
\r
141 // add object to the list
\r
142 fCompList->AddLast(pObj);
\r
147 //_____________________________________________________________________________
\r
148 void AliPerformanceTask::UserCreateOutputObjects()
\r
150 // Create histograms
\r
153 // create output list
\r
154 fOutput = new TList;
\r
155 fOutput->SetOwner();
\r
156 fPitList = fOutput->MakeIterator();
\r
158 // create output list
\r
159 //fOutputSummary = new TTree;
\r
161 // add comparison objects to the output
\r
162 AliPerformanceObject *pObj=0;
\r
164 TIterator *pitCompList = fCompList->MakeIterator();
\r
165 pitCompList->Reset();
\r
166 while(( pObj = (AliPerformanceObject *)pitCompList->Next()) != NULL) {
\r
167 fOutput->Add(pObj);
\r
170 Printf("UserCreateOutputObjects(): Number of output comparison objects: %d \n", count);
\r
172 PostData(1, fOutput);
\r
173 PostData(0, fOutputSummary);
\r
176 //_____________________________________________________________________________
\r
177 void AliPerformanceTask::UserExec(Option_t *)
\r
180 // Called for each event
\r
183 // Decide whether to use HLT or Offline ESD
\r
186 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*>
\r
187 (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
\r
190 printf("ERROR: Could not get ESDInputHandler");
\r
193 fESD = esdH->GetHLTEvent();
\r
196 fESD = (AliESDEvent*) (InputEvent());
\r
200 fESDfriend = static_cast<AliESDfriend*>(fESD->FindListObject("AliESDfriend"));
\r
202 Printf("ERROR: ESD friends not available");
\r
212 Printf("ERROR: ESD event not available");
\r
216 if (fUseMCInfo && !fMC) {
\r
217 Printf("ERROR: MC event not available");
\r
224 Printf("ERROR: ESD friends not available");
\r
228 // Process analysis
\r
229 Bool_t process = kTRUE;
\r
231 // Check for centrality
\r
232 if (fUseCentrality) {
\r
233 if ( CalculateCentralityBin() != fUseCentralityBin ) {
\r
238 // Process comparison
\r
240 AliPerformanceObject *pObj=0;
\r
242 while(( pObj = (AliPerformanceObject *)fPitList->Next()) != NULL) {
\r
243 pObj->Exec(fMC,fESD,fESDfriend,fUseMCInfo,fUseESDfriend);
\r
247 // Post output data.
\r
248 PostData(1, fOutput);
\r
251 //_____________________________________________________________________________
\r
252 void AliPerformanceTask::Terminate(Option_t *)
\r
254 // Called once at the end
\r
256 if ( !fUseTerminate )
\r
259 // check output data
\r
260 fOutputSummary = dynamic_cast<TTree*> (GetOutputData(0));
\r
261 fOutput = dynamic_cast<TList*> (GetOutputData(1));
\r
263 Printf("ERROR: AliPerformanceTask::Terminate(): fOutput data not available ..." );
\r
266 if (fOutputSummary) { delete fOutputSummary; fOutputSummary=0; }
\r
267 AliPerformanceObject* pObj=0;
\r
268 AliPerformanceTPC* pTPC = 0;
\r
269 AliPerformanceDEdx* pDEdx = 0;
\r
270 AliPerformanceMatch* pMatch = 0;
\r
271 AliPerformanceMatch* pPull = 0;
\r
272 AliPerformanceMatch* pConstrain = 0;
\r
273 TIterator* itOut = fOutput->MakeIterator();
\r
275 while(( pObj = dynamic_cast<AliPerformanceObject*>(itOut->Next())) != NULL) {
\r
276 pObj->AnalyseFinal();
\r
277 /* if (! pTPC) { pTPC = dynamic_cast<AliPerformanceTPC*>(pObj); }
\r
278 if (! pDEdx) { pDEdx = dynamic_cast<AliPerformanceDEdx*>(pObj); }
\r
279 if (! pMatch) { pMatch = dynamic_cast<AliPerformanceMatch*>(pObj); }
\r
280 if ((! pPull) && pMatch ) { pPull = dynamic_cast<AliPerformanceMatch*>(pObj);*/
\r
282 if (!strcmp(pObj->GetName(),"AliPerformanceTPC")) { pTPC = dynamic_cast<AliPerformanceTPC*>(pObj); }
\r
283 if (!strcmp(pObj->GetName(),"AliPerformanceDEdxTPCInner")) { pDEdx = dynamic_cast<AliPerformanceDEdx*>(pObj); }
\r
284 if (!strcmp(pObj->GetName(),"AliPerformanceMatchTPCITS")) { pMatch = dynamic_cast<AliPerformanceMatch*>(pObj); }
\r
285 if (!strcmp(pObj->GetName(),"AliPerformanceMatchITSTPC")) { pPull = dynamic_cast<AliPerformanceMatch*>(pObj);}
\r
286 if (!strcmp(pObj->GetName(),"AliPerformanceMatchTPCConstrain")) { pConstrain = dynamic_cast<AliPerformanceMatch*>(pObj);}
\r
291 printf("DO NOT USE OCDB \n");
\r
295 if (! AliCDBManager::Instance()->GetDefaultStorage()) { AliCDBManager::Instance()->SetDefaultStorage("raw://"); }
\r
297 TString tmpFile = gSystem->TempDirectory() + TString("/TPCQASummary.") + uuid.AsString() + TString(".root");
\r
298 AliTPCPerformanceSummary::WriteToFile(pTPC, pDEdx, pMatch, pPull, pConstrain, tmpFile.Data());
\r
299 TChain* chain = new TChain("tpcQA");
\r
301 chain->Add(tmpFile.Data());
\r
302 TTree *tree = chain->CopyTree("1");
\r
303 if (chain) { delete chain; chain=0; }
\r
304 fOutputSummary = tree;
\r
306 // Post output data.
\r
307 PostData(0, fOutputSummary);
\r
311 //_____________________________________________________________________________
\r
312 void AliPerformanceTask::FinishTaskOutput()
\r
314 // called once at the end of each job (on the workernode)
\r
316 // projects THnSparse to TH1,2,3
\r
318 fOutput = dynamic_cast<TList*> (GetOutputData(1));
\r
320 Printf("ERROR: AliPerformanceTask::FinishTaskOutput(): fOutput data not available ..." );
\r
324 AliPerformanceObject* pObj=0;
\r
325 TIterator* itOut = fOutput->MakeIterator();
\r
327 while(( pObj = dynamic_cast<AliPerformanceObject*>(itOut->Next())) != NULL) {
\r
328 pObj->SetRunNumber(fCurrentRunNumber);
\r
332 // Post output data.
\r
333 PostData(1, fOutput);
\r
336 //_____________________________________________________________________________
\r
337 Bool_t AliPerformanceTask::Notify()
\r
339 static Int_t count = 0;
\r
341 Printf("Processing %d. file", count);
\r
346 //________________________________________________________________________
\r
347 Int_t AliPerformanceTask::CalculateCentralityBin(){
\r
348 // Get Centrality bin
\r
350 Int_t centrality = -1;
\r
351 Float_t centralityF = -1;
\r
353 if (fUseCentrality == 0)
\r
356 AliCentrality *esdCentrality = fESD->GetCentrality();
\r
358 // New : 2010-11-18 JMT
\r
359 if ( fUseCentrality == 1 )
\r
360 centralityF = esdCentrality->GetCentralityPercentile("V0M");
\r
361 else if ( fUseCentrality == 2 )
\r
362 centralityF = esdCentrality->GetCentralityPercentile("CL1");
\r
363 else if ( fUseCentrality == 3 )
\r
364 centralityF = esdCentrality->GetCentralityPercentile("TRK");
\r
365 if (centralityF == 0.)
\r
366 centralityF = 100.;
\r
368 if ( centralityF >= 0. && centralityF < 5.) centrality = 0;
\r
369 else if ( centralityF >= 5. && centralityF < 10.) centrality = 5;
\r
370 else if ( centralityF >= 10. && centralityF < 20.) centrality = 10;
\r
371 else if ( centralityF >= 20. && centralityF < 30.) centrality = 20;
\r
372 else if ( centralityF >= 30. && centralityF < 40.) centrality = 30;
\r
373 else if ( centralityF >= 40. && centralityF < 50.) centrality = 40;
\r
374 else if ( centralityF >= 50. && centralityF < 60.) centrality = 50;
\r
375 else if ( centralityF >= 60. && centralityF < 70.) centrality = 60;
\r
376 else if ( centralityF >= 70. && centralityF < 80.) centrality = 70;
\r
377 else if ( centralityF >= 80. && centralityF < 90.) centrality = 80;
\r
378 else if ( centralityF >= 90. && centralityF <=100.) centrality = 90;
\r