1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
23 #include "THnSparse.h"
40 #include "AliAnalysisManager.h"
42 #include "AliAODHandler.h"
43 #include "AliAODInputHandler.h"
44 #include "AliInputEventHandler.h"
46 #include "AliESDEvent.h"
47 #include "AliESDInputHandler.h"
48 #include "AliMultiplicity.h"
49 #include "AliCentrality.h"
50 #include "AliAnalysisTaskDptDptCorrelations.h"
53 #include "AliPIDResponse.h"
55 #include "AliESDVertex.h"
56 #include "AliESDEvent.h"
57 #include "AliESDInputHandler.h"
58 #include "AliAODEvent.h"
59 #include "AliAODTrack.h"
60 #include "AliAODInputHandler.h"
62 #include "AliESDEvent.h"
63 #include "AliAODEvent.h"
65 #include "AliESDtrackCuts.h"
66 #include "AliAODMCHeader.h"
69 #include "AliGenHijingEventHeader.h"
70 #include "AliGenEventHeader.h"
72 #include "AliAODPid.h"
73 #include "AliPIDResponse.h"
74 #include "AliAODpidUtil.h"
75 #include "AliPIDCombined.h"
79 ClassImp(AliAnalysisTaskDptDptCorrelations)
81 AliAnalysisTaskDptDptCorrelations::AliAnalysisTaskDptDptCorrelations()
82 : AliAnalysisTaskSE(),
84 fESDEvent(0), //! ESD Event
88 _twoPi ( 2.0 * 3.1415927),
95 _rejectPairConversion ( 0),
100 _centralityMethod ( 4),
101 _centralityMin ( 0.),
102 _centralityMax ( 0.),
103 _requestedCharge_1 ( 1),
104 _requestedCharge_2 ( -1),
145 _correctionWeight_1(0),
146 _correctionWeight_2(0),
147 _nBins_M0(500), _min_M0(0), _max_M0(10000), _width_M0(20),
148 _nBins_M1(500), _min_M1(0), _max_M1(10000), _width_M1(20),
149 _nBins_M2(500), _min_M2(0), _max_M2(10000), _width_M2(20),
150 _nBins_M3(500), _min_M3(0), _max_M3(10000), _width_M3(20),
151 _nBins_M4(100), _min_M4(0), _max_M4(1), _width_M4(0.01),
152 _nBins_M5(100), _min_M5(0), _max_M5(1), _width_M5(0.01),
153 _nBins_M6(100), _min_M6(0), _max_M6(1), _width_M6(0.01),
154 _nBins_vertexZ(40), _min_vertexZ(-10), _max_vertexZ(10), _width_vertexZ(0.5),
156 _nBins_pt_1(18), _min_pt_1(0.2), _max_pt_1(2.0), _width_pt_1(0.1),
157 _nBins_phi_1(72), _min_phi_1(0), _max_phi_1(2.*3.1415927),_width_phi_1(2.*3.1415927/72.),
158 _nBins_eta_1(0), _min_eta_1(0), _max_eta_1(0), _width_eta_1(0.1),
161 _nBins_etaPhiPt_1(0),
162 _nBins_zEtaPhiPt_1(0),
164 _nBins_pt_2(18), _min_pt_2(0.2), _max_pt_2(2.0), _width_pt_2(0.1),
165 _nBins_phi_2(72), _min_phi_2(0), _max_phi_2(2.*3.1415927),_width_phi_2(2.*3.1415927/72),
166 _nBins_eta_2(0), _min_eta_2(0), _max_eta_2(0), _width_eta_2(0.1),
169 _nBins_etaPhiPt_2(0),
170 _nBins_zEtaPhiPt_2(0),
190 __s1pt_1_vsEtaPhi(0),
191 __n1_1_vsZEtaPhiPt(0),
194 __s1pt_2_vsEtaPhi(0),
195 __n1_2_vsZEtaPhiPt(0),
198 __s2ptpt_12_vsEtaPhi(0),
199 __s2PtN_12_vsEtaPhi(0),
200 __s2NPt_12_vsEtaPhi(0),
203 _eventAccounting ( 0),
220 _n1_1_vsEtaVsPhi ( 0),
221 _s1pt_1_vsEtaVsPhi ( 0),
222 _n1_1_vsZVsEtaVsPhiVsPt ( 0),
231 _n1_2_vsEtaVsPhi ( 0),
232 _s1pt_2_vsEtaVsPhi ( 0),
233 _n1_2_vsZVsEtaVsPhiVsPt ( 0),
241 _n2_12_vsEtaPhi ( 0),
242 _n2_12_vsPtVsPt ( 0),
243 _s2PtPt_12_vsEtaPhi( 0),
244 _s2PtN_12_vsEtaPhi ( 0),
245 _s2NPt_12_vsEtaPhi ( 0),
251 _s2PtPtNw_12_vsM ( 0),
252 _s2PtNNw_12_vsM ( 0),
253 _s2NPtNw_12_vsM ( 0),
276 intBinCorrName("NA"),
326 _title_etaPhi_1("NA"),
328 _title_SumPt_1("NA"),
329 _title_AvgPt_1("NA"),
331 _title_AvgSumPt_1("NA"),
336 _title_etaPhi_2("NA"),
338 _title_SumPt_2("NA"),
339 _title_AvgPt_2("NA"),
341 _title_AvgSumPt_2("NA"),
343 _title_etaPhi_12("NA"),
345 _title_AvgN2_12("NA"),
346 _title_AvgSumPtPt_12("NA"),
347 _title_AvgSumPtN_12("NA"),
348 _title_AvgNSumPt_12("NA"),
358 printf("Default constructor called \n");
360 printf("passed \n ");
364 AliAnalysisTaskDptDptCorrelations::AliAnalysisTaskDptDptCorrelations(const TString & name)
365 : AliAnalysisTaskSE(name),
371 _twoPi ( 2.0 * 3.1415927),
376 _sameFilter ( false),
378 _rejectPairConversion ( 0),
381 _vertexXYMin ( -10.),
383 _centralityMethod ( 4),
384 _centralityMin ( 0.),
385 _centralityMax ( 1.),
386 _requestedCharge_1 ( 1),
387 _requestedCharge_2 ( -1),
395 _trackFilterBit ( 0),
430 _correctionWeight_1(0),
431 _correctionWeight_2(0),
432 _nBins_M0(500), _min_M0(0), _max_M0(10000), _width_M0(20),
433 _nBins_M1(500), _min_M1(0), _max_M1(10000), _width_M1(20),
434 _nBins_M2(500), _min_M2(0), _max_M2(10000), _width_M2(20),
435 _nBins_M3(500), _min_M3(0), _max_M3(10000), _width_M3(20),
436 _nBins_M4(100), _min_M4(0), _max_M4(1), _width_M4(0.01),
437 _nBins_M5(100), _min_M5(0), _max_M5(1), _width_M5(0.01),
438 _nBins_M6(100), _min_M6(0), _max_M6(1), _width_M6(0.01),
439 _nBins_vertexZ(40), _min_vertexZ(-10), _max_vertexZ(10), _width_vertexZ(0.5),
441 _nBins_pt_1(18), _min_pt_1(0.2), _max_pt_1(2.0), _width_pt_1(0.1),
442 _nBins_phi_1(72), _min_phi_1(0), _max_phi_1(2.*3.1415927),_width_phi_1(2.*3.1415927/72.),
443 _nBins_eta_1(0), _min_eta_1(0), _max_eta_1(0), _width_eta_1(0.1),
446 _nBins_etaPhiPt_1(0),
447 _nBins_zEtaPhiPt_1(0),
449 _nBins_pt_2(18), _min_pt_2(0.2), _max_pt_2(2.0), _width_pt_2(0.1),
450 _nBins_phi_2(72), _min_phi_2(0), _max_phi_2(2.*3.1415927),_width_phi_2(2.*3.1415927/72),
451 _nBins_eta_2(0), _min_eta_2(0), _max_eta_2(0), _width_eta_2(0.1),
454 _nBins_etaPhiPt_2(0),
455 _nBins_zEtaPhiPt_2(0),
475 __s1pt_1_vsEtaPhi(0),
476 __n1_1_vsZEtaPhiPt(0),
479 __s1pt_2_vsEtaPhi(0),
480 __n1_2_vsZEtaPhiPt(0),
483 __s2ptpt_12_vsEtaPhi(0),
484 __s2PtN_12_vsEtaPhi(0),
485 __s2NPt_12_vsEtaPhi(0),
488 _eventAccounting ( 0),
505 _n1_1_vsEtaVsPhi ( 0),
506 _s1pt_1_vsEtaVsPhi ( 0),
507 _n1_1_vsZVsEtaVsPhiVsPt ( 0),
516 _n1_2_vsEtaVsPhi ( 0),
517 _s1pt_2_vsEtaVsPhi ( 0),
518 _n1_2_vsZVsEtaVsPhiVsPt ( 0),
526 _n2_12_vsEtaPhi ( 0),
527 _n2_12_vsPtVsPt ( 0),
528 _s2PtPt_12_vsEtaPhi( 0),
529 _s2PtN_12_vsEtaPhi ( 0),
530 _s2NPt_12_vsEtaPhi ( 0),
536 _s2PtPtNw_12_vsM ( 0),
537 _s2PtNNw_12_vsM ( 0),
538 _s2NPtNw_12_vsM ( 0),
561 intBinCorrName("NA"),
611 _title_etaPhi_1("NA"),
613 _title_SumPt_1("NA"),
614 _title_AvgPt_1("NA"),
616 _title_AvgSumPt_1("NA"),
621 _title_etaPhi_2("NA"),
623 _title_SumPt_2("NA"),
624 _title_AvgPt_2("NA"),
626 _title_AvgSumPt_2("NA"),
628 _title_etaPhi_12("NA"),
630 _title_AvgN2_12("NA"),
631 _title_AvgSumPtPt_12("NA"),
632 _title_AvgSumPtN_12("NA"),
633 _title_AvgNSumPt_12("NA"),
643 printf("2nd constructor called ");
645 DefineOutput(0, TList::Class());
651 AliAnalysisTaskDptDptCorrelations::~AliAnalysisTaskDptDptCorrelations()
656 void AliAnalysisTaskDptDptCorrelations::UserCreateOutputObjects()
659 _outputHistoList = new TList();
660 _outputHistoList->SetOwner();
662 _nBins_M0 = 500; _min_M0 = 0.; _max_M0 = 5000.; _width_M0 = (_max_M0-_min_M0)/_nBins_M0;
663 _nBins_M1 = 500; _min_M1 = 0.; _max_M1 = 5000.; _width_M1 = (_max_M1-_min_M1)/_nBins_M1;
664 _nBins_M2 = 500; _min_M2 = 0.; _max_M2 = 5000.; _width_M2 = (_max_M2-_min_M2)/_nBins_M2;
665 _nBins_M3 = 500; _min_M3 = 0.; _max_M3 = 5000.; _width_M3 = (_max_M3-_min_M3)/_nBins_M3;
666 _nBins_M4 = 100; _min_M4 = 0.; _max_M4 = 100.; _width_M4 = (_max_M4-_min_M4)/_nBins_M4;
667 _nBins_M5 = 100; _min_M5 = 0.; _max_M5 = 100.; _width_M5 = (_max_M5-_min_M5)/_nBins_M5;
668 _nBins_M6 = 100; _min_M6 = 0.; _max_M6 = 100.; _width_M6 = (_max_M6-_min_M6)/_nBins_M6;
670 _min_vertexZ = _vertexZMin;
671 _max_vertexZ = _vertexZMax;
672 _width_vertexZ = 0.5;
673 _nBins_vertexZ = int(0.5+ (_max_vertexZ - _min_vertexZ)/_width_vertexZ);
674 _nBins_pt_1 = int(0.5+ (_max_pt_1 -_min_pt_1 )/_width_pt_1);
675 _nBins_eta_1 = int(0.5+ (_max_eta_1-_min_eta_1)/_width_eta_1);
676 _width_phi_1 = (_max_phi_1 - _min_phi_1) /_nBins_phi_1;
677 _nBins_etaPhi_1 = _nBins_phi_1 * _nBins_eta_1;
678 _nBins_etaPhiPt_1 = _nBins_etaPhi_1 * _nBins_pt_1;
679 _nBins_zEtaPhiPt_1 = _nBins_vertexZ * _nBins_etaPhiPt_1;
681 _nBins_pt_2 = int(0.5+ (_max_pt_2 -_min_pt_2 )/_width_pt_2);
682 _nBins_eta_2 = int(0.5+ (_max_eta_2-_min_eta_2)/_width_eta_2);
683 _width_phi_2 = (_max_phi_2 - _min_phi_2) /_nBins_phi_2;
684 _nBins_etaPhi_2 = _nBins_phi_2 * _nBins_eta_2;
685 _nBins_etaPhiPt_2 = _nBins_etaPhi_2 * _nBins_pt_2;
686 _nBins_zEtaPhiPt_2 = _nBins_vertexZ * _nBins_etaPhiPt_2;
687 _nBins_etaPhi_12 = _nBins_etaPhi_1 * _nBins_etaPhi_2;
689 _id_1 = new int[arraySize];
690 _charge_1 = new int[arraySize];
691 _iEtaPhi_1 = new int[arraySize];
692 _iPt_1 = new int[arraySize];
693 _pt_1 = new float[arraySize];
694 _px_1 = new float[arraySize];
695 _py_1 = new float[arraySize];
696 _pz_1 = new float[arraySize];
697 _correction_1 = new float[arraySize];
699 __n1_1_vsPt = getDoubleArray(_nBins_pt_1, 0.);
700 __n1_1_vsEtaPhi = getDoubleArray(_nBins_etaPhi_1, 0.);
701 __s1pt_1_vsEtaPhi = getDoubleArray(_nBins_etaPhi_1, 0.);
702 __n1_1_vsZEtaPhiPt = getFloatArray(_nBins_zEtaPhiPt_1, 0.);
705 if (_requestedCharge_2!=_requestedCharge_1)
709 _id_2 = new int[arraySize];
710 _charge_2 = new int[arraySize];
711 _iEtaPhi_2 = new int[arraySize];
712 _iPt_2 = new int[arraySize];
713 _pt_2 = new float[arraySize];
714 _px_2 = new float[arraySize];
715 _py_2 = new float[arraySize];
716 _pz_2 = new float[arraySize];
717 _correction_2 = new float[arraySize];
719 __n1_2_vsPt = getDoubleArray(_nBins_pt_2, 0.);
720 __n1_2_vsEtaPhi = getDoubleArray(_nBins_etaPhi_2, 0.);
721 __s1pt_2_vsEtaPhi = getDoubleArray(_nBins_etaPhi_2, 0.);
722 __n1_2_vsZEtaPhiPt = getFloatArray(_nBins_zEtaPhiPt_2, 0.);
726 __n2_12_vsPtPt = getDoubleArray(_nBins_pt_1*_nBins_pt_2,0.);
727 __n2_12_vsEtaPhi = getFloatArray(_nBins_etaPhi_12, 0.);
728 __s2ptpt_12_vsEtaPhi = getFloatArray(_nBins_etaPhi_12, 0.);
729 __s2PtN_12_vsEtaPhi = getFloatArray(_nBins_etaPhi_12, 0.);
730 __s2NPt_12_vsEtaPhi = getFloatArray(_nBins_etaPhi_12, 0.);
732 // Setup all the labels needed.
736 pair_12_Name = "_12";
753 binCorrName = "binCorr";
754 intBinCorrName = "intBinCorr";
759 s1ptNwName = "sumPtNw";
760 s1DptName = "sumDpt";
761 s2PtPtName = "sumPtPt";
762 s2PtPtNwName = "sumPtPtNw";
763 s2DptDptName = "sumDptDpt";
764 s2NPtName = "sumNPt";
765 s2NPtNwName = "sumNPtNw";
766 s2PtNName = "sumPtN";
767 s2NPtNwName = "sumNPtNw";
768 s2PtNNwName = "sumPtNNw";
770 ptptName = "avgPtPt";
771 pt1pt1Name = "avgPtavgPt";
773 DptDptName = "avgDptDpt";
774 RDptDptName = "relDptDpt"; // ratio of avgDptDpt by avgPt*avgPt
779 _title_counts = "yield";
785 _title_m4 = "V0Centrality";
786 _title_m5 = "TrkCentrality";
787 _title_m6 = "SpdCentrality";
789 _title_eta_1 = "#eta_{1}";
790 _title_phi_1 = "#varphi_{1} (radian)";
791 _title_etaPhi_1 = "#eta_{1}#times#varphi_{1}";
792 _title_pt_1 = "p_{t,1} (GeV/c)";
793 _title_n_1 = "n_{1}";
794 _title_SumPt_1 = "#Sigma p_{t,1} (GeV/c)";
795 _title_AvgPt_1 = "#LT p_{t,1} #GT (GeV/c)";
796 _title_AvgN_1 = "#LT n_{1} #GT";
797 _title_AvgSumPt_1 = "#LT #Sigma p_{t,1} #GT (GeV/c)";
799 _title_eta_2 = "#eta_{2}";
800 _title_phi_2 = "#varphi_{2} (radian)";
801 _title_etaPhi_2 = "#eta_{2}#times#varphi_{2}";
802 _title_pt_2 = "p_{t,2} (GeV/c)";
803 _title_n_2 = "n_{2}";
804 _title_SumPt_2 = "#Sigma p_{t,1} (GeV/c)";
805 _title_AvgPt_2 = "#LT p_{t,2} #GT (GeV/c)";
806 _title_AvgN_2 = "#LT n_{2} #GT";
807 _title_AvgSumPt_2 = "#LT #Sigma p_{t,2} #GT (GeV/c)";
809 _title_etaPhi_12 = "#eta_{1}#times#varphi_{1}#times#eta_{2}#times#varphi_{2}";
811 _title_AvgN2_12 = "#LT n_{2} #GT";;
812 _title_AvgSumPtPt_12 = "#LT #Sigma p_{t,1}p_{t,2} #GT";;
813 _title_AvgSumPtN_12 = "#LT #Sigma p_{t,1}N #GT";;
814 _title_AvgNSumPt_12 = "#LT N#Sigma p_{t,2} #GT";;
822 vsEtaPhi = "_vsEtaPhi";
823 vsPtVsPt = "_vsPtVsPt";
828 int iZ, iEtaPhi, iPt;
829 int iZ1,iEtaPhi1,iPt1;
833 _correctionWeight_1 = new float[_nBins_vertexZ*_nBins_etaPhi_1*_nBins_pt_1];
835 b = _nBins_etaPhi_1*_nBins_pt_1;
836 for (iZ=0,iZ1=1; iZ<_nBins_vertexZ; iZ++, iZ1++)
838 for (iEtaPhi=0,iEtaPhi1=1; iEtaPhi<_nBins_etaPhi_1; iEtaPhi++, iEtaPhi1++)
840 for (iPt=0,iPt1=1; iPt<_nBins_pt_1; iPt++, iPt1++)
842 _correctionWeight_1[iZ*b+iEtaPhi*a+iPt] = _weight_1->GetBinContent(iZ1,iEtaPhi1,iPt1);
849 AliError("AliAnalysisTaskDptDptCorrelations:: _weight_1 is a null pointer.");
856 _correctionWeight_2 = new float[_nBins_vertexZ*_nBins_etaPhi_2*_nBins_pt_2];
858 b = _nBins_etaPhi_2*_nBins_pt_2;
859 for (iZ=0,iZ1=1; iZ<_nBins_vertexZ; iZ++, iZ1++)
861 for (iEtaPhi=0,iEtaPhi1=1; iEtaPhi<_nBins_etaPhi_2; iEtaPhi++, iEtaPhi1++)
863 for (iPt=0,iPt1=1; iPt<_nBins_pt_2; iPt++, iPt1++)
865 _correctionWeight_2[iZ*b+iEtaPhi*a+iPt] = _weight_2->GetBinContent(iZ1,iEtaPhi1,iPt1);
872 AliError("AliAnalysisTaskDptDptCorrelations:: _weight_1 is a null pointer.");
879 PostData(0,_outputHistoList);
881 //cout<< "AliAnalysisTaskDptDptCorrelations::CreateOutputObjects() DONE " << endl;
885 void AliAnalysisTaskDptDptCorrelations::createHistograms()
887 AliInfo(" AliAnalysisTaskDptDptCorrelations::createHistoHistograms() Creating Event Histos");
890 name = "eventAccounting";
892 _eventAccounting = createHisto1D(name,name,10, -0.5, 9.5, "event Code", _title_counts);
894 name = "m0"; _m0 = createHisto1D(name,name,_nBins_M1, _min_M1, _max_M1, _title_m0, _title_counts);
895 name = "m1"; _m1 = createHisto1D(name,name,_nBins_M1, _min_M1, _max_M1, _title_m1, _title_counts);
896 name = "m2"; _m2 = createHisto1D(name,name,_nBins_M2, _min_M2, _max_M2, _title_m2, _title_counts);
897 name = "m3"; _m3 = createHisto1D(name,name,_nBins_M3, _min_M3, _max_M3, _title_m3, _title_counts);
898 name = "m4"; _m4 = createHisto1D(name,name,_nBins_M4, _min_M4, _max_M4, _title_m4, _title_counts);
899 name = "m5"; _m5 = createHisto1D(name,name,_nBins_M5, _min_M5, _max_M5, _title_m5, _title_counts);
900 name = "m6"; _m6 = createHisto1D(name,name,_nBins_M6, _min_M6, _max_M6, _title_m6, _title_counts);
901 name = "zV"; _vertexZ = createHisto1D(name,name,100, -10, 10, "z-Vertex (cm)", _title_counts);
904 name = "Eta"; _etadis = createHisto1F(name,name, 200, -1.0, 1.0, "#eta","counts");
905 name = "Phi"; _phidis = createHisto1F(name,name, 360, 0.0, 6.4, "#phi","counts");
906 name = "DCAz"; _dcaz = createHisto1F(name,name, 500, -5.0, 5.0, "dcaZ","counts");
907 name = "DCAxy"; _dcaxy = createHisto1F(name,name, 500, -5.0, 5.0, "dcaXY","counts");
909 //name = "Nclus1"; _Ncluster1 = createHisto1F(name,name, 200, 0, 200, "Ncluster1","counts");
910 //name = "Nclus2"; _Ncluster2 = createHisto1F(name,name, 200, 0, 200, "Ncluster2","counts");
914 name = n1Name+part_1_Name+vsPt; _n1_1_vsPt = createHisto1F(name,name, _nBins_pt_1, _min_pt_1, _max_pt_1, _title_pt_1, _title_AvgN_1);
915 name = n1Name+part_1_Name+vsZ+vsEtaPhi+vsPt; _n1_1_vsZVsEtaVsPhiVsPt = createHisto3F(name,name, _nBins_vertexZ,_min_vertexZ,_max_vertexZ, _nBins_etaPhi_1, 0., double(_nBins_etaPhi_1), _nBins_pt_1, _min_pt_1, _max_pt_1, "zVertex", _title_etaPhi_1, _title_pt_1);
917 name = n1Name+part_2_Name+vsPt; _n1_2_vsPt = createHisto1F(name,name, _nBins_pt_2, _min_pt_2, _max_pt_2, _title_pt_2, _title_AvgN_2);
918 name = n1Name+part_2_Name+vsZ+vsEtaPhi+vsPt; _n1_2_vsZVsEtaVsPhiVsPt = createHisto3F(name,name, _nBins_vertexZ,_min_vertexZ,_max_vertexZ, _nBins_etaPhi_2, 0., double(_nBins_etaPhi_2), _nBins_pt_2, _min_pt_2, _max_pt_2, "zVertex", _title_etaPhi_2, _title_pt_2);
923 name = n1Name+part_1_Name+vsEtaPhi; _n1_1_vsEtaVsPhi = createHisto2F(name,name, _nBins_eta_1, _min_eta_1, _max_eta_1, _nBins_phi_1, _min_phi_1, _max_phi_1, _title_eta_1, _title_phi_1, _title_AvgN_1);
924 name = s1ptName+part_1_Name+vsEtaPhi; _s1pt_1_vsEtaVsPhi = createHisto2F(name,name, _nBins_eta_1, _min_eta_1, _max_eta_1, _nBins_phi_1, _min_phi_1, _max_phi_1, _title_eta_1, _title_phi_1, _title_AvgSumPt_1);
925 name = n1Name+part_1_Name+vsM; _n1_1_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgN_1);
926 name = s1ptName+part_1_Name+vsM; _s1pt_1_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgSumPt_1);
927 name = n1NwName+part_1_Name+vsM; _n1Nw_1_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgN_1);
928 name = s1ptNwName+part_1_Name+vsM; _s1ptNw_1_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgSumPt_1);
930 name = n1Name+part_2_Name+vsEtaPhi; _n1_2_vsEtaVsPhi = createHisto2F(name,name, _nBins_eta_2, _min_eta_2, _max_eta_2, _nBins_phi_2, _min_phi_2, _max_phi_2, _title_eta_2, _title_phi_2, _title_AvgN_2);
931 name = s1ptName+part_2_Name+vsEtaPhi; _s1pt_2_vsEtaVsPhi = createHisto2F(name,name, _nBins_eta_2, _min_eta_2, _max_eta_2, _nBins_phi_2, _min_phi_2, _max_phi_2, _title_eta_2, _title_phi_2, _title_AvgSumPt_2);
932 name = n1Name+part_2_Name + vsM; _n1_2_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgN_2);
933 name = s1ptName+part_2_Name + vsM; _s1pt_2_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgSumPt_2);
934 name = n1NwName+part_2_Name+vsM; _n1Nw_2_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgN_1);
935 name = s1ptNwName+part_2_Name+vsM; _s1ptNw_2_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgSumPt_1);
937 name = n2Name+pair_12_Name+vsEtaPhi; _n2_12_vsEtaPhi = createHisto1F(name,name, _nBins_etaPhi_12, 0., double(_nBins_etaPhi_12), _title_etaPhi_12, _title_AvgN2_12);
938 name = s2PtPtName+pair_12_Name + vsEtaPhi;_s2PtPt_12_vsEtaPhi = createHisto1F(name,name, _nBins_etaPhi_12, 0., double(_nBins_etaPhi_12), _title_etaPhi_12, _title_AvgSumPtPt_12);
939 name = s2PtNName+pair_12_Name + vsEtaPhi; _s2PtN_12_vsEtaPhi = createHisto1F(name,name, _nBins_etaPhi_12, 0., double(_nBins_etaPhi_12), _title_etaPhi_12, _title_AvgSumPtN_12);
940 name = s2NPtName+pair_12_Name + vsEtaPhi; _s2NPt_12_vsEtaPhi = createHisto1F(name,name, _nBins_etaPhi_12, 0., double(_nBins_etaPhi_12), _title_etaPhi_12, _title_AvgNSumPt_12);
941 name = n2Name+pair_12_Name+vsPtVsPt; _n2_12_vsPtVsPt = createHisto2F(name,name, _nBins_pt_1, _min_pt_1, _max_pt_1, _nBins_pt_2, _min_pt_2, _max_pt_2, _title_pt_1, _title_pt_2, _title_AvgN2_12);
943 name = n2Name+pair_12_Name + vsM; _n2_12_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgN2_12);
944 name = s2PtPtName+pair_12_Name + vsM; _s2PtPt_12_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgSumPtPt_12);
945 name = s2PtNName+pair_12_Name + vsM; _s2PtN_12_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgSumPtN_12);
946 name = s2NPtName+pair_12_Name + vsM; _s2NPt_12_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgNSumPt_12);
948 name = n2NwName+pair_12_Name + vsM; _n2Nw_12_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgN2_12);
949 name = s2PtPtNwName+pair_12_Name + vsM; _s2PtPtNw_12_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgSumPtPt_12);
950 name = s2PtNNwName+pair_12_Name + vsM; _s2PtNNw_12_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgSumPtN_12);
951 name = s2NPtNwName+pair_12_Name + vsM; _s2NPtNw_12_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgNSumPt_12);
953 name = "mInv"; _invMass = createHisto1F(name,name, 50, 0.41, 0.55, "M_{inv}","counts");
954 name = "mInvElec"; _invMassElec = createHisto1F(name,name, 500, 0., 1.000, "M_{inv}","counts");
957 AliInfo(" AliAnalysisTaskDptDptCorrelations::createHistoHistograms() All Done");
959 //-----------------------//
961 void AliAnalysisTaskDptDptCorrelations::finalizeHistograms()
964 AliInfo("AliAnalysisTaskDptDptCorrelations::finalizeHistograms() starting");
965 AliInfo(Form("CorrelationAnalyzers::finalizeHistograms() _eventCount : %d",int(_eventCount)));
970 fillHistoWithArray(_n1_1_vsPt, __n1_1_vsPt, _nBins_pt_1);
971 fillHistoWithArray(_n1_1_vsZVsEtaVsPhiVsPt, __n1_1_vsZEtaPhiPt, _nBins_vertexZ, _nBins_etaPhi_1, _nBins_pt_1);
972 fillHistoWithArray(_n1_2_vsPt, __n1_1_vsPt, _nBins_pt_1);
973 fillHistoWithArray(_n1_2_vsZVsEtaVsPhiVsPt, __n1_1_vsZEtaPhiPt, _nBins_vertexZ, _nBins_etaPhi_1, _nBins_pt_1);
977 fillHistoWithArray(_n1_1_vsPt, __n1_1_vsPt, _nBins_pt_1);
978 fillHistoWithArray(_n1_1_vsZVsEtaVsPhiVsPt, __n1_1_vsZEtaPhiPt, _nBins_vertexZ, _nBins_etaPhi_1, _nBins_pt_1);
979 fillHistoWithArray(_n1_2_vsPt, __n1_2_vsPt, _nBins_pt_2);
980 fillHistoWithArray(_n1_2_vsZVsEtaVsPhiVsPt, __n1_2_vsZEtaPhiPt, _nBins_vertexZ, _nBins_etaPhi_2, _nBins_pt_2);
987 fillHistoWithArray(_n1_1_vsEtaVsPhi, __n1_1_vsEtaPhi, _nBins_eta_1, _nBins_phi_1);
988 fillHistoWithArray(_s1pt_1_vsEtaVsPhi, __s1pt_1_vsEtaPhi, _nBins_eta_1, _nBins_phi_1);
989 fillHistoWithArray(_n1_2_vsEtaVsPhi, __n1_1_vsEtaPhi, _nBins_eta_1, _nBins_phi_1);
990 fillHistoWithArray(_s1pt_2_vsEtaVsPhi, __s1pt_1_vsEtaPhi, _nBins_eta_1, _nBins_phi_1);
994 fillHistoWithArray(_n1_1_vsEtaVsPhi, __n1_1_vsEtaPhi, _nBins_eta_1, _nBins_phi_1);
995 fillHistoWithArray(_s1pt_1_vsEtaVsPhi, __s1pt_1_vsEtaPhi, _nBins_eta_1, _nBins_phi_1);
996 fillHistoWithArray(_n1_2_vsEtaVsPhi, __n1_2_vsEtaPhi, _nBins_eta_2, _nBins_phi_2);
997 fillHistoWithArray(_s1pt_2_vsEtaVsPhi, __s1pt_2_vsEtaPhi, _nBins_eta_2, _nBins_phi_2);
999 fillHistoWithArray(_n2_12_vsEtaPhi, __n2_12_vsEtaPhi, _nBins_etaPhi_12);
1000 fillHistoWithArray(_s2PtPt_12_vsEtaPhi, __s2ptpt_12_vsEtaPhi, _nBins_etaPhi_12);
1001 fillHistoWithArray(_s2PtN_12_vsEtaPhi, __s2PtN_12_vsEtaPhi, _nBins_etaPhi_12);
1002 fillHistoWithArray(_s2NPt_12_vsEtaPhi, __s2NPt_12_vsEtaPhi, _nBins_etaPhi_12);
1003 fillHistoWithArray(_n2_12_vsPtVsPt, __n2_12_vsPtPt, _nBins_pt_1, _nBins_pt_2);
1006 AliInfo("AliAnalysisTaskDptDptCorrelations::finalizeHistograms() Done ");
1011 void AliAnalysisTaskDptDptCorrelations::UserExec(Option_t */*option*/)
1015 int iPhi, iEta, iEtaPhi, iPt, charge;
1016 float q, phi, pt, eta, corr, corrPt, px, py, pz;
1018 int id_1, q_1, iEtaPhi_1, iPt_1;
1019 float pt_1, px_1, py_1, pz_1, corr_1;
1020 int id_2, q_2, iEtaPhi_2, iPt_2;
1021 float pt_2, px_2, py_2, pz_2, corr_2;
1023 int iVertex, iVertexP1, iVertexP2;
1025 float massElecSq = 1.94797849000000016e-02;
1028 const AliAODVertex* vertex;
1031 AliAnalysisManager* manager = AliAnalysisManager::GetAnalysisManager();
1035 AliAODInputHandler* inputHandler = dynamic_cast<AliAODInputHandler*> (manager->GetInputEventHandler());
1036 if (!inputHandler) {
1040 fAODEvent = dynamic_cast<AliAODEvent*>(InputEvent());
1041 //AliAODEvent* fAODEvent = dynamic_cast<AliAODEvent*>(InputEvent());
1047 fPIDResponse =inputHandler->GetPIDResponse();
1049 AliFatal("This Task needs the PID response attached to the inputHandler");
1053 // count all events looked at here
1056 if (_eventAccounting)
1058 _eventAccounting->Fill(0);// count all calls to this function
1066 _eventAccounting->Fill(1);// count all calls to this function with a valid pointer
1067 //reset single particle counters
1069 __n1_1 = __n1_2 = __s1pt_1 = __s1pt_2 = __n1Nw_1 = __n1Nw_2 = __s1ptNw_1 = __s1ptNw_2 = 0;
1071 float v0Centr = -999.;
1072 float v0ACentr = -999.;
1073 float trkCentr = -999.;
1074 float spdCentr = -999.;
1076 float vertexX = -999;
1077 float vertexY = -999;
1078 float vertexZ = -999;
1079 //float vertexXY = -999;
1080 //float dcaZ = -999;
1081 //float dcaXY = -999;
1082 float centrality = -999;
1087 AliCentrality* centralityObject = fAODEvent->GetHeader()->GetCentralityP();
1088 if (centralityObject)
1090 //cout << "AliAnalysisTaskDptDptCorrelations::UserExec(Option_t *option) - 6" << endl;
1092 v0Centr = centralityObject->GetCentralityPercentile("V0M");
1093 v0ACentr = centralityObject->GetCentralityPercentile("V0A");
1094 trkCentr = centralityObject->GetCentralityPercentile("TRK");
1095 spdCentr = centralityObject->GetCentralityPercentile("CL1");
1099 _nTracks =fAODEvent->GetNumberOfTracks();//NEW Test
1106 _field = fAODEvent->GetMagneticField();
1109 switch (_centralityMethod)
1111 case 0: centrality = _mult0; break;
1112 case 1: centrality = _mult1; break;
1113 case 2: centrality = _mult2; break;
1114 case 3: centrality = _mult3; break;
1115 case 4: centrality = _mult4; break;
1116 case 5: centrality = _mult5; break;
1117 case 6: centrality = _mult6; break;
1118 case 7: centrality = _mult4a; break;
1122 if ( centrality < _centralityMin || centrality > _centralityMax )
1126 _eventAccounting->Fill(2);// count all events with right centrality
1128 // filter on z and xy vertex
1129 vertex = (AliAODVertex*) fAODEvent->GetPrimaryVertex();
1131 //vertex->GetXYZ(V);
1136 vertex->GetCovarianceMatrix(fCov);
1137 if(vertex->GetNContributors() > 0)
1141 vertexX = vertex->GetX();
1142 vertexY = vertex->GetY();
1143 vertexZ = vertex->GetZ();
1145 if(TMath::Abs(vertexZ) > 10)
1153 _vertexZ->Fill(vertexZ);
1155 iVertex = int((vertexZ-_min_vertexZ)/_width_vertexZ);
1156 iVertexP1 = iVertex*_nBins_etaPhiPt_1;
1157 iVertexP2 = iVertex*_nBins_etaPhiPt_2;
1158 if (iVertex<0 || iVertex>=_nBins_vertexZ)
1160 AliError("AliAnalysisTaskDptDptCorrelations::Exec(Option_t *option) iVertex<0 || iVertex>=_nBins_vertexZ ");
1163 _eventAccounting->Fill(3);// count all calls to this function with a valid pointer
1164 //======================
1166 //*********************************************************
1167 TExMap *trackMap = new TExMap();//Mapping matrix----
1169 for(Int_t i = 0; i < _nTracks; i++)
1171 AliAODTrack* aodTrack = dynamic_cast<AliAODTrack *>(fAODEvent->GetTrack(i));
1173 AliError(Form("ERROR: Could not retrieve AODtrack %d",i));
1176 Int_t gID = aodTrack->GetID();
1177 if (aodTrack->TestFilterBit(1)) trackMap->Add(gID, i);//Global tracks
1180 AliAODTrack* newAodTrack;
1182 //Track Loop starts here
1183 for (int iTrack=0; iTrack< _nTracks; iTrack++)
1185 AliAODTrack* t = dynamic_cast<AliAODTrack *>(fAODEvent->GetTrack(iTrack));
1187 AliError(Form("Could not receive track %d", iTrack));
1191 bitOK = t->TestFilterBit(_trackFilterBit);
1192 if (!bitOK) continue; //128bit or 272bit
1194 Int_t gID = t->GetID();
1195 newAodTrack = gID >= 0 ?t : fAODEvent->GetTrack(trackMap->GetValue(-1-gID));
1206 //dcaZ = t->ZAtDCA();
1208 Double_t nsigmaelectron = TMath::Abs(fPIDResponse->NumberOfSigmasTPC(newAodTrack,(AliPID::EParticleType)AliPID::kElectron));
1209 Double_t nsigmapion = TMath::Abs(fPIDResponse->NumberOfSigmasTPC(newAodTrack,(AliPID::EParticleType)AliPID::kPion));
1210 Double_t nsigmakaon = TMath::Abs(fPIDResponse->NumberOfSigmasTPC(newAodTrack,(AliPID::EParticleType)AliPID::kKaon));
1211 Double_t nsigmaproton = TMath::Abs(fPIDResponse->NumberOfSigmasTPC(newAodTrack,(AliPID::EParticleType)AliPID::kProton));
1213 //nsigma cut to reject electron
1215 if(nsigmaelectron < fNSigmaCut
1216 && nsigmapion > fNSigmaCut
1217 && nsigmakaon > fNSigmaCut
1218 && nsigmaproton > fNSigmaCut ) continue;
1220 if(charge == 0) continue;
1222 // Kinematics cuts used
1223 if( pt < _min_pt_1 || pt > _max_pt_1) continue;
1224 if( eta < _min_eta_1 || eta > _max_eta_1) continue;
1227 newAodTrack->GetXYZ(pos);
1229 Double_t DCAX = pos[0] - vertexX;
1230 Double_t DCAY = pos[1] - vertexY;
1231 Double_t DCAZ = pos[2] - vertexZ;
1233 Double_t DCAXY = TMath::Sqrt((DCAX*DCAX) + (DCAY*DCAY));
1235 if (DCAZ < _dcaZMin ||
1237 DCAXY > _dcaXYMax ) continue;
1238 //==== QA ===========================
1240 _dcaxy->Fill(DCAXY);
1243 //===================================
1244 //*************************************************
1247 if (_requestedCharge_1 == charge)
1250 iPhi = int( phi/_width_phi_1);
1252 if (iPhi<0 || iPhi>=_nBins_phi_1 )
1254 AliWarning("AliAnalysisTaskDptDptCorrelations::analyze() iPhi<0 || iPhi>=_nBins_phi_1");
1258 iEta = int((eta-_min_eta_1)/_width_eta_1);
1259 if (iEta<0 || iEta>=_nBins_eta_1)
1261 AliWarning(Form("AliAnalysisTaskDptDptCorrelations::analyze(AliceEvent * event) Mismatched iEta: %d", iEta));
1264 iPt = int((pt -_min_pt_1 )/_width_pt_1 );
1265 if (iPt<0 || iPt >=_nBins_pt_1)
1267 AliWarning(Form("AliAnalysisTaskDptDptCorrelations::analyze(AliceEvent * event) Mismatched iPt: %d",iPt));
1270 iEtaPhi = iEta*_nBins_phi_1+iPhi;
1271 iZEtaPhiPt = iVertexP1 + iEtaPhi*_nBins_pt_1 + iPt;
1273 if (_correctionWeight_1)
1274 corr = _correctionWeight_1[iZEtaPhiPt];
1277 if (iZEtaPhiPt<0 || iZEtaPhiPt>=_nBins_zEtaPhiPt_1)
1279 AliWarning("AliAnalysisTaskDptDptCorrelations::analyze(AliceEvent * event) iZEtaPhiPt<0 || iZEtaPhiPt>=_nBins_zEtaPhiPt_1");
1287 __n1_1_vsPt[iPt] += corr; //cout << "step 15" << endl;
1288 __n1_1_vsZEtaPhiPt[iZEtaPhiPt] += corr; //cout << "step 12" << endl;
1295 _charge_1[k1] = charge;
1296 _iEtaPhi_1[k1] = iEtaPhi;
1302 _correction_1[k1] = corr;
1304 __n1_1_vsEtaPhi[iEtaPhi] += corr;
1306 __s1pt_1_vsEtaPhi[iEtaPhi] += corrPt;
1312 AliError(Form("AliAnalysisTaskDptDptCorrelations::analyze(AliceEvent * event) k1 >=arraySize; arraySize: %d",arraySize));
1318 if (!_sameFilter && _requestedCharge_2 == charge)
1321 iPhi = int( phi/_width_phi_2);
1323 if (iPhi<0 || iPhi>=_nBins_phi_2 )
1325 AliWarning("AliAnalysisTaskDptDptCorrelations::analyze() iPhi<0 || iPhi>=_nBins_phi_1");
1329 iEta = int((eta-_min_eta_2)/_width_eta_2);
1330 if (iEta<0 || iEta>=_nBins_eta_2)
1332 AliWarning(Form("AliAnalysisTaskDptDptCorrelations::analyze(AliceEvent * event) Mismatched iEta: %d", iEta));
1335 iPt = int((pt -_min_pt_2 )/_width_pt_2 );
1336 if (iPt<0 || iPt >=_nBins_pt_2)
1338 AliWarning(Form("AliAnalysisTaskDptDptCorrelations::analyze(AliceEvent * event) Mismatched iPt: %d",iPt));
1342 iEtaPhi = iEta*_nBins_phi_2+iPhi;
1343 iZEtaPhiPt = iVertexP2 + iEtaPhi*_nBins_pt_2 + iPt;
1344 if (iZEtaPhiPt<0 || iZEtaPhiPt>=_nBins_zEtaPhiPt_2)
1346 AliWarning("AliAnalysisTaskDptDptCorrelations::analyze(AliceEvent * event) iZEtaPhiPt<0 || iZEtaPhiPt>=_nBins_zEtaPhiPt_2");
1351 if (_correctionWeight_2)
1352 corr = _correctionWeight_2[iZEtaPhiPt];
1358 __n1_2_vsPt[iPt] += corr; //cout << "step 15" << endl;
1359 __n1_2_vsZEtaPhiPt[iZEtaPhiPt] += corr; //cout << "step 12" << endl;
1364 _id_2[k2] = iTrack; //cout << "step 1" << endl;
1365 _charge_2[k2] = charge; //cout << "step 2" << endl;
1366 _iEtaPhi_2[k2] = iEtaPhi; //cout << "step 3" << endl;
1367 _iPt_2[k2] = iPt; //cout << "step 4" << endl;
1368 _pt_2[k2] = pt; //cout << "step 5" << endl;
1369 _px_2[k2] = px; //cout << "step 6" << endl;
1370 _py_2[k2] = py; //cout << "step 7" << endl;
1371 _pz_2[k2] = pz; //cout << "step 8" << endl;
1372 _correction_2[k2] = corr; //cout << "step 9" << endl;
1373 __n1_2 += corr; //cout << "step 10" << endl;
1374 __s1pt_2 += corrPt; //cout << "step 13" << endl;
1376 __n1_2_vsEtaPhi[iEtaPhi] += corr; //cout << "step 11" << endl;
1377 __s1pt_2_vsEtaPhi[iEtaPhi] += corrPt; //cout << "step 14" << endl;
1382 AliWarning(Form("-W- k2 >=arraySize; arraySize: %d",arraySize));
1387 //cout << "done with track" << endl;
1393 //cout << "Filling histograms now" << endl;
1401 //_vertexZ->Fill(vertexZ);
1405 // nothing to do here.
1411 _n1_1_vsM->Fill(centrality, __n1_1);
1412 _s1pt_1_vsM->Fill(centrality, __s1pt_1);
1413 _n1Nw_1_vsM->Fill(centrality, __n1Nw_1);
1414 _s1ptNw_1_vsM->Fill(centrality, __s1ptNw_1);
1415 _n1_2_vsM->Fill(centrality, __n1_1);
1416 _s1pt_2_vsM->Fill(centrality, __s1pt_1);
1417 _n1Nw_2_vsM->Fill(centrality, __n1Nw_1);
1418 _s1ptNw_2_vsM->Fill(centrality, __s1ptNw_1);
1419 // reset pair counters
1420 __n2_12 = __s2ptpt_12 = __s2NPt_12 = __s2PtN_12 = 0;
1421 __n2Nw_12 = __s2ptptNw_12 = __s2NPtNw_12 = __s2PtNNw_12 = 0;
1424 for (int i1=0; i1<k1; i1++)
1426 ////cout << " i1:" << i1 << endl;
1427 id_1 = _id_1[i1]; ////cout << " id_1:" << id_1 << endl;
1428 q_1 = _charge_1[i1]; ////cout << " q_1:" << q_1 << endl;
1429 iEtaPhi_1 = _iEtaPhi_1[i1]; ////cout << " iEtaPhi_1:" << iEtaPhi_1 << endl;
1430 iPt_1 = _iPt_1[i1]; ////cout << " iPt_1:" << iPt_1 << endl;
1431 corr_1 = _correction_1[i1]; ////cout << " corr_1:" << corr_1 << endl;
1432 pt_1 = _pt_1[i1]; ////cout << " pt_1:" << pt_1 << endl;
1433 //dedx_1 = _dedx_1[i1]; ////cout << " dedx_1:" << dedx_1 << endl;
1435 for (int i2=i1+1; i2<k1; i2++)
1437 ////cout << " i2:" << i2 << endl;
1438 id_2 = _id_1[i2]; ////cout << " id_2:" << id_2 << endl;
1441 q_2 = _charge_1[i2]; ////cout << " q_1:" << q_1 << endl;
1442 iEtaPhi_2 = _iEtaPhi_1[i2]; ////cout << " iEtaPhi_1:" << iEtaPhi_1 << endl;
1443 iPt_2 = _iPt_1[i2]; ////cout << " iPt_1:" << iPt_1 << endl;
1444 corr_2 = _correction_1[i2]; ////cout << " corr_1:" << corr_1 << endl;
1445 pt_2 = _pt_1[i2]; ////cout << " pt_1:" << pt_1 << endl;
1446 //dedx_2 = _dedx_1[i2]; ////cout << " dedx_2:" << dedx_2 << endl;
1447 corr = corr_1*corr_2;
1448 if (q_2>q_1 || (q_1>0 && q_2>0 && pt_2<=pt_1) || (q_1<0 && q_2<0 && pt_2>=pt_1))
1450 ij = iEtaPhi_1*_nBins_etaPhi_1 + iEtaPhi_2; ////cout << " ij:" << ij<< endl;
1452 else // swap particles
1454 ij = iEtaPhi_2*_nBins_etaPhi_1 + iEtaPhi_1; ////cout << " ij:" << ij<< endl;
1458 __n2_12_vsEtaPhi[ij] += corr;
1460 __s2ptpt_12 += corr*ptpt;
1461 __s2PtN_12 += corr*pt_1;
1462 __s2NPt_12 += corr*pt_2;
1463 __s2ptpt_12_vsEtaPhi[ij] += corr*ptpt;
1464 __s2PtN_12_vsEtaPhi[ij] += corr*pt_1;
1465 __s2NPt_12_vsEtaPhi[ij] += corr*pt_2;
1466 __n2_12_vsPtPt[iPt_1*_nBins_pt_2 + iPt_2] += corr;
1469 __s2ptptNw_12 += ptpt;
1470 __s2PtNNw_12 += pt_1;
1471 __s2NPtNw_12 += pt_2;
1479 for (int i1=0; i1<k1; i1++)
1481 ////cout << " i1:" << i1 << endl;
1482 id_1 = _id_1[i1]; ////cout << " id_1:" << id_1 << endl;
1483 q_1 = _charge_1[i1]; ////cout << " q_1:" << q_1 << endl;
1484 iEtaPhi_1 = _iEtaPhi_1[i1]; ////cout << " iEtaPhi_1:" << iEtaPhi_1 << endl;
1485 iPt_1 = _iPt_1[i1]; ////cout << " iPt_1:" << iPt_1 << endl;
1486 corr_1 = _correction_1[i1]; ////cout << " corr_1:" << corr_1 << endl;
1487 pt_1 = _pt_1[i1]; ////cout << " pt_1:" << pt_1 << endl;
1488 //dedx_1 = _dedx_1[i1]; ////cout << " dedx_1:" << dedx_1 << endl;
1490 for (int i2=i1+1; i2<k1; i2++)
1492 ////cout << " i2:" << i2 << endl;
1493 id_2 = _id_1[i2]; ////cout << " id_2:" << id_2 << endl;
1496 q_2 = _charge_1[i2]; ////cout << " q_2:" << q_2 << endl;
1497 iEtaPhi_2 = _iEtaPhi_1[i2]; ////cout << " iEtaPhi_2:" << iEtaPhi_2 << endl;
1498 iPt_2 = _iPt_1[i2]; ////cout << " iPt_2:" << iPt_2 << endl;
1499 corr_2 = _correction_1[i2]; ////cout << " corr_2:" << corr_2 << endl;
1500 pt_2 = _pt_1[i2]; ////cout << " pt_2:" << pt_2 << endl;
1501 //dedx_2 = _dedx_1[i2]; ////cout << " dedx_2:" << dedx_2 << endl;
1502 corr = corr_1*corr_2;
1503 if ( q_2<q_1 || (q_1>0 && q_2>0 && pt_2>=pt_1) || (q_1<0 && q_2<0 && pt_2<=pt_1))
1505 ij = iEtaPhi_1*_nBins_etaPhi_1 + iEtaPhi_2; ////cout << " ij:" << ij<< endl;
1507 else // swap particles
1509 ij = iEtaPhi_2*_nBins_etaPhi_1 + iEtaPhi_1; ////cout << " ij:" << ij<< endl;
1513 __n2_12_vsEtaPhi[ij] += corr;
1515 __s2ptpt_12 += corr*ptpt;
1516 __s2PtN_12 += corr*pt_1;
1517 __s2NPt_12 += corr*pt_2;
1518 __s2ptpt_12_vsEtaPhi[ij] += corr*ptpt;
1519 __s2PtN_12_vsEtaPhi[ij] += corr*pt_1;
1520 __s2NPt_12_vsEtaPhi[ij] += corr*pt_2;
1521 __n2_12_vsPtPt[iPt_1*_nBins_pt_2 + iPt_2] += corr;
1524 __s2ptptNw_12 += ptpt;
1525 __s2PtNNw_12 += pt_1;
1526 __s2NPtNw_12 += pt_2;
1533 else // filter 1 and 2 are different -- must do all particle pairs...
1535 _n1_1_vsM->Fill(centrality, __n1_1);
1536 _s1pt_1_vsM->Fill(centrality, __s1pt_1);
1537 _n1Nw_1_vsM->Fill(centrality, __n1Nw_1);
1538 _s1ptNw_1_vsM->Fill(centrality, __s1ptNw_1);
1539 _n1_2_vsM->Fill(centrality, __n1_2);
1540 _s1pt_2_vsM->Fill(centrality, __s1pt_2);
1541 _n1Nw_2_vsM->Fill(centrality, __n1Nw_2);
1542 _s1ptNw_2_vsM->Fill(centrality, __s1ptNw_2);
1543 // reset pair counters
1544 __n2_12 = __s2ptpt_12 = __s2NPt_12 = __s2PtN_12 = 0;
1545 __n2Nw_12 = __s2ptptNw_12 = __s2NPtNw_12 = __s2PtNNw_12 = 0;
1546 for (int i1=0; i1<k1; i1++)
1548 ////cout << " i1:" << i1 << endl;
1549 id_1 = _id_1[i1]; ////cout << " id_1:" << id_1 << endl;
1550 q_1 = _charge_1[i1]; ////cout << " q_1:" << q_1 << endl;
1551 iEtaPhi_1 = _iEtaPhi_1[i1]; ////cout << " iEtaPhi_1:" << iEtaPhi_1 << endl;
1552 iPt_1 = _iPt_1[i1]; ////cout << " iPt_1:" << iPt_1 << endl;
1553 corr_1 = _correction_1[i1]; ////cout << " corr_1:" << corr_1 << endl;
1554 pt_1 = _pt_1[i1]; ////cout << " pt_1:" << pt_1 << endl;
1555 px_1 = _px_1[i1]; ////cout << " px_1:" << px_1 << endl;
1556 py_1 = _py_1[i1]; ////cout << " py_1:" << py_1 << endl;
1557 pz_1 = _pz_1[i1]; ////cout << " pz_1:" << pz_1 << endl;
1558 //dedx_1 = _dedx_1[i1]; ////cout << " dedx_1:" << dedx_1 << endl;
1561 for (int i2=0; i2<k2; i2++)
1563 ////cout << " i2:" << i2 << endl;
1564 id_2 = _id_2[i2]; ////cout << " id_2:" << id_2 << endl;
1565 if (id_1!=id_2) // exclude auto correlation
1567 q_2 = _charge_2[i2]; ////cout << " q_2:" << q_2 << endl;
1568 iEtaPhi_2 = _iEtaPhi_2[i2]; ////cout << " iEtaPhi_2:" << iEtaPhi_2 << endl;
1569 iPt_2 = _iPt_2[i2]; ////cout << " iPt_2:" << iPt_2 << endl;
1570 corr_2 = _correction_2[i2]; ////cout << " corr_2:" << corr_2 << endl;
1571 pt_2 = _pt_2[i2]; ////cout << " pt_2:" << pt_2 << endl;
1572 px_2 = _px_2[i2]; ////cout << " px_2:" << px_2 << endl;
1573 py_2 = _py_2[i2]; ////cout << " py_2:" << py_2 << endl;
1574 pz_2 = _pz_2[i2]; ////cout << " pz_2:" << pz_2 << endl;
1575 //dedx_2 = _dedx_2[i2]; ////cout << " dedx_2:" << dedx_2 << endl;
1578 if (_rejectPairConversion)
1580 float e1Sq = massElecSq + pt_1*pt_1 + pz_1*pz_1;
1581 float e2Sq = massElecSq + pt_2*pt_2 + pz_2*pz_2;
1582 float mInvSq = 2*(massElecSq + sqrt(e1Sq*e2Sq) - px_1*px_2 - py_1*py_2 - pz_1*pz_2 );
1583 float mInv = sqrt(mInvSq);
1584 _invMass->Fill(mInv);
1587 corr = corr_1*corr_2;
1588 ij = iEtaPhi_1*_nBins_etaPhi_1 + iEtaPhi_2; ////cout << " ij:" << ij<< endl;
1590 __n2_12_vsEtaPhi[ij] += corr;
1592 __s2ptpt_12 += corr*ptpt;
1593 __s2PtN_12 += corr*pt_1;
1594 __s2NPt_12 += corr*pt_2;
1595 __s2ptpt_12_vsEtaPhi[ij] += corr*ptpt;
1596 __s2PtN_12_vsEtaPhi[ij] += corr*pt_1;
1597 __s2NPt_12_vsEtaPhi[ij] += corr*pt_2;
1598 __n2_12_vsPtPt[iPt_1*_nBins_pt_2 + iPt_2] += corr;
1600 __s2ptptNw_12 += ptpt;
1601 __s2PtNNw_12 += pt_1;
1602 __s2NPtNw_12 += pt_2;
1609 _n2_12_vsM->Fill(centrality, __n2_12);
1610 _s2PtPt_12_vsM->Fill(centrality, __s2ptpt_12);
1611 _s2PtN_12_vsM->Fill(centrality, __s2NPt_12);
1612 _s2NPt_12_vsM->Fill(centrality, __s2PtN_12);
1614 _n2Nw_12_vsM->Fill(centrality, __n2Nw_12);
1615 _s2PtPtNw_12_vsM->Fill(centrality, __s2ptptNw_12);
1616 _s2PtNNw_12_vsM->Fill(centrality, __s2NPtNw_12);
1617 _s2NPtNw_12_vsM->Fill(centrality, __s2PtNNw_12);
1622 AliInfo("AliAnalysisTaskDptDptCorrelations::UserExec() -----------------Event Done ");
1623 PostData(0,_outputHistoList);
1627 void AliAnalysisTaskDptDptCorrelations::FinishTaskOutput()
1629 AliInfo("AliAnalysisTaskDptDptCorrelations::FinishTaskOutput() Starting.");
1630 Printf("= 0 ====================================================================");
1631 finalizeHistograms();
1632 AliInfo("= 1 ====================================================================");
1633 PostData(0,_outputHistoList);
1634 AliInfo("= 2 ====================================================================");
1635 AliInfo("AliAnalysisTaskDptDptCorrelations::FinishTaskOutput() Done.");
1638 void AliAnalysisTaskDptDptCorrelations::Terminate(Option_t* /*option*/)
1640 AliInfo("AliAnalysisTaskDptDptCorrelations::Terminate() Starting/Done.");
1645 //===================================================================================================
1646 void AliAnalysisTaskDptDptCorrelations::fillHistoWithArray(TH1 * h, float * array, int size)
1649 float v1, ev1, v2, ev2, sum, esum;
1650 for (i=0, i1=1; i<size; ++i,++i1)
1652 v1 = array[i]; ev1 = sqrt(v1);
1653 v2 = h->GetBinContent(i1);
1654 ev2 = h->GetBinError(i1);
1656 esum = sqrt(ev1*ev1+ev2*ev2);
1657 h->SetBinContent(i1,sum);
1658 h->SetBinError(i1,esum);
1662 void AliAnalysisTaskDptDptCorrelations::fillHistoWithArray(TH2 * h, float * array, int size1, int size2)
1666 float v1, ev1, v2, ev2, sum, esum;
1667 for (i=0, i1=1; i<size1; ++i,++i1)
1669 for (j=0, j1=1; j<size2; ++j,++j1)
1671 v1 = array[i*size2+j]; ev1 = sqrt(v1);
1672 v2 = h->GetBinContent(i1,j1);
1673 ev2 = h->GetBinError(i1,j1);
1675 esum = sqrt(ev1*ev1+ev2*ev2);
1676 h->SetBinContent(i1,j1,sum);
1677 h->SetBinError(i1,j1,esum);
1682 void AliAnalysisTaskDptDptCorrelations::fillHistoWithArray(TH3 * h, float * array, int size1, int size2, int size3)
1687 float v1, ev1, v2, ev2, sum, esum;
1688 int size23 = size2*size3;
1689 for (i=0, i1=1; i<size1; ++i,++i1)
1691 for (j=0, j1=1; j<size2; ++j,++j1)
1693 for (k=0, k1=1; k<size3; ++k,++k1)
1695 v1 = array[i*size23+j*size3+k]; ev1 = sqrt(v1);
1696 v2 = h->GetBinContent(i1,j1,k1);
1697 ev2 = h->GetBinError(i1,j1,k1);
1699 esum = sqrt(ev1*ev1+ev2*ev2);
1700 h->SetBinContent(i1,j1,k1,sum);
1701 h->SetBinError(i1,j1,k1,esum);
1707 void AliAnalysisTaskDptDptCorrelations::fillHistoWithArray(TH1 * h, double * array, int size)
1710 double v1, ev1, v2, ev2, sum, esum;
1711 for (i=0, i1=1; i<size; ++i,++i1)
1713 v1 = array[i]; ev1 = sqrt(v1);
1714 v2 = h->GetBinContent(i1);
1715 ev2 = h->GetBinError(i1);
1717 esum = sqrt(ev1*ev1+ev2*ev2);
1718 h->SetBinContent(i1,sum);
1719 h->SetBinError(i1,esum);
1723 void AliAnalysisTaskDptDptCorrelations::fillHistoWithArray(TH2 * h, double * array, int size1, int size2)
1727 double v1, ev1, v2, ev2, sum, esum;
1728 for (i=0, i1=1; i<size1; ++i,++i1)
1730 for (j=0, j1=1; j<size2; ++j,++j1)
1732 v1 = array[i*size2+j]; ev1 = sqrt(v1);
1733 v2 = h->GetBinContent(i1,j1);
1734 ev2 = h->GetBinError(i1,j1);
1736 esum = sqrt(ev1*ev1+ev2*ev2);
1737 h->SetBinContent(i1,j1,sum);
1738 h->SetBinError(i1,j1,esum);
1743 void AliAnalysisTaskDptDptCorrelations::fillHistoWithArray(TH3 * h, double * array, int size1, int size2, int size3)
1748 double v1, ev1, v2, ev2, sum, esum;
1749 int size23 = size2*size3;
1750 for (i=0, i1=1; i<size1; ++i,++i1)
1752 for (j=0, j1=1; j<size2; ++j,++j1)
1754 for (k=0, k1=1; k<size3; ++k,++k1)
1756 v1 = array[i*size23+j*size3+k]; ev1 = sqrt(v1);
1757 v2 = h->GetBinContent(i1,j1,k1);
1758 ev2 = h->GetBinError(i1,j1,k1);
1760 esum = sqrt(ev1*ev1+ev2*ev2);
1761 h->SetBinContent(i1,j1,k1,sum);
1762 h->SetBinError(i1,j1,k1,esum);
1768 //________________________________________________________________________
1769 double * AliAnalysisTaskDptDptCorrelations::getDoubleArray(int size, double v)
1771 /// Allocate an array of type double with n values
1772 /// Initialize the array to the given value
1773 double * array = new double [size];
1774 for (int i=0;i<size;++i) array[i]=v;
1778 //________________________________________________________________________
1779 float * AliAnalysisTaskDptDptCorrelations::getFloatArray(int size, float v)
1781 /// Allocate an array of type float with n values
1782 /// Initialize the array to the given value
1783 float * array = new float [size];
1784 for (int i=0;i<size;++i) array[i]=v;
1789 //________________________________________________________________________
1790 TH1D * AliAnalysisTaskDptDptCorrelations::createHisto1D(const TString & name, const TString & title,
1791 int n, double xMin, double xMax,
1792 const TString & xTitle, const TString & yTitle)
1794 //CreateHisto new 1D historgram
1795 AliInfo(Form("createHisto 1D histo %s nBins: %d xMin: %f xMax: %f",name.Data(),n,xMin,xMax));
1796 TH1D * h = new TH1D(name,title,n,xMin,xMax);
1797 h->GetXaxis()->SetTitle(xTitle);
1798 h->GetYaxis()->SetTitle(yTitle);
1804 //________________________________________________________________________
1805 TH1D * AliAnalysisTaskDptDptCorrelations::createHisto1D(const TString & name, const TString & title,
1806 int n, double * bins,
1807 const TString & xTitle, const TString & yTitle)
1809 AliInfo(Form("createHisto 1D histo %s with %d non uniform nBins",name.Data(),n));
1810 TH1D * h = new TH1D(name,title,n,bins);
1811 h->GetXaxis()->SetTitle(xTitle);
1812 h->GetYaxis()->SetTitle(yTitle);
1818 //________________________________________________________________________
1819 TH2D * AliAnalysisTaskDptDptCorrelations::createHisto2D(const TString & name, const TString & title,
1820 int nx, double xMin, double xMax, int ny, double yMin, double yMax,
1821 const TString & xTitle, const TString & yTitle, const TString & zTitle)
1823 AliInfo(Form("createHisto 2D histo %s nx: %d xMin: %f10.4 xMax: %f10.4 ny: %d yMin: %f10.4 yMax: %f10.4",name.Data(),nx,xMin,xMax,ny,yMin,yMax));
1824 TH2D * h = new TH2D(name,title,nx,xMin,xMax,ny,yMin,yMax);
1825 h->GetXaxis()->SetTitle(xTitle);
1826 h->GetYaxis()->SetTitle(yTitle);
1827 h->GetZaxis()->SetTitle(zTitle);
1832 //________________________________________________________________________
1833 TH2D * AliAnalysisTaskDptDptCorrelations::createHisto2D(const TString & name, const TString & title,
1834 int nx, double* xbins, int ny, double yMin, double yMax,
1835 const TString & xTitle, const TString & yTitle, const TString & zTitle)
1837 AliInfo(Form("createHisto 2D histo %s with %d non uniform nBins",name.Data(),nx));
1839 h = new TH2D(name,title,nx,xbins,ny,yMin,yMax);
1840 h->GetXaxis()->SetTitle(xTitle);
1841 h->GetYaxis()->SetTitle(yTitle);
1842 h->GetZaxis()->SetTitle(zTitle);
1848 //________________________________________________________________________
1849 TH1F * AliAnalysisTaskDptDptCorrelations::createHisto1F(const TString & name, const TString & title,
1850 int n, double xMin, double xMax,
1851 const TString & xTitle, const TString & yTitle)
1853 //CreateHisto new 1D historgram
1854 AliInfo(Form("createHisto 1D histo %s nBins: %d xMin: %f xMax: %f",name.Data(),n,xMin,xMax));
1855 TH1F * h = new TH1F(name,title,n,xMin,xMax);
1856 h->GetXaxis()->SetTitle(xTitle);
1857 h->GetYaxis()->SetTitle(yTitle);
1863 //________________________________________________________________________
1864 TH1F * AliAnalysisTaskDptDptCorrelations::createHisto1F(const TString & name, const TString & title,
1865 int n, double * bins,
1866 const TString & xTitle, const TString & yTitle)
1868 AliInfo(Form("createHisto 1D histo %s with %d non uniform nBins",name.Data(),n));
1869 TH1F * h = new TH1F(name,title,n,bins);
1870 h->GetXaxis()->SetTitle(xTitle);
1871 h->GetYaxis()->SetTitle(yTitle);
1877 //________________________________________________________________________
1878 TH2F * AliAnalysisTaskDptDptCorrelations::createHisto2F(const TString & name, const TString & title,
1879 int nx, double xMin, double xMax, int ny, double yMin, double yMax,
1880 const TString & xTitle, const TString & yTitle, const TString & zTitle)
1882 AliInfo(Form("createHisto 2D histo %s nx: %d xMin: %f10.4 xMax: %f10.4 ny: %d yMin: %f10.4 yMax: %f10.4",name.Data(),nx,xMin,xMax,ny,yMin,yMax));
1883 TH2F * h = new TH2F(name,title,nx,xMin,xMax,ny,yMin,yMax);
1884 h->GetXaxis()->SetTitle(xTitle);
1885 h->GetYaxis()->SetTitle(yTitle);
1886 h->GetZaxis()->SetTitle(zTitle);
1891 //________________________________________________________________________
1892 TH2F * AliAnalysisTaskDptDptCorrelations::createHisto2F(const TString & name, const TString & title,
1893 int nx, double* xbins, int ny, double yMin, double yMax,
1894 const TString & xTitle, const TString & yTitle, const TString & zTitle)
1896 AliInfo(Form("createHisto 2D histo %s with %d non uniform nBins",name.Data(),nx));
1898 h = new TH2F(name,title,nx,xbins,ny,yMin,yMax);
1899 h->GetXaxis()->SetTitle(xTitle);
1900 h->GetYaxis()->SetTitle(yTitle);
1901 h->GetZaxis()->SetTitle(zTitle);
1907 //________________________________________________________________________
1908 TH3F * AliAnalysisTaskDptDptCorrelations::createHisto3F(const TString & name, const TString & title,
1909 int nx, double xMin, double xMax,
1910 int ny, double yMin, double yMax,
1911 int nz, double zMin, double zMax,
1912 const TString & xTitle, const TString & yTitle, const TString & zTitle)
1914 AliInfo(Form("createHisto 2D histo %s nx: %d xMin: %f10.4 xMax: %f10.4 ny: %d yMin: %f10.4 yMax: %f10.4 nz: %d zMin: %f10.4 zMax: %f10.4",name.Data(),nx,xMin,xMax,ny,yMin,yMax,nz,zMin,zMax));
1915 TH3F * h = new TH3F(name,title,nx,xMin,xMax,ny,yMin,yMax,nz,zMin,zMax);
1916 h->GetXaxis()->SetTitle(xTitle);
1917 h->GetYaxis()->SetTitle(yTitle);
1918 h->GetZaxis()->SetTitle(zTitle);
1924 //________________________________________________________________________
1925 TProfile * AliAnalysisTaskDptDptCorrelations::createProfile(const TString & name, const TString & description,
1926 int nx,double xMin,double xMax,
1927 const TString & xTitle, const TString & yTitle)
1929 AliInfo(Form("createHisto 1D profile %s nBins: %d xMin: %f10.4 xMax: %f10.4",name.Data(),nx,xMin,xMax));
1930 TProfile * h = new TProfile(name,description,nx,xMin,xMax);
1931 h->GetXaxis()->SetTitle(xTitle);
1932 h->GetYaxis()->SetTitle(yTitle);
1937 //________________________________________________________________________
1938 TProfile * AliAnalysisTaskDptDptCorrelations::createProfile(const TString & name,const TString & description,
1939 int nx, double* bins,
1940 const TString & xTitle, const TString & yTitle)
1942 AliInfo(Form("createHisto 1D profile %s with %d non uniform bins",name.Data(),nx));
1943 TProfile * h = new TProfile(name,description,nx,bins);
1944 h->GetXaxis()->SetTitle(xTitle);
1945 h->GetYaxis()->SetTitle(yTitle);
1951 void AliAnalysisTaskDptDptCorrelations::addToList(TH1 *h)
1953 if (_outputHistoList)
1955 _outputHistoList->Add(h);
1958 AliInfo("addToList(TH1 *h) _outputHistoList is null!!!!! Should abort ship");