]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/RESONANCES/macros/mini/ConfigPhiPbPbTPC.C
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGLF / RESONANCES / macros / mini / ConfigPhiPbPbTPC.C
1 //
2 // *** Configuration script for phi->KK analysis with 2010 runs ***
3 // 
4 // A configuration script for RSN package needs to define the followings:
5 //
6 // (1) decay tree of each resonance to be studied, which is needed to select
7 //     true pairs and to assign the right mass to all candidate daughters
8 // (2) cuts at all levels: single daughters, tracks, events
9 // (3) output objects: histograms or trees
10 //
11 Bool_t ConfigPhiPbPbTPC
12 (  
13    AliRsnMiniAnalysisTask *task,
14    Bool_t                  isMC,
15    Bool_t                  isESD,
16    const char             *suffix,
17    AliRsnCutSet           *cutsPair
18 )
19 {
20    // manage suffix
21    if (strlen(suffix) > 0) suffix = Form("_%s", suffix);
22    
23    // 
24    // -- Define track cuts -------------------------------------------------------------------------
25    //
26    
27    // BB parameterization depends on data sample (MC, data)
28    Double_t bbPar[5];
29    if (isMC) {
30       bbPar[0] = 2.15898 / 50.0;
31       bbPar[1] = 1.75295E1;
32       bbPar[2] = 3.40030E-9;
33       bbPar[3] = 1.96178;
34       bbPar[4] = 3.91720;
35    } else {
36       bbPar[0] = 1.41543 / 50.0;
37       bbPar[1] = 2.63394E1;
38       bbPar[2] = 5.0411E-11;
39       bbPar[3] = 2.12543;
40       bbPar[4] = 4.88663;
41    }
42    
43    // standard kaon cut
44    AliRsnCutKaonForPhi2010 *cut = new AliRsnCutKaonForPhi2010(Form("cut%s", suffix), 3.0, 3.0, 0.8);
45    
46    // setup (set manually the TPC PID)
47    cut->SetMode(AliRsnCutKaonForPhi2010::kOnlyTPC);
48    cut->InitMyPID(isMC, isESD);
49    cut->MyPID()->GetTPCResponse().SetBetheBlochParameters(bbPar[0], bbPar[1], bbPar[2], bbPar[3], bbPar[4]);
50    
51    // cut set
52    AliRsnCutSet *cutSet = new AliRsnCutSet(Form("set%s", suffix), AliRsnTarget::kDaughter);
53    cutSet->AddCut(cut);
54    cutSet->SetCutScheme(cut->GetName());
55    
56    // add to task
57    Int_t icut = task->AddTrackCuts(cutSet);
58    ::Info("Config", "Cut ID = %d", icut);
59    
60    //
61    // -- Values ------------------------------------------------------------------------------------
62    //
63    
64    /* invariant mass   */ Int_t imID   = task->CreateValue(AliRsnMiniValue::kInvMass, kFALSE);
65    /* IM resolution    */ Int_t resID  = task->CreateValue(AliRsnMiniValue::kInvMassRes, kTRUE);
66    /* transv. momentum */ Int_t ptID   = task->CreateValue(AliRsnMiniValue::kPt, kFALSE);
67    /* centrality       */ Int_t centID = task->CreateValue(AliRsnMiniValue::kMult, kFALSE);
68    
69    //
70    // -- Create all needed outputs -----------------------------------------------------------------
71    //
72    
73    // use an array for more compact writing, which are different on mixing and charges
74    // [0] = unlike
75    // [1] = mixing
76    // [2] = like ++
77    // [3] = like --
78    Bool_t  use     [6] = { 1      ,  0      ,  0      ,  0      ,  isMC   ,  isMC     };
79    Bool_t  useIM   [6] = { 1      ,  1      ,  1      ,  1      ,  1      ,  0        };
80    TString name    [6] = {"Unlike", "Mixing", "LikePP", "LikeMM", "Trues" , "Res"     };
81    TString comp    [6] = {"PAIR"  , "MIX"   , "PAIR"  , "PAIR"  , "TRUE"  , "TRUE"    };
82    TString output  [6] = {"SPARSE", "SPARSE", "SPARSE", "SPARSE", "SPARSE", "SPARSE"  };
83    Char_t  charge1 [6] = {'+'     , '+'     , '+'     , '-'     , '+'     , '+'       };
84    Char_t  charge2 [6] = {'-'     , '-'     , '+'     , '-'     , '-'     , '-'       };
85    Int_t   cutID   [6] = { icut   ,  icut   ,  icut   ,  icut   ,  icut   ,  icut     };
86    
87    for (Int_t i = 0; i < 6; i++) {
88       if (!use[i]) continue;
89       // create output
90       AliRsnMiniOutput *out = task->CreateOutput(Form("phi_%s%s", name[i].Data(), suffix), output[i].Data(), comp[i].Data());
91       // selection settings
92       out->SetCutID(0, cutID[i]);
93       out->SetCutID(1, cutID[i]);
94       out->SetDaughter(0, AliRsnDaughter::kKaon);
95       out->SetDaughter(1, AliRsnDaughter::kKaon);
96       out->SetCharge(0, charge1[i]);
97       out->SetCharge(1, charge2[i]);
98       out->SetMotherPDG(333);
99       out->SetMotherMass(1.019455);
100       // pair cuts
101       out->SetPairCuts(cutsPair);
102       // axis X: invmass (or resolution)
103       if (useIM) 
104          out->AddAxis(imID, 500, 0.9,  1.4);
105       else
106          out->AddAxis(resID, 200, -0.02, 0.02);
107       // axis Y: transverse momentum
108       out->AddAxis(ptID, 100, 0.0, 10.0);
109       // axis Z: centrality 
110       out->AddAxis(centID, 100, 0.0, 100.0);
111    }
112    
113    return kTRUE;
114 }