]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/RESONANCES/macros/mini/ConfigDeltaPP7TeV.C
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGLF / RESONANCES / macros / mini / ConfigDeltaPP7TeV.C
1 //
2 // *** Configuration script for delta(++), delta(--) and delta(0) 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 ConfigDeltaPP7TeV
12 (  
13    AliRsnMiniAnalysisTask *task, 
14    Bool_t                  isMC, 
15    const char             *suffix,
16    AliRsnCutSet           *cutsPair
17 )
18 {
19    // manage suffix
20    if (strlen(suffix) > 0) suffix = Form("_%s", suffix);
21
22
23    // integrated pion cut TOF
24    AliRsnCutDelta *cutPiTOF = new AliRsnCutDelta("cutPionTOF",AliPID::kPion,kFALSE);
25    // cut set TOF
26    AliRsnCutSet *cutSetPiTOF = new AliRsnCutSet("setPionForDeltaTOF", AliRsnTarget::kDaughter);
27    cutSetPiTOF->AddCut(cutPiTOF);
28    cutSetPiTOF->SetCutScheme(cutPiTOF->GetName());   
29    // add to task TOF
30    Int_t iCutPiTOF = task->AddTrackCuts(cutSetPiTOF);
31    
32    // integrated proton cut TOF
33    AliRsnCutDelta *cutPTOF = new AliRsnCutDelta("cutProtonTOF",AliPID::kProton,kFALSE);
34    // cut set TOF
35    AliRsnCutSet *cutSetPTOF = new AliRsnCutSet("setProtonForDeltaTOF", AliRsnTarget::kDaughter);
36    cutSetPTOF->AddCut(cutPTOF);
37    cutSetPTOF->SetCutScheme(cutPTOF->GetName());
38    // add to task TOF
39    Int_t iCutPTOF = task->AddTrackCuts(cutSetPTOF);
40
41 /////////////////////////////////////////////////////SYSTEMATICS
42  
43    AliRsnCutDelta *cutPiTOF1 = new AliRsnCutDelta("cutPionTOF1",AliPID::kPion,kFALSE);
44    cutPiTOF1->SetTPCNSigmaProton(2.0);
45    cutPiTOF1->SetTPCNSigmaPion(3.0);
46    cutPiTOF1->SetTOFNSigmaProton(2.0);
47    cutPiTOF1->SetTOFNSigmaPion(3.0);
48    AliRsnCutSet *cutSetPiTOF1 = new AliRsnCutSet("setPionForDeltaTOF1", AliRsnTarget::kDaughter);
49    cutSetPiTOF1->AddCut(cutPiTOF1);
50    cutSetPiTOF1->SetCutScheme(cutPiTOF1->GetName());
51    Int_t iCutPiTOF1 = task->AddTrackCuts(cutSetPiTOF1);
52
53    AliRsnCutDelta *cutPTOF1 = new AliRsnCutDelta("cutProtonTOF1",AliPID::kProton,kFALSE);
54    cutPTOF1->SetTPCNSigmaProton(2.0);
55    cutPTOF1->SetTPCNSigmaPion(3.0);
56    cutPTOF1->SetTOFNSigmaProton(2.0);
57    cutPTOF1->SetTOFNSigmaPion(3.0);
58    AliRsnCutSet *cutSetPTOF1 = new AliRsnCutSet("setProtonForDeltaTOF1", AliRsnTarget::kDaughter);
59    cutSetPTOF1->AddCut(cutPTOF1);
60    cutSetPTOF1->SetCutScheme(cutPTOF1->GetName());
61    Int_t iCutPTOF1 = task->AddTrackCuts(cutSetPTOF1);
62
63    AliRsnCutDelta *cutPiTOF2 = new AliRsnCutDelta("cutPionTOF2",AliPID::kPion,kFALSE);
64    cutPiTOF2->SetTPCNSigmaProton(3.0);
65    cutPiTOF2->SetTPCNSigmaPion(2.0);
66    cutPiTOF2->SetTOFNSigmaProton(3.0);
67    cutPiTOF2->SetTOFNSigmaPion(2.0);
68    AliRsnCutSet *cutSetPiTOF2 = new AliRsnCutSet("setPionForDeltaTOF2", AliRsnTarget::kDaughter);
69    cutSetPiTOF2->AddCut(cutPiTOF2);
70    cutSetPiTOF2->SetCutScheme(cutPiTOF2->GetName());
71    Int_t iCutPiTOF2 = task->AddTrackCuts(cutSetPiTOF2);
72
73    AliRsnCutDelta *cutPTOF2 = new AliRsnCutDelta("cutProtonTOF2",AliPID::kProton,kFALSE);
74    cutPTOF2->SetTPCNSigmaProton(3.0);
75    cutPTOF2->SetTPCNSigmaPion(2.0);
76    cutPTOF2->SetTOFNSigmaProton(3.0);
77    cutPTOF2->SetTOFNSigmaPion(2.0);
78    AliRsnCutSet *cutSetPTOF2 = new AliRsnCutSet("setProtonForDeltaTOF2", AliRsnTarget::kDaughter);
79    cutSetPTOF2->AddCut(cutPTOF2);
80    cutSetPTOF2->SetCutScheme(cutPTOF2->GetName());
81    Int_t iCutPTOF2 = task->AddTrackCuts(cutSetPTOF2); 
82
83    AliRsnCutDelta *cutPiTOF3 = new AliRsnCutDelta("cutPionTOF3",AliPID::kPion,kFALSE);
84    cutPiTOF3->SetTPCNSigmaProton(4.0);
85    cutPiTOF3->SetTPCNSigmaPion(3.0);
86    cutPiTOF3->SetTOFNSigmaProton(4.0);
87    cutPiTOF3->SetTOFNSigmaPion(3.0);
88    AliRsnCutSet *cutSetPiTOF3 = new AliRsnCutSet("setPionForDeltaTOF3", AliRsnTarget::kDaughter);
89    cutSetPiTOF3->AddCut(cutPiTOF3);
90    cutSetPiTOF3->SetCutScheme(cutPiTOF3->GetName());
91    Int_t iCutPiTOF3 = task->AddTrackCuts(cutSetPiTOF3);
92
93    AliRsnCutDelta *cutPTOF3 = new AliRsnCutDelta("cutProtonTOF3",AliPID::kProton,kFALSE);
94    cutPiTOF3->SetTPCNSigmaProton(4.0);
95    cutPiTOF3->SetTPCNSigmaPion(3.0);
96    cutPiTOF3->SetTOFNSigmaProton(4.0);
97    cutPiTOF3->SetTOFNSigmaPion(3.0);
98    AliRsnCutSet *cutSetPTOF3 = new AliRsnCutSet("setProtonForDeltaTOF3", AliRsnTarget::kDaughter);
99    cutSetPTOF3->AddCut(cutPTOF3);
100    cutSetPTOF3->SetCutScheme(cutPTOF3->GetName());
101    Int_t iCutPTOF3 = task->AddTrackCuts(cutSetPTOF3); 
102
103    AliRsnCutDelta *cutPiTOF4 = new AliRsnCutDelta("cutPionTOF4",AliPID::kPion,kFALSE);
104    cutPiTOF4->SetTPCNSigmaProton(3.0);
105    cutPiTOF4->SetTPCNSigmaPion(4.0);
106    cutPiTOF4->SetTOFNSigmaProton(3.0);
107    cutPiTOF4->SetTOFNSigmaPion(4.0);
108    AliRsnCutSet *cutSetPiTOF4 = new AliRsnCutSet("setPionForDeltaTOF4", AliRsnTarget::kDaughter);
109    cutSetPiTOF4->AddCut(cutPiTOF4);
110    cutSetPiTOF4->SetCutScheme(cutPiTOF4->GetName());
111    Int_t iCutPiTOF4 = task->AddTrackCuts(cutSetPiTOF4);
112
113    AliRsnCutDelta *cutPTOF4 = new AliRsnCutDelta("cutProtonTOF4",AliPID::kProton,kFALSE);
114    cutPTOF4->SetTPCNSigmaProton(3.0);
115    cutPTOF4->SetTPCNSigmaPion(4.0);
116    cutPTOF4->SetTOFNSigmaProton(3.0);
117    cutPTOF4->SetTOFNSigmaPion(4.0);
118    AliRsnCutSet *cutSetPTOF4 = new AliRsnCutSet("setProtonForDeltaTOF4", AliRsnTarget::kDaughter);
119    cutSetPTOF4->AddCut(cutPTOF4);
120    cutSetPTOF4->SetCutScheme(cutPTOF4->GetName());
121    Int_t iCutPTOF4 = task->AddTrackCuts(cutSetPTOF4);
122
123    AliRsnCutDelta *cutPiTOF5 = new AliRsnCutDelta("cutPionTOF5",AliPID::kPion,kFALSE);
124    cutPiTOF5->SetTOFMomProton(2.3);
125    AliRsnCutSet *cutSetPiTOF5 = new AliRsnCutSet("setPionForDeltaTOF5", AliRsnTarget::kDaughter);
126    cutSetPiTOF5->AddCut(cutPiTOF5);
127    cutSetPiTOF5->SetCutScheme(cutPiTOF5->GetName());
128    Int_t iCutPiTOF5 = task->AddTrackCuts(cutSetPiTOF5);
129
130    AliRsnCutDelta *cutPTOF5 = new AliRsnCutDelta("cutProtonTOF5",AliPID::kProton,kFALSE);
131    cutPTOF5->SetTOFMomProton(2.3);
132    AliRsnCutSet *cutSetPTOF5 = new AliRsnCutSet("setProtonForDeltaTOF5", AliRsnTarget::kDaughter);
133    cutSetPTOF5->AddCut(cutPTOF5);
134    cutSetPTOF5->SetCutScheme(cutPTOF5->GetName());
135    Int_t iCutPTOF5 = task->AddTrackCuts(cutSetPTOF5);
136
137    AliRsnCutDelta *cutPiTOF6 = new AliRsnCutDelta("cutPionTOF6",AliPID::kPion,kFALSE);
138    cutPiTOF6->SetTOFMomProton(2.7);
139    AliRsnCutSet *cutSetPiTOF6 = new AliRsnCutSet("setPionForDeltaTOF6", AliRsnTarget::kDaughter);
140    cutSetPiTOF6->AddCut(cutPiTOF6);
141    cutSetPiTOF6->SetCutScheme(cutPiTOF6->GetName());
142    Int_t iCutPiTOF6 = task->AddTrackCuts(cutSetPiTOF6);
143
144    AliRsnCutDelta *cutPTOF6 = new AliRsnCutDelta("cutProtonTOF6",AliPID::kProton,kFALSE);
145    cutPTOF6->SetTOFMomProton(2.7);
146    AliRsnCutSet *cutSetPTOF6 = new AliRsnCutSet("setProtonForDeltaTOF6", AliRsnTarget::kDaughter);
147    cutSetPTOF6->AddCut(cutPTOF6);
148    cutSetPTOF6->SetCutScheme(cutPTOF6->GetName());
149    Int_t iCutPTOF6 = task->AddTrackCuts(cutSetPTOF6);
150
151 /////////////////////////////////////////////////////SYSTEMATICS END
152
153
154    
155    //
156    // -- Values ------------------------------------------------------------------------------------
157    //
158 //CreateValue(AliRsnMiniValue::EType type, Bool_t useMC = kFALSE); 
159
160    /* IM resolution    */ Int_t resID  = task->CreateValue(AliRsnMiniValue::kInvMassRes, kTRUE);
161    /* invariant mass   */ Int_t imID   = task->CreateValue(AliRsnMiniValue::kInvMass, kFALSE);
162    /* transv. momentum */ Int_t ptID   = task->CreateValue(AliRsnMiniValue::kPt, kFALSE);
163    /* centrality       */ Int_t centID = task->CreateValue(AliRsnMiniValue::kMult, kFALSE);
164
165
166    /* 1st daughter pt  */ Int_t fdpt   = task->CreateValue(AliRsnMiniValue::kFirstDaughterPt  , kFALSE);
167    /* 2nd daughter pt  */ Int_t sdpt   = task->CreateValue(AliRsnMiniValue::kSecondDaughterPt , kFALSE);
168    /* 1st daughter p   */ Int_t fdp    = task->CreateValue(AliRsnMiniValue::kFirstDaughterP   , kFALSE);
169    /* 2nd daughter p   */ Int_t sdp    = task->CreateValue(AliRsnMiniValue::kSecondDaughterP  , kFALSE);
170    
171    /* pseudorapidity   */ Int_t etaID  = task->CreateValue(AliRsnMiniValue::kEta , kFALSE);
172    /* rapidity         */ Int_t yID    = task->CreateValue(AliRsnMiniValue::kY   , kFALSE);
173
174
175    
176    //
177    // -- Create all needed outputs -----------------------------------------------------------------
178    //
179    
180    // use an array for more compact writing, which are different on mixing and charges
181    // [0] = unlike
182    // [1] = mixing
183    // [2] = like ++
184    // [3] = like --
185       
186    Bool_t  use      [12] = { 1          ,  1          ,  1         ,  1         ,  1         ,  1         ,  isMC       ,   isMC       ,  1          ,   1         ,  isMC       ,   isMC         };
187    Bool_t  useIM    [12] = { 1          ,  1          ,  1         ,  1         ,  1         ,  1         ,  1          ,   1          ,  1          ,   1         ,  1          ,   1            };
188    TString name     [12] = {"UnlikePM"  , "UnlikeMP"  , "MixingPM" , "MixingMP" , "LikePP"   , "LikeMM"   , "TruesPP"   ,  "TruesMM"   , "MixingPP"  ,  "MixingMM" , "TRUESPM"   ,  "TRUESMP"     };
189    TString comp     [12] = {"PAIR"      , "PAIR"      , "MIX"      , "MIX"      , "PAIR"     , "PAIR"     , "TRUE"      ,  "TRUE"      , "MIX"       ,  "MIX"      , "TRUE"      ,  "TRUE"        };
190    TString output   [12] = {"SPARSE"    , "SPARSE"    , "SPARSE"   , "SPARSE"   , "SPARSE"   , "SPARSE"   , "SPARSE"    ,  "SPARSE"    , "SPARSE"    ,  "SPARSE"   , "SPARSE"    ,  "SPARSE"      };
191    Char_t  charge1  [12] = {'+'         , '-'         , '+'        , '-'        , '+'        , '-'        , '+'         ,  '-'         , '+'         ,  '-'        , '+'         ,  '-'           };
192    Char_t  charge2  [12] = {'-'         , '+'         , '-'        , '+'        , '+'        , '-'        , '+'         ,  '-'         , '+'         ,  '-'        , '-'         ,  '+'           };
193    Int_t   cutID1   [12] = {iCutPTOF    , iCutPTOF    , iCutPTOF   , iCutPTOF   , iCutPTOF   , iCutPTOF   , iCutPTOF    ,  iCutPTOF    , iCutPTOF    ,  iCutPTOF   , iCutPTOF    ,  iCutPTOF      };
194    Int_t   cutID2   [12] = {iCutPiTOF   , iCutPiTOF   , iCutPiTOF  , iCutPiTOF  , iCutPiTOF  , iCutPiTOF  , iCutPiTOF   ,  iCutPiTOF   , iCutPiTOF   ,  iCutPiTOF  , iCutPiTOF   ,  iCutPiTOF     };
195    Int_t   cutID3   [12] = {iCutPTOF1   , iCutPTOF1   , iCutPTOF1  , iCutPTOF1  , iCutPTOF1  , iCutPTOF1  , iCutPTOF1   ,  iCutPTOF1   , iCutPTOF1   ,  iCutPTOF1  , iCutPTOF1   ,  iCutPTOF1     };
196    Int_t   cutID4   [12] = {iCutPiTOF1  , iCutPiTOF1  , iCutPiTOF1 , iCutPiTOF1 , iCutPiTOF1 , iCutPiTOF1 , iCutPiTOF1  ,  iCutPiTOF1  , iCutPiTOF1  ,  iCutPiTOF1 , iCutPiTOF1  ,  iCutPiTOF1    };
197    Int_t   cutID5   [12] = {iCutPTOF2   , iCutPTOF2   , iCutPTOF2  , iCutPTOF2  , iCutPTOF2  , iCutPTOF2  , iCutPTOF2   ,  iCutPTOF2   , iCutPTOF2   ,  iCutPTOF2  , iCutPTOF2   ,  iCutPTOF2     };
198    Int_t   cutID6   [12] = {iCutPiTOF2  , iCutPiTOF2  , iCutPiTOF2 , iCutPiTOF2 , iCutPiTOF2 , iCutPiTOF2 , iCutPiTOF2  ,  iCutPiTOF2  , iCutPiTOF2  ,  iCutPiTOF2 , iCutPiTOF2  ,  iCutPiTOF2    };
199    Int_t   cutID7   [12] = {iCutPTOF3   , iCutPTOF3   , iCutPTOF3  , iCutPTOF3  , iCutPTOF3  , iCutPTOF3  , iCutPTOF3   ,  iCutPTOF3   , iCutPTOF3   ,  iCutPTOF3  , iCutPTOF3   ,  iCutPTOF3     };
200    Int_t   cutID8   [12] = {iCutPiTOF3  , iCutPiTOF3  , iCutPiTOF3 , iCutPiTOF3 , iCutPiTOF3 , iCutPiTOF3 , iCutPiTOF3  ,  iCutPiTOF3  , iCutPiTOF3  ,  iCutPiTOF3 , iCutPiTOF3  ,  iCutPiTOF3    };
201    Int_t   cutID9   [12] = {iCutPTOF4   , iCutPTOF4   , iCutPTOF4  , iCutPTOF4  , iCutPTOF4  , iCutPTOF4  , iCutPTOF4   ,  iCutPTOF4   , iCutPTOF4   ,  iCutPTOF4  , iCutPTOF4   ,  iCutPTOF4     };
202    Int_t   cutID10  [12] = {iCutPiTOF4  , iCutPiTOF4  , iCutPiTOF4 , iCutPiTOF4 , iCutPiTOF4 , iCutPiTOF4 , iCutPiTOF4  ,  iCutPiTOF4  , iCutPiTOF4  ,  iCutPiTOF4 , iCutPiTOF4  ,  iCutPiTOF4    };
203    Int_t   cutID11  [12] = {iCutPTOF5   , iCutPTOF5   , iCutPTOF5  , iCutPTOF5  , iCutPTOF5  , iCutPTOF5  , iCutPTOF5   ,  iCutPTOF5   , iCutPTOF5   ,  iCutPTOF5  , iCutPTOF5   ,  iCutPTOF5     };
204    Int_t   cutID12  [12] = {iCutPiTOF5  , iCutPiTOF5  , iCutPiTOF5 , iCutPiTOF5 , iCutPiTOF5 , iCutPiTOF5 , iCutPiTOF5  ,  iCutPiTOF5  , iCutPiTOF5  ,  iCutPiTOF5 , iCutPiTOF5  ,  iCutPiTOF5    };
205    Int_t   cutID13  [12] = {iCutPTOF6   , iCutPTOF6   , iCutPTOF6  , iCutPTOF6  , iCutPTOF6  , iCutPTOF6  , iCutPTOF6   ,  iCutPTOF6   , iCutPTOF6   ,  iCutPTOF6  , iCutPTOF6   ,  iCutPTOF6     };
206    Int_t   cutID14  [12] = {iCutPiTOF6  , iCutPiTOF6  , iCutPiTOF6 , iCutPiTOF6 , iCutPiTOF6 , iCutPiTOF6 , iCutPiTOF6  ,  iCutPiTOF6  , iCutPiTOF6  ,  iCutPiTOF6 , iCutPiTOF6  ,  iCutPiTOF6    };
207    Int_t   ipdg     [12] = {2114        , -2114       , 2114       , -2114      , 2224       , -2224      , 2224        ,  -2224       , 2224        ,   -2224     , 2114        ,  -2114         };
208    
209    
210
211    for (Int_t i = 0; i < 12; i++) {
212       if (!use[i]) continue;
213  
214        AliRsnMiniOutput *outTOF = task->CreateOutput(Form("deltatof_%s%s", name[i].Data(), suffix), output[i].Data(), comp[i].Data());
215       // selection settings
216       outTOF->SetCutID(0, cutID1[i]);
217       outTOF->SetCutID(1, cutID2[i]);
218       outTOF->SetDaughter(0, AliRsnDaughter::kProton);
219       outTOF->SetDaughter(1, AliRsnDaughter::kPion);
220       outTOF->SetCharge(0, charge1[i]);
221       outTOF->SetCharge(1, charge2[i]);
222       outTOF->SetMotherPDG(ipdg[i]);
223       outTOF->SetMotherMass(1.232);
224       // pair cuts
225       outTOF->SetPairCuts(cutsPair);
226       // axis X: invmass (or resolution)
227       if (useIM[i])
228          outTOF->AddAxis(imID, 100, 1.0, 2.0);
229       else
230          outTOF->AddAxis(resID, 200, -0.02, 0.02);
231       // axis Y: transverse momentum
232       outTOF->AddAxis(ptID, 100, 0.0, 10.0);
233      // axis Z: Multiplicity
234       outTOF->AddAxis( centID , 120 , 0.0, 120  ); 
235       
236       outTOF->AddAxis( fdpt   , 100 , 0.0, 10.0 ); 
237       outTOF->AddAxis( sdpt   , 100 , 0.0, 10.0 ); 
238       outTOF->AddAxis( fdp    , 100 , 0.0, 10.0 ); 
239       outTOF->AddAxis( sdp    , 100 , 0.0, 10.0 ); 
240       outTOF->AddAxis( etaID  , 20  ,-1.0, 1.0  );
241       outTOF->AddAxis( yID    , 10  ,-0.5, 0.5 ); 
242  
243  
244  
245       AliRsnMiniOutput *outTOF1 = task->CreateOutput(Form("config1deltatof_%s%s", name[i].Data(), suffix), output[i].Data(), comp[i].Data());
246       outTOF1->SetCutID(0, cutID3[i]);
247       outTOF1->SetCutID(1, cutID4[i]);
248       outTOF1->SetDaughter(0, AliRsnDaughter::kProton);
249       outTOF1->SetDaughter(1, AliRsnDaughter::kPion);
250       outTOF1->SetCharge(0, charge1[i]);
251       outTOF1->SetCharge(1, charge2[i]);
252       outTOF1->SetMotherPDG(ipdg[i]);
253       outTOF1->SetMotherMass(1.232);
254       outTOF1->SetPairCuts(cutsPair);
255       if (useIM[i])
256          outTOF1->AddAxis(imID, 100, 1.0, 2.0);
257       else
258          outTOF1->AddAxis(resID, 200, -0.02, 0.02);
259       outTOF1->AddAxis(ptID, 100, 0.0, 10.0);
260       outTOF1->AddAxis(centID, 120, 0.0, 120); 
261
262
263       AliRsnMiniOutput *outTOF2 = task->CreateOutput(Form("config2deltatof_%s%s", name[i].Data(), suffix), output[i].Data(), comp[i].Data());
264       outTOF2->SetCutID(0, cutID5[i]);
265       outTOF2->SetCutID(1, cutID6[i]);
266       outTOF2->SetDaughter(0, AliRsnDaughter::kProton);
267       outTOF2->SetDaughter(1, AliRsnDaughter::kPion);
268       outTOF2->SetCharge(0, charge1[i]);
269       outTOF2->SetCharge(1, charge2[i]);
270       outTOF2->SetMotherPDG(ipdg[i]);
271       outTOF2->SetMotherMass(1.232);
272       outTOF2->SetPairCuts(cutsPair);
273       if (useIM[i])
274          outTOF2->AddAxis(imID, 100, 1.0, 2.0);
275       else
276          outTOF2->AddAxis(resID, 200, -0.02, 0.02);
277       outTOF2->AddAxis(ptID, 100, 0.0, 10.0);
278       outTOF2->AddAxis(centID, 120, 0.0, 120);
279
280
281       AliRsnMiniOutput *outTOF3 = task->CreateOutput(Form("config3deltatof_%s%s", name[i].Data(), suffix), output[i].Data(), comp[i].Data());
282       outTOF3->SetCutID(0, cutID7[i]);
283       outTOF3->SetCutID(1, cutID8[i]);
284       outTOF3->SetDaughter(0, AliRsnDaughter::kProton);
285       outTOF3->SetDaughter(1, AliRsnDaughter::kPion);
286       outTOF3->SetCharge(0, charge1[i]);
287       outTOF3->SetCharge(1, charge2[i]);
288       outTOF3->SetMotherPDG(ipdg[i]);
289       outTOF3->SetMotherMass(1.232);
290       outTOF3->SetPairCuts(cutsPair);
291       if (useIM[i])
292          outTOF3->AddAxis(imID, 100, 1.0, 2.0);
293       else
294          outTOF3->AddAxis(resID, 200, -0.02, 0.02);
295       outTOF3->AddAxis(ptID, 100, 0.0, 10.0);
296       outTOF3->AddAxis(centID, 120, 0.0, 120);
297
298
299       AliRsnMiniOutput *outTOF4 = task->CreateOutput(Form("config4deltatof_%s%s", name[i].Data(), suffix), output[i].Data(), comp[i].Data());
300       outTOF4->SetCutID(0, cutID9[i]);
301       outTOF4->SetCutID(1, cutID10[i]);
302       outTOF4->SetDaughter(0, AliRsnDaughter::kProton);
303       outTOF4->SetDaughter(1, AliRsnDaughter::kPion);
304       outTOF4->SetCharge(0, charge1[i]);
305       outTOF4->SetCharge(1, charge2[i]);
306       outTOF4->SetMotherPDG(ipdg[i]);
307       outTOF4->SetMotherMass(1.232);
308       outTOF4->SetPairCuts(cutsPair);
309       if (useIM[i])
310          outTOF4->AddAxis(imID, 100, 1.0, 2.0);
311       else
312          outTOF4->AddAxis(resID, 200, -0.02, 0.02);
313       outTOF4->AddAxis(ptID, 100, 0.0, 10.0);
314       outTOF4->AddAxis(centID, 120, 0.0, 120);
315
316
317       AliRsnMiniOutput *outTOF5 = task->CreateOutput(Form("config5deltatof_%s%s", name[i].Data(), suffix), output[i].Data(), comp[i].Data());
318       outTOF5->SetCutID(0, cutID11[i]);
319       outTOF5->SetCutID(1, cutID12[i]);
320       outTOF5->SetDaughter(0, AliRsnDaughter::kProton);
321       outTOF5->SetDaughter(1, AliRsnDaughter::kPion);
322       outTOF5->SetCharge(0, charge1[i]);
323       outTOF5->SetCharge(1, charge2[i]);
324       outTOF5->SetMotherPDG(ipdg[i]);
325       outTOF5->SetMotherMass(1.232);
326       outTOF5->SetPairCuts(cutsPair);
327       if (useIM[i])
328          outTOF5->AddAxis(imID, 100, 1.0, 2.0);
329       else
330          outTOF5->AddAxis(resID, 200, -0.02, 0.02);
331       outTOF5->AddAxis(ptID, 100, 0.0, 10.0);
332       outTOF5->AddAxis(centID, 120, 0.0, 120);
333
334       AliRsnMiniOutput *outTOF6 = task->CreateOutput(Form("config6deltatof_%s%s", name[i].Data(), suffix), output[i].Data(), comp[i].Data());
335       outTOF6->SetCutID(0, cutID13[i]);
336       outTOF6->SetCutID(1, cutID14[i]);
337       outTOF6->SetDaughter(0, AliRsnDaughter::kProton);
338       outTOF6->SetDaughter(1, AliRsnDaughter::kPion);
339       outTOF6->SetCharge(0, charge1[i]);
340       outTOF6->SetCharge(1, charge2[i]);
341       outTOF6->SetMotherPDG(ipdg[i]);
342       outTOF6->SetMotherMass(1.232);
343       outTOF6->SetPairCuts(cutsPair);
344       if (useIM[i])
345          outTOF6->AddAxis(imID, 100, 1.0, 2.0);
346       else
347          outTOF6->AddAxis(resID, 200, -0.02, 0.02);
348       outTOF6->AddAxis(ptID, 100, 0.0, 10.0);
349       outTOF6->AddAxis(centID, 120, 0.0, 120);
350
351    }
352    
353    return kTRUE;
354 }