]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/macros/AnalysisTrainCAF.C
Completion of previous checkin
[u/mrichter/AliRoot.git] / PWG4 / macros / AnalysisTrainCAF.C
1 //______________________________________________________________________________
2 void AnalysisTrainCAF(Int_t nEvents = 10000, Int_t nOffset = 0, char *ds = "/PWG4/kleinb/LHC08v_jetjet15-50")
3 {
4 // Example of running analysis train in CAF. To run in debug mode:
5 //  - export ROOTSYS=debug  on your local client
6 //  - un-comment gProof->ClearPackages()
7 //  - un-comment lines with debugging info
8
9 // WHY AOD is not a exchange container when running from ESD->AOD
10
11     Bool_t debug         = kTRUE;
12     Bool_t useMC         = kTRUE;
13     Bool_t readTR        = kFALSE;
14     Bool_t bPROOF        = kTRUE;
15
16
17     
18     Int_t iAODanalysis   = 0;
19     Int_t iAODhandler    = 1;
20     Int_t iESDfilter     = 1;  // Only active if iAODanalysis=0
21     Int_t iJETAN         = 1;
22     Int_t iJETANAOD      = 0;
23     Int_t iJETANMC       = 1;
24     Int_t iDIJETAN       = 0;
25     Int_t iPWG4SPECTRUM  = 1;
26     Int_t iPWG4UE        = 0;
27
28     if (iAODanalysis) {
29        useMC = kFALSE;
30        readTR = kFALSE;
31        iESDfilter = 0;
32        iMUONfilter = 0;
33     }    
34     if (iJETANAOD) iESDfilter=1;
35     if (iESDfilter) iAODhandler=1;
36
37     // Dataset from CAF
38     TString dataset(ds);
39     // CKB quick hack for local analysis
40     //    gROOT->LoadMacro("CreateESDChain.C");
41     // TChain *chain = CreateESDChain("jetjet15-50.txt",1000);
42     TChain *chain = 0;
43
44     printf("==================================================================\n");
45     printf("===========    RUNNING ANALYSIS TRAIN IN CAF MODE    =============\n");
46     printf("==================================================================\n");
47     if (iAODanalysis) printf("=  AOD analysis on dataset: %s\n", dataset.Data());
48     else              printf("=  ESD analysis on dataset: %s\n", dataset.Data());
49     if (iESDfilter)   printf("=  ESD filter                                                    =\n");
50     if (iJETAN)       printf("=  Jet analysis                                                  =\n");
51     if (iJETANAOD)    printf("=  Jet analysis from AOD                                        =\n");
52     if (iJETANMC)     printf("=  Jet analysis from Kinematics                                  =\n");
53     if (iDIJETAN)     printf("=  DiJet analysis                                                =\n");
54     if (iPWG4SPECTRUM)printf("=  PWG4 Jet spectrum analysis                                    =\n");
55     if (iPWG4UE)      printf("=  PWG4 UE                                                       =\n");
56     printf("==================================================================\n");
57     if (useMC) printf(":: use MC    TRUE\n");
58     else       printf(":: use MC    FALSE\n");
59     if (readTR) printf(":: read TR   TRUE\n");
60     else        printf(":: read TR   FALSE\n");
61     if (debug) printf(":: debugging TRUE\n");
62     else       printf(":: debugging FALSE\n");
63     
64     // Load common libraries
65     gSystem->Load("libTree.so");
66     gSystem->Load("libGeom.so");
67     gSystem->Load("libVMC.so");
68     gSystem->Load("libPhysics.so");
69
70
71     // Reset user processes if CAF if not responding anymore
72     // TProof::Reset("alicecaf"); 
73     // One may enable a different ROOT version on CAF
74
75     //    const char* proofNode = "alicecaf";
76     const char* proofNode = "localhost";
77
78
79
80
81     // Connect to proof
82     if(bPROOF){
83       TProof::Mgr(proofNode)->ShowROOTVersions();
84       //      TProof::Mgr(proofNode)->SetROOTVersion("v5-21-01-alice_dbg");
85       TProof::Open(proofNode); 
86
87       // Clear packages if changing ROOT version on CAF or local
88       //      gProof->ClearPackages();
89       // Enable proof debugging if needed
90       //    gProof->SetLogLevel(5);
91       // To debug the train in PROOF mode, type in a root session:
92       // root[0] TProof::Mgr("lxb6064")->GetSessionLogs()->Display("*",0,10000);
93       // Common packages
94       // --- Enable the STEERBase Package
95       gProof->UploadPackage("${ALICE_ROOT}/STEERBase.par");
96       gProof->EnablePackage("STEERBase");
97       // --- Enable the ESD Package
98       gProof->UploadPackage("${ALICE_ROOT}/ESD.par");
99       gProof->EnablePackage("ESD");
100       // --- Enable the AOD Package
101       gProof->UploadPackage("${ALICE_ROOT}/AOD.par");
102       gProof->EnablePackage("AOD");
103       // --- Enable the ANALYSIS Package
104       gProof->UploadPackage("${ALICE_ROOT}/ANALYSIS.par");
105       gProof->EnablePackage("ANALYSIS");
106       // --- Enable the ANALYSISalice Package
107       gProof->UploadPackage("${ALICE_ROOT}/ANALYSISalice.par");
108       gProof->EnablePackage("ANALYSISalice");
109
110
111       // --- Enable the JETAN Package
112       if (iJETAN||iJETANMC) {
113         gProof->UploadPackage("${ALICE_ROOT}/JETAN.par");
114         gProof->EnablePackage("JETAN");
115       }   
116       // --- Enable particle correlation analysis
117       if (iPWG4UE||iPWG4SPECTRUM) {
118         gProof->UploadPackage("${ALICE_ROOT}/PWG4JetTasks.par");
119         gProof->EnablePackage("PWG4JetTasks");
120       }   
121
122     }
123     else{
124       gSystem->Load("libSTEERBase");
125       gSystem->Load("libESD");
126       gSystem->Load("libAOD");
127       gSystem->Load("libANALYSIS");
128       gSystem->Load("libANALYSISalice");
129
130
131       // --- Enable the JETAN Package
132       if (iJETAN||iJETANMC||iJETANAOD) gSystem->Load("libJETAN");
133       // --- Enable particle correlation analysis
134       if (iPWG4UE||iPWG4SPECTRUM)gSystem->Load("libPWG4JetTasks");  
135     }
136
137
138     // Make the analysis manager
139     AliAnalysisManager *mgr  = new AliAnalysisManager("Analysis Train", "A test setup for the analysis train");
140     if (iAODanalysis) {
141     // AOD input handler
142        AliAODInputHandler *aodH = new AliAODInputHandler();
143        mgr->SetInputEventHandler(aodH);
144     } else {   
145     // ESD input handler
146        AliESDInputHandler *esdHandler = new AliESDInputHandler();
147        mgr->SetInputEventHandler(esdHandler);
148 //       esdHandler->SetInactiveBranches("FMD CaloCluster");
149     }
150     // Monte Carlo handler
151     if (useMC && !iAODanalysis) {
152        AliMCEventHandler* mcHandler = new AliMCEventHandler();
153        mgr->SetMCtruthEventHandler(mcHandler);
154        mcHandler->SetReadTR(readTR); 
155     }   
156     // Top container for input 
157     AliAnalysisDataContainer *cinput = mgr->GetCommonInputContainer();
158     
159     // This container is managed by the AOD handler
160     AliAnalysisDataContainer *cout_aod = 0;
161     if (iAODhandler) {
162        // AOD output handler
163        AliAODHandler* aodHandler   = new AliAODHandler();
164        aodHandler->SetFillAOD(kFALSE);
165        mgr->SetOutputEventHandler(aodHandler);       
166        aodHandler->SetOutputFileName(Form("AliAODs_pwg4_%07d-%07d.root",nOffset,nOffset+nEvents));
167        cout_aod = mgr->GetCommonOutputContainer();
168        cout_aod->SetSpecialOutput();
169     }   
170
171     // Debugging if needed
172     if (debug) mgr->SetDebugLevel(0);
173 //    AliLog::EnableDebug(kTRUE);
174     AliLog::SetGlobalLogLevel(0);
175
176
177     if (iESDfilter && !iAODanalysis) {
178        // Set of cuts plugged into the ESD filter
179        // 
180        // standard
181        AliESDtrackCuts* esdTrackCutsL = new AliESDtrackCuts("AliESDtrackCuts", "Loose");
182        esdTrackCutsL->SetMinNClustersTPC(50);   
183        esdTrackCutsL->SetMaxChi2PerClusterTPC(3.5);
184        esdTrackCutsL->SetMaxCovDiagonalElements(2,2,0.5,0.5,2);
185        esdTrackCutsL->SetRequireTPCRefit(kTRUE);
186        esdTrackCutsL->SetMinNsigmaToVertex(3);
187        esdTrackCutsL->SetRequireSigmaToVertex(kTRUE);
188        esdTrackCutsL->SetAcceptKingDaughters(kFALSE);
189        //
190        AliAnalysisFilter* trackFilter = new AliAnalysisFilter("trackFilter");
191        trackFilter->AddCuts(esdTrackCutsL);
192        //
193        // ESD filter task putting standard info to output AOD (with cuts)
194        AliAnalysisTaskESDfilter *esdfilter = new AliAnalysisTaskESDfilter("ESD Filter");
195        esdfilter->SetTrackFilter(trackFilter);
196        esdfilter->SetDebugLevel(10);
197        mgr->AddTask(esdfilter);
198        // Connect to data containers
199        mgr->ConnectInput  (esdfilter,  0, cinput  );
200        mgr->ConnectOutput (esdfilter,  0, cout_aod );
201     }   
202     // Jet analysis
203     AliAnalysisDataContainer *c_aodjet = 0;
204     if (iJETAN && !iAODanalysis) {
205        AliAnalysisTaskJets *jetana = new AliAnalysisTaskJets("JetAnalysis");
206        jetana->SetDebugLevel(10);
207        jetana->SetConfigFile("${ALICE_ROOT}/JETAN/ConfigJetAnalysis.C");
208        mgr->AddTask(jetana);
209        // Output histograms list for jet analysis                       
210        AliAnalysisDataContainer *cout_jet = mgr->CreateContainer("jethist", TList::Class(),
211                                                                  AliAnalysisManager::kOutputContainer,Form("jethist_%07d-%07d.root",nOffset,nOffset+nEvents));
212        // Dummy AOD output container for jet analysis (no client yet)
213        c_aodjet = mgr->CreateContainer("cAODjet", TTree::Class(),
214                            AliAnalysisManager::kExchangeContainer);
215        // Connect to data containers
216        mgr->ConnectInput  (jetana,     0, cinput  );
217        mgr->ConnectOutput (jetana,     0, c_aodjet );
218        // mgr->ConnectOutput (jetana,     0, cout_aod );
219        mgr->ConnectOutput (jetana,     1, cout_jet );
220     }   
221     // JETANALYSIS from the AOD
222     if (iJETANAOD) {
223        AliAnalysisTaskJets *jetanaAOD = new AliAnalysisTaskJets("AODJetAnalysis");
224        jetanaAOD->SetDebugLevel(10);
225        jetanaAOD->SetNonStdBranch("jetsAOD");    
226        jetanaAOD->SetConfigFile("${ALICE_ROOT}/JETAN/ConfigJetAnalysisAOD.C");
227        mgr->AddTask(jetanaAOD);
228        // Output histograms list for jet analysis                       
229        AliAnalysisDataContainer *cout_jetAOD = mgr->CreateContainer("jethistAOD", TList::Class(),
230                                                                     AliAnalysisManager::kOutputContainer, Form("jethistAOD_%07d-%07d.root",nOffset,nOffset+nEvents));
231        // Dummy AOD output container for jet analysis (no client yet)
232        c_aodjet0 = mgr->CreateContainer("cAODjet0", TTree::Class(),
233                            AliAnalysisManager::kExchangeContainer);
234        // Connect to data containers
235        mgr->ConnectInput  (jetanaAOD,     0, cout_aod  );
236        mgr->ConnectOutput (jetanaAOD,     0, c_aodjet0 );
237        // mgr->ConnectOutput (jetana,     0, cout_aod );
238        mgr->ConnectOutput (jetanaAOD,     1, cout_jetAOD );
239     }   
240     // Jet analysisMC
241     AliAnalysisDataContainer *c_aodjetMC = 0;
242     if (iJETANMC && useMC) {
243        AliAnalysisTaskJets *jetanaMC = new AliAnalysisTaskJets("JetAnalysisMC");
244        jetanaMC->SetDebugLevel(10);
245        jetanaMC->SetConfigFile("${ALICE_ROOT}/JETAN/ConfigJetAnalysisMC.C");
246        jetanaMC->SetNonStdBranch("jetsMC");
247        mgr->AddTask(jetanaMC);
248        // Output histograms list for jet analysis                       
249        AliAnalysisDataContainer *cout_jetMC = mgr->CreateContainer("jethistMC", TList::Class(),
250                                                                    AliAnalysisManager::kOutputContainer,Form("jethistMC_%07d-%07d.root",nOffset,nOffset+nEvents));
251        // Dummy AOD output container for jet analysis (no client yet)
252        c_aodjetMC = mgr->CreateContainer("cAODjetMC", TTree::Class(),
253                            AliAnalysisManager::kExchangeContainer);
254        // Connect to data containers
255        mgr->ConnectInput  (jetanaMC,     0, cinput  );
256        mgr->ConnectOutput (jetanaMC,     0, c_aodjetMC );
257        // mgr->ConnectOutput (jetanaMC,     0, cout_aod );
258        mgr->ConnectOutput (jetanaMC,     1, cout_jetMC );
259     }   
260     // Dijet analysis
261     if(iDIJETAN){
262       AliAnalysisTaskDiJets *dijetana = new AliAnalysisTaskDiJets("DiJetAnalysis");
263       dijetana->SetDebugLevel(2);
264       
265       mgr->AddTask(dijetana);
266
267       //
268       // Create containers for input/output
269       AliAnalysisDataContainer *c_aod_dijet = mgr->CreateContainer("cAODdijet", TTree::Class(),
270                                                                    AliAnalysisManager::kExchangeContainer);
271       mgr->ConnectInput  (dijetana,  0, cinput  );
272       mgr->ConnectOutput (dijetana,  0, c_aod_dijet);
273     }
274
275
276     if (iPWG4SPECTRUM) {
277       AliAnalysisTaskJetSpectrum* pwg4spec = new  AliAnalysisTaskJetSpectrum("Jet Spectrum");
278       
279       // default parameters use a switch via iPWGSPECTRUM
280       // or a config file
281       pwg4spec->SetAnalysisType(AliAnalysisTaskJetSpectrum::kAnaMC);
282       //      if(iAODanalysis)pwg4spec->SetAODInput(kTRUE);
283       pwg4spec->SetDebugLevel(11); 
284       //      pwg4spec->SetBranchRec("jetsMC"); 
285       //      pwg4spec->SetBranchGen("jetsMC"); 
286       mgr->AddTask(pwg4spec);
287
288       AliAnalysisDataContainer *coutput1_Spec = mgr->CreateContainer("pwg4spec", TList::Class(),AliAnalysisManager::kOutputContainer,Form( "pwg4spec_%07d-%07d.root",nOffset,nOffset+nEvents));
289       // coutput1_Spec->SetSpecialOutput();
290       // Dummy AOD output container for jet analysis (no client yet)
291       c_aodSpec = mgr->CreateContainer("cAODjetSpec", TTree::Class(),
292                                         AliAnalysisManager::kExchangeContainer);
293       mgr->ConnectInput  (pwg4spec,  0, cinput);    
294       // mgr->ConnectInput  (pwg4spec,  0, c_aodjet);    
295       mgr->ConnectOutput (pwg4spec,  0, c_aodSpec );
296       mgr->ConnectOutput (pwg4spec,  1, coutput1_Spec );
297     }   
298
299     
300     // Particle correlation analysis
301     if (iPWG4UE) {
302       AliAnalysisTaskUE* ueana = new  AliAnalysisTaskUE("Underlying Event");
303       
304
305       // default parameters use a switch via iPWGUE
306       // or a config file
307       Int_t anaType =1; 
308       Int_t regType =1;
309       Double_t jetEtaCut=0.2;
310       Double_t trackPtCut=0.5; 
311       Double_t trackEtaCut= 0.9; 
312       Double_t rad=0.7; 
313       Double_t deltaPhiCut = 2.616;
314
315       ueana->SetDebugLevel(10); 
316       ueana->SetPtRangeInHist(25, 0., 250.);
317       ueana->SetAnaTopology(anaType);      
318       ueana->SetRegionType(regType);        
319       ueana->SetJet1EtaCut(jetEtaCut);     
320       ueana->SetTrackPtCut(trackPtCut); 
321       ueana->SetPtSumOrdering(2);
322       ueana->SetConeRadius(rad);   
323       ueana->SetTrackEtaCut(trackEtaCut);
324       ueana->SetJet2DeltaPhiCut(deltaPhiCut);
325       mgr->AddTask(ueana);
326
327
328       AliAnalysisDataContainer *coutput1_UE = mgr->CreateContainer("histosUE", TList::Class(),AliAnalysisManager::kOutputContainer, Form("pwg4UE_%07d-%07d.root",nOffset,nOffset+nEvents));
329
330       mgr->ConnectInput  (ueana,  0, cinput);    
331       //      mgr->ConnectInput  (ueana,  0, c_aodjet);    
332       mgr->ConnectOutput (ueana,     0, coutput1_UE );
333     }   
334
335     // Run the analysis
336     //    
337     if (mgr->InitAnalysis()) {
338        mgr->PrintStatus();
339        if(bPROOF)mgr->StartAnalysis("proof",dataset.Data(), nEvents,nOffset);
340        else mgr->StartAnalysis("local",chain);
341     }   
342 }