]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/vertexingHF/RunAnalysisAODVertexingHF.C
Added merging via jdl (Renu)
[u/mrichter/AliRoot.git] / PWG3 / vertexingHF / RunAnalysisAODVertexingHF.C
1 class AliAnalysisGrid;
2 class AliAnalysisAlien;
3
4 void RunAnalysisAODVertexingHF()
5 {
6   //
7   // Test macro for AliAnalysisTaskSE's for heavy-flavour candidates
8   // It has the structure of a Analysis Train:
9   // - in this macro, change things related to running mode
10   //   and input preparation 
11   // - add your task using a AddTaskXXX macro 
12   //
13   // A.Dainese, andrea.dainese@lnl.infn.it
14   // "grid" mode added by R.Bala, bala@to.infn.it
15   //
16
17
18   gSystem->SetIncludePath("-I. -I$ROOTSYS/include -I$ALICE_ROOT -I$ALICE_ROOT/include -I$ALICE_ROOT/ITS -I$ALICE_ROOT/TPC -I$ALICE_ROOT/CONTAINERS -I$ALICE_ROOT/STEER -I$ALICE_ROOT/TRD -I$ALICE_ROOT/macros -I$ALICE_ROOT/ANALYSIS -I$ALICE_ROOT/PWG3 -I$ALICE_ROOT/PWG3/base -I$ALICE_ROOT/PWG3/vertexingHF -g"); 
19   //
20   TString trainName = "D2H";
21   TString analysisMode = "grid"; // "local", "grid", or "proof"
22   TString inputMode    = "list"; // "list", "xml", or "dataset"
23   Long64_t nentries=123567890,firstentry=0;
24   Bool_t useParFiles=kFALSE;
25   Bool_t useAlienPlugin=kTRUE;
26   TString pluginmode="full";
27   Bool_t saveProofToAlien=kFALSE;
28   TString proofOutdir = "";
29   TString loadMacroPath="$ALICE_ROOT/PWG3/vertexingHF/macros/";
30   //TString loadMacroPath="./"; // this is normally needed for CAF
31   //
32
33   if(analysisMode=="grid") {
34     // Connect to AliEn
35     TGrid::Connect("alien://");
36   } else if(analysisMode=="proof") {
37     // Connect to the PROOF cluster
38     if(inputMode!="dataset") {printf("Input mode must be dataset, for proof analysis\n"); return;}
39     gEnv->SetValue("XSec.GSI.DelegProxy","2");
40     TProof::Open("alicecaf");
41     //TProof::Reset("alicecaf");
42     if(saveProofToAlien) {
43       TGrid::Connect("alien://");
44       if(gGrid) {
45         TString homedir = gGrid->GetHomeDirectory();
46         TString workdir = homedir + trainName;
47         if(!gGrid->Cd(workdir)) {
48           gGrid->Cd(homedir);
49           if(gGrid->Mkdir(workdir)) {
50             gGrid->Cd(trainName);
51             ::Info("VertexingTrain::Connect()", "Directory %s created", gGrid->Pwd());
52           }
53         }          
54         gGrid->Mkdir("proof_output");
55         gGrid->Cd("proof_output");
56         proofOutdir = Form("alien://%s", gGrid->Pwd());
57       } 
58     }
59   }
60
61
62   // AliRoot libraries
63   if(analysisMode=="local" || analysisMode=="grid") {
64     TString loadLibraries="LoadLibraries.C"; loadLibraries.Prepend(loadMacroPath.Data());
65     gROOT->LoadMacro(loadLibraries.Data());
66     LoadLibraries(useParFiles);
67   } else if (analysisMode=="proof") {
68     gSystem->Load("libTree.so");
69     gSystem->Load("libGeom.so");
70     gSystem->Load("libPhysics.so");
71     gSystem->Load("libVMC.so");    
72     gSystem->Load("libMinuit.so");    
73     // Enable the needed packages
74     //gProof->ClearPackages();
75     TString parDir="/afs/cern.ch/user/d/dainesea/code/";
76     TString parFile;
77     if(!useParFiles) {
78       gProof->UploadPackage("AF-v4-17");
79       gProof->EnablePackage("AF-v4-17");
80       // --- Enable the PWG3vertexingHF Package
81       parFile="PWG3vertexingHF.par"; parFile.Prepend(parDir.Data());
82       gProof->UploadPackage(parFile.Data());
83       gProof->EnablePackage("PWG3vertexingHF");
84     } else {
85       // --- Enable the STEERBase Package
86       parFile="STEERBase.par"; parFile.Prepend(parDir.Data());
87       gProof->UploadPackage(parFile.Data());
88       gProof->EnablePackage("STEERBase");
89       // --- Enable the ESD Package
90       parFile="ESD.par"; parFile.Prepend(parDir.Data());
91       gProof->UploadPackage(parFile.Data());
92       gProof->EnablePackage("ESD");
93       // --- Enable the AOD Package
94       parFile="AOD.par"; parFile.Prepend(parDir.Data());
95       gProof->UploadPackage(parFile.Data());
96       gProof->EnablePackage("AOD");
97       // --- Enable the ANALYSIS Package
98       parFile="ANALYSIS.par"; parFile.Prepend(parDir.Data());
99       gProof->UploadPackage(parFile.Data());
100       gProof->EnablePackage("ANALYSIS");
101       // --- Enable the ANALYSISalice Package
102       parFile="ANALYSISalice.par"; parFile.Prepend(parDir.Data());
103       gProof->UploadPackage(parFile.Data());
104       gProof->EnablePackage("ANALYSISalice");
105       // --- Enable the CORRFW Package
106       parFile="CORRFW.par"; parFile.Prepend(parDir.Data());
107       gProof->UploadPackage(parFile.Data());
108       gProof->EnablePackage("CORRFW");
109       // --- Enable the PWG3base Package
110       parFile="PWG3base.par"; parFile.Prepend(parDir.Data());
111       gProof->UploadPackage(parFile.Data());
112       gProof->EnablePackage("PWG3base");
113       // --- Enable the PWG3vertexingHF Package
114       parFile="PWG3vertexingHF.par"; parFile.Prepend(parDir.Data());
115       gProof->UploadPackage(parFile.Data());
116       gProof->EnablePackage("PWG3vertexingHF");
117       // --- Enable the PWG3muon Package
118       parFile="PWG3muon.par"; parFile.Prepend(parDir.Data());
119       gProof->UploadPackage(parFile.Data());
120       gProof->EnablePackage("PWG3muon");
121     }
122     gProof->ShowEnabledPackages(); // show a list of enabled packages
123   }
124
125
126   // Create Alien plugin, if requested
127   if(useAlienPlugin) {  
128     if(analysisMode!="grid") {printf("Analysis mode must be grid, to use alien plugin\n"); return;}
129     AliAnalysisGrid *alienHandler = CreateAlienHandler(pluginmode,useParFiles);  
130     if(!alienHandler) return;
131   }
132
133
134   //-------------------------------------------------------------------
135   // Prepare input
136   TChain *chainAOD = 0;
137   TString dataset; // for proof
138
139   if(!useAlienPlugin) {
140     TString makeAODInputChain="../MakeAODInputChain.C"; makeAODInputChain.Prepend(loadMacroPath.Data());
141     if(inputMode=="list") {
142       // Local files
143       gROOT->LoadMacro(makeAODInputChain.Data());
144       chainAOD = MakeAODInputChain();// with this it reads ./AliAOD.root and ./AliAOD.VertexingHF.root
145       //chainAOD = MakeAODInputChain("alien:///alice/cern.ch/user/r/rbala/newtrain/out_lhc08x/180100/",1,1);
146       printf("ENTRIES %d\n",chainAOD->GetEntries());
147     } else if(inputMode=="xml") {
148       // xml
149       gROOT->LoadMacro(makeAODInputChain.Data());
150       chainAOD = MakeAODInputChain("collection_aod.xml","collection_aodHF.xml");
151     } else if(inputMode=="dataset") {
152       // CAF dataset
153       //gProof->ShowDataSets();
154       dataset="/ITS/dainesea/AODVertexingHF_LHC08x_180100";
155     }
156   }
157
158   // Create the analysis manager
159   AliAnalysisManager *mgr  = new AliAnalysisManager("My Manager","My Manager");
160   mgr->SetDebugLevel(10);
161   // Connect plug-in to the analysis manager
162   if(useAlienPlugin) mgr->SetGridHandler(alienHandler);
163
164   // Input
165   AliAODInputHandler *inputHandler = new AliAODInputHandler();
166   if(analysisMode=="proof" ) {
167     inputHandler->AddFriend("./AliAOD.VertexingHF.root");
168     //inputHandler->AddFriend("deltas/AliAOD.VertexingHF.root");
169     if(saveProofToAlien) mgr->SetSpecialOutputLocation(proofOutdir);
170   }
171   mgr->SetInputEventHandler(inputHandler);
172   //-------------------------------------------------------------------
173
174   
175   //-------------------------------------------------------------------
176   // Analysis tasks (wagons of the train)   
177   //
178   TString taskName;
179   
180   ////// ADD THE FULL D2H TRAIN
181   /*taskName="../AddD2HTrain.C"; taskName.Prepend(loadMacroPath.Data());
182   gROOT->LoadMacro(taskName.Data());
183   Bool_t readMC=kFALSE;
184   AddD2HTrain(readMC);//,1,0,0,0,0,0,0,0,0,0,0);*/
185   
186   ////// OR ADD INDIVIDUAL TASKS
187   
188   
189   /*  taskName="AddTaskCompareHF.C"; taskName.Prepend(loadMacroPath.Data());
190     gROOT->LoadMacro(taskName.Data());
191     AliAnalysisTaskSECompareHF *cmpTask = AddTaskCompareHF();
192   */ 
193     taskName="AddTaskD0Mass.C"; taskName.Prepend(loadMacroPath.Data());
194     gROOT->LoadMacro(taskName.Data());
195     AliAnalysisTaskSED0Mass *d0massTask = AddTaskD0Mass();
196     AliAnalysisTaskSED0Mass *d0massLikeSignTask = AddTaskD0Mass(1); 
197     /*
198     taskName="AddTaskDplus.C"; taskName.Prepend(loadMacroPath.Data());
199     gROOT->LoadMacro(taskName.Data());
200     AliAnalysisTaskSEDplus *dplusTask = AddTaskDplus();
201   
202     taskName="AddTaskDs.C"; taskName.Prepend(loadMacroPath.Data());
203     gROOT->LoadMacro(taskName.Data());
204     AliAnalysisTaskSEDs *dsTask = AddTaskDs();
205     
206     //taskName="AddTaskSelectHF.C"; taskName.Prepend(loadMacroPath.Data());
207     //gROOT->LoadMacro(taskName.Data());
208     //AliAnalysisTaskSESelectHF *seleTask = AddTaskSelectHF();
209     
210     taskName="AddTaskBkgLikeSignD0.C"; taskName.Prepend(loadMacroPath.Data());
211     gROOT->LoadMacro(taskName.Data());
212     AliAnalysisTaskSEBkgLikeSignD0 *lsD0Task = AddTaskBkgLikeSignD0();
213     
214     taskName="AddTaskBkgLikeSignJPSI.C"; taskName.Prepend(loadMacroPath.Data());
215     gROOT->LoadMacro(taskName.Data());
216     AliAnalysisTaskSEBkgLikeSignJPSI *lsJPSITask = AddTaskBkgLikeSignJPSI();
217     
218     //taskName="AddTaskBtoJPSItoEle.C"; taskName.Prepend(loadMacroPath.Data());
219     //gROOT->LoadMacro(taskName.Data());
220     //AliAnalysisTaskSEBtoJPSItoEle *jpsiTask = AddTaskBtoJPSItoEle();
221     
222     taskName="AddTaskCFMultiVarMultiStep.C"; taskName.Prepend(loadMacroPath.Data());
223     gROOT->LoadMacro(taskName.Data());
224     AliCFHeavyFlavourTaskMultiVarMultiStep *cfmvmsTask = AddTaskCFMultiVarMultiStep();
225     
226
227     taskName="AddTaskSECharmFraction.C"; 
228     taskName.Prepend(loadMacroPath.Data());
229     gROOT->LoadMacro(taskName.Data());
230     Int_t switchMC[5]={0,0,0,0,0};
231     Int_t ppPbPb=1;// 0 for pp, 1 for PbPb, used to siwtch on/off the removal of daughters from the primary vertex
232     AliAnalysisTaskSECharmFraction *cFractTask = AddTaskSECharmFraction("standard",switchMC,readMC,kTRUE,kFALSE,"D0toKpiCharmFractCuts.root","c",ppPbPb);
233     // arguments: filename,switchMC,readmc,usepid,likesign,cutfilename,containerprefix
234
235     
236     // attach a private task (not committed)
237     // (the files MyTask.h MyTask.cxx AddMyTask.C have to be declared in plugin
238     // configuration, see below)
239     
240     if(analysisMode.Data()=="proof") {
241     gProof->LoadMacro("MyTask.cxx++g");
242     } else {
243     gROOT->LoadMacro("MyTask.cxx++g");
244     }
245     gROOT->LoadMacro("AddMyTask.C");
246     MyTask *myTask = AddMyTask();
247     
248     
249     if(analysisMode.Data()=="proof") {
250     gProof->LoadMacro("AliDStarJets.cxx++g");
251     } else {
252     gROOT->LoadMacro("AliDStarJets.cxx++g");
253     }
254     gROOT->LoadMacro("AddTaskDStarJets.C");
255     AliDStarJets *myTask = AddTaskDStarJets();
256   */
257   //-------------------------------------------------------------------
258   
259   //
260   // Run the analysis
261   //    
262   if(chainAOD) printf("CHAIN HAS %d ENTRIES\n",(Int_t)chainAOD->GetEntries());
263   
264   if(!mgr->InitAnalysis()) return;
265   mgr->PrintStatus();
266   if(analysisMode=="grid" && !useAlienPlugin) analysisMode="local";
267   if(analysisMode!="proof") {
268     mgr->StartAnalysis(analysisMode.Data(),chainAOD,nentries,firstentry);
269   } else {
270     // proof
271     mgr->StartAnalysis(analysisMode.Data(),dataset.Data(),nentries,firstentry);
272   }
273   
274   return;
275 }
276 //_____________________________________________________________________________
277 //
278 AliAnalysisGrid* CreateAlienHandler(TString pluginmode="test",Bool_t useParFiles=kFALSE)
279 {
280   // Check if user has a valid token, otherwise make one. This has limitations.
281   // One can always follow the standard procedure of calling alien-token-init then
282   //   source /tmp/gclient_env_$UID in the current shell.
283    AliAnalysisAlien *plugin = new AliAnalysisAlien();
284    // Set the run mode (can be "full", "test", "offline", "submit" or "terminate")
285    plugin->SetRunMode(pluginmode.Data());
286    plugin->SetUser();
287    // Set versions of used packages
288    plugin->SetAPIVersion("V1.1x");
289    plugin->SetROOTVersion();
290    plugin->SetAliROOTVersion();
291
292    // Declare input data to be processed.
293    //************************************************
294    // Set data search pattern for DATA
295    //************************************************   
296    plugin->SetGridDataDir("/alice/data/2010/LHC10d"); // specify LHC period
297    plugin->SetDataPattern("pass2/AOD018/*AliAOD.root"); // specify reco pass and AOD set
298    plugin->SetFriendChainName("./AliAOD.VertexingHF.root");
299    // OR plugin->SetFriendChainName("deltas/AliAOD.VertexingHF.root");
300    // Adds only the good runs from the Monalisa Run Condition Table
301    // More than one period can be added but the period name has to be removed from GridDataDir (to be tested)
302    Int_t totruns=0;
303    //totruns += AddGoodRuns(plugin,"LHC10b"); // specify LHC period
304    //totruns += AddGoodRuns(plugin,"LHC10c"); // specify LHC period
305    totruns += AddGoodRuns(plugin,"LHC10d"); // specify LHC period
306    plugin->SetNrunsPerMaster(totruns);
307    
308    //************************************************
309    // Set data search pattern for MONTECARLO
310    //************************************************
311    /* 
312    plugin->SetGridDataDir("/alice/sim/LHC10d3"); // specify MC sample
313    plugin->SetDataPattern("AOD005/*AliAOD.root"); // specify AOD set
314    plugin->SetFriendChainName("./AliAOD.VertexingHF.root");
315    // OR plugin->SetFriendChainName("deltas/AliAOD.VertexingHF.root");
316    // Adds only the good runs from the Monalisa Run Condition Table 
317    // More than one period can be added!
318    Int_t totruns=0;
319    totruns += AddGoodRuns(plugin,"LHC10b","LHC10d3"); // specify LHC period for anchor runs; and the name of the MC production
320    //totruns += AddGoodRuns(plugin,"LHC10c","LHC10f7"); // specify LHC period for anchor runs;  and the name of the MC production
321    //totruns += AddGoodRuns(plugin,"LHC10d","LHC10f7"); // specify LHC period for anchor runs;  and the name of the MC production
322    plugin->SetNrunsPerMaster(totruns);
323    */
324    //
325    // Define alien work directory where all files will be copied. Relative to alien $HOME.
326    plugin->SetGridWorkingDir("myHFanalysis");
327    // Declare alien output directory. Relative to working directory.
328    plugin->SetGridOutputDir("output"); // In this case will be $HOME/work/output
329    // Declare the analysis source files names separated by blancs. To be compiled runtime
330    // using ACLiC on the worker nodes.
331    //plugin->SetAnalysisSource("AliDStarJets.cxx");
332    // Declare all libraries (other than the default ones for the framework. These will be
333    // loaded by the generated analysis macro. Add all extra files (task .cxx/.h) here.
334    plugin->SetAdditionalLibs("libPWG3base.so libPWG3vertexingHF.so");
335    // use par files
336    if(useParFiles) {
337      plugin->EnablePackage("STEERBase.par");
338      plugin->EnablePackage("ESD.par");
339      plugin->EnablePackage("AOD.par");
340      plugin->EnablePackage("ANALYSIS.par");
341      plugin->EnablePackage("ANALYSISalice.par");
342      plugin->EnablePackage("CORRFW.par");
343      plugin->EnablePackage("PWG3base.par");
344      plugin->EnablePackage("PWG3vertexingHF.par");
345    }
346    plugin->AddIncludePath("-I. -I$ROOTSYS/include -I$ALICE_ROOT -I$ALICE_ROOT/include -I$ALICE_ROOT/ITS -I$ALICE_ROOT/TPC -I$ALICE_ROOT/CONTAINERS -I$ALICE_ROOT/STEER -I$ALICE_ROOT/TRD -I$ALICE_ROOT/macros -I$ALICE_ROOT/ANALYSIS -I$ALICE_ROOT/PWG3 -I$ALICE_ROOT/PWG3/base -I$ALICE_ROOT/PWG3/vertexingHF -g");
347
348    plugin->SetDefaultOutputs(kTRUE);
349    // merging via jdl
350    plugin->SetMergeViaJDL(kTRUE);
351    plugin->SetOneStageMerging(kFALSE);
352    plugin->SetMaxMergeStages(2);
353
354    // Optionally set a name for the generated analysis macro (default MyAnalysis.C)
355    plugin->SetAnalysisMacro("AnalysisHF.C");
356    // Optionally set maximum number of input files/subjob (default 100, put 0 to ignore)
357    // Optionally modify the name of the generated JDL (default analysis.jdl)
358    plugin->SetJDLName("TaskHF.jdl");
359
360    return plugin;
361 }
362 //----------------------------------------------------------------------------
363 Int_t AddGoodRuns(AliAnalysisAlien* plugin,TString lhcPeriod,TString mcprod="") {
364   //
365   // Adds good runs from the Monalisa Run Condition Table
366   //
367   if(mcprod=="") plugin->SetRunPrefix("000"); // DATA
368
369   Int_t firstrun=0,lastrun=9999999;
370   Int_t nruns=0,ngoodruns=0;
371
372   if(mcprod=="LHC10d3") {firstrun=117054;lastrun=117222;}
373   if(mcprod=="LHC10d5") {firstrun=117086;lastrun=117222;}
374
375
376   if(lhcPeriod=="LHC10b") {
377     nruns=31;
378     Int_t runlist[31]={117222, 117220, 117116, 117112, 117109, 117099, 117092, 117086, 117077, 117065, 117063, 117060, 117059, 117054, 117053, 117052, 117050, 117048, 116645, 116643, 116574, 116571, 116562, 116403, 116288, 116102, 115401, 115393, 115193, 115186, 114931};
379    
380     for(Int_t k=0;k<nruns;k++){
381       if(runlist[k]<firstrun || runlist[k]>lastrun) continue;
382       plugin->AddRunNumber(runlist[k]);
383       ngoodruns++;
384     }
385     plugin->SetNrunsPerMaster(ngoodruns);
386   }
387
388   if(lhcPeriod=="LHC10c") { 
389     nruns=36;
390     Int_t runlist[36]={120829, 120825, 120824, 120823, 120822, 120821, 120820, 120758, 120750, 120741, 120671, 120617, 120616, 120505, 120504, 120503, 120244, 120079, 120076, 120073, 120072, 120069, 120067, 119862, 119859, 119856, 119853, 119849, 119846, 119845, 119844, 119842, 119841, 119163, 119161, 119159};
391    
392     for(Int_t k=0;k<nruns;k++){
393       if(runlist[k]<firstrun || runlist[k]>lastrun) continue;
394       plugin->AddRunNumber(runlist[k]);
395       ngoodruns++;
396     }
397     plugin->SetNrunsPerMaster(ngoodruns);
398   }
399
400   if(lhcPeriod=="LHC10dhighmu") { // only runs with high mu
401     nruns=17;
402     Int_t runlist[17]={124750, 124746, 124702, 124608, 124607, 124606, 124605, 124604, 124381, 124380, 124378, 124367, 124362, 124358, 124355, 124191, 124187};
403    
404     for(Int_t k=0;k<nruns;k++){
405       if(runlist[k]<firstrun || runlist[k]>lastrun) continue;
406       plugin->AddRunNumber(runlist[k]);
407       ngoodruns++;
408     }
409     plugin->SetNrunsPerMaster(ngoodruns);
410   }
411
412   if(lhcPeriod=="LHC10d") { // runs with high mu excluded
413     nruns=55;
414     Int_t runlist[55]={126437, 126432, 126425, 126424, 126422, 126409, 126408, 126407, 126406, 126405, 126404, 126403, 126359, 126352, 126351, 126350, 126285, 126284, 126283, 126168, 126167, 126160, 126158, 126097, 126090, 126088, 126082, 126081, 126078, 126073, 126008, 126007, 126004, 125855, 125851, 125850, 125849, 125848, 125847, 125844, 125843, 125842, 125633, 125632, 125630, 125296, 125134, 125101, 125100, 125097, 125085, 125023, 124751, 122375, 122374};
415    
416     for(Int_t k=0;k<nruns;k++){
417       if(runlist[k]<firstrun || runlist[k]>lastrun) continue;
418       plugin->AddRunNumber(runlist[k]);
419       ngoodruns++;
420     }
421     plugin->SetNrunsPerMaster(ngoodruns);
422   }
423
424
425   return ngoodruns;
426 }