]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/RESONANCES/macros/mini/AddMonitorOutput.C
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGLF / RESONANCES / macros / mini / AddMonitorOutput.C
1 /***************************************************************************
2               fbellini@cern.ch - last modified on 09/04/2013
3
4  Macro to add monitoring histograms for track and PID cuts to the rsn task 
5  Tuned for monitoring TOF KStar analysis of PbPb 2010 data
6
7 Options ("opt" argument):
8 - dim1 --> use TH1 only (no p or pt dependence)
9 - dim3 --> use TH3 (p dependence) for TOF vs TPC pid histos
10 - NoTOFSIGMA  --> disable the nsigma_TOF vs p histos
11 - NoTPCSIGMA  --> disable the nsigma_TPC vs p_TPC histos
12 - NoSIGN --> disable splitting for single track charge
13 - NoTrackQ  --> disable track quality monitoring (DCA, nCls) histos
14 /***************************************************************************/
15
16 void AddMonitorOutput(Bool_t useMCMon = 0, TObjArray *mon=0,TString opt="",AliRsnLoopDaughter *lm=0)
17 {
18   //Set options
19   Bool_t useTH1 = opt.Contains("dim1");
20   Bool_t useTH3 = opt.Contains("dim3");
21   Bool_t monitorSign = !opt.Contains("NoSIGN");
22   Bool_t monitorTrackQ = !opt.Contains("NoTrackQ");
23   Bool_t monitorTPCpid = !opt.Contains("NoTPCSIGMA");
24   Bool_t monitorTOFpid = !opt.Contains("NoTPCSIGMA");
25   
26   // Multiplicity/centrality
27   AliRsnValueEvent *multi = new AliRsnValueEvent("multi",AliRsnValueEvent::kMult);
28   multi->SetBins(0.0, 400.0, 1);
29   // Momentum
30   AliRsnValueDaughter *axisMomTPC = new AliRsnValueDaughter("pTPC", AliRsnValueDaughter::kPtpc);
31   axisMomTPC->SetBins(0.0, 10.0, 0.02);
32   AliRsnValueDaughter *axisMomP = new AliRsnValueDaughter("p", AliRsnValueDaughter::kP);
33   axisMomP->SetBins(0.0, 10.0, 0.02);
34   // Momentum Pt
35   AliRsnValueDaughter *axisMomPt = new AliRsnValueDaughter("pt", AliRsnValueDaughter::kPt);
36   axisMomPt->SetBins(0.0,10.0,0.02);
37   // Eta
38   AliRsnValueDaughter *axisMomEta = new AliRsnValueDaughter("eta", AliRsnValueDaughter::kEta);
39   axisMomEta->SetBins(-1.0, 1.0, 0.1);
40   //ITS clusters
41   AliRsnValueDaughter *axisITScls = new AliRsnValueDaughter("ITScls", AliRsnValueDaughter::kNITSclusters);
42   axisITScls->SetBins(0.0, 12.0, 1.0);
43   //TPC clusters
44   AliRsnValueDaughter *axisTPCcls = new AliRsnValueDaughter("TPCcls", AliRsnValueDaughter::kNTPCclusters);
45   axisTPCcls->SetBins(0.0, 160.0, 1.0);
46   //TPC crossed rows
47   AliRsnValueDaughter *axisTPCcrossedRows = new AliRsnValueDaughter("TPCcrossedRows", AliRsnValueDaughter::kNTPCcrossedRows);
48   axisTPCcrossedRows->SetBins(0.0, 160.0, 1.0);
49   //TPC crossed rows / findable clusters
50   AliRsnValueDaughter *axisTPCcrossedRows2Fcls = new AliRsnValueDaughter("TPCcrossedRows2Fcls", AliRsnValueDaughter::kNTPCcrossedRowsFclusters);
51   axisTPCcrossedRows2Fcls->SetBins(0.0, 1.0, 0.1);
52   //ITS chi2
53   AliRsnValueDaughter *axisITSchi2 = new AliRsnValueDaughter("ITSchi2", AliRsnValueDaughter::kITSchi2);
54   axisITSchi2->SetBins(0.0, 12.0, 0.1);
55   //TPC chi2
56   AliRsnValueDaughter *axisTPCchi2 = new AliRsnValueDaughter("TPCchi2", AliRsnValueDaughter::kTPCchi2);
57   axisTPCchi2->SetBins(0.0, 10.0, 0.1);
58   //DCA xy
59   AliRsnValueDaughter *axisDCAxy = new AliRsnValueDaughter("DCAxy", AliRsnValueDaughter::kDCAXY);
60   axisDCAxy->SetBins(-2.5, 2.5, 0.001);
61   //DCA z
62   AliRsnValueDaughter *axisDCAz = new AliRsnValueDaughter("DCAz", AliRsnValueDaughter::kDCAZ);
63   axisDCAz->SetBins(-5.0, 5.0, 0.01);
64   //Charge
65   AliRsnValueDaughter *axisCharge = new AliRsnValueDaughter("charge",AliRsnValueDaughter::kCharge);
66   axisCharge->SetBins(-1.5, 1.5, 1.0);
67   //Phi
68   AliRsnValueDaughter *axisPhi = new AliRsnValueDaughter("phi", AliRsnValueDaughter::kPhi);
69   axisPhi->SetBins(0.0, 360.0, 1.0);
70   AliRsnValueDaughter *axisPhiOuterTPC = new AliRsnValueDaughter("phiOuterTPC", AliRsnValueDaughter::kPhiOuterTPC);
71   axisPhiOuterTPC->SetBins(0.0, 360.0, 1.0);
72   
73   /* dEdx tpc */
74   AliRsnValueDaughter *axisSigTPC = new AliRsnValueDaughter("sTPC", AliRsnValueDaughter::kTPCsignal);
75   axisSigTPC->SetBins(0.0, 500.0, 2.0);
76   // kTPCnsigmaPi
77   AliRsnValueDaughter *axisTPCnsigmaPi = new AliRsnValueDaughter("pi", AliRsnValueDaughter::kTPCnsigmaPi);
78   axisTPCnsigmaPi->SetBins(-10.,10., 0.1);
79   // kTPCnsigmaK
80   AliRsnValueDaughter *axisTPCnsigmaK = new AliRsnValueDaughter("K", AliRsnValueDaughter::kTPCnsigmaK);
81   axisTPCnsigmaK->SetBins(-10.,10., 0.1);
82   // kTPCnsigmaP
83   AliRsnValueDaughter *axisTPCnsigmaP = new AliRsnValueDaughter("p", AliRsnValueDaughter::kTPCnsigmaP);
84   axisTPCnsigmaP->SetBins(-10.,10., 0.1);
85
86   /* tof signal - time */
87   AliRsnValueDaughter *axisSigTOF = new AliRsnValueDaughter("sTOF", AliRsnValueDaughter::kTOFsignal);
88   axisSigTOF->SetBins(10000.0, 50000.0, 250);//in ps
89   
90   // kTOFdeltaPi
91   AliRsnValueDaughter *axisTOFdeltaPi = new AliRsnValueDaughter("Dpi", AliRsnValueDaughter::kTOFdeltaPi);
92   axisTOFdeltaPi->SetBins(-3000.,3000., 10.);
93   // kTOFdeltaK
94   AliRsnValueDaughter *axisTOFdeltaK = new AliRsnValueDaughter("DK", AliRsnValueDaughter::kTOFdeltaK);
95   axisTOFdeltaK->SetBins(-3000.,3000., 10.);
96   // kTOFdeltaP
97   AliRsnValueDaughter *axisTOFdeltaP = new AliRsnValueDaughter("Dp", AliRsnValueDaughter::kTOFdeltaP);
98   axisTOFdeltaP->SetBins(-3000.,3000., 10.);
99
100   // kTOFnsigmaPi
101   AliRsnValueDaughter *axisTOFnsigmaPi = new AliRsnValueDaughter("pi", AliRsnValueDaughter::kTOFnsigmaPi);
102   axisTOFnsigmaPi->SetBins(-10.,10., 0.1);
103   // kTOFnsigmaK
104   AliRsnValueDaughter *axisTOFnsigmaK = new AliRsnValueDaughter("K", AliRsnValueDaughter::kTOFnsigmaK);
105   axisTOFnsigmaK->SetBins(-10.,10., 0.1);
106   // kTOFnsigmaP
107   AliRsnValueDaughter *axisTOFnsigmaP = new AliRsnValueDaughter("p", AliRsnValueDaughter::kTOFnsigmaP);
108   axisTOFnsigmaP->SetBins(-10.,10., 0.1);
109    
110   /****************************************************************/
111   /***************       MONITOR AOB           ********************/
112   /****************************************************************/
113   AliRsnListOutput *outMonitorPTvsMult = new AliRsnListOutput("PtVsMult",AliRsnListOutput::kHistoDefault);
114   outMonitorPTvsMult->AddValue(axisMomPt);
115   outMonitorPTvsMult->AddValue(multi);
116   if (mon) mon->Add(outMonitorPTvsMult);
117   if (lm) lm->AddOutput(outMonitorPTvsMult);
118
119   //    if (lm) lm->SetTrueMC(kTRUE);
120   /****************************************************************/
121   /***************       MONITOR kinematics    ********************/
122   /****************************************************************/
123   // output: TH1D for momentum
124   AliRsnListOutput *outMonitorP = new AliRsnListOutput("P", AliRsnListOutput::kHistoDefault);
125   outMonitorP->AddValue(axisMomP);
126   if (monitorSign) outMonitorP->AddValue(axisCharge);
127   if (mon) mon->Add(outMonitorP);
128   if (lm) lm->AddOutput(outMonitorP);
129   
130   // output:  TH1D for pt
131   AliRsnListOutput *outMonitorPt = new AliRsnListOutput("Pt", AliRsnListOutput::kHistoDefault);
132   outMonitorPt->AddValue(axisMomPt);
133   if (monitorSign) outMonitorPt->AddValue(axisCharge);
134   if (mon) mon->Add(outMonitorPt);
135   if (lm) lm->AddOutput(outMonitorPt);
136   
137   // output: TH1D for pseudorapidity
138   AliRsnListOutput *outMonitorEta = new AliRsnListOutput("Eta", AliRsnListOutput::kHistoDefault);
139   outMonitorEta->AddValue(axisMomEta);
140   if (monitorSign) outMonitorEta->AddValue(axisCharge);
141   if (mon) mon->Add(outMonitorEta);
142   if (lm) lm->AddOutput(outMonitorEta);
143   
144   // output:  TH1D for phi at vertex
145   AliRsnListOutput *outMonitorPhi = new AliRsnListOutput("Phi", AliRsnListOutput::kHistoDefault);
146   outMonitorPhi->AddValue(axisPhi);
147   if (monitorSign) outMonitorPhi->AddValue(axisCharge);
148   if (mon) mon->Add(outMonitorPhi);
149   if (lm) lm->AddOutput(outMonitorPhi);
150   
151   // output:  TH1D for phiOuterTPC at TPC outer radius
152   AliRsnListOutput *outMonitorPhiOuterTPC = new AliRsnListOutput("PhiOuterTPC", AliRsnListOutput::kHistoDefault);
153   outMonitorPhiOuterTPC->AddValue(axisPhiOuterTPC);
154   if (monitorSign) outMonitorPhiOuterTPC->AddValue(axisCharge);
155   if (mon) mon->Add(outMonitorPhiOuterTPC);
156   if (lm) lm->AddOutput(outMonitorPhiOuterTPC);
157   
158   // output:  TH2D for phi vs pt
159   AliRsnListOutput *outMonitorPhiVsPt = new AliRsnListOutput("PhiVsPt", AliRsnListOutput::kHistoDefault);
160   outMonitorPhiVsPt->AddValue(axisMomPt);
161   outMonitorPhiVsPt->AddValue(axisPhi);
162   if (monitorSign) outMonitorPhiVsPt->AddValue(axisCharge);
163   if (mon) mon->Add(outMonitorPhiVsPt);
164   if (lm) lm->AddOutput(outMonitorPhiVsPt);
165
166   /****************************************************************/
167   /***************      MONITOR TRACK QUALITY  ********************/
168   /****************************************************************/
169   if (monitorTrackQ) {
170     // output: 2D histogram of DCAxy vs pt
171     AliRsnListOutput *outMonitorDCAxy = new AliRsnListOutput("DCAxyVsPt", AliRsnListOutput::kHistoDefault);
172     if (!useTH1) outMonitorDCAxy->AddValue(axisMomPt);
173     outMonitorDCAxy->AddValue(axisDCAxy);
174     if (mon) mon->Add(outMonitorDCAxy);
175     if (lm) lm->AddOutput(outMonitorDCAxy);
176
177     // output: 2D histogram of DCAz vs P
178     AliRsnListOutput *outMonitorDCAz = new AliRsnListOutput("DCAzVsP", AliRsnListOutput::kHistoDefault);
179     if (!useTH1) outMonitorDCAz->AddValue(axisMomP);
180     outMonitorDCAz->AddValue(axisDCAz);
181     if (mon) mon->Add(outMonitorDCAz);
182     if (lm) lm->AddOutput(outMonitorDCAz);
183
184     // output: 2D histogram of ITS cls vs pt
185     AliRsnListOutput *outMonitorITScls = new AliRsnListOutput("ITSclsVsPt", AliRsnListOutput::kHistoDefault);
186     if (!useTH1) outMonitorITScls->AddValue(axisMomPt);
187     outMonitorITScls->AddValue(axisITScls);
188     if (mon) mon->Add(outMonitorITScls);
189     if (lm) lm->AddOutput(outMonitorITScls);
190
191     // output: 2D histogram of TPC cls vs. pt
192     AliRsnListOutput *outMonitorTPCcls = new AliRsnListOutput("TPCclsVsPt", AliRsnListOutput::kHistoDefault);
193     if (!useTH1) outMonitorTPCcls->AddValue(axisMomPt);
194     outMonitorTPCcls->AddValue(axisTPCcls);
195     if (mon) mon->Add(outMonitorTPCcls);
196     if (lm) lm->AddOutput(outMonitorTPCcls);
197
198     // output: 2D histogram of TPC cls vs. TPC momentum
199     AliRsnListOutput *outMonitorTPCclsVsPtpc = new AliRsnListOutput("TPCclsVsPtpc", AliRsnListOutput::kHistoDefault);
200     if (!useTH1) outMonitorTPCclsVsPtpc->AddValue(axisMomTPC);
201     outMonitorTPCclsVsPtpc->AddValue(axisTPCcls);
202     if (mon) mon->Add(outMonitorTPCclsVsPtpc);
203     if (lm) lm->AddOutput(outMonitorTPCclsVsPtpc);
204
205     // output: 2D histogram of TPC crossed rows vs. TPC momentum
206     AliRsnListOutput *outMonitorTPCcrossedRowsVsPtpc = new AliRsnListOutput("TPCcrossedRowsVsPtpc", AliRsnListOutput::kHistoDefault);
207     if (!useTH1) outMonitorTPCcrossedRowsVsPtpc->AddValue(axisMomTPC);
208     outMonitorTPCcrossedRowsVsPtpc->AddValue(axisTPCcrossedRows);
209     if (mon) mon->Add(outMonitorTPCcrossedRowsVsPtpc);
210     if (lm) lm->AddOutput(outMonitorTPCcrossedRowsVsPtpc);
211
212     // output: 2D histogram of TPC crossed rows/Findable cls vs. TPC momentum
213     AliRsnListOutput *outMonitorTPCcrossedRows2FclsVsPtpc = new AliRsnListOutput("TPCcrossedRows2FclsVsPtpc", AliRsnListOutput::kHistoDefault);
214     if (!useTH1) outMonitorTPCcrossedRows2FclsVsPtpc->AddValue(axisMomTPC);
215     outMonitorTPCcrossedRows2FclsVsPtpc->AddValue(axisTPCcrossedRows2Fcls);
216     if (mon) mon->Add(outMonitorTPCcrossedRows2FclsVsPtpc);
217     if (lm) lm->AddOutput(outMonitorTPCcrossedRows2FclsVsPtpc);
218
219     // output: 2D histogram of ITS chi2 vs pt
220     AliRsnListOutput *outMonitorITSchi2 = new AliRsnListOutput("ITSchi2VsPt", AliRsnListOutput::kHistoDefault);
221     if (!useTH1) outMonitorITSchi2->AddValue(axisMomPt);
222     outMonitorITSchi2->AddValue(axisITSchi2);
223     if (mon) mon->Add(outMonitorITSchi2);
224     if (lm) lm->AddOutput(outMonitorITSchi2);
225
226     // output: 2D histogram of TPC chi2 vs. pt
227     AliRsnListOutput *outMonitorTPCchi2 = new AliRsnListOutput("TPCchi2VsPt", AliRsnListOutput::kHistoDefault);
228     if (!useTH1) outMonitorTPCchi2->AddValue(axisMomPt);
229     outMonitorTPCchi2->AddValue(axisTPCchi2);
230     if (mon) mon->Add(outMonitorTPCchi2);
231     if (lm) lm->AddOutput(outMonitorTPCchi2);
232
233     // output: 2D histogram of TPC chi2 vs. TPC momentum
234     AliRsnListOutput *outMonitorTPCchi2VsPtpc = new AliRsnListOutput("TPCchi2VsPtpc", AliRsnListOutput::kHistoDefault);
235    if (!useTH1) outMonitorTPCchi2VsPtpc->AddValue(axisMomTPC);
236     outMonitorTPCchi2VsPtpc->AddValue(axisTPCchi2);
237     if (mon) mon->Add(outMonitorTPCchi2VsPtpc);
238     if (lm) lm->AddOutput(outMonitorTPCchi2VsPtpc);
239   }
240   /****************************************************************/
241   /***************       MONITOR TPC           ********************/
242   /****************************************************************/
243   if (monitorTPCpid) {
244     // output: 2D histogram of TPC signal vs. TPC momentum
245     AliRsnListOutput *outMonitordEdxTPC = new AliRsnListOutput("dEdx_VsPtpc", AliRsnListOutput::kHistoDefault);
246     outMonitordEdxTPC->AddValue(axisMomTPC);
247     outMonitordEdxTPC->AddValue(axisSigTPC);
248     if (mon) mon->Add(outMonitordEdxTPC);
249     if (lm) lm->AddOutput(outMonitordEdxTPC);
250
251     // output: 2D histogram of TPC nsigma pi vs. TPC momentum
252     AliRsnListOutput *outMonitorTPCnsigmaPi = new AliRsnListOutput("TPC_nsigmaPi_VsPtpc", AliRsnListOutput::kHistoDefault);
253     outMonitorTPCnsigmaPi->AddValue(axisMomTPC);
254     outMonitorTPCnsigmaPi->AddValue(axisTPCnsigmaPi);
255     if (mon) mon->Add(outMonitorTPCnsigmaPi);
256     if (lm) lm->AddOutput(outMonitorTPCnsigmaPi);
257
258     // output: 2D histogram of TPC nsigma K vs. TPC momentum
259     AliRsnListOutput *outMonitorTPCnsigmaK = new AliRsnListOutput("TPC_nsigmaK_VsPtpc", AliRsnListOutput::kHistoDefault);
260     outMonitorTPCnsigmaK->AddValue(axisMomTPC);
261     outMonitorTPCnsigmaK->AddValue(axisTPCnsigmaK);
262     if (mon) mon->Add(outMonitorTPCnsigmaK);
263     if (lm) lm->AddOutput(outMonitorTPCnsigmaK);
264
265     // output: 2D histogram of TPC nsigma pro vs. TPC momentum
266     AliRsnListOutput *outMonitorTPCnsigmaP = new AliRsnListOutput("TPC_nsigmaPro_VsPtpc", AliRsnListOutput::kHistoDefault);
267     outMonitorTPCnsigmaP->AddValue(axisMomTPC);
268     outMonitorTPCnsigmaP->AddValue(axisTPCnsigmaP);
269     if (mon) mon->Add(outMonitorTPCnsigmaP);
270     if (lm) lm->AddOutput(outMonitorTPCnsigmaP);
271   }
272   /****************************************************************/
273   /***************       MONITOR TOF           ********************/
274   /****************************************************************/
275   // output:2D histogram of TOF Nsigma pi vs. TPC Nsigma pi vs momentum
276   AliRsnListOutput *outMonitorTOFvsTPCnsigmaPi = new AliRsnListOutput("TOFnsigmaPi_TPCnsigmaPi", AliRsnListOutput::kHistoDefault);
277   outMonitorTOFvsTPCnsigmaPi->AddValue(axisTOFnsigmaPi);
278   outMonitorTOFvsTPCnsigmaPi->AddValue(axisTPCnsigmaPi);
279   if (useTH3) outMonitorTOFvsTPCnsigmaPi->AddValue(axisMomP);
280   if (mon) mon->Add(outMonitorTOFvsTPCnsigmaPi);
281   if (lm) lm->AddOutput(outMonitorTOFvsTPCnsigmaPi);
282   
283   // output:2D histogram of TOF Nsigma pi vs. TPC Nsigma pi vs momentum
284   AliRsnListOutput *outMonitorTOFvsTPCnsigmaK = new AliRsnListOutput("TOFnsigmaK_TPCnsigmaK", AliRsnListOutput::kHistoDefault);
285   outMonitorTOFvsTPCnsigmaK->AddValue(axisTOFnsigmaK);
286   outMonitorTOFvsTPCnsigmaK->AddValue(axisTPCnsigmaK);
287   if (useTH3) outMonitorTOFvsTPCnsigmaK->AddValue(axisMomP);
288   if (mon) mon->Add(outMonitorTOFvsTPCnsigmaK);
289   if (lm) lm->AddOutput(outMonitorTOFvsTPCnsigmaK);
290
291   AliRsnListOutput *outMonitorTOFvsTPCnsigmaP = new AliRsnListOutput("TOFnsigmaP_TPCnsigmaP", AliRsnListOutput::kHistoDefault);
292   outMonitorTOFvsTPCnsigmaP->AddValue(axisTOFnsigmaP);
293   outMonitorTOFvsTPCnsigmaP->AddValue(axisTPCnsigmaP);
294   if (useTH3) outMonitorTOFvsTPCnsigmaP->AddValue(axisMomP);
295   if (mon) mon->Add(outMonitorTOFvsTPCnsigmaP);
296   if (lm) lm->AddOutput(outMonitorTOFvsTPCnsigmaP);
297
298   // // output: 2D histogram of TOF signal vs. momentum
299   // AliRsnListOutput *outMonitorTimeTOF = new AliRsnListOutput("time_VsP", AliRsnListOutput::kHistoDefault);
300   // outMonitorTimeTOF->AddValue(axisMomP);
301   // outMonitorTimeTOF->AddValue(axisSigTOF);
302   // if (mon) mon->Add(outMonitorTimeTOF);
303   // if (lm) lm->AddOutput(outMonitorTimeTOF);
304
305   // // output: 2D histogram of TOF signal vs. pt
306   // AliRsnListOutput *outMonitorTimeTOFPt = new AliRsnListOutput("time_VsPt", AliRsnListOutput::kHistoDefault);
307   // outMonitorTimeTOFPt->AddValue(axisMomPt);
308   // outMonitorTimeTOFPt->AddValue(axisSigTOF);
309   // if (mon) mon->Add(outMonitorTimeTOFPt);
310   // if (lm) lm->AddOutput(outMonitorTimeTOFPt);
311
312   if (monitorTOFpid) {
313     // output: 2D histogram of TOF Nsigma pi vs. TPC momentum
314     AliRsnListOutput *outMonitorTOFnsigmaPi = new AliRsnListOutput("TOF_nsigmaPi_vsP", AliRsnListOutput::kHistoDefault);
315     outMonitorTOFnsigmaPi->AddValue(axisMomP);
316     outMonitorTOFnsigmaPi->AddValue(axisTOFnsigmaPi);
317     if (mon) mon->Add(outMonitorTOFnsigmaPi);
318     if (lm) lm->AddOutput(outMonitorTOFnsigmaPi);
319      
320     // output: 2D histogram of TOF signal vs. TOF momentum
321     AliRsnListOutput *outMonitorTOFnsigmaK = new AliRsnListOutput("TOF_nsigmaK_vsP", AliRsnListOutput::kHistoDefault);
322     outMonitorTOFnsigmaK->AddValue(axisMomP);
323     outMonitorTOFnsigmaK->AddValue(axisTOFnsigmaK);
324     if (mon) mon->Add(outMonitorTOFnsigmaK);
325     if (lm) lm->AddOutput(outMonitorTOFnsigmaK);
326       
327     // output: 2D histogram of TOF signal vs. TOF momentum
328     AliRsnListOutput *outMonitorTOFnsigmaP = new AliRsnListOutput("TOF_nsigmaPro_vsP", AliRsnListOutput::kHistoDefault);
329     outMonitorTOFnsigmaP->AddValue(axisMomP);
330     outMonitorTOFnsigmaP->AddValue(axisTOFnsigmaP);
331     if (mon) mon->Add(outMonitorTOFnsigmaP);
332     if (lm) lm->AddOutput(outMonitorTOFnsigmaP);
333
334    // output: 2D histogram of TOF Delta pi vs. TPC momentum
335     AliRsnListOutput *outMonitorTOFdeltaPi = new AliRsnListOutput("TOF_deltaPi_vsP", AliRsnListOutput::kHistoDefault);
336     outMonitorTOFdeltaPi->AddValue(axisMomP);
337     outMonitorTOFdeltaPi->AddValue(axisTOFdeltaPi);
338     if (mon) mon->Add(outMonitorTOFdeltaPi);
339     if (lm) lm->AddOutput(outMonitorTOFdeltaPi);
340      
341     // output: 2D histogram of TOF signal vs. TOF momentum
342     AliRsnListOutput *outMonitorTOFdeltaK = new AliRsnListOutput("TOF_deltaK_vsP", AliRsnListOutput::kHistoDefault);
343     outMonitorTOFdeltaK->AddValue(axisMomP);
344     outMonitorTOFdeltaK->AddValue(axisTOFdeltaK);
345     if (mon) mon->Add(outMonitorTOFdeltaK);
346     if (lm) lm->AddOutput(outMonitorTOFdeltaK);
347       
348     // output: 2D histogram of TOF signal vs. TOF momentum
349     AliRsnListOutput *outMonitorTOFdeltaP = new AliRsnListOutput("TOF_deltaPro_vsP", AliRsnListOutput::kHistoDefault);
350     outMonitorTOFdeltaP->AddValue(axisMomP);
351     outMonitorTOFdeltaP->AddValue(axisTOFdeltaP);
352     if (mon) mon->Add(outMonitorTOFdeltaP);
353     if (lm) lm->AddOutput(outMonitorTOFdeltaP);
354   }
355
356   /****************************************************************/
357   /***************       MONITOR MC            ********************/
358   /****************************************************************/
359   if (useMCMon) {
360     //Momentum 
361     AliRsnValueDaughter *axisMomPMC = new AliRsnValueDaughter("pMC", AliRsnValueDaughter::kP);
362     axisMomPMC->SetUseMCInfo(kTRUE);
363     axisMomPMC->SetBins(0.0,10.0,0.01);
364     // Momentum Pt
365     AliRsnValueDaughter *axisMomPtMC = new AliRsnValueDaughter("ptMC", AliRsnValueDaughter::kPt);
366     axisMomPtMC->SetUseMCInfo(kTRUE);
367     axisMomPtMC->SetBins(0.0,10.0,0.01);
368     // Eta
369     AliRsnValueDaughter *axisEtaMC = new AliRsnValueDaughter("etaMC", AliRsnValueDaughter::kEta);
370     axisEtaMC->SetUseMCInfo(kTRUE);
371     axisEtaMC->SetBins(-1.0,1.0,0.1);
372
373     // output: 2D histo for kine acceptance
374     AliRsnListOutput *outMonitorPtVsEtaMC = new AliRsnListOutput("Pt_VsEtaMC", AliRsnListOutput::kHistoDefault);
375     outMonitorPtVsEtaMC->AddValue(axisMomPtMC);
376     outMonitorPtVsEtaMC->AddValue(axisEtaMC);
377     if (mon) mon->Add(outMonitorPtVsEtaMC);
378     if (lm) lm->AddOutput(outMonitorPtVsEtaMC);
379
380     // output: 1D histo pt from MC
381     AliRsnListOutput *outMonitorPtMC = new AliRsnListOutput("PtMC", AliRsnListOutput::kHistoDefault);
382     outMonitorPtMC->AddValue(axisMomPtMC);
383     if (mon) mon->Add(outMonitorPtMC);
384     if (lm) lm->AddOutput(outMonitorPtMC);
385  
386     // output: 1D histo eta from MC
387     AliRsnListOutput *outMonitorEtaMC = new AliRsnListOutput("EtaMC", AliRsnListOutput::kHistoDefault);
388     outMonitorEtaMC->AddValue(axisEtaMC);
389     if (mon) mon->Add(outMonitorEtaMC);
390     if (lm) lm->AddOutput(outMonitorEtaMC);
391
392     // output: 1D histo pt from MC
393     AliRsnListOutput *outMonitorPtMC = new AliRsnListOutput("PtMC", AliRsnListOutput::kHistoDefault);
394     outMonitorPtMC->AddValue(axisMomPtMC);
395     if (mon) mon->Add(outMonitorPtMC);
396     if (lm) lm->AddOutput(outMonitorPtMC);
397   }
398  
399 }