]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/JetTasks/AliAnalysisTaskCorrectionsUE.cxx
New analysis devoted to shower shape studies
[u/mrichter/AliRoot.git] / PWG4 / JetTasks / AliAnalysisTaskCorrectionsUE.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 /* $Id:$ */
17
18 ////////////////////////////////////////////////////////////////////////
19 //
20 // Analysis class to Correct Underlying Event studies
21 // (requires input AOD)
22 //
23 // Different analysis are performed according to input.
24 //
25 // Run on MC:
26 // - fraction of diffractive events after different cuts
27 // - tracking efficiency and contamination
28 // 
29 // Run on DATA:
30 // - fraction of events after different cuts
31 // - vertex reconstruction efficiency
32 //
33 // vallero@physi.uni-heidelberg.de
34 // 
35 ////////////////////////////////////////////////////////////////////////
36
37 #include <TROOT.h>
38 #include <TChain.h>
39 #include <TCanvas.h>
40 #include <TFile.h>
41 #include <TList.h>
42 #include <TMath.h>
43 #include <TProfile.h>
44 #include <TTree.h>
45 #include <TVector3.h>
46
47 #include "AliAnalyseUE.h"
48 #include "AliAnalysisHelperJetTasks.h"
49 #include "AliAnalysisManager.h"
50 #include "AliAnalysisTaskCorrectionsUE.h"
51 #include "AliAODHandler.h"
52 #include "AliESDHandler.h"
53 #include "AliESDtrack.h"
54 #include "AliESDEvent.h"
55 #include "AliAODInputHandler.h"
56 #include "AliESDInputHandler.h"
57 #include "AliCFContainer.h"
58 #include "AliCFManager.h"
59 #include "AliGenDPMjetEventHeader.h"
60 #include "AliGenPythiaEventHeader.h"
61 #include "AliHistogramsUE.h"
62 #include "AliInputEventHandler.h"
63 #include "AliLog.h"
64 #include "AliMCEventHandler.h"
65 #include "AliAnalysisHelperJetTasks.h"
66
67 ClassImp( AliAnalysisTaskCorrectionsUE)
68
69 // Define global pointer
70 AliAnalysisTaskCorrectionsUE* AliAnalysisTaskCorrectionsUE::fgTaskCorrectionsUE=NULL;
71
72 //____________________________________________________________________
73 AliAnalysisTaskCorrectionsUE:: AliAnalysisTaskCorrectionsUE(const char* name):
74 AliAnalysisTask(name,""),
75 fAnaUE(0x0),
76 fAOD(0x0), 
77 fESDHandler(0x0),
78 fAODBranch("jets"),
79 fCFManager(0x0),
80 fDebug(0),
81 fHistosUE(0x0),
82 fInputHandler(0x0),
83 fListOfHistos(0x0),  
84 fMcHandler(0x0),
85 fMcEvent(0x0),
86 fBinsPtInHist(30),     
87 fIsNorm2Area(kTRUE),
88 fMaxJetPtInHist(300.), 
89 fMinJetPtInHist(0.),
90 fConstrainDistance(kTRUE),
91 fMinDistance(0.2),
92 fSimulateChJetPt(kFALSE),
93 fUseAliStack(kTRUE),
94 fUseMCParticleBranch(kFALSE),
95 fnTracksVertex(3),  // QA tracks pointing to principal vertex (= 3 default) 
96 fZVertex(5.),
97 fAnaType(1),         
98 fConePosition(1),
99 fConeRadius(0.7),
100 fFilterBit(0xFF),
101 fJetsOnFly(kFALSE),
102 fRegionType(1),
103 fUseChargeHadrons(kFALSE),
104 fUseChPartJet(kFALSE),
105 fUsePositiveCharge(kTRUE),
106 fUseSingleCharge(kFALSE),
107 fOrdering(1),
108 fChJetPtMin(5.0),
109 fJet1EtaCut(0.2),
110 fJet2DeltaPhiCut(2.616),    // 150 degrees
111 fJet2RatioPtCut(0.8),
112 fJet3PtCut(15.),
113 fTrackEtaCut(0.9),
114 fTrackPtCut(0.),
115 fAvgTrials(1)
116 {
117   // Default constructor
118   // Define input and output slots here
119   // Input slot #0 works with a TChain
120   DefineInput(0, TChain::Class());
121   // Output slot #0 writes into a TList container
122   DefineOutput(0, TList::Class());
123   DefineOutput(1, AliCFContainer::Class());
124
125 }
126
127 //____________________________________________________________________
128 AliAnalysisTaskCorrectionsUE:: AliAnalysisTaskCorrectionsUE(const AliAnalysisTaskCorrectionsUE & original):
129 AliAnalysisTask(),
130 fAnaUE(original.fAnaUE),
131 fAOD(original.fAOD),            
132 fESDHandler(original.fESDHandler),            
133 fAODBranch(original.fAODBranch),
134 fCFManager(original.fCFManager),
135 fDebug(original.fDebug),
136 fHistosUE(original.fHistosUE),
137 fInputHandler(original.fInputHandler),
138 fListOfHistos(original.fListOfHistos),  
139 fMcHandler(original.fMcHandler),
140 fMcEvent(original.fMcEvent),
141 fBinsPtInHist(original.fBinsPtInHist),     
142 fIsNorm2Area(original.fIsNorm2Area),
143 fMaxJetPtInHist(original.fMaxJetPtInHist), 
144 fMinJetPtInHist(original.fMinJetPtInHist),
145 fConstrainDistance(original.fConstrainDistance),
146 fMinDistance(original.fMinDistance),
147 fSimulateChJetPt(original.fSimulateChJetPt),
148 fUseAliStack(original.fUseAliStack),
149 fUseMCParticleBranch(original.fUseMCParticleBranch),
150 fnTracksVertex(original.fnTracksVertex),  // QA tracks pointing to principal vertex (= 3 default) 
151 fZVertex(original.fZVertex),
152 fAnaType(original.fAnaType),         
153 fConePosition(original.fConePosition),
154 fConeRadius(original.fConeRadius),
155 fFilterBit(original.fFilterBit),
156 fJetsOnFly(original.fJetsOnFly),
157 fRegionType(original.fRegionType),
158 fUseChargeHadrons(original.fUseChargeHadrons),
159 fUseChPartJet(original.fUseChPartJet),
160 fUsePositiveCharge(original.fUsePositiveCharge),
161 fUseSingleCharge(original.fUseSingleCharge),
162 fOrdering(original.fOrdering),
163 fChJetPtMin(original.fChJetPtMin),
164 fJet1EtaCut(original.fJet1EtaCut),
165 fJet2DeltaPhiCut(original.fJet2DeltaPhiCut),    // 150 degrees
166 fJet2RatioPtCut(original.fJet2RatioPtCut),
167 fJet3PtCut(original.fJet3PtCut),
168 fTrackEtaCut(original.fTrackEtaCut),
169 fTrackPtCut(original.fTrackPtCut),
170 fAvgTrials(original.fAvgTrials)
171 {
172   // Copy constructor
173 }
174
175
176 //______________________________________________________________
177 AliAnalysisTaskCorrectionsUE & AliAnalysisTaskCorrectionsUE::operator = (const AliAnalysisTaskCorrectionsUE & /*source*/)
178 {
179   // assignment operator
180   return *this;
181 }
182
183 //______________________________________________________________
184 Bool_t AliAnalysisTaskCorrectionsUE::Notify()
185 {
186   //
187   // Implemented Notify() to read the cross sections
188   // and number of trials from pyxsec.root
189   // Copy from AliAnalysisTaskJFSystematics
190   fAvgTrials = 1;
191   TTree *tree = AliAnalysisManager::GetAnalysisManager()->GetTree();
192   Float_t xsection = 0;
193   Float_t trials  = 1;
194   if(tree){
195         TFile *curfile = tree->GetCurrentFile();
196         if (!curfile) {
197                 Error("Notify","No current file");
198                 return kFALSE;
199                 }
200         
201         AliAnalysisHelperJetTasks::PythiaInfoFromFile(curfile->GetName(),xsection,trials); 
202         fHistosUE->GetXsec()->Fill("<#sigma>",xsection);
203
204         // construct average trials
205         Float_t nEntries = (Float_t)tree->GetTree()->GetEntries();
206         if(trials>=nEntries && nEntries>0.)fAvgTrials = trials/nEntries;
207         }
208   
209   return kTRUE;
210
211 }
212
213 //____________________________________________________________________
214 void AliAnalysisTaskCorrectionsUE::ConnectInputData(Option_t* /*option*/)
215 {
216   // Connect the input data  
217   
218   // We need AODs with tracks and jets.
219   // Since AODs can either be connected to the InputEventHandler 
220   // or to the OutputEventHandler ( the AOD is created by a previous task in the train )
221   // we need to check where it is and get the pointer to AODEvent in the right way
222   
223   // Delta AODs are also accepted
224   
225  
226   if (fDebug > 1) AliInfo("ConnectInputData() ");
227   
228   //Get the input handler
229   TObject* handler = AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
230
231   if( handler && handler->InheritsFrom("AliAODInputHandler") ) { // input AOD
232         fAOD = ((AliAODInputHandler*)handler)->GetEvent();
233         if(!fJetsOnFly){
234                 if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODInputHandler");
235                 }else{
236                 if (fDebug > 1) AliInfo(" ==== Tracks from AliAODInputHandler / Jets on-the-fly");
237                 }
238         } else {  //output AOD
239         handler = AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler();
240         if( handler && handler->InheritsFrom("AliAODHandler") ) {
241                 fAOD = ((AliAODHandler*)handler)->GetAOD();
242                 if (!fJetsOnFly){
243                         if (fDebug > 1) AliInfo(" ==== Tracks and Jets from AliAODHandler");
244                         } else {
245                         if (fDebug > 1) AliInfo(" ==== Tracks from AliAODHandler / Jets on-the-fly");
246                         }
247                 }else {
248                 AliFatal("I can't get any AOD Event Handler");
249                 return;
250                 }
251         }
252
253   //Get ESD input handler
254   AliInputEventHandler *inputHandler = dynamic_cast<AliInputEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
255    AliESDEvent *esdEvent = (AliESDEvent*)inputHandler->GetEvent();
256    if (!esdEvent && fDebug > 1) {
257         AliInfo("********************** No ESD event: cannot retrive DCA values from AOD !!! ");
258   }else fAnaUE->SetESDEvent(esdEvent);
259
260   //Get MC handler 
261   fMcHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
262  
263   //Initialize AliAnalysisUE class
264   fAnaUE->Initialize(fAnaType, fAOD, fConeRadius, fDebug, fFilterBit, fJet1EtaCut, fJet2DeltaPhiCut, fJet2RatioPtCut, fJet3PtCut, fOrdering, fRegionType, fSimulateChJetPt, fTrackEtaCut, fTrackPtCut, fUseChargeHadrons, fUseChPartJet, fUsePositiveCharge, fUseSingleCharge, fHistosUE);
265  
266 }
267
268 //____________________________________________________________________
269 void  AliAnalysisTaskCorrectionsUE::CreateOutputObjects()
270 {
271   // Create the output container
272   
273   if (fDebug > 1) AliInfo("CreateOutPutData()");
274    
275   // Create pointer to AliAnalysisUE, a class implementing the main analysis algorithms
276   fAnaUE = new AliAnalyseUE();
277   if (!fAnaUE){
278      AliError("UE analysis class not initialized!!!");
279      return;
280   }
281   // Create pointer to AliHistogramsUE, a class handling histogram creation/filling
282   fHistosUE = new AliHistogramsUE();
283   if (!fHistosUE){
284      AliError("UE histograms class not initialized!!!");
285      return;
286   }
287
288   // Create list of output histograms
289   if (fListOfHistos != NULL){
290         delete fListOfHistos;
291         fListOfHistos = NULL;
292         }
293   if (!fListOfHistos){
294         fListOfHistos = new TList();
295         fListOfHistos->SetOwner(kTRUE); 
296         }
297   
298   
299   //Initialize output histograms
300   fHistosUE->CreateHistogramsCorrections(fListOfHistos,fBinsPtInHist, fMinJetPtInHist, fMaxJetPtInHist, fTrackEtaCut);
301   AddSettingsTree();
302
303   //Configure the CF manager
304   if (fCFManager != NULL){
305         delete fCFManager;
306         fCFManager = NULL;
307         }
308   if (!fCFManager){
309         fCFManager = new AliCFManager();
310         }
311   fHistosUE->CreateCorrectionsContainer(fCFManager,fBinsPtInHist, fMinJetPtInHist, fMaxJetPtInHist, fTrackEtaCut,fJet1EtaCut);
312
313   //Post outputs
314   PostData(0,fListOfHistos);
315   PostData(1,fCFManager->GetEventContainer());
316
317 }
318
319 //____________________________________________________________________
320 void  AliAnalysisTaskCorrectionsUE::Exec(Option_t */*option*/)
321 {
322
323   Bool_t flag=kTRUE;
324
325   //Determine the corrections
326   if (fMcHandler){
327         fMcEvent = fMcHandler->MCEvent();
328         if ( fDebug > 3 ) AliInfo( " Processing MC event..." );
329   }else{
330         if ( fDebug > 3 ) AliInfo( " Processing DATA event..." );
331         }
332   
333
334   flag = EvaluateCorrections(); 
335
336   
337   if (flag){ //executed if event passes trigger+vertex selection
338         //Fetch the pythia header info and get the trials
339         Float_t nTrials = 1;
340         if (fMcHandler && fMcEvent) {
341                 AliGenPythiaEventHeader*  pythiaGenHeader = AliAnalysisHelperJetTasks::GetPythiaEventHeader(fMcEvent);
342                 if(pythiaGenHeader) nTrials = pythiaGenHeader->Trials();
343                 }
344         fHistosUE->GetTrials()->Fill("#sum{ntrials}",fAvgTrials);
345         PostData(0,fListOfHistos);
346         PostData(1,fCFManager->GetEventContainer());
347         }
348 }
349
350 //____________________________________________________________________
351 void  AliAnalysisTaskCorrectionsUE::AddSettingsTree()
352 {
353   //Write settings to output list
354   TTree *settingsTree   = new TTree("UECorrectionsSettings","Analysis Settings in UE corrections");
355   settingsTree->Branch("fFilterBit", &fFilterBit,"FilterBit/I");
356   settingsTree->Branch("fConeRadius", &fConeRadius,"Rad/D");
357   settingsTree->Branch("fJet1EtaCut", &fJet1EtaCut, "LeadJetEtaCut/D");
358   settingsTree->Branch("fJet2DeltaPhiCut", &fJet2DeltaPhiCut, "DeltaPhi/D");
359   settingsTree->Branch("fJet2RatioPtCut", &fJet2RatioPtCut, "Jet2Ratio/D");
360   settingsTree->Branch("fJet3PtCut", &fJet3PtCut, "Jet3PtCut/D");
361   settingsTree->Branch("fTrackPtCut", &fTrackPtCut, "TrackPtCut/D");
362   settingsTree->Branch("fTrackEtaCut", &fTrackEtaCut, "TrackEtaCut/D");
363   settingsTree->Branch("fAnaType", &fAnaType, "Ana/I");        
364   settingsTree->Branch("fRegionType", &fRegionType,"Reg/I");
365   settingsTree->Branch("fOrdering", &fOrdering,"OrderMeth/I");
366   settingsTree->Branch("fUseChPartJet", &fUseChPartJet,"UseChPart/O");
367   settingsTree->Branch("fUseChargeHadrons", &fUseChargeHadrons,"UseChHadrons/O");
368   settingsTree->Branch("fUseSingleCharge", &fUseSingleCharge,"UseSingleCh/O");
369   settingsTree->Branch("fUsePositiveCharge", &fUsePositiveCharge,"UsePositiveCh/O");
370   settingsTree->Fill();
371   fListOfHistos->Add(settingsTree);
372
373
374
375 //____________________________________________________________________
376 Bool_t  AliAnalysisTaskCorrectionsUE::EvaluateCorrections()
377 {
378   
379   /////////////////////////////////////////////////////////////////////////
380   // 
381   // EVENT SELECTION:
382   // CF containers are filled to get the number of entries after every cut.
383   // *** Cuts: ***
384   // 0 - triggered
385   // 1 - physics selection
386   // 2 - vertex selection
387   // 3 - event topology
388   // 4 - leading track pT cut
389   // 5 - leading track correctly identified
390   // *** Variables: ***
391   // 0 - leading track pT (reco)
392   // 1 - leading track eta (reco)
393   // 2 - process type:
394   //    1: non diffractive
395   //    2: double diffractive
396   //    4: single diffractive
397   // 3 - leading track pT (true)
398   // 4 - leading track eta (true)
399   // 5 - delta eta (reco-true)
400   // 6 - delta phi (reco-true)
401   // 7 - radius (reco-true)
402   //
403   // TRACK-LEVEL CORRECTIONS:
404   // Fill histograms similar to AliAnalysisTaskUE.
405   //
406   /////////////////////////////////////////////////////////////////////////
407
408   Double_t containerInput[8];// relevant variables (see above)  
409   
410   //PROCESS TYPE (ND,SD,DD)
411   AliAnalysisHelperJetTasks::MCProcessType eventId = AliAnalysisHelperJetTasks::kInvalidProcess;
412   if (fMcHandler && fMcEvent) {
413                 AliGenEventHeader* genHeader = fMcEvent->GenEventHeader();
414                 eventId = AliAnalysisHelperJetTasks::GetPythiaEventProcessType(genHeader,kFALSE);
415                 if (eventId<0){
416                 eventId = AliAnalysisHelperJetTasks::GetDPMjetEventProcessType(genHeader,kFALSE);
417                 }
418                 if (eventId<0 && fDebug>1)AliInfo("No Pythia or Phojet Header retrived!"); 
419   }else if (fDebug>1) AliInfo("No MC handler or Event!"); 
420
421   //Initialize container inputs 
422   for (Int_t i =0; i<8; i++){  
423         containerInput[i]=-999.;
424         }
425   
426   //Assign process type
427   if (eventId == 1 ) containerInput[2]=1.; //Non diffractive
428   if (eventId == 2 ) containerInput[2]=2.; //Double diffractive
429   if (eventId == 4 ) containerInput[2]=4.; //Single diffractive
430   
431   // Execute analysis for current event ******************************
432   
433   // Get jets and order by pT
434   TVector3 jetVect[3];
435   *jetVect = fAnaUE->GetOrderedClusters(fAODBranch, fUseChPartJet, fChJetPtMin );
436   
437   //now define leading track pT and eta
438   containerInput[0]=jetVect[0].Pt();
439   containerInput[1]=jetVect[0].Eta();
440
441   fCFManager->GetEventContainer()->Fill(containerInput,kCFStepTriggered);  //fill CF container 
442   
443   // Physics selection ************************************************
444   fInputHandler = (AliInputEventHandler*)AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler();
445   if (!fAnaUE->TriggerSelection(fInputHandler)) return kFALSE;
446   
447   fCFManager->GetEventContainer()->Fill(containerInput, kCFStepPhysSelect); //fill CF container
448   
449   // Event selection (vertex) *****************************************
450   
451   if(!fAnaUE->VertexSelection(fAOD,fnTracksVertex,fZVertex)) return kFALSE;
452   //if(!fAnaUE->VertexSelectionOld(fAOD)) return kFALSE; // temporary to compare with old task and to have same cuts for MC !!!
453
454   fCFManager->GetEventContainer()->Fill(containerInput, kCFStepVertexSelect); //fill CF container
455    
456    // Select events according to analysis type *************************
457    // (in the leading track analysis it should not happen that there are no "jets")
458    if( jetVect[0].Pt() < 0. ) {
459         if( fDebug > 1 ) AliInfo("\n   Skipping Event, not jet found...");
460         return kTRUE;
461     } else {
462         if (fDebug >1 ) AliInfo(Form("\n   Pt Leading Jet = %6.1f eta=%5.3f ",  jetVect[0].Pt(), jetVect[0].Eta() ));
463         }
464   
465   
466   if ( ! (fAnaUE->AnaTypeSelection(jetVect))) return kTRUE;
467    
468    
469    // For the events selected check the real MC leading particle
470    // (only when running on MC)
471    TVector3 jetVectMC[3];
472    if (fMcEvent){
473         *jetVectMC = fAnaUE->GetLeadingTracksMC(fMcEvent); 
474         if (fAnaUE->AnaTypeSelection(jetVectMC) ){
475                 //now define leading track pT and eta
476                 containerInput[3]=jetVectMC[0].Pt();
477                 containerInput[4]=jetVectMC[0].Eta();
478                 //Check distance between real and reco leading-particle
479                 containerInput[5] = jetVect[0].Eta()-jetVectMC[0].Eta();
480                 containerInput[6] = TMath::Abs(jetVect[0].Phi()-jetVectMC[0].Phi());
481                 if (containerInput[6] >= 2.*TMath::Pi())containerInput[6] -= 2.*TMath::Pi(); 
482                 containerInput[7] = sqrt(TMath::Power(containerInput[5],2.) + TMath::Power(containerInput[6],2.));
483                 } 
484         }  
485    
486    fCFManager->GetEventContainer()->Fill(containerInput, kCFStepAnaTopology); //fill CF container
487    
488    //Cut on the leading track pT 
489    if (!(jetVect[0].Pt() >= 1.)) return kTRUE;
490    // (only when running on MC)   
491    if (fMcEvent){ 
492         if (!(jetVectMC[0].Pt()>= 1.)){
493                 containerInput[3]=containerInput[4]=-999.;
494                 containerInput[5]=containerInput[6]=containerInput[7]=-999.;
495                 }
496         }       
497    
498    fCFManager->GetEventContainer()->Fill(containerInput, kCFStepLtPtCut1); //fill CF container
499
500    // Check if leading track is correctly identified
501    // (only when running on MC)
502    if (fMcEvent){
503         Int_t labelLt = fAnaUE->GetLtLabel();
504         Int_t labelLtMC = fAnaUE->GetLtMCLabel();
505         if (labelLt == labelLtMC){
506                 fCFManager->GetEventContainer()->Fill(containerInput, kCFStepLtCorrect); //fill CF container
507                 } 
508         }
509
510    // For track efficiency and contamination
511    // (only when running on MC)
512    if (fMcEvent){
513         //Run once on reco...
514         fAnaUE->FindMaxMinRegions( jetVect,  fConePosition, 0, 0 );//normal track cut 
515         fAnaUE->FindMaxMinRegions( jetVect,  fConePosition, 0, 1 );//for efficiency: cut on pT and eta are set on MC true  
516         //and once on MC true
517         fAnaUE->FindMaxMinRegions( jetVect,  fConePosition, 1, 0 );
518   }else{
519         // For d0 distribution, runs only on real data
520         fAnaUE->FindMaxMinRegions( jetVect,  fConePosition, 0, 0 );//run on real data 
521         }
522
523
524    return kTRUE;
525 }
526
527 //____________________________________________________________________
528 AliAnalysisTaskCorrectionsUE* AliAnalysisTaskCorrectionsUE::Instance()
529
530   
531   //Create instance
532   if (fgTaskCorrectionsUE) {
533         return fgTaskCorrectionsUE;
534   } else {
535         fgTaskCorrectionsUE = new AliAnalysisTaskCorrectionsUE();
536         return fgTaskCorrectionsUE;
537         }
538 }
539
540 //____________________________________________________________________
541 void  AliAnalysisTaskCorrectionsUE::Terminate(Option_t */*option*/)
542 {
543   // Terminate analysis
544   
545   if (fDebug >1) AliAnalysisHelperJetTasks::PrintDirectorySize("PWG4_JetTasksOutput.root");
546   
547   if (!gROOT->IsBatch()){
548         fListOfHistos = dynamic_cast<TList*> (GetOutputData(0));
549         if (!fListOfHistos){
550                 AliError("Histogram List is not available");
551                 return;
552                 }
553
554   } else {
555         AliInfo(" Batch mode, not histograms will be shown...");
556   }
557
558   
559 }
560
561 //____________________________________________________________________
562 void  AliAnalysisTaskCorrectionsUE::WriteSettings()
563
564
565 }