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