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),
428 _correctionWeight_1(0),
429 _correctionWeight_2(0),
430 _nBins_M0(500), _min_M0(0), _max_M0(10000), _width_M0(20),
431 _nBins_M1(500), _min_M1(0), _max_M1(10000), _width_M1(20),
432 _nBins_M2(500), _min_M2(0), _max_M2(10000), _width_M2(20),
433 _nBins_M3(500), _min_M3(0), _max_M3(10000), _width_M3(20),
434 _nBins_M4(100), _min_M4(0), _max_M4(1), _width_M4(0.01),
435 _nBins_M5(100), _min_M5(0), _max_M5(1), _width_M5(0.01),
436 _nBins_M6(100), _min_M6(0), _max_M6(1), _width_M6(0.01),
437 _nBins_vertexZ(40), _min_vertexZ(-10), _max_vertexZ(10), _width_vertexZ(0.5),
439 _nBins_pt_1(18), _min_pt_1(0.2), _max_pt_1(2.0), _width_pt_1(0.1),
440 _nBins_phi_1(72), _min_phi_1(0), _max_phi_1(2.*3.1415927),_width_phi_1(2.*3.1415927/72.),
441 _nBins_eta_1(0), _min_eta_1(0), _max_eta_1(0), _width_eta_1(0.1),
444 _nBins_etaPhiPt_1(0),
445 _nBins_zEtaPhiPt_1(0),
447 _nBins_pt_2(18), _min_pt_2(0.2), _max_pt_2(2.0), _width_pt_2(0.1),
448 _nBins_phi_2(72), _min_phi_2(0), _max_phi_2(2.*3.1415927),_width_phi_2(2.*3.1415927/72),
449 _nBins_eta_2(0), _min_eta_2(0), _max_eta_2(0), _width_eta_2(0.1),
452 _nBins_etaPhiPt_2(0),
453 _nBins_zEtaPhiPt_2(0),
473 __s1pt_1_vsEtaPhi(0),
474 __n1_1_vsZEtaPhiPt(0),
477 __s1pt_2_vsEtaPhi(0),
478 __n1_2_vsZEtaPhiPt(0),
481 __s2ptpt_12_vsEtaPhi(0),
482 __s2PtN_12_vsEtaPhi(0),
483 __s2NPt_12_vsEtaPhi(0),
486 _eventAccounting ( 0),
503 _n1_1_vsEtaVsPhi ( 0),
504 _s1pt_1_vsEtaVsPhi ( 0),
505 _n1_1_vsZVsEtaVsPhiVsPt ( 0),
514 _n1_2_vsEtaVsPhi ( 0),
515 _s1pt_2_vsEtaVsPhi ( 0),
516 _n1_2_vsZVsEtaVsPhiVsPt ( 0),
524 _n2_12_vsEtaPhi ( 0),
525 _n2_12_vsPtVsPt ( 0),
526 _s2PtPt_12_vsEtaPhi( 0),
527 _s2PtN_12_vsEtaPhi ( 0),
528 _s2NPt_12_vsEtaPhi ( 0),
534 _s2PtPtNw_12_vsM ( 0),
535 _s2PtNNw_12_vsM ( 0),
536 _s2NPtNw_12_vsM ( 0),
559 intBinCorrName("NA"),
609 _title_etaPhi_1("NA"),
611 _title_SumPt_1("NA"),
612 _title_AvgPt_1("NA"),
614 _title_AvgSumPt_1("NA"),
619 _title_etaPhi_2("NA"),
621 _title_SumPt_2("NA"),
622 _title_AvgPt_2("NA"),
624 _title_AvgSumPt_2("NA"),
626 _title_etaPhi_12("NA"),
628 _title_AvgN2_12("NA"),
629 _title_AvgSumPtPt_12("NA"),
630 _title_AvgSumPtN_12("NA"),
631 _title_AvgNSumPt_12("NA"),
641 printf("2nd constructor called ");
643 DefineOutput(0, TList::Class());
649 AliAnalysisTaskDptDptCorrelations::~AliAnalysisTaskDptDptCorrelations()
654 void AliAnalysisTaskDptDptCorrelations::UserCreateOutputObjects()
657 _outputHistoList = new TList();
658 _outputHistoList->SetOwner();
660 _nBins_M0 = 500; _min_M0 = 0.; _max_M0 = 5000.; _width_M0 = (_max_M0-_min_M0)/_nBins_M0;
661 _nBins_M1 = 500; _min_M1 = 0.; _max_M1 = 5000.; _width_M1 = (_max_M1-_min_M1)/_nBins_M1;
662 _nBins_M2 = 500; _min_M2 = 0.; _max_M2 = 5000.; _width_M2 = (_max_M2-_min_M2)/_nBins_M2;
663 _nBins_M3 = 500; _min_M3 = 0.; _max_M3 = 5000.; _width_M3 = (_max_M3-_min_M3)/_nBins_M3;
664 _nBins_M4 = 100; _min_M4 = 0.; _max_M4 = 100.; _width_M4 = (_max_M4-_min_M4)/_nBins_M4;
665 _nBins_M5 = 100; _min_M5 = 0.; _max_M5 = 100.; _width_M5 = (_max_M5-_min_M5)/_nBins_M5;
666 _nBins_M6 = 100; _min_M6 = 0.; _max_M6 = 100.; _width_M6 = (_max_M6-_min_M6)/_nBins_M6;
668 _min_vertexZ = _vertexZMin;
669 _max_vertexZ = _vertexZMax;
670 _width_vertexZ = 0.5;
671 _nBins_vertexZ = int(0.5+ (_max_vertexZ - _min_vertexZ)/_width_vertexZ);
672 _nBins_pt_1 = int(0.5+ (_max_pt_1 -_min_pt_1 )/_width_pt_1);
673 _nBins_eta_1 = int(0.5+ (_max_eta_1-_min_eta_1)/_width_eta_1);
674 _width_phi_1 = (_max_phi_1 - _min_phi_1) /_nBins_phi_1;
675 _nBins_etaPhi_1 = _nBins_phi_1 * _nBins_eta_1;
676 _nBins_etaPhiPt_1 = _nBins_etaPhi_1 * _nBins_pt_1;
677 _nBins_zEtaPhiPt_1 = _nBins_vertexZ * _nBins_etaPhiPt_1;
679 _nBins_pt_2 = int(0.5+ (_max_pt_2 -_min_pt_2 )/_width_pt_2);
680 _nBins_eta_2 = int(0.5+ (_max_eta_2-_min_eta_2)/_width_eta_2);
681 _width_phi_2 = (_max_phi_2 - _min_phi_2) /_nBins_phi_2;
682 _nBins_etaPhi_2 = _nBins_phi_2 * _nBins_eta_2;
683 _nBins_etaPhiPt_2 = _nBins_etaPhi_2 * _nBins_pt_2;
684 _nBins_zEtaPhiPt_2 = _nBins_vertexZ * _nBins_etaPhiPt_2;
685 _nBins_etaPhi_12 = _nBins_etaPhi_1 * _nBins_etaPhi_2;
687 _id_1 = new int[arraySize];
688 _charge_1 = new int[arraySize];
689 _iEtaPhi_1 = new int[arraySize];
690 _iPt_1 = new int[arraySize];
691 _pt_1 = new float[arraySize];
692 _px_1 = new float[arraySize];
693 _py_1 = new float[arraySize];
694 _pz_1 = new float[arraySize];
695 _correction_1 = new float[arraySize];
697 __n1_1_vsPt = getDoubleArray(_nBins_pt_1, 0.);
698 __n1_1_vsEtaPhi = getDoubleArray(_nBins_etaPhi_1, 0.);
699 __s1pt_1_vsEtaPhi = getDoubleArray(_nBins_etaPhi_1, 0.);
700 __n1_1_vsZEtaPhiPt = getFloatArray(_nBins_zEtaPhiPt_1, 0.);
703 if (_requestedCharge_2!=_requestedCharge_1)
707 _id_2 = new int[arraySize];
708 _charge_2 = new int[arraySize];
709 _iEtaPhi_2 = new int[arraySize];
710 _iPt_2 = new int[arraySize];
711 _pt_2 = new float[arraySize];
712 _px_2 = new float[arraySize];
713 _py_2 = new float[arraySize];
714 _pz_2 = new float[arraySize];
715 _correction_2 = new float[arraySize];
717 __n1_2_vsPt = getDoubleArray(_nBins_pt_2, 0.);
718 __n1_2_vsEtaPhi = getDoubleArray(_nBins_etaPhi_2, 0.);
719 __s1pt_2_vsEtaPhi = getDoubleArray(_nBins_etaPhi_2, 0.);
720 __n1_2_vsZEtaPhiPt = getFloatArray(_nBins_zEtaPhiPt_2, 0.);
724 __n2_12_vsPtPt = getDoubleArray(_nBins_pt_1*_nBins_pt_2,0.);
725 __n2_12_vsEtaPhi = getFloatArray(_nBins_etaPhi_12, 0.);
726 __s2ptpt_12_vsEtaPhi = getFloatArray(_nBins_etaPhi_12, 0.);
727 __s2PtN_12_vsEtaPhi = getFloatArray(_nBins_etaPhi_12, 0.);
728 __s2NPt_12_vsEtaPhi = getFloatArray(_nBins_etaPhi_12, 0.);
730 // Setup all the labels needed.
734 pair_12_Name = "_12";
751 binCorrName = "binCorr";
752 intBinCorrName = "intBinCorr";
757 s1ptNwName = "sumPtNw";
758 s1DptName = "sumDpt";
759 s2PtPtName = "sumPtPt";
760 s2PtPtNwName = "sumPtPtNw";
761 s2DptDptName = "sumDptDpt";
762 s2NPtName = "sumNPt";
763 s2NPtNwName = "sumNPtNw";
764 s2PtNName = "sumPtN";
765 s2NPtNwName = "sumNPtNw";
766 s2PtNNwName = "sumPtNNw";
768 ptptName = "avgPtPt";
769 pt1pt1Name = "avgPtavgPt";
771 DptDptName = "avgDptDpt";
772 RDptDptName = "relDptDpt"; // ratio of avgDptDpt by avgPt*avgPt
777 _title_counts = "yield";
783 _title_m4 = "V0Centrality";
784 _title_m5 = "TrkCentrality";
785 _title_m6 = "SpdCentrality";
787 _title_eta_1 = "#eta_{1}";
788 _title_phi_1 = "#varphi_{1} (radian)";
789 _title_etaPhi_1 = "#eta_{1}#times#varphi_{1}";
790 _title_pt_1 = "p_{t,1} (GeV/c)";
791 _title_n_1 = "n_{1}";
792 _title_SumPt_1 = "#Sigma p_{t,1} (GeV/c)";
793 _title_AvgPt_1 = "#LT p_{t,1} #GT (GeV/c)";
794 _title_AvgN_1 = "#LT n_{1} #GT";
795 _title_AvgSumPt_1 = "#LT #Sigma p_{t,1} #GT (GeV/c)";
797 _title_eta_2 = "#eta_{2}";
798 _title_phi_2 = "#varphi_{2} (radian)";
799 _title_etaPhi_2 = "#eta_{2}#times#varphi_{2}";
800 _title_pt_2 = "p_{t,2} (GeV/c)";
801 _title_n_2 = "n_{2}";
802 _title_SumPt_2 = "#Sigma p_{t,1} (GeV/c)";
803 _title_AvgPt_2 = "#LT p_{t,2} #GT (GeV/c)";
804 _title_AvgN_2 = "#LT n_{2} #GT";
805 _title_AvgSumPt_2 = "#LT #Sigma p_{t,2} #GT (GeV/c)";
807 _title_etaPhi_12 = "#eta_{1}#times#varphi_{1}#times#eta_{2}#times#varphi_{2}";
809 _title_AvgN2_12 = "#LT n_{2} #GT";;
810 _title_AvgSumPtPt_12 = "#LT #Sigma p_{t,1}p_{t,2} #GT";;
811 _title_AvgSumPtN_12 = "#LT #Sigma p_{t,1}N #GT";;
812 _title_AvgNSumPt_12 = "#LT N#Sigma p_{t,2} #GT";;
820 vsEtaPhi = "_vsEtaPhi";
821 vsPtVsPt = "_vsPtVsPt";
826 int iZ, iEtaPhi, iPt;
827 int iZ1,iEtaPhi1,iPt1;
831 _correctionWeight_1 = new float[_nBins_vertexZ*_nBins_etaPhi_1*_nBins_pt_1];
833 b = _nBins_etaPhi_1*_nBins_pt_1;
834 for (iZ=0,iZ1=1; iZ<_nBins_vertexZ; iZ++, iZ1++)
836 for (iEtaPhi=0,iEtaPhi1=1; iEtaPhi<_nBins_etaPhi_1; iEtaPhi++, iEtaPhi1++)
838 for (iPt=0,iPt1=1; iPt<_nBins_pt_1; iPt++, iPt1++)
840 _correctionWeight_1[iZ*b+iEtaPhi*a+iPt] = _weight_1->GetBinContent(iZ1,iEtaPhi1,iPt1);
847 AliError("AliAnalysisTaskDptDptCorrelations:: _weight_1 is a null pointer.");
854 _correctionWeight_2 = new float[_nBins_vertexZ*_nBins_etaPhi_2*_nBins_pt_2];
856 b = _nBins_etaPhi_2*_nBins_pt_2;
857 for (iZ=0,iZ1=1; iZ<_nBins_vertexZ; iZ++, iZ1++)
859 for (iEtaPhi=0,iEtaPhi1=1; iEtaPhi<_nBins_etaPhi_2; iEtaPhi++, iEtaPhi1++)
861 for (iPt=0,iPt1=1; iPt<_nBins_pt_2; iPt++, iPt1++)
863 _correctionWeight_2[iZ*b+iEtaPhi*a+iPt] = _weight_2->GetBinContent(iZ1,iEtaPhi1,iPt1);
870 AliError("AliAnalysisTaskDptDptCorrelations:: _weight_1 is a null pointer.");
877 PostData(0,_outputHistoList);
879 //cout<< "AliAnalysisTaskDptDptCorrelations::CreateOutputObjects() DONE " << endl;
883 void AliAnalysisTaskDptDptCorrelations::createHistograms()
885 AliInfo(" AliAnalysisTaskDptDptCorrelations::createHistoHistograms() Creating Event Histos");
888 name = "eventAccounting";
890 _eventAccounting = createHisto1D(name,name,10, -0.5, 9.5, "event Code", _title_counts);
892 name = "m0"; _m0 = createHisto1D(name,name,_nBins_M1, _min_M1, _max_M1, _title_m0, _title_counts);
893 name = "m1"; _m1 = createHisto1D(name,name,_nBins_M1, _min_M1, _max_M1, _title_m1, _title_counts);
894 name = "m2"; _m2 = createHisto1D(name,name,_nBins_M2, _min_M2, _max_M2, _title_m2, _title_counts);
895 name = "m3"; _m3 = createHisto1D(name,name,_nBins_M3, _min_M3, _max_M3, _title_m3, _title_counts);
896 name = "m4"; _m4 = createHisto1D(name,name,_nBins_M4, _min_M4, _max_M4, _title_m4, _title_counts);
897 name = "m5"; _m5 = createHisto1D(name,name,_nBins_M5, _min_M5, _max_M5, _title_m5, _title_counts);
898 name = "m6"; _m6 = createHisto1D(name,name,_nBins_M6, _min_M6, _max_M6, _title_m6, _title_counts);
899 name = "zV"; _vertexZ = createHisto1D(name,name,100, -10, 10, "z-Vertex (cm)", _title_counts);
902 name = "Eta"; _etadis = createHisto1F(name,name, 200, -1.0, 1.0, "#eta","counts");
903 name = "Phi"; _phidis = createHisto1F(name,name, 360, 0.0, 6.4, "#phi","counts");
904 name = "DCAz"; _dcaz = createHisto1F(name,name, 500, -5.0, 5.0, "dcaZ","counts");
905 name = "DCAxy"; _dcaxy = createHisto1F(name,name, 500, -5.0, 5.0, "dcaXY","counts");
907 //name = "Nclus1"; _Ncluster1 = createHisto1F(name,name, 200, 0, 200, "Ncluster1","counts");
908 //name = "Nclus2"; _Ncluster2 = createHisto1F(name,name, 200, 0, 200, "Ncluster2","counts");
912 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);
913 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);
915 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);
916 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);
921 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);
922 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);
923 name = n1Name+part_1_Name+vsM; _n1_1_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgN_1);
924 name = s1ptName+part_1_Name+vsM; _s1pt_1_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgSumPt_1);
925 name = n1NwName+part_1_Name+vsM; _n1Nw_1_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgN_1);
926 name = s1ptNwName+part_1_Name+vsM; _s1ptNw_1_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgSumPt_1);
928 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);
929 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);
930 name = n1Name+part_2_Name + vsM; _n1_2_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgN_2);
931 name = s1ptName+part_2_Name + vsM; _s1pt_2_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgSumPt_2);
932 name = n1NwName+part_2_Name+vsM; _n1Nw_2_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgN_1);
933 name = s1ptNwName+part_2_Name+vsM; _s1ptNw_2_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgSumPt_1);
935 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);
936 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);
937 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);
938 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);
939 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);
941 name = n2Name+pair_12_Name + vsM; _n2_12_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgN2_12);
942 name = s2PtPtName+pair_12_Name + vsM; _s2PtPt_12_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgSumPtPt_12);
943 name = s2PtNName+pair_12_Name + vsM; _s2PtN_12_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgSumPtN_12);
944 name = s2NPtName+pair_12_Name + vsM; _s2NPt_12_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgNSumPt_12);
946 name = n2NwName+pair_12_Name + vsM; _n2Nw_12_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgN2_12);
947 name = s2PtPtNwName+pair_12_Name + vsM; _s2PtPtNw_12_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgSumPtPt_12);
948 name = s2PtNNwName+pair_12_Name + vsM; _s2PtNNw_12_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgSumPtN_12);
949 name = s2NPtNwName+pair_12_Name + vsM; _s2NPtNw_12_vsM = createProfile(name,name, _nBins_M4, _min_M4, _max_M4, _title_m4, _title_AvgNSumPt_12);
951 name = "mInv"; _invMass = createHisto1F(name,name, 50, 0.41, 0.55, "M_{inv}","counts");
952 name = "mInvElec"; _invMassElec = createHisto1F(name,name, 500, 0., 1.000, "M_{inv}","counts");
955 AliInfo(" AliAnalysisTaskDptDptCorrelations::createHistoHistograms() All Done");
957 //-----------------------//
959 void AliAnalysisTaskDptDptCorrelations::finalizeHistograms()
962 AliInfo("AliAnalysisTaskDptDptCorrelations::finalizeHistograms() starting");
963 AliInfo(Form("CorrelationAnalyzers::finalizeHistograms() _eventCount : %d",int(_eventCount)));
968 fillHistoWithArray(_n1_1_vsPt, __n1_1_vsPt, _nBins_pt_1);
969 fillHistoWithArray(_n1_1_vsZVsEtaVsPhiVsPt, __n1_1_vsZEtaPhiPt, _nBins_vertexZ, _nBins_etaPhi_1, _nBins_pt_1);
970 fillHistoWithArray(_n1_2_vsPt, __n1_1_vsPt, _nBins_pt_1);
971 fillHistoWithArray(_n1_2_vsZVsEtaVsPhiVsPt, __n1_1_vsZEtaPhiPt, _nBins_vertexZ, _nBins_etaPhi_1, _nBins_pt_1);
975 fillHistoWithArray(_n1_1_vsPt, __n1_1_vsPt, _nBins_pt_1);
976 fillHistoWithArray(_n1_1_vsZVsEtaVsPhiVsPt, __n1_1_vsZEtaPhiPt, _nBins_vertexZ, _nBins_etaPhi_1, _nBins_pt_1);
977 fillHistoWithArray(_n1_2_vsPt, __n1_2_vsPt, _nBins_pt_2);
978 fillHistoWithArray(_n1_2_vsZVsEtaVsPhiVsPt, __n1_2_vsZEtaPhiPt, _nBins_vertexZ, _nBins_etaPhi_2, _nBins_pt_2);
985 fillHistoWithArray(_n1_1_vsEtaVsPhi, __n1_1_vsEtaPhi, _nBins_eta_1, _nBins_phi_1);
986 fillHistoWithArray(_s1pt_1_vsEtaVsPhi, __s1pt_1_vsEtaPhi, _nBins_eta_1, _nBins_phi_1);
987 fillHistoWithArray(_n1_2_vsEtaVsPhi, __n1_1_vsEtaPhi, _nBins_eta_1, _nBins_phi_1);
988 fillHistoWithArray(_s1pt_2_vsEtaVsPhi, __s1pt_1_vsEtaPhi, _nBins_eta_1, _nBins_phi_1);
992 fillHistoWithArray(_n1_1_vsEtaVsPhi, __n1_1_vsEtaPhi, _nBins_eta_1, _nBins_phi_1);
993 fillHistoWithArray(_s1pt_1_vsEtaVsPhi, __s1pt_1_vsEtaPhi, _nBins_eta_1, _nBins_phi_1);
994 fillHistoWithArray(_n1_2_vsEtaVsPhi, __n1_2_vsEtaPhi, _nBins_eta_2, _nBins_phi_2);
995 fillHistoWithArray(_s1pt_2_vsEtaVsPhi, __s1pt_2_vsEtaPhi, _nBins_eta_2, _nBins_phi_2);
997 fillHistoWithArray(_n2_12_vsEtaPhi, __n2_12_vsEtaPhi, _nBins_etaPhi_12);
998 fillHistoWithArray(_s2PtPt_12_vsEtaPhi, __s2ptpt_12_vsEtaPhi, _nBins_etaPhi_12);
999 fillHistoWithArray(_s2PtN_12_vsEtaPhi, __s2PtN_12_vsEtaPhi, _nBins_etaPhi_12);
1000 fillHistoWithArray(_s2NPt_12_vsEtaPhi, __s2NPt_12_vsEtaPhi, _nBins_etaPhi_12);
1001 fillHistoWithArray(_n2_12_vsPtVsPt, __n2_12_vsPtPt, _nBins_pt_1, _nBins_pt_2);
1004 AliInfo("AliAnalysisTaskDptDptCorrelations::finalizeHistograms() Done ");
1009 void AliAnalysisTaskDptDptCorrelations::UserExec(Option_t */*option*/)
1013 int iPhi, iEta, iEtaPhi, iPt, charge;
1014 float q, phi, pt, eta, corr, corrPt, px, py, pz;
1016 int id_1, q_1, iEtaPhi_1, iPt_1;
1017 float pt_1, px_1, py_1, pz_1, corr_1;
1018 int id_2, q_2, iEtaPhi_2, iPt_2;
1019 float pt_2, px_2, py_2, pz_2, corr_2;
1021 int iVertex, iVertexP1, iVertexP2;
1023 float massElecSq = 1.94797849000000016e-02;
1026 const AliAODVertex* vertex;
1029 AliAnalysisManager* manager = AliAnalysisManager::GetAnalysisManager();
1033 AliAODInputHandler* inputHandler = dynamic_cast<AliAODInputHandler*> (manager->GetInputEventHandler());
1034 if (!inputHandler) {
1038 fAODEvent = dynamic_cast<AliAODEvent*>(InputEvent());
1039 //AliAODEvent* fAODEvent = dynamic_cast<AliAODEvent*>(InputEvent());
1045 fPIDResponse =inputHandler->GetPIDResponse();
1047 AliFatal("This Task needs the PID response attached to the inputHandler");
1051 // count all events looked at here
1054 if (_eventAccounting)
1056 _eventAccounting->Fill(0);// count all calls to this function
1064 _eventAccounting->Fill(1);// count all calls to this function with a valid pointer
1065 //reset single particle counters
1067 __n1_1 = __n1_2 = __s1pt_1 = __s1pt_2 = __n1Nw_1 = __n1Nw_2 = __s1ptNw_1 = __s1ptNw_2 = 0;
1069 float v0Centr = -999.;
1070 float v0ACentr = -999.;
1071 float trkCentr = -999.;
1072 float spdCentr = -999.;
1074 float vertexX = -999;
1075 float vertexY = -999;
1076 float vertexZ = -999;
1077 //float vertexXY = -999;
1078 //float dcaZ = -999;
1079 //float dcaXY = -999;
1080 float centrality = -999;
1085 AliCentrality* centralityObject = fAODEvent->GetHeader()->GetCentralityP();
1086 if (centralityObject)
1088 //cout << "AliAnalysisTaskDptDptCorrelations::UserExec(Option_t *option) - 6" << endl;
1090 v0Centr = centralityObject->GetCentralityPercentile("V0M");
1091 v0ACentr = centralityObject->GetCentralityPercentile("V0A");
1092 trkCentr = centralityObject->GetCentralityPercentile("TRK");
1093 spdCentr = centralityObject->GetCentralityPercentile("CL1");
1097 _nTracks =fAODEvent->GetNumberOfTracks();//NEW Test
1104 _field = fAODEvent->GetMagneticField();
1107 switch (_centralityMethod)
1109 case 0: centrality = _mult0; break;
1110 case 1: centrality = _mult1; break;
1111 case 2: centrality = _mult2; break;
1112 case 3: centrality = _mult3; break;
1113 case 4: centrality = _mult4; break;
1114 case 5: centrality = _mult5; break;
1115 case 6: centrality = _mult6; break;
1116 case 7: centrality = _mult4a; break;
1120 if ( centrality < _centralityMin || centrality > _centralityMax )
1124 _eventAccounting->Fill(2);// count all events with right centrality
1126 // filter on z and xy vertex
1127 vertex = (AliAODVertex*) fAODEvent->GetPrimaryVertex();
1129 //vertex->GetXYZ(V);
1134 vertex->GetCovarianceMatrix(fCov);
1135 if(vertex->GetNContributors() > 0)
1139 vertexX = vertex->GetX();
1140 vertexY = vertex->GetY();
1141 vertexZ = vertex->GetZ();
1143 if(TMath::Abs(vertexZ) > 10)
1151 _vertexZ->Fill(vertexZ);
1153 iVertex = int((vertexZ-_min_vertexZ)/_width_vertexZ);
1154 iVertexP1 = iVertex*_nBins_etaPhiPt_1;
1155 iVertexP2 = iVertex*_nBins_etaPhiPt_2;
1156 if (iVertex<0 || iVertex>=_nBins_vertexZ)
1158 AliError("AliAnalysisTaskDptDptCorrelations::Exec(Option_t *option) iVertex<0 || iVertex>=_nBins_vertexZ ");
1161 _eventAccounting->Fill(3);// count all calls to this function with a valid pointer
1162 //======================
1164 //*********************************************************
1165 TExMap *trackMap = new TExMap();//Mapping matrix----
1167 for(Int_t i = 0; i < _nTracks; i++)
1169 AliAODTrack* aodTrack = dynamic_cast<AliAODTrack *>(fAODEvent->GetTrack(i));
1171 AliError(Form("ERROR: Could not retrieve AODtrack %d",i));
1174 Int_t gID = aodTrack->GetID();
1175 if (aodTrack->TestFilterBit(1)) trackMap->Add(gID, i);//Global tracks
1178 AliAODTrack* newAodTrack;
1180 //Track Loop starts here
1181 for (int iTrack=0; iTrack< _nTracks; iTrack++)
1183 AliAODTrack* t = dynamic_cast<AliAODTrack *>(fAODEvent->GetTrack(iTrack));
1185 AliError(Form("Could not receive track %d", iTrack));
1189 bitOK = t->TestFilterBit(_trackFilterBit);
1190 if (!bitOK) continue; //128bit or 272bit
1192 Int_t gID = t->GetID();
1193 newAodTrack = gID >= 0 ?t : fAODEvent->GetTrack(trackMap->GetValue(-1-gID));
1204 //dcaZ = t->ZAtDCA();
1206 Double_t nsigmaelectron = TMath::Abs(fPIDResponse->NumberOfSigmasTPC(newAodTrack,(AliPID::EParticleType)AliPID::kElectron));
1207 Double_t nsigmapion = TMath::Abs(fPIDResponse->NumberOfSigmasTPC(newAodTrack,(AliPID::EParticleType)AliPID::kPion));
1208 Double_t nsigmakaon = TMath::Abs(fPIDResponse->NumberOfSigmasTPC(newAodTrack,(AliPID::EParticleType)AliPID::kKaon));
1209 Double_t nsigmaproton = TMath::Abs(fPIDResponse->NumberOfSigmasTPC(newAodTrack,(AliPID::EParticleType)AliPID::kProton));
1211 //nsigma cut to reject electron
1213 if(nsigmaelectron < fNSigmaCut
1214 && nsigmapion > fNSigmaCut
1215 && nsigmakaon > fNSigmaCut
1216 && nsigmaproton > fNSigmaCut ) continue;
1218 if(charge == 0) continue;
1220 // Kinematics cuts used
1221 if( pt < _min_pt_1 || pt > _max_pt_1) continue;
1222 if( eta < _min_eta_1 || eta > _max_eta_1) continue;
1225 newAodTrack->GetXYZ(pos);
1227 Double_t DCAX = pos[0] - vertexX;
1228 Double_t DCAY = pos[1] - vertexY;
1229 Double_t DCAZ = pos[2] - vertexZ;
1231 Double_t DCAXY = TMath::Sqrt((DCAX*DCAX) + (DCAY*DCAY));
1233 if (DCAZ < _dcaZMin ||
1235 DCAXY > _dcaXYMax ) continue;
1236 //==== QA ===========================
1238 _dcaxy->Fill(DCAXY);
1241 //===================================
1242 //*************************************************
1245 if (_requestedCharge_1 == charge)
1248 iPhi = int( phi/_width_phi_1);
1250 if (iPhi<0 || iPhi>=_nBins_phi_1 )
1252 AliWarning("AliAnalysisTaskDptDptCorrelations::analyze() iPhi<0 || iPhi>=_nBins_phi_1");
1256 iEta = int((eta-_min_eta_1)/_width_eta_1);
1257 if (iEta<0 || iEta>=_nBins_eta_1)
1259 AliWarning(Form("AliAnalysisTaskDptDptCorrelations::analyze(AliceEvent * event) Mismatched iEta: %d", iEta));
1262 iPt = int((pt -_min_pt_1 )/_width_pt_1 );
1263 if (iPt<0 || iPt >=_nBins_pt_1)
1265 AliWarning(Form("AliAnalysisTaskDptDptCorrelations::analyze(AliceEvent * event) Mismatched iPt: %d",iPt));
1268 iEtaPhi = iEta*_nBins_phi_1+iPhi;
1269 iZEtaPhiPt = iVertexP1 + iEtaPhi*_nBins_pt_1 + iPt;
1271 if (_correctionWeight_1)
1272 corr = _correctionWeight_1[iZEtaPhiPt];
1275 if (iZEtaPhiPt<0 || iZEtaPhiPt>=_nBins_zEtaPhiPt_1)
1277 AliWarning("AliAnalysisTaskDptDptCorrelations::analyze(AliceEvent * event) iZEtaPhiPt<0 || iZEtaPhiPt>=_nBins_zEtaPhiPt_1");
1285 __n1_1_vsPt[iPt] += corr; //cout << "step 15" << endl;
1286 __n1_1_vsZEtaPhiPt[iZEtaPhiPt] += corr; //cout << "step 12" << endl;
1293 _charge_1[k1] = charge;
1294 _iEtaPhi_1[k1] = iEtaPhi;
1300 _correction_1[k1] = corr;
1302 __n1_1_vsEtaPhi[iEtaPhi] += corr;
1304 __s1pt_1_vsEtaPhi[iEtaPhi] += corrPt;
1310 AliError(Form("AliAnalysisTaskDptDptCorrelations::analyze(AliceEvent * event) k1 >=arraySize; arraySize: %d",arraySize));
1316 if (!_sameFilter && _requestedCharge_2 == charge)
1319 iPhi = int( phi/_width_phi_2);
1321 if (iPhi<0 || iPhi>=_nBins_phi_2 )
1323 AliWarning("AliAnalysisTaskDptDptCorrelations::analyze() iPhi<0 || iPhi>=_nBins_phi_1");
1327 iEta = int((eta-_min_eta_2)/_width_eta_2);
1328 if (iEta<0 || iEta>=_nBins_eta_2)
1330 AliWarning(Form("AliAnalysisTaskDptDptCorrelations::analyze(AliceEvent * event) Mismatched iEta: %d", iEta));
1333 iPt = int((pt -_min_pt_2 )/_width_pt_2 );
1334 if (iPt<0 || iPt >=_nBins_pt_2)
1336 AliWarning(Form("AliAnalysisTaskDptDptCorrelations::analyze(AliceEvent * event) Mismatched iPt: %d",iPt));
1340 iEtaPhi = iEta*_nBins_phi_2+iPhi;
1341 iZEtaPhiPt = iVertexP2 + iEtaPhi*_nBins_pt_2 + iPt;
1342 if (iZEtaPhiPt<0 || iZEtaPhiPt>=_nBins_zEtaPhiPt_2)
1344 AliWarning("AliAnalysisTaskDptDptCorrelations::analyze(AliceEvent * event) iZEtaPhiPt<0 || iZEtaPhiPt>=_nBins_zEtaPhiPt_2");
1349 if (_correctionWeight_2)
1350 corr = _correctionWeight_2[iZEtaPhiPt];
1356 __n1_2_vsPt[iPt] += corr; //cout << "step 15" << endl;
1357 __n1_2_vsZEtaPhiPt[iZEtaPhiPt] += corr; //cout << "step 12" << endl;
1362 _id_2[k2] = iTrack; //cout << "step 1" << endl;
1363 _charge_2[k2] = charge; //cout << "step 2" << endl;
1364 _iEtaPhi_2[k2] = iEtaPhi; //cout << "step 3" << endl;
1365 _iPt_2[k2] = iPt; //cout << "step 4" << endl;
1366 _pt_2[k2] = pt; //cout << "step 5" << endl;
1367 _px_2[k2] = px; //cout << "step 6" << endl;
1368 _py_2[k2] = py; //cout << "step 7" << endl;
1369 _pz_2[k2] = pz; //cout << "step 8" << endl;
1370 _correction_2[k2] = corr; //cout << "step 9" << endl;
1371 __n1_2 += corr; //cout << "step 10" << endl;
1372 __s1pt_2 += corrPt; //cout << "step 13" << endl;
1374 __n1_2_vsEtaPhi[iEtaPhi] += corr; //cout << "step 11" << endl;
1375 __s1pt_2_vsEtaPhi[iEtaPhi] += corrPt; //cout << "step 14" << endl;
1380 AliWarning(Form("-W- k2 >=arraySize; arraySize: %d",arraySize));
1385 //cout << "done with track" << endl;
1391 //cout << "Filling histograms now" << endl;
1399 //_vertexZ->Fill(vertexZ);
1403 // nothing to do here.
1409 _n1_1_vsM->Fill(centrality, __n1_1);
1410 _s1pt_1_vsM->Fill(centrality, __s1pt_1);
1411 _n1Nw_1_vsM->Fill(centrality, __n1Nw_1);
1412 _s1ptNw_1_vsM->Fill(centrality, __s1ptNw_1);
1413 _n1_2_vsM->Fill(centrality, __n1_1);
1414 _s1pt_2_vsM->Fill(centrality, __s1pt_1);
1415 _n1Nw_2_vsM->Fill(centrality, __n1Nw_1);
1416 _s1ptNw_2_vsM->Fill(centrality, __s1ptNw_1);
1417 // reset pair counters
1418 __n2_12 = __s2ptpt_12 = __s2NPt_12 = __s2PtN_12 = 0;
1419 __n2Nw_12 = __s2ptptNw_12 = __s2NPtNw_12 = __s2PtNNw_12 = 0;
1422 for (int i1=0; i1<k1; i1++)
1424 ////cout << " i1:" << i1 << endl;
1425 id_1 = _id_1[i1]; ////cout << " id_1:" << id_1 << endl;
1426 q_1 = _charge_1[i1]; ////cout << " q_1:" << q_1 << endl;
1427 iEtaPhi_1 = _iEtaPhi_1[i1]; ////cout << " iEtaPhi_1:" << iEtaPhi_1 << endl;
1428 iPt_1 = _iPt_1[i1]; ////cout << " iPt_1:" << iPt_1 << endl;
1429 corr_1 = _correction_1[i1]; ////cout << " corr_1:" << corr_1 << endl;
1430 pt_1 = _pt_1[i1]; ////cout << " pt_1:" << pt_1 << endl;
1431 //dedx_1 = _dedx_1[i1]; ////cout << " dedx_1:" << dedx_1 << endl;
1433 for (int i2=i1+1; i2<k1; i2++)
1435 ////cout << " i2:" << i2 << endl;
1436 id_2 = _id_1[i2]; ////cout << " id_2:" << id_2 << endl;
1439 q_2 = _charge_1[i2]; ////cout << " q_1:" << q_1 << endl;
1440 iEtaPhi_2 = _iEtaPhi_1[i2]; ////cout << " iEtaPhi_1:" << iEtaPhi_1 << endl;
1441 iPt_2 = _iPt_1[i2]; ////cout << " iPt_1:" << iPt_1 << endl;
1442 corr_2 = _correction_1[i2]; ////cout << " corr_1:" << corr_1 << endl;
1443 pt_2 = _pt_1[i2]; ////cout << " pt_1:" << pt_1 << endl;
1444 //dedx_2 = _dedx_1[i2]; ////cout << " dedx_2:" << dedx_2 << endl;
1445 corr = corr_1*corr_2;
1446 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))
1448 ij = iEtaPhi_1*_nBins_etaPhi_1 + iEtaPhi_2; ////cout << " ij:" << ij<< endl;
1450 else // swap particles
1452 ij = iEtaPhi_2*_nBins_etaPhi_1 + iEtaPhi_1; ////cout << " ij:" << ij<< endl;
1456 __n2_12_vsEtaPhi[ij] += corr;
1458 __s2ptpt_12 += corr*ptpt;
1459 __s2PtN_12 += corr*pt_1;
1460 __s2NPt_12 += corr*pt_2;
1461 __s2ptpt_12_vsEtaPhi[ij] += corr*ptpt;
1462 __s2PtN_12_vsEtaPhi[ij] += corr*pt_1;
1463 __s2NPt_12_vsEtaPhi[ij] += corr*pt_2;
1464 __n2_12_vsPtPt[iPt_1*_nBins_pt_2 + iPt_2] += corr;
1467 __s2ptptNw_12 += ptpt;
1468 __s2PtNNw_12 += pt_1;
1469 __s2NPtNw_12 += pt_2;
1477 for (int i1=0; i1<k1; i1++)
1479 ////cout << " i1:" << i1 << endl;
1480 id_1 = _id_1[i1]; ////cout << " id_1:" << id_1 << endl;
1481 q_1 = _charge_1[i1]; ////cout << " q_1:" << q_1 << endl;
1482 iEtaPhi_1 = _iEtaPhi_1[i1]; ////cout << " iEtaPhi_1:" << iEtaPhi_1 << endl;
1483 iPt_1 = _iPt_1[i1]; ////cout << " iPt_1:" << iPt_1 << endl;
1484 corr_1 = _correction_1[i1]; ////cout << " corr_1:" << corr_1 << endl;
1485 pt_1 = _pt_1[i1]; ////cout << " pt_1:" << pt_1 << endl;
1486 //dedx_1 = _dedx_1[i1]; ////cout << " dedx_1:" << dedx_1 << endl;
1488 for (int i2=i1+1; i2<k1; i2++)
1490 ////cout << " i2:" << i2 << endl;
1491 id_2 = _id_1[i2]; ////cout << " id_2:" << id_2 << endl;
1494 q_2 = _charge_1[i2]; ////cout << " q_2:" << q_2 << endl;
1495 iEtaPhi_2 = _iEtaPhi_1[i2]; ////cout << " iEtaPhi_2:" << iEtaPhi_2 << endl;
1496 iPt_2 = _iPt_1[i2]; ////cout << " iPt_2:" << iPt_2 << endl;
1497 corr_2 = _correction_1[i2]; ////cout << " corr_2:" << corr_2 << endl;
1498 pt_2 = _pt_1[i2]; ////cout << " pt_2:" << pt_2 << endl;
1499 //dedx_2 = _dedx_1[i2]; ////cout << " dedx_2:" << dedx_2 << endl;
1500 corr = corr_1*corr_2;
1501 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))
1503 ij = iEtaPhi_1*_nBins_etaPhi_1 + iEtaPhi_2; ////cout << " ij:" << ij<< endl;
1505 else // swap particles
1507 ij = iEtaPhi_2*_nBins_etaPhi_1 + iEtaPhi_1; ////cout << " ij:" << ij<< endl;
1511 __n2_12_vsEtaPhi[ij] += corr;
1513 __s2ptpt_12 += corr*ptpt;
1514 __s2PtN_12 += corr*pt_1;
1515 __s2NPt_12 += corr*pt_2;
1516 __s2ptpt_12_vsEtaPhi[ij] += corr*ptpt;
1517 __s2PtN_12_vsEtaPhi[ij] += corr*pt_1;
1518 __s2NPt_12_vsEtaPhi[ij] += corr*pt_2;
1519 __n2_12_vsPtPt[iPt_1*_nBins_pt_2 + iPt_2] += corr;
1522 __s2ptptNw_12 += ptpt;
1523 __s2PtNNw_12 += pt_1;
1524 __s2NPtNw_12 += pt_2;
1531 else // filter 1 and 2 are different -- must do all particle pairs...
1533 _n1_1_vsM->Fill(centrality, __n1_1);
1534 _s1pt_1_vsM->Fill(centrality, __s1pt_1);
1535 _n1Nw_1_vsM->Fill(centrality, __n1Nw_1);
1536 _s1ptNw_1_vsM->Fill(centrality, __s1ptNw_1);
1537 _n1_2_vsM->Fill(centrality, __n1_2);
1538 _s1pt_2_vsM->Fill(centrality, __s1pt_2);
1539 _n1Nw_2_vsM->Fill(centrality, __n1Nw_2);
1540 _s1ptNw_2_vsM->Fill(centrality, __s1ptNw_2);
1541 // reset pair counters
1542 __n2_12 = __s2ptpt_12 = __s2NPt_12 = __s2PtN_12 = 0;
1543 __n2Nw_12 = __s2ptptNw_12 = __s2NPtNw_12 = __s2PtNNw_12 = 0;
1544 for (int i1=0; i1<k1; i1++)
1546 ////cout << " i1:" << i1 << endl;
1547 id_1 = _id_1[i1]; ////cout << " id_1:" << id_1 << endl;
1548 q_1 = _charge_1[i1]; ////cout << " q_1:" << q_1 << endl;
1549 iEtaPhi_1 = _iEtaPhi_1[i1]; ////cout << " iEtaPhi_1:" << iEtaPhi_1 << endl;
1550 iPt_1 = _iPt_1[i1]; ////cout << " iPt_1:" << iPt_1 << endl;
1551 corr_1 = _correction_1[i1]; ////cout << " corr_1:" << corr_1 << endl;
1552 pt_1 = _pt_1[i1]; ////cout << " pt_1:" << pt_1 << endl;
1553 px_1 = _px_1[i1]; ////cout << " px_1:" << px_1 << endl;
1554 py_1 = _py_1[i1]; ////cout << " py_1:" << py_1 << endl;
1555 pz_1 = _pz_1[i1]; ////cout << " pz_1:" << pz_1 << endl;
1556 //dedx_1 = _dedx_1[i1]; ////cout << " dedx_1:" << dedx_1 << endl;
1559 for (int i2=0; i2<k2; i2++)
1561 ////cout << " i2:" << i2 << endl;
1562 id_2 = _id_2[i2]; ////cout << " id_2:" << id_2 << endl;
1563 if (id_1!=id_2) // exclude auto correlation
1565 q_2 = _charge_2[i2]; ////cout << " q_2:" << q_2 << endl;
1566 iEtaPhi_2 = _iEtaPhi_2[i2]; ////cout << " iEtaPhi_2:" << iEtaPhi_2 << endl;
1567 iPt_2 = _iPt_2[i2]; ////cout << " iPt_2:" << iPt_2 << endl;
1568 corr_2 = _correction_2[i2]; ////cout << " corr_2:" << corr_2 << endl;
1569 pt_2 = _pt_2[i2]; ////cout << " pt_2:" << pt_2 << endl;
1570 px_2 = _px_2[i2]; ////cout << " px_2:" << px_2 << endl;
1571 py_2 = _py_2[i2]; ////cout << " py_2:" << py_2 << endl;
1572 pz_2 = _pz_2[i2]; ////cout << " pz_2:" << pz_2 << endl;
1573 //dedx_2 = _dedx_2[i2]; ////cout << " dedx_2:" << dedx_2 << endl;
1576 if (_rejectPairConversion)
1578 float e1Sq = massElecSq + pt_1*pt_1 + pz_1*pz_1;
1579 float e2Sq = massElecSq + pt_2*pt_2 + pz_2*pz_2;
1580 float mInvSq = 2*(massElecSq + sqrt(e1Sq*e2Sq) - px_1*px_2 - py_1*py_2 - pz_1*pz_2 );
1581 float mInv = sqrt(mInvSq);
1582 _invMass->Fill(mInv);
1585 corr = corr_1*corr_2;
1586 ij = iEtaPhi_1*_nBins_etaPhi_1 + iEtaPhi_2; ////cout << " ij:" << ij<< endl;
1588 __n2_12_vsEtaPhi[ij] += corr;
1590 __s2ptpt_12 += corr*ptpt;
1591 __s2PtN_12 += corr*pt_1;
1592 __s2NPt_12 += corr*pt_2;
1593 __s2ptpt_12_vsEtaPhi[ij] += corr*ptpt;
1594 __s2PtN_12_vsEtaPhi[ij] += corr*pt_1;
1595 __s2NPt_12_vsEtaPhi[ij] += corr*pt_2;
1596 __n2_12_vsPtPt[iPt_1*_nBins_pt_2 + iPt_2] += corr;
1598 __s2ptptNw_12 += ptpt;
1599 __s2PtNNw_12 += pt_1;
1600 __s2NPtNw_12 += pt_2;
1607 _n2_12_vsM->Fill(centrality, __n2_12);
1608 _s2PtPt_12_vsM->Fill(centrality, __s2ptpt_12);
1609 _s2PtN_12_vsM->Fill(centrality, __s2NPt_12);
1610 _s2NPt_12_vsM->Fill(centrality, __s2PtN_12);
1612 _n2Nw_12_vsM->Fill(centrality, __n2Nw_12);
1613 _s2PtPtNw_12_vsM->Fill(centrality, __s2ptptNw_12);
1614 _s2PtNNw_12_vsM->Fill(centrality, __s2NPtNw_12);
1615 _s2NPtNw_12_vsM->Fill(centrality, __s2PtNNw_12);
1620 AliInfo("AliAnalysisTaskDptDptCorrelations::UserExec() -----------------Event Done ");
1621 PostData(0,_outputHistoList);
1625 void AliAnalysisTaskDptDptCorrelations::FinishTaskOutput()
1627 AliInfo("AliAnalysisTaskDptDptCorrelations::FinishTaskOutput() Starting.");
1628 Printf("= 0 ====================================================================");
1629 finalizeHistograms();
1630 AliInfo("= 1 ====================================================================");
1631 PostData(0,_outputHistoList);
1632 AliInfo("= 2 ====================================================================");
1633 AliInfo("AliAnalysisTaskDptDptCorrelations::FinishTaskOutput() Done.");
1636 void AliAnalysisTaskDptDptCorrelations::Terminate(Option_t* /*option*/)
1638 AliInfo("AliAnalysisTaskDptDptCorrelations::Terminate() Starting/Done.");
1643 //===================================================================================================
1644 void AliAnalysisTaskDptDptCorrelations::fillHistoWithArray(TH1 * h, float * array, int size)
1647 float v1, ev1, v2, ev2, sum, esum;
1648 for (i=0, i1=1; i<size; ++i,++i1)
1650 v1 = array[i]; ev1 = sqrt(v1);
1651 v2 = h->GetBinContent(i1);
1652 ev2 = h->GetBinError(i1);
1654 esum = sqrt(ev1*ev1+ev2*ev2);
1655 h->SetBinContent(i1,sum);
1656 h->SetBinError(i1,esum);
1660 void AliAnalysisTaskDptDptCorrelations::fillHistoWithArray(TH2 * h, float * array, int size1, int size2)
1664 float v1, ev1, v2, ev2, sum, esum;
1665 for (i=0, i1=1; i<size1; ++i,++i1)
1667 for (j=0, j1=1; j<size2; ++j,++j1)
1669 v1 = array[i*size2+j]; ev1 = sqrt(v1);
1670 v2 = h->GetBinContent(i1,j1);
1671 ev2 = h->GetBinError(i1,j1);
1673 esum = sqrt(ev1*ev1+ev2*ev2);
1674 h->SetBinContent(i1,j1,sum);
1675 h->SetBinError(i1,j1,esum);
1680 void AliAnalysisTaskDptDptCorrelations::fillHistoWithArray(TH3 * h, float * array, int size1, int size2, int size3)
1685 float v1, ev1, v2, ev2, sum, esum;
1686 int size23 = size2*size3;
1687 for (i=0, i1=1; i<size1; ++i,++i1)
1689 for (j=0, j1=1; j<size2; ++j,++j1)
1691 for (k=0, k1=1; k<size3; ++k,++k1)
1693 v1 = array[i*size23+j*size3+k]; ev1 = sqrt(v1);
1694 v2 = h->GetBinContent(i1,j1,k1);
1695 ev2 = h->GetBinError(i1,j1,k1);
1697 esum = sqrt(ev1*ev1+ev2*ev2);
1698 h->SetBinContent(i1,j1,k1,sum);
1699 h->SetBinError(i1,j1,k1,esum);
1705 void AliAnalysisTaskDptDptCorrelations::fillHistoWithArray(TH1 * h, double * array, int size)
1708 double v1, ev1, v2, ev2, sum, esum;
1709 for (i=0, i1=1; i<size; ++i,++i1)
1711 v1 = array[i]; ev1 = sqrt(v1);
1712 v2 = h->GetBinContent(i1);
1713 ev2 = h->GetBinError(i1);
1715 esum = sqrt(ev1*ev1+ev2*ev2);
1716 h->SetBinContent(i1,sum);
1717 h->SetBinError(i1,esum);
1721 void AliAnalysisTaskDptDptCorrelations::fillHistoWithArray(TH2 * h, double * array, int size1, int size2)
1725 double v1, ev1, v2, ev2, sum, esum;
1726 for (i=0, i1=1; i<size1; ++i,++i1)
1728 for (j=0, j1=1; j<size2; ++j,++j1)
1730 v1 = array[i*size2+j]; ev1 = sqrt(v1);
1731 v2 = h->GetBinContent(i1,j1);
1732 ev2 = h->GetBinError(i1,j1);
1734 esum = sqrt(ev1*ev1+ev2*ev2);
1735 h->SetBinContent(i1,j1,sum);
1736 h->SetBinError(i1,j1,esum);
1741 void AliAnalysisTaskDptDptCorrelations::fillHistoWithArray(TH3 * h, double * array, int size1, int size2, int size3)
1746 double v1, ev1, v2, ev2, sum, esum;
1747 int size23 = size2*size3;
1748 for (i=0, i1=1; i<size1; ++i,++i1)
1750 for (j=0, j1=1; j<size2; ++j,++j1)
1752 for (k=0, k1=1; k<size3; ++k,++k1)
1754 v1 = array[i*size23+j*size3+k]; ev1 = sqrt(v1);
1755 v2 = h->GetBinContent(i1,j1,k1);
1756 ev2 = h->GetBinError(i1,j1,k1);
1758 esum = sqrt(ev1*ev1+ev2*ev2);
1759 h->SetBinContent(i1,j1,k1,sum);
1760 h->SetBinError(i1,j1,k1,esum);
1766 //________________________________________________________________________
1767 double * AliAnalysisTaskDptDptCorrelations::getDoubleArray(int size, double v)
1769 /// Allocate an array of type double with n values
1770 /// Initialize the array to the given value
1771 double * array = new double [size];
1772 for (int i=0;i<size;++i) array[i]=v;
1776 //________________________________________________________________________
1777 float * AliAnalysisTaskDptDptCorrelations::getFloatArray(int size, float v)
1779 /// Allocate an array of type float with n values
1780 /// Initialize the array to the given value
1781 float * array = new float [size];
1782 for (int i=0;i<size;++i) array[i]=v;
1787 //________________________________________________________________________
1788 TH1D * AliAnalysisTaskDptDptCorrelations::createHisto1D(const TString & name, const TString & title,
1789 int n, double xMin, double xMax,
1790 const TString & xTitle, const TString & yTitle)
1792 //CreateHisto new 1D historgram
1793 AliInfo(Form("createHisto 1D histo %s nBins: %d xMin: %f xMax: %f",name.Data(),n,xMin,xMax));
1794 TH1D * h = new TH1D(name,title,n,xMin,xMax);
1795 h->GetXaxis()->SetTitle(xTitle);
1796 h->GetYaxis()->SetTitle(yTitle);
1802 //________________________________________________________________________
1803 TH1D * AliAnalysisTaskDptDptCorrelations::createHisto1D(const TString & name, const TString & title,
1804 int n, double * bins,
1805 const TString & xTitle, const TString & yTitle)
1807 AliInfo(Form("createHisto 1D histo %s with %d non uniform nBins",name.Data(),n));
1808 TH1D * h = new TH1D(name,title,n,bins);
1809 h->GetXaxis()->SetTitle(xTitle);
1810 h->GetYaxis()->SetTitle(yTitle);
1816 //________________________________________________________________________
1817 TH2D * AliAnalysisTaskDptDptCorrelations::createHisto2D(const TString & name, const TString & title,
1818 int nx, double xMin, double xMax, int ny, double yMin, double yMax,
1819 const TString & xTitle, const TString & yTitle, const TString & zTitle)
1821 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));
1822 TH2D * h = new TH2D(name,title,nx,xMin,xMax,ny,yMin,yMax);
1823 h->GetXaxis()->SetTitle(xTitle);
1824 h->GetYaxis()->SetTitle(yTitle);
1825 h->GetZaxis()->SetTitle(zTitle);
1830 //________________________________________________________________________
1831 TH2D * AliAnalysisTaskDptDptCorrelations::createHisto2D(const TString & name, const TString & title,
1832 int nx, double* xbins, int ny, double yMin, double yMax,
1833 const TString & xTitle, const TString & yTitle, const TString & zTitle)
1835 AliInfo(Form("createHisto 2D histo %s with %d non uniform nBins",name.Data(),nx));
1837 h = new TH2D(name,title,nx,xbins,ny,yMin,yMax);
1838 h->GetXaxis()->SetTitle(xTitle);
1839 h->GetYaxis()->SetTitle(yTitle);
1840 h->GetZaxis()->SetTitle(zTitle);
1846 //________________________________________________________________________
1847 TH1F * AliAnalysisTaskDptDptCorrelations::createHisto1F(const TString & name, const TString & title,
1848 int n, double xMin, double xMax,
1849 const TString & xTitle, const TString & yTitle)
1851 //CreateHisto new 1D historgram
1852 AliInfo(Form("createHisto 1D histo %s nBins: %d xMin: %f xMax: %f",name.Data(),n,xMin,xMax));
1853 TH1F * h = new TH1F(name,title,n,xMin,xMax);
1854 h->GetXaxis()->SetTitle(xTitle);
1855 h->GetYaxis()->SetTitle(yTitle);
1861 //________________________________________________________________________
1862 TH1F * AliAnalysisTaskDptDptCorrelations::createHisto1F(const TString & name, const TString & title,
1863 int n, double * bins,
1864 const TString & xTitle, const TString & yTitle)
1866 AliInfo(Form("createHisto 1D histo %s with %d non uniform nBins",name.Data(),n));
1867 TH1F * h = new TH1F(name,title,n,bins);
1868 h->GetXaxis()->SetTitle(xTitle);
1869 h->GetYaxis()->SetTitle(yTitle);
1875 //________________________________________________________________________
1876 TH2F * AliAnalysisTaskDptDptCorrelations::createHisto2F(const TString & name, const TString & title,
1877 int nx, double xMin, double xMax, int ny, double yMin, double yMax,
1878 const TString & xTitle, const TString & yTitle, const TString & zTitle)
1880 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));
1881 TH2F * h = new TH2F(name,title,nx,xMin,xMax,ny,yMin,yMax);
1882 h->GetXaxis()->SetTitle(xTitle);
1883 h->GetYaxis()->SetTitle(yTitle);
1884 h->GetZaxis()->SetTitle(zTitle);
1889 //________________________________________________________________________
1890 TH2F * AliAnalysisTaskDptDptCorrelations::createHisto2F(const TString & name, const TString & title,
1891 int nx, double* xbins, int ny, double yMin, double yMax,
1892 const TString & xTitle, const TString & yTitle, const TString & zTitle)
1894 AliInfo(Form("createHisto 2D histo %s with %d non uniform nBins",name.Data(),nx));
1896 h = new TH2F(name,title,nx,xbins,ny,yMin,yMax);
1897 h->GetXaxis()->SetTitle(xTitle);
1898 h->GetYaxis()->SetTitle(yTitle);
1899 h->GetZaxis()->SetTitle(zTitle);
1905 //________________________________________________________________________
1906 TH3F * AliAnalysisTaskDptDptCorrelations::createHisto3F(const TString & name, const TString & title,
1907 int nx, double xMin, double xMax,
1908 int ny, double yMin, double yMax,
1909 int nz, double zMin, double zMax,
1910 const TString & xTitle, const TString & yTitle, const TString & zTitle)
1912 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));
1913 TH3F * h = new TH3F(name,title,nx,xMin,xMax,ny,yMin,yMax,nz,zMin,zMax);
1914 h->GetXaxis()->SetTitle(xTitle);
1915 h->GetYaxis()->SetTitle(yTitle);
1916 h->GetZaxis()->SetTitle(zTitle);
1922 //________________________________________________________________________
1923 TProfile * AliAnalysisTaskDptDptCorrelations::createProfile(const TString & name, const TString & description,
1924 int nx,double xMin,double xMax,
1925 const TString & xTitle, const TString & yTitle)
1927 AliInfo(Form("createHisto 1D profile %s nBins: %d xMin: %f10.4 xMax: %f10.4",name.Data(),nx,xMin,xMax));
1928 TProfile * h = new TProfile(name,description,nx,xMin,xMax);
1929 h->GetXaxis()->SetTitle(xTitle);
1930 h->GetYaxis()->SetTitle(yTitle);
1935 //________________________________________________________________________
1936 TProfile * AliAnalysisTaskDptDptCorrelations::createProfile(const TString & name,const TString & description,
1937 int nx, double* bins,
1938 const TString & xTitle, const TString & yTitle)
1940 AliInfo(Form("createHisto 1D profile %s with %d non uniform bins",name.Data(),nx));
1941 TProfile * h = new TProfile(name,description,nx,bins);
1942 h->GetXaxis()->SetTitle(xTitle);
1943 h->GetYaxis()->SetTitle(yTitle);
1949 void AliAnalysisTaskDptDptCorrelations::addToList(TH1 *h)
1951 if (_outputHistoList)
1953 _outputHistoList->Add(h);
1956 AliInfo("addToList(TH1 *h) _outputHistoList is null!!!!! Should abort ship");