]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/SPECTRA/PiKaPr/TestAOD/runAODProof.C
PostData() fix, runAODProof updated
[u/mrichter/AliRoot.git] / PWGLF / SPECTRA / PiKaPr / TestAOD / runAODProof.C
1 #include "AliSpectraAODPID.h" 
2
3
4 void runAODProof(Int_t c=1, const char * proofMode = "full")
5 { //1 data AOD049
6   //2 MC AOD048
7   //3 data AOD086
8   //4 MC AOD090
9   
10   gEnv->SetValue("XSec.GSI.DelegProxy", "2");
11   
12   gSystem->Load("libTree.so");
13    gSystem->Load("libGeom.so");
14    gSystem->Load("libVMC.so");
15    gSystem->Load("libPhysics.so");
16    gSystem->Load("libSTEERBase.so");
17    gSystem->Load("libESD.so");
18    gSystem->Load("libAOD.so");
19    gSystem->Load("libANALYSIS.so");
20    gSystem->Load("libOADB.so");
21    gSystem->Load("libANALYSISalice.so");
22    gSystem->AddIncludePath("-I$ALICE_ROOT/include");
23
24
25    AliAnalysisAlien * handler = new AliAnalysisAlien("test");
26    handler->SetOverwriteMode();
27    handler->SetRunMode(proofMode);
28    handler->SetProofReset(0);
29    //handler->SetROOTVersion("v5-33-02a");
30    //handler->SetAliROOTVersion("v5-03-11-AN");
31    handler->SetAliROOTVersion("v5-03-17-AN");
32    
33    //handler->SetNproofWorkers(2);
34    handler->SetNproofWorkersPerSlave(1);
35    handler->SetProofCluster(Form("%s@alice-caf.cern.ch", gSystem->Getenv("CAFUSER")));
36    //handler->SetProofCluster(Form("%s@skaf.saske.sk",gSystem->Getenv("CAFUSER")));
37    
38    // Set handler for Real DATA:
39    if (c == 1){
40      handler->SetProofDataSet("/alice/data/LHC10h_000138653_p2_AOD049#aodTree");
41      //handler->SetProofDataSet("/alice/data/LHC10h_000138653_p2_AOD049#aodTree|/alice/data/LHC10h_000138662_p2_AOD049#aodTree|/alice/data/LHC10h_000138666_p2_AOD049#aodTree|/alice/data/LHC10h_000138795_p2_AOD049#aodTree|/alice/data/LHC10h_000139107_p2_AOD049#aodTree|/alice/data/LHC10h_000139110_p2_AOD049#aodTree");
42    }
43    if (c == 2){
44      handler->SetProofDataSet("/alice/sim/LHC11a10a_000138653_AOD048#aodTree|/alice/sim/LHC11a10a_000138662_AOD048#aodTree|/alice/sim/LHC11a10a_000138666_AOD048#aodTree|/alice/sim/LHC11a10a_000138795_AOD048#aodTree|/alice/sim/LHC11a10a_000139107_AOD048#aodTree|/alice/sim/LHC11a10a_000139110_AOD048#aodTree");      
45    }
46    if (c == 3){
47      handler->SetProofDataSet("/default/lmilano/LHC10h_000138653_AOD086_p2#aodTree|/default/lmilano/LHC10h_000138662_AOD086_p2#aodTree|/default/lmilano/LHC10h_000138666_AOD086_p2#aodTree|/default/lmilano/LHC10h_000139107_AOD086_p2#aodTree");      
48    }
49    if (c == 4){
50      handler->SetProofDataSet("/default/lmilano/LHC11a10a_bis_138653_AOD090#aodTree|/default/lmilano/LHC11a10a_bis_138662_AOD090#aodTree|/default/lmilano/LHC11a10a_bis_138666_AOD090#aodTree|/default/lmilano/LHC11a10a_bis_139107_AOD090#aodTree");      
51    }
52    
53    gROOT->LoadMacro("AliSpectraAODTrackCuts.cxx+g");
54    gROOT->LoadMacro("AliSpectraAODEventCuts.cxx+g");
55    gROOT->LoadMacro("AliSpectraAODHistoManager.cxx+g");
56    gROOT->LoadMacro("AliSpectraAODPID.cxx+g");
57    gROOT->LoadMacro("AliAnalysisTaskSpectraAOD.cxx+g");
58    
59    handler->SetAliRootMode("default");
60    handler->SetAdditionalLibs("AliSpectraAODHistoManager.cxx AliSpectraAODHistoManager.h AliSpectraAODPID.cxx AliSpectraAODPID.h AliSpectraAODEventCuts.cxx AliSpectraAODEventCuts.h AliSpectraAODTrackCuts.cxx AliSpectraAODTrackCuts.h AliAnalysisTaskSpectraAOD.cxx AliAnalysisTaskSpectraAOD.h");
61    handler->SetAnalysisSource("Histograms.h HistogramNames.h AliSpectraAODHistoManager.cxx+ AliSpectraAODEventCuts.cxx+ AliSpectraAODTrackCuts.cxx+ AliSpectraAODPID.cxx+ AliAnalysisTaskSpectraAOD.cxx+");
62    handler->SetFileForTestMode("filelist.txt"); // list of local files for testing
63    //  handler->SetAliRootMode("");
64    handler->SetClearPackages();
65    
66    AliAnalysisManager *mgr = new AliAnalysisManager("testAnalysis");
67    mgr->SetGridHandler(handler);
68    AliAODInputHandler* aodH = new AliAODInputHandler();
69    mgr->SetInputEventHandler(aodH);
70    
71    // Add PID task
72    gROOT->LoadMacro("$ALICE_ROOT/ANALYSIS/macros/AddTaskPIDResponse.C");
73    Bool_t isMC = kFALSE;
74    if (c == 2 || c == 4) isMC = kTRUE;   
75    AliAnalysisTask * taskPID = AddTaskPIDResponse(isMC);
76    mgr->AddTask(taskPID);
77    
78    //LOOP OVER SELECTION
79    Double_t CentCutMin[4]={0,0,20,20};
80    Double_t CentCutMax[4]={100,5,40,40};
81    Double_t QvecCutMin[4]={0,0,0,3};
82    Double_t QvecCutMax[4]={100,100,2,100};
83    // Double_t CentCutMin[4]={5,10,20,30};
84    // Double_t CentCutMax[4]={10,20,30,40};
85    // Double_t QvecCutMin[4]={0,0,0,0};
86    // Double_t QvecCutMax[4]={100,100,100,100};
87    AliSpectraAODPID * pid = new AliSpectraAODPID(AODPIDType_t::kNSigmaTPCTOF);
88    pid->SetNSigmaCut(5.);
89    
90    for(Int_t iCut=1;iCut<2;iCut++){
91      AliAnalysisTaskSpectraAOD *task = new AliAnalysisTaskSpectraAOD("TaskAODExercise");
92      mgr->AddTask(task);
93      //physics selection
94      task->SelectCollisionCandidates();     
95      // set pid object
96      task->SetPID(pid);
97      // Set the cuts
98      AliSpectraAODEventCuts * vcuts = new AliSpectraAODEventCuts("Event Cuts");
99      AliSpectraAODTrackCuts  * tcuts = new AliSpectraAODTrackCuts("Track Cuts");
100      //if(c==1 || c==2)tcuts->SetTrackType(5); //AOD 046 & 047. Standard Cut with loose DCA
101      if(c==1 || c==2)tcuts->SetTrackType(5); //AOD 046 & 047. Standard Cut with tight DCA
102      if(c==3 || c==4)tcuts->SetTrackType(10); //AOD 086 & 090. Standard Raa cut
103      tcuts->SetEta(.8);
104      //tcuts->SetDCA(.1);
105      tcuts->SetPt(5);
106      tcuts->SetY(.5);
107      tcuts->SetPtTOFMatching(0.6);   
108      tcuts->SetQvecMin(QvecCutMin[iCut]);   
109      tcuts->SetQvecMax(QvecCutMax[iCut]);    
110      vcuts->SetCentralityCutMax(CentCutMax[iCut]);  
111      vcuts->SetCentralityCutMin(CentCutMin[iCut]);
112      task->SetEventCuts(vcuts);
113      task->SetTrackCuts(tcuts);
114      vcuts->PrintCuts();
115      tcuts->PrintCuts();
116      
117      // check for MC or real data
118      if (c == 2 || c == 4)
119        {
120          task->SetIsMC(kTRUE);
121          AliAnalysisDataContainer *cinput = mgr->GetCommonInputContainer();
122          AliAnalysisDataContainer *coutputpt1 = mgr->CreateContainer(Form("chistpt%d",iCut), AliSpectraAODHistoManager::Class(),  AliAnalysisManager::kOutputContainer, 
123                                                                      Form("Pt.AOD.1._MC_Cent%.0fto%.0f_QVec%.1fto%.1f.root",CentCutMin[iCut],CentCutMax[iCut],QvecCutMin[iCut],QvecCutMax[iCut]));
124          AliAnalysisDataContainer *coutputpt2 = mgr->CreateContainer(Form("cvcutpt%d",iCut), AliSpectraAODEventCuts::Class(),    AliAnalysisManager::kOutputContainer, 
125                                                                      Form("Pt.AOD.1._MC_Cent%.0fto%.0f_QVec%.1fto%.1f.root",CentCutMin[iCut],CentCutMax[iCut],QvecCutMin[iCut],QvecCutMax[iCut]));
126          AliAnalysisDataContainer *coutputpt3 = mgr->CreateContainer(Form("ctcutpt%d",iCut), AliSpectraAODTrackCuts::Class(),     AliAnalysisManager::kOutputContainer, 
127                                                                      Form("Pt.AOD.1._MC_Cent%.0fto%.0f_QVec%.1fto%.1f.root",CentCutMin[iCut],CentCutMax[iCut],QvecCutMin[iCut],QvecCutMax[iCut]));
128          AliAnalysisDataContainer *coutputpt4 = mgr->CreateContainer(Form("cpidpt%d",iCut),  AliSpectraAODPID::Class(),     AliAnalysisManager::kOutputContainer, 
129                                                                      Form("Pt.AOD.1._MC_Cent%.0fto%.0f_QVec%.1fto%.1f.root",CentCutMin[iCut],CentCutMax[iCut],QvecCutMin[iCut],QvecCutMax[iCut]));
130        }
131      if (c == 1 || c==3)
132        {
133          AliAnalysisDataContainer *cinput = mgr->GetCommonInputContainer();
134          AliAnalysisDataContainer *coutputpt1 = mgr->CreateContainer(Form("chistpt%d",iCut), AliSpectraAODHistoManager::Class(),  AliAnalysisManager::kOutputContainer, 
135                                                                      Form("Pt.AOD.1._data_ptcut_Cent%.0fto%.0f_QVec%.1fto%.1f.root",CentCutMin[iCut],CentCutMax[iCut],QvecCutMin[iCut],QvecCutMax[iCut]));
136          AliAnalysisDataContainer *coutputpt2 = mgr->CreateContainer(Form("cvcutpt%d",iCut), AliSpectraAODEventCuts::Class(),    AliAnalysisManager::kOutputContainer, 
137                                                                      Form("Pt.AOD.1._data_ptcut_Cent%.0fto%.0f_QVec%.1fto%.1f.root",CentCutMin[iCut],CentCutMax[iCut],QvecCutMin[iCut],QvecCutMax[iCut]));
138          AliAnalysisDataContainer *coutputpt3 = mgr->CreateContainer(Form("ctcutpt%d",iCut), AliSpectraAODTrackCuts::Class(),     AliAnalysisManager::kOutputContainer, 
139                                                                      Form("Pt.AOD.1._data_ptcut_Cent%.0fto%.0f_QVec%.1fto%.1f.root",CentCutMin[iCut],CentCutMax[iCut],QvecCutMin[iCut],QvecCutMax[iCut]));
140          AliAnalysisDataContainer *coutputpt4 = mgr->CreateContainer(Form("cpidpt%d",iCut),  AliSpectraAODPID::Class(),     AliAnalysisManager::kOutputContainer, 
141                                                                      Form("Pt.AOD.1._data_Cent%.0fto%.0f_QVec%.1fto%.1f.root",CentCutMin[iCut],CentCutMax[iCut],QvecCutMin[iCut],QvecCutMax[iCut]));
142          
143        }
144      mgr->ConnectInput(task, 0, cinput);
145      mgr->ConnectOutput(task, 1, coutputpt1);
146      mgr->ConnectOutput(task, 2, coutputpt2);
147      mgr->ConnectOutput(task, 3, coutputpt3);
148      mgr->ConnectOutput(task, 4, coutputpt4);
149    }
150    mgr->SetDebugLevel(2);
151    
152    if (!mgr->InitAnalysis()) return;
153    mgr->PrintStatus();
154    mgr->StartAnalysis("proof");
155 }