1 /**************************************************************************
2 * Author: Panos Christakoglou. *
3 * Contributors are mentioned in the code where appropriate. *
5 * Permission to use, copy, modify and distribute this software and its *
6 * documentation strictly for non-commercial purposes is hereby granted *
7 * without fee, provided that the above copyright notice appears in all *
8 * copies and that both the copyright notice and this permission notice *
9 * appear in the supporting documentation. The authors make no claims *
10 * about the suitability of this software for any purpose. It is *
11 * provided "as is" without express or implied warranty. *
12 **************************************************************************/
16 //-----------------------------------------------------------------
17 // Balance Function class
18 // This is the class to deal with the Balance Function analysis
19 // Origin: Panos Christakoglou, UOA-CERN, Panos.Christakoglou@cern.ch
20 //-----------------------------------------------------------------
24 #include <Riostream.h>
28 #include <TLorentzVector.h>
29 #include <TObjArray.h>
30 #include <TGraphErrors.h>
33 #include "AliVParticle.h"
34 #include "AliMCParticle.h"
35 #include "AliESDtrack.h"
36 #include "AliAODTrack.h"
38 #include "AliBalance.h"
42 //____________________________________________________________________//
43 AliBalance::AliBalance() :
47 bConversionCut(kFALSE),
48 fAnalysisLevel("ESD"),
54 // Default constructor
56 for(Int_t i = 0; i < ANALYSIS_TYPES; i++){
58 fNumberOfBins[i] = 180;
66 fNumberOfBins[i] = 20;
72 fP2Step[i] = TMath::Abs(fP2Start - fP2Stop) / (Double_t)fNumberOfBins[i];
79 for(Int_t j = 0; j < MAXIMUM_NUMBER_OF_STEPS; j++) {
99 //____________________________________________________________________//
100 AliBalance::AliBalance(const AliBalance& balance):
102 bShuffle(balance.bShuffle),
103 bHBTcut(balance.bHBTcut),
104 bConversionCut(balance.bConversionCut),
105 fAnalysisLevel(balance.fAnalysisLevel),
106 fAnalyzedEvents(balance.fAnalyzedEvents),
107 fCentralityId(balance.fCentralityId),
108 fCentStart(balance.fCentStart),
109 fCentStop(balance.fCentStop) {
111 for(Int_t i = 0; i < ANALYSIS_TYPES; i++){
112 fNn[i] = balance.fNn[i];
113 fNp[i] = balance.fNp[i];
115 fP1Start[i] = balance.fP1Start[i];
116 fP1Stop[i] = balance.fP1Stop[i];
117 fNumberOfBins[i] = balance.fNumberOfBins[i];
118 fP2Start[i] = balance.fP2Start[i];
119 fP2Stop[i] = balance.fP2Stop[i];
120 fP2Step[i] = balance.fP2Step[i];
121 fCentStart = balance.fCentStart;
122 fCentStop = balance.fCentStop;
124 fHistP[i] = balance.fHistP[i];
125 fHistN[i] = balance.fHistN[i];
126 fHistPN[i] = balance.fHistPN[i];
127 fHistNP[i] = balance.fHistNP[i];
128 fHistPP[i] = balance.fHistPP[i];
129 fHistNN[i] = balance.fHistNN[i];
131 for(Int_t j = 0; j < MAXIMUM_NUMBER_OF_STEPS; j++) {
143 //____________________________________________________________________//
144 AliBalance::~AliBalance() {
147 for(Int_t i = 0; i < ANALYSIS_TYPES; i++){
159 //____________________________________________________________________//
160 void AliBalance::SetInterval(Int_t iAnalysisType,
161 Double_t p1Start, Double_t p1Stop,
162 Int_t ibins, Double_t p2Start, Double_t p2Stop) {
163 // Sets the analyzed interval.
164 // Set the same Information for all analyses
166 if(iAnalysisType == -1){
167 for(Int_t i = 0; i < ANALYSIS_TYPES; i++){
168 fP1Start[i] = p1Start;
170 fNumberOfBins[i] = ibins;
171 fP2Start[i] = p2Start;
173 fP2Step[i] = TMath::Abs(p2Start - p2Stop) / (Double_t)fNumberOfBins[i];
176 // Set the Information for one analysis
177 else if((iAnalysisType > -1) && (iAnalysisType < ANALYSIS_TYPES)) {
178 fP1Start[iAnalysisType] = p1Start;
179 fP1Stop[iAnalysisType] = p1Stop;
180 fNumberOfBins[iAnalysisType] = ibins;
181 fP2Start[iAnalysisType] = p2Start;
182 fP2Stop[iAnalysisType] = p2Stop;
183 fP2Step[iAnalysisType] = TMath::Abs(p2Start - p2Stop) / (Double_t)fNumberOfBins[iAnalysisType];
186 AliError("Wrong ANALYSIS number!");
190 //____________________________________________________________________//
191 void AliBalance::InitHistograms() {
192 //Initialize the histograms
194 for(Int_t iAnalysisType = 0; iAnalysisType < ANALYSIS_TYPES; iAnalysisType++) {
195 histName = "fHistP"; histName += gBFAnalysisType[iAnalysisType];
196 if(bShuffle) histName.Append("_shuffle");
197 if(fCentralityId) histName += fCentralityId.Data();
198 fHistP[iAnalysisType] = new TH2D(histName.Data(),"",fCentStop-fCentStart,fCentStart,fCentStop,100,fP1Start[iAnalysisType],fP1Stop[iAnalysisType]);
200 histName = "fHistN"; histName += gBFAnalysisType[iAnalysisType];
201 if(bShuffle) histName.Append("_shuffle");
202 if(fCentralityId) histName += fCentralityId.Data();
203 fHistN[iAnalysisType] = new TH2D(histName.Data(),"",fCentStop-fCentStart,fCentStart,fCentStop,100,fP1Start[iAnalysisType],fP1Stop[iAnalysisType]);
205 histName = "fHistPN"; histName += gBFAnalysisType[iAnalysisType];
206 if(bShuffle) histName.Append("_shuffle");
207 if(fCentralityId) histName += fCentralityId.Data();
208 fHistPN[iAnalysisType] = new TH2D(histName.Data(),"",fCentStop-fCentStart,fCentStart,fCentStop,fNumberOfBins[iAnalysisType],fP2Start[iAnalysisType],fP2Stop[iAnalysisType]);
210 histName = "fHistNP"; histName += gBFAnalysisType[iAnalysisType];
211 if(bShuffle) histName.Append("_shuffle");
212 if(fCentralityId) histName += fCentralityId.Data();
213 fHistNP[iAnalysisType] = new TH2D(histName.Data(),"",fCentStop-fCentStart,fCentStart,fCentStop,fNumberOfBins[iAnalysisType],fP2Start[iAnalysisType],fP2Stop[iAnalysisType]);
215 histName = "fHistPP"; histName += gBFAnalysisType[iAnalysisType];
216 if(bShuffle) histName.Append("_shuffle");
217 if(fCentralityId) histName += fCentralityId.Data();
218 fHistPP[iAnalysisType] = new TH2D(histName.Data(),"",fCentStop-fCentStart,fCentStart,fCentStop,fNumberOfBins[iAnalysisType],fP2Start[iAnalysisType],fP2Stop[iAnalysisType]);
220 histName = "fHistNN"; histName += gBFAnalysisType[iAnalysisType];
221 if(bShuffle) histName.Append("_shuffle");
222 if(fCentralityId) histName += fCentralityId.Data();
223 fHistNN[iAnalysisType] = new TH2D(histName.Data(),"",fCentStop-fCentStart,fCentStart,fCentStop,fNumberOfBins[iAnalysisType],fP2Start[iAnalysisType],fP2Stop[iAnalysisType]);
227 //____________________________________________________________________//
228 void AliBalance::PrintAnalysisSettings() {
230 Printf("======================================");
231 Printf("Analysis level: %s",fAnalysisLevel.Data());
232 Printf("======================================");
233 for(Int_t ibin = 0; ibin < ANALYSIS_TYPES; ibin++){
234 Printf("Interval info for variable %d",ibin);
235 Printf("Analyzed interval (min.): %lf",fP2Start[ibin]);
236 Printf("Analyzed interval (max.): %lf",fP2Stop[ibin]);
237 Printf("Number of bins: %d",fNumberOfBins[ibin]);
238 Printf("Step: %lf",fP2Step[ibin]);
241 Printf("======================================");
244 //____________________________________________________________________//
245 void AliBalance::CalculateBalance(Float_t fCentrality,vector<Double_t> **chargeVector,Float_t bSign) {
246 // Calculates the balance function
251 // Initialize histograms if not done yet
253 AliWarning("Histograms not yet initialized! --> Will be done now");
254 AliWarning("This works only in local mode --> Add 'gBalance->InitHistograms()' in your configBalanceFunction");
258 Int_t gNtrack = chargeVector[0]->size();
259 //Printf("(AliBalance) Number of tracks: %d",gNtrack);
261 for(i = 0; i < gNtrack;i++){
262 Short_t charge = chargeVector[0]->at(i);
263 Double_t rapidity = chargeVector[1]->at(i);
264 Double_t pseudorapidity = chargeVector[2]->at(i);
265 Double_t phi = chargeVector[3]->at(i);
267 //0:y - 1:eta - 2:Qlong - 3:Qout - 4:Qside - 5:Qinv - 6:phi
268 for(Int_t iAnalysisType = 0; iAnalysisType < ANALYSIS_TYPES; iAnalysisType++) {
269 if(iAnalysisType == kEta) {
270 if((pseudorapidity >= fP1Start[iAnalysisType]) && (pseudorapidity <= fP1Stop[iAnalysisType])) {
272 fNp[iAnalysisType] += 1.;
273 fHistP[iAnalysisType]->Fill(fCentrality,pseudorapidity);
276 fNn[iAnalysisType] += 1.;
277 fHistN[iAnalysisType]->Fill(fCentrality,pseudorapidity);
280 }//analysis type: eta
281 else if(iAnalysisType == kPhi) {
282 if((phi >= fP1Start[iAnalysisType]) && (phi <= fP1Stop[iAnalysisType])) {
284 fNp[iAnalysisType] += 1.;
285 fHistP[iAnalysisType]->Fill(fCentrality,phi);
288 fNn[iAnalysisType] += 1.;
289 fHistN[iAnalysisType]->Fill(fCentrality,phi);
292 }//analysis type: phi
294 if((rapidity >= fP1Start[iAnalysisType]) && (rapidity <= fP1Stop[iAnalysisType])) {
296 fNp[iAnalysisType] += 1.;
297 fHistP[iAnalysisType]->Fill(fCentrality,rapidity);
300 fNn[iAnalysisType] += 1.;
301 fHistN[iAnalysisType]->Fill(fCentrality,rapidity);
304 }//analysis type: y, qside, qout, qlong, qinv
305 }//analysis type loop
308 //Printf("Np: %lf - Nn: %lf",fNp[0],fNn[0]);
310 Double_t dy = 0., deta = 0.;
311 Double_t qLong = 0., qOut = 0., qSide = 0., qInv = 0.;
315 Double_t eta1 = 0., rap1 = 0.;
316 Double_t px1 = 0., py1 = 0., pz1 = 0.;
318 Double_t energy1 = 0.;
322 Double_t eta2 = 0., rap2 = 0.;
323 Double_t px2 = 0., py2 = 0., pz2 = 0.;
325 Double_t energy2 = 0.;
327 //0:y - 1:eta - 2:Qlong - 3:Qout - 4:Qside - 5:Qinv - 6:phi
328 for(i = 1; i < gNtrack; i++) {
330 charge1 = chargeVector[0]->at(i);
331 rap1 = chargeVector[1]->at(i);
332 eta1 = chargeVector[2]->at(i);
333 phi1 = chargeVector[3]->at(i);
334 px1 = chargeVector[4]->at(i);
335 py1 = chargeVector[5]->at(i);
336 pz1 = chargeVector[6]->at(i);
337 pt1 = chargeVector[7]->at(i);
338 energy1 = chargeVector[8]->at(i);
340 for(j = 0; j < i; j++) {
342 charge2 = chargeVector[0]->at(j);
343 rap2 = chargeVector[1]->at(j);
344 eta2 = chargeVector[2]->at(j);
345 phi2 = chargeVector[3]->at(j);
346 px2 = chargeVector[4]->at(j);
347 py2 = chargeVector[5]->at(j);
348 pz2 = chargeVector[6]->at(j);
349 pt2 = chargeVector[7]->at(i);
350 energy2 = chargeVector[8]->at(j);
352 // filling the arrays
355 dy = TMath::Abs(rap1 - rap2);
358 deta = TMath::Abs(eta1 - eta2);
361 Double_t eTot = energy1 + energy2;
362 Double_t pxTot = px1 + px2;
363 Double_t pyTot = py1 + py2;
364 Double_t pzTot = pz1 + pz2;
365 Double_t q0Tot = energy1 - energy2;
366 Double_t qxTot = px1 - px2;
367 Double_t qyTot = py1 - py2;
368 Double_t qzTot = pz1 - pz2;
370 Double_t eTot2 = eTot*eTot;
371 Double_t pTot2 = pxTot*pxTot + pyTot*pyTot + pzTot*pzTot;
372 Double_t pzTot2 = pzTot*pzTot;
374 Double_t q0Tot2 = q0Tot*q0Tot;
375 Double_t qTot2 = qxTot*qxTot + qyTot*qyTot + qzTot*qzTot;
377 Double_t snn = eTot2 - pTot2;
378 Double_t ptTot2 = pTot2 - pzTot2 ;
379 Double_t ptTot = TMath::Sqrt( ptTot2 );
381 qLong = TMath::Abs(eTot*qzTot - pzTot*q0Tot)/TMath::Sqrt(snn + ptTot2);
384 qOut = TMath::Sqrt(snn/(snn + ptTot2)) * TMath::Abs(pxTot*qxTot + pyTot*qyTot)/ptTot;
387 qSide = TMath::Abs(pxTot*qyTot - pyTot*qxTot)/ptTot;
390 qInv = TMath::Sqrt(TMath::Abs(-q0Tot2 + qTot2 ));
393 dphi = TMath::Abs(phi1 - phi2);
394 if(dphi>180) dphi = 360 - dphi; //dphi should be between 0 and 180!
398 //if( dphi < 3 || deta < 0.01 ){ // VERSION 1
401 // VERSION 2 (Taken from DPhiCorrelations)
402 // the variables & cuthave been developed by the HBT group
403 // see e.g. https://indico.cern.ch/materialDisplay.py?contribId=36&sessionId=6&materialId=slides&confId=142700
406 if (TMath::Abs(deta) < 0.02 * 2.5 * 3) //twoTrackEfficiencyCutValue = 0.02 [default for dphicorrelations]
408 // check first boundaries to see if is worth to loop and find the minimum
409 Float_t dphistar1 = GetDPhiStar(phi1*TMath::DegToRad(), pt1, charge1, phi2*TMath::DegToRad(), pt2, charge2, 0.8, bSign);
410 Float_t dphistar2 = GetDPhiStar(phi1*TMath::DegToRad(), pt1, charge1, phi2*TMath::DegToRad(), pt2, charge2, 2.5, bSign);
412 const Float_t kLimit = 0.02 * 3;
414 Float_t dphistarminabs = 1e5;
415 Float_t dphistarmin = 1e5;
416 if (TMath::Abs(dphistar1) < kLimit || TMath::Abs(dphistar2) < kLimit || dphistar1 * dphistar2 < 0)
418 for (Double_t rad=0.8; rad<2.51; rad+=0.01)
420 Float_t dphistar = GetDPhiStar(phi1*TMath::DegToRad(), pt1, charge1, phi2*TMath::DegToRad(), pt2, charge2, rad, bSign);
422 Float_t dphistarabs = TMath::Abs(dphistar);
424 if (dphistarabs < dphistarminabs)
426 dphistarmin = dphistar;
427 dphistarminabs = dphistarabs;
431 if (dphistarminabs < 0.02 && TMath::Abs(deta) < 0.02)
433 //Printf("Removed track pair %d %d with %f %f %f %f %d %f %f %d %f", i, j, deta, dphistarminabs, phi1, pt1, charge1, phi2, pt2, charge2, bSign);
442 if (charge1 * charge2 < 0)
444 Float_t m0 = 0.510e-3;
445 Float_t tantheta1 = 1e10;
447 if (eta1 < -1e-10 || eta1 > 1e-10)
448 tantheta1 = 2 * TMath::Exp(-eta1) / ( 1 - TMath::Exp(-2*eta1));
450 Float_t tantheta2 = 1e10;
451 if (eta2 < -1e-10 || eta2 > 1e-10)
452 tantheta2 = 2 * TMath::Exp(-eta2) / ( 1 - TMath::Exp(-2*eta2));
454 Float_t e1squ = m0 * m0 + pt1 * pt1 * (1.0 + 1.0 / tantheta1 / tantheta1);
455 Float_t e2squ = m0 * m0 + pt2 * pt2 * (1.0 + 1.0 / tantheta2 / tantheta2);
457 Float_t mass = 2 * m0 * m0 + 2 * ( TMath::Sqrt(e1squ * e2squ) - ( pt1 * pt2 * ( TMath::Cos(phi1 - phi2) + 1.0 / tantheta1 / tantheta2 ) ) );
459 if (mass < 0.04*0.04){
460 //Printf("Removed track pair %d %d with %f %f %d %d ", i, j, deta, mass, charge1, charge2);
467 //0:y - 1:eta - 2:Qlong - 3:Qout - 4:Qside - 5:Qinv - 6:phi
468 if((rap1 >= fP1Start[kRapidity]) && (rap1 <= fP1Stop[kRapidity]) && (rap2 >= fP1Start[kRapidity]) && (rap2 <= fP1Stop[kRapidity])) {
471 if( dy > fP2Start[kRapidity] && dy < fP2Stop[kRapidity]){
472 iBin = Int_t((dy-fP2Start[kRapidity])/fP2Step[kRapidity]);
473 if(iBin >=0 && iBin < MAXIMUM_NUMBER_OF_STEPS){
475 if((charge1 > 0)&&(charge2 > 0)) {
476 fNpp[kRapidity][iBin] += 1.;
477 fHistPP[kRapidity]->Fill(fCentrality,dy);
479 else if((charge1 < 0)&&(charge2 < 0)) {
480 fNnn[kRapidity][iBin] += 1.;
481 fHistNN[kRapidity]->Fill(fCentrality,dy);
483 else if((charge1 > 0)&&(charge2 < 0)) {
484 fNpn[kRapidity][iBin] += 1.;
485 fHistPN[kRapidity]->Fill(fCentrality,dy);
487 else if((charge1 < 0)&&(charge2 > 0)) {
488 fNpn[kRapidity][iBin] += 1.;
489 fHistPN[kRapidity]->Fill(fCentrality,dy);
496 if((eta1 >= fP1Start[kEta]) && (eta1 <= fP1Stop[kEta]) && (eta2 >= fP1Start[kEta]) && (eta2 <= fP1Stop[kEta])) {
497 if( deta > fP2Start[kEta] && deta < fP2Stop[kEta]){
498 iBin = Int_t((deta-fP2Start[kEta])/fP2Step[kEta]);
499 if(iBin >=0 && iBin < MAXIMUM_NUMBER_OF_STEPS){
500 if((charge1 > 0)&&(charge2 > 0)) {
501 fNpp[kEta][iBin] += 1.;
502 fHistPP[kEta]->Fill(fCentrality,deta);
504 if((charge1 < 0)&&(charge2 < 0)) {
505 fNnn[kEta][iBin] += 1.;
506 fHistNN[kEta]->Fill(fCentrality,deta);
508 if((charge1 > 0)&&(charge2 < 0)) {
509 fNpn[kEta][iBin] += 1.;
510 fHistPN[kEta]->Fill(fCentrality,deta);
512 if((charge1 < 0)&&(charge2 > 0)) {
513 fNpn[kEta][iBin] += 1.;
514 fHistPN[kEta]->Fill(fCentrality,deta);
520 // Qlong, out, side, inv
521 // Check the p1 intervall for rapidity here (like for single tracks above)
522 if((rap1 >= fP1Start[kRapidity]) && (rap1 <= fP1Stop[kRapidity]) && (rap2 >= fP1Start[kRapidity]) && (rap2 <= fP1Stop[kRapidity])) {
523 if( qLong > fP2Start[kQlong] && qLong < fP2Stop[kQlong]){
524 iBin = Int_t((qLong-fP2Start[kQlong])/fP2Step[kQlong]);
525 if(iBin >=0 && iBin < MAXIMUM_NUMBER_OF_STEPS){
526 if((charge1 > 0)&&(charge2 > 0)) {
527 fNpp[kQlong][iBin] += 1.;
528 fHistPP[kQlong]->Fill(fCentrality,qLong);
530 if((charge1 < 0)&&(charge2 < 0)) {
531 fNnn[kQlong][iBin] += 1.;
532 fHistNN[kQlong]->Fill(fCentrality,qLong);
534 if((charge1 > 0)&&(charge2 < 0)) {
535 fNpn[kQlong][iBin] += 1.;
536 fHistPN[kQlong]->Fill(fCentrality,qLong);
538 if((charge1 < 0)&&(charge2 > 0)) {
539 fNpn[kQlong][iBin] += 1.;
540 fHistPN[kQlong]->Fill(fCentrality,qLong);
546 if((rap1 >= fP1Start[kRapidity]) && (rap1 <= fP1Stop[kRapidity]) && (rap2 >= fP1Start[kRapidity]) && (rap2 <= fP1Stop[kRapidity])) {
547 if( qOut > fP2Start[kQout] && qOut < fP2Stop[kQout]){
548 iBin = Int_t((qOut-fP2Start[kQout])/fP2Step[kQout]);
549 if(iBin >=0 && iBin < MAXIMUM_NUMBER_OF_STEPS){
550 if((charge1 > 0)&&(charge2 > 0)) {
551 fNpp[kQout][iBin] += 1.;
552 fHistPP[kQout]->Fill(fCentrality,qOut);
554 if((charge1 < 0)&&(charge2 < 0)) {
555 fNnn[kQout][iBin] += 1.;
556 fHistNN[kQout]->Fill(fCentrality,qOut);
558 if((charge1 > 0)&&(charge2 < 0)) {
559 fNpn[kQout][iBin] += 1.;
560 fHistPN[kQout]->Fill(fCentrality,qOut);
562 if((charge1 < 0)&&(charge2 > 0)) {
563 fNpn[kQout][iBin] += 1.;
564 fHistPN[kQout]->Fill(fCentrality,qOut);
570 if((rap1 >= fP1Start[kRapidity]) && (rap1 <= fP1Stop[kRapidity]) && (rap2 >= fP1Start[kRapidity]) && (rap2 <= fP1Stop[kRapidity])) {
571 if( qSide > fP2Start[kQside] && qSide < fP2Stop[kQside]){
572 iBin = Int_t((qSide-fP2Start[kQside])/fP2Step[kQside]);
573 if(iBin >=0 && iBin < MAXIMUM_NUMBER_OF_STEPS){
574 if((charge1 > 0)&&(charge2 > 0)) {
575 fNpp[kQside][iBin] += 1.;
576 fHistPP[kQside]->Fill(fCentrality,qSide);
578 if((charge1 < 0)&&(charge2 < 0)) {
579 fNnn[kQside][iBin] += 1.;
580 fHistNN[kQside]->Fill(fCentrality,qSide);
582 if((charge1 > 0)&&(charge2 < 0)) {
583 fNpn[kQside][iBin] += 1.;
584 fHistPN[kQside]->Fill(fCentrality,qSide);
586 if((charge1 < 0)&&(charge2 > 0)) {
587 fNpn[kQside][iBin] += 1.;
588 fHistPN[kQside]->Fill(fCentrality,qSide);
594 if((rap1 >= fP1Start[kRapidity]) && (rap1 <= fP1Stop[kRapidity]) && (rap2 >= fP1Start[kRapidity]) && (rap2 <= fP1Stop[kRapidity])) {
595 if( qInv > fP2Start[kQinv] && qInv < fP2Stop[kQinv]){
596 iBin = Int_t((qInv-fP2Start[kQinv])/fP2Step[kQinv]);
597 if(iBin >=0 && iBin < MAXIMUM_NUMBER_OF_STEPS){
598 if((charge1 > 0)&&(charge2 > 0)) {
599 fNpp[kQinv][iBin] += 1.;
600 fHistPP[kQinv]->Fill(fCentrality,qInv);
602 if((charge1 < 0)&&(charge2 < 0)) {
603 fNnn[kQinv][iBin] += 1.;
604 fHistNN[kQinv]->Fill(fCentrality,qInv);
606 if((charge1 > 0)&&(charge2 < 0)) {
607 fNpn[kQinv][iBin] += 1.;
608 fHistPN[kQinv]->Fill(fCentrality,qInv);
610 if((charge1 < 0)&&(charge2 > 0)) {
611 fNpn[kQinv][iBin] += 1.;
612 fHistPN[kQinv]->Fill(fCentrality,qInv);
619 if((phi1 >= fP1Start[kPhi]) && (phi1 <= fP1Stop[kPhi]) && (phi2 >= fP1Start[kPhi]) && (phi2 <= fP1Stop[kPhi])) {
620 if( dphi > fP2Start[kPhi] && dphi < fP2Stop[kPhi]){
621 iBin = Int_t((dphi-fP2Start[kPhi])/fP2Step[kPhi]);
622 if(iBin >=0 && iBin < MAXIMUM_NUMBER_OF_STEPS){
623 if((charge1 > 0)&&(charge2 > 0)) {
624 fNpp[kPhi][iBin] += 1.;
625 fHistPP[kPhi]->Fill(fCentrality,dphi);
627 if((charge1 < 0)&&(charge2 < 0)) {
628 fNnn[kPhi][iBin] += 1.;
629 fHistNN[kPhi]->Fill(fCentrality,dphi);
631 if((charge1 > 0)&&(charge2 < 0)) {
632 fNpn[kPhi][iBin] += 1.;
633 fHistPN[kPhi]->Fill(fCentrality,dphi);
635 if((charge1 < 0)&&(charge2 > 0)) {
636 fNpn[kPhi][iBin] += 1.;
637 fHistPN[kPhi]->Fill(fCentrality,dphi);
642 }//end of 2nd particle loop
643 }//end of 1st particle loop
644 //Printf("Number of analyzed events: %i",fAnalyzedEvents);
645 //Printf("DeltaEta NN[0] = %.0f, PP[0] = %.0f, NP[0] = %.0f, PN[0] = %.0f",fNnn[kEta][0],fNpp[kEta][0],fNnp[kEta][0],fNpn[kEta][0]);
649 //____________________________________________________________________//
650 Double_t AliBalance::GetBalance(Int_t iAnalysisType, Int_t p2) {
651 // Returns the value of the balance function in bin p2
652 fB[iAnalysisType][p2] = 0.5*(((fNpn[iAnalysisType][p2] - 2.*fNnn[iAnalysisType][p2])/fNn[iAnalysisType]) + ((fNpn[iAnalysisType][p2] - 2.*fNpp[iAnalysisType][p2])/fNp[iAnalysisType]))/fP2Step[iAnalysisType];
654 return fB[iAnalysisType][p2];
657 //____________________________________________________________________//
658 Double_t AliBalance::GetError(Int_t iAnalysisType, Int_t p2) {
659 // Returns the error on the BF value for bin p2
660 // The errors for fNn and fNp are neglected here (0.1 % of total error)
661 /*ferror[iAnalysisType][p2] = TMath::Sqrt(Double_t(fNpp[iAnalysisType][p2])/(Double_t(fNp[iAnalysisType])*Double_t(fNp[iAnalysisType]))
662 + Double_t(fNnn[iAnalysisType][p2])/(Double_t(fNn[iAnalysisType])*Double_t(fNn[iAnalysisType]))
663 + Double_t(fNpn[iAnalysisType][p2])/(Double_t(fNp[iAnalysisType])*Double_t(fNp[iAnalysisType]))
664 + Double_t(fNnp[iAnalysisType][p2])/(Double_t(fNp[iAnalysisType])*Double_t(fNp[iAnalysisType]))
665 //+ TMath::Power(fNpn[iAnalysisType][p2]-fNpp[iAnalysisType][p2],2)/TMath::Power(Double_t(fNp[iAnalysisType]),3)
666 //+ TMath::Power(fNnp[iAnalysisType][p2]-fNnn[iAnalysisType][p2],2)/TMath::Power(Double_t(fNn[iAnalysisType]),3)
667 ) /fP2Step[iAnalysisType];*/
669 ferror[iAnalysisType][p2] = TMath::Sqrt( Double_t(fNpp[iAnalysisType][p2])/(Double_t(fNp[iAnalysisType])*Double_t(fNp[iAnalysisType])) +
670 Double_t(fNnn[iAnalysisType][p2])/(Double_t(fNn[iAnalysisType])*Double_t(fNn[iAnalysisType])) +
671 Double_t(fNpn[iAnalysisType][p2])*TMath::Power((0.5/Double_t(fNp[iAnalysisType]) + 0.5/Double_t(fNn[iAnalysisType])),2))/fP2Step[iAnalysisType];
673 return ferror[iAnalysisType][p2];
675 //____________________________________________________________________//
676 TGraphErrors *AliBalance::DrawBalance(Int_t iAnalysisType) {
679 Double_t x[MAXIMUM_NUMBER_OF_STEPS];
680 Double_t xer[MAXIMUM_NUMBER_OF_STEPS];
681 Double_t b[MAXIMUM_NUMBER_OF_STEPS];
682 Double_t ber[MAXIMUM_NUMBER_OF_STEPS];
684 if((fNp[iAnalysisType] == 0)||(fNn[iAnalysisType] == 0)) {
685 cerr<<"Couldn't find any particles in the analyzed interval!!!"<<endl;
689 for(Int_t i = 0; i < fNumberOfBins[iAnalysisType]; i++) {
690 b[i] = GetBalance(iAnalysisType,i);
691 ber[i] = GetError(iAnalysisType,i);
692 x[i] = fP2Start[iAnalysisType] + fP2Step[iAnalysisType]*i + fP2Step[iAnalysisType]/2;
696 TGraphErrors *gr = new TGraphErrors(fNumberOfBins[iAnalysisType],x,b,xer,ber);
697 gr->GetXaxis()->SetTitleColor(1);
698 if(iAnalysisType==0) {
699 gr->SetTitle("Balance function B(#Delta y)");
700 gr->GetXaxis()->SetTitle("#Delta y");
701 gr->GetYaxis()->SetTitle("B(#Delta y)");
703 if(iAnalysisType==1) {
704 gr->SetTitle("Balance function B(#Delta #eta)");
705 gr->GetXaxis()->SetTitle("#Delta #eta");
706 gr->GetYaxis()->SetTitle("B(#Delta #eta)");
708 if(iAnalysisType==2) {
709 gr->SetTitle("Balance function B(q_{long})");
710 gr->GetXaxis()->SetTitle("q_{long} (GeV/c)");
711 gr->GetYaxis()->SetTitle("B(q_{long}) ((GeV/c)^{-1})");
713 if(iAnalysisType==3) {
714 gr->SetTitle("Balance function B(q_{out})");
715 gr->GetXaxis()->SetTitle("q_{out} (GeV/c)");
716 gr->GetYaxis()->SetTitle("B(q_{out}) ((GeV/c)^{-1})");
718 if(iAnalysisType==4) {
719 gr->SetTitle("Balance function B(q_{side})");
720 gr->GetXaxis()->SetTitle("q_{side} (GeV/c)");
721 gr->GetYaxis()->SetTitle("B(q_{side}) ((GeV/c)^{-1})");
723 if(iAnalysisType==5) {
724 gr->SetTitle("Balance function B(q_{inv})");
725 gr->GetXaxis()->SetTitle("q_{inv} (GeV/c)");
726 gr->GetYaxis()->SetTitle("B(q_{inv}) ((GeV/c)^{-1})");
728 if(iAnalysisType==6) {
729 gr->SetTitle("Balance function B(#Delta #phi)");
730 gr->GetXaxis()->SetTitle("#Delta #phi");
731 gr->GetYaxis()->SetTitle("B(#Delta #phi)");
737 //____________________________________________________________________//
738 void AliBalance::PrintResults(Int_t iAnalysisType, TH1D *gHistBalance) {
739 //Prints the calculated width of the BF and its error
740 Double_t x[MAXIMUM_NUMBER_OF_STEPS];
741 Double_t gSumXi = 0.0, gSumBi = 0.0, gSumBiXi = 0.0;
742 Double_t gSumBiXi2 = 0.0, gSumBi2Xi2 = 0.0;
743 Double_t gSumDeltaBi2 = 0.0, gSumXi2DeltaBi2 = 0.0;
744 Double_t deltaBalP2 = 0.0, integral = 0.0;
745 Double_t deltaErrorNew = 0.0;
747 cout<<"=================================================="<<endl;
748 for(Int_t i = 1; i <= fNumberOfBins[iAnalysisType]; i++) {
749 x[i-1] = fP2Start[iAnalysisType] + fP2Step[iAnalysisType]*i + fP2Step[iAnalysisType]/2;
750 //cout<<"B: "<<gHistBalance->GetBinContent(i)<<"\t Error: "<<gHistBalance->GetBinError(i)<<"\t bin: "<<gHistBalance->GetBinCenter(i)<<endl;
752 //cout<<"=================================================="<<endl;
753 for(Int_t i = 2; i <= fNumberOfBins[iAnalysisType]; i++) {
754 gSumXi += gHistBalance->GetBinCenter(i);
755 gSumBi += gHistBalance->GetBinContent(i);
756 gSumBiXi += gHistBalance->GetBinContent(i)*gHistBalance->GetBinCenter(i);
757 gSumBiXi2 += gHistBalance->GetBinContent(i)*TMath::Power(gHistBalance->GetBinCenter(i),2);
758 gSumBi2Xi2 += TMath::Power(gHistBalance->GetBinContent(i),2)*TMath::Power(gHistBalance->GetBinCenter(i),2);
759 gSumDeltaBi2 += TMath::Power(gHistBalance->GetBinError(i),2);
760 gSumXi2DeltaBi2 += TMath::Power(gHistBalance->GetBinCenter(i),2) * TMath::Power(gHistBalance->GetBinError(i),2);
762 deltaBalP2 += fP2Step[iAnalysisType]*TMath::Power(gHistBalance->GetBinError(i),2);
763 integral += fP2Step[iAnalysisType]*gHistBalance->GetBinContent(i);
765 for(Int_t i = 1; i < fNumberOfBins[iAnalysisType]; i++)
766 deltaErrorNew += gHistBalance->GetBinError(i)*(gHistBalance->GetBinCenter(i)*gSumBi - gSumBiXi)/TMath::Power(gSumBi,2);
768 Double_t integralError = TMath::Sqrt(deltaBalP2);
770 Double_t delta = gSumBiXi / gSumBi;
771 Double_t deltaError = (gSumBiXi / gSumBi) * TMath::Sqrt(TMath::Power((TMath::Sqrt(gSumXi2DeltaBi2)/gSumBiXi),2) + TMath::Power((gSumDeltaBi2/gSumBi),2) );
772 cout<<"Analysis type: "<<gBFAnalysisType[iAnalysisType].Data()<<endl;
773 cout<<"Width: "<<delta<<"\t Error: "<<deltaError<<endl;
774 cout<<"New error: "<<deltaErrorNew<<endl;
775 cout<<"Integral: "<<integral<<"\t Error: "<<integralError<<endl;
776 cout<<"=================================================="<<endl;
779 //____________________________________________________________________//
780 TH1D *AliBalance::GetBalanceFunctionHistogram(Int_t iAnalysisType,Double_t centrMin, Double_t centrMax, Double_t etaWindow) {
781 //Returns the BF histogram, extracted from the 6 TH2D objects
782 //(private members) of the AliBalance class.
784 // Acceptance correction:
785 // - only for analysis type = kEta
786 // - only if etaWindow > 0 (default = -1.)
787 // - calculated as proposed by STAR
789 TString gAnalysisType[ANALYSIS_TYPES] = {"y","eta","qlong","qout","qside","qinv","phi"};
790 TString histName = "gHistBalanceFunctionHistogram";
791 histName += gAnalysisType[iAnalysisType];
793 SetInterval(iAnalysisType, fHistP[iAnalysisType]->GetYaxis()->GetXmin(),
794 fHistP[iAnalysisType]->GetYaxis()->GetXmin(),
795 fHistPP[iAnalysisType]->GetNbinsY(),
796 fHistPP[iAnalysisType]->GetYaxis()->GetXmin(),
797 fHistPP[iAnalysisType]->GetYaxis()->GetXmax());
799 // determine the projection thresholds
800 Int_t binMinX, binMinY, binMinZ;
801 Int_t binMaxX, binMaxY, binMaxZ;
803 fHistPP[iAnalysisType]->GetBinXYZ(fHistPP[iAnalysisType]->FindBin(centrMin),binMinX,binMinY,binMinZ);
804 fHistPP[iAnalysisType]->GetBinXYZ(fHistPP[iAnalysisType]->FindBin(centrMax),binMaxX,binMaxY,binMaxZ);
806 TH1D *gHistBalanceFunctionHistogram = new TH1D(histName.Data(),"",fHistPP[iAnalysisType]->GetNbinsY(),fHistPP[iAnalysisType]->GetYaxis()->GetXmin(),fHistPP[iAnalysisType]->GetYaxis()->GetXmax());
807 switch(iAnalysisType) {
809 gHistBalanceFunctionHistogram->GetXaxis()->SetTitle("#Delta y");
810 gHistBalanceFunctionHistogram->GetYaxis()->SetTitle("B(#Delta y)");
813 gHistBalanceFunctionHistogram->GetXaxis()->SetTitle("#Delta #eta");
814 gHistBalanceFunctionHistogram->GetYaxis()->SetTitle("B(#Delta #eta)");
817 gHistBalanceFunctionHistogram->GetXaxis()->SetTitle("q_{long} (GeV/c)");
818 gHistBalanceFunctionHistogram->GetYaxis()->SetTitle("B(q_{long})");
821 gHistBalanceFunctionHistogram->GetXaxis()->SetTitle("q_{out} (GeV/c)");
822 gHistBalanceFunctionHistogram->GetYaxis()->SetTitle("B(q_{out})");
825 gHistBalanceFunctionHistogram->GetXaxis()->SetTitle("q_{side} (GeV/c)");
826 gHistBalanceFunctionHistogram->GetYaxis()->SetTitle("B(q_{side})");
829 gHistBalanceFunctionHistogram->GetXaxis()->SetTitle("q_{inv} (GeV/c)");
830 gHistBalanceFunctionHistogram->GetYaxis()->SetTitle("B(q_{inv})");
833 gHistBalanceFunctionHistogram->GetXaxis()->SetTitle("#Delta #phi (deg.)");
834 gHistBalanceFunctionHistogram->GetYaxis()->SetTitle("B(#Delta #phi)");
840 TH1D *hTemp1 = dynamic_cast<TH1D *>(fHistPN[iAnalysisType]->ProjectionY(Form("%s_Cent_%.0f_%.0f",fHistPN[iAnalysisType]->GetName(),centrMin,centrMax),binMinX,binMaxX));
841 TH1D *hTemp2 = dynamic_cast<TH1D *>(fHistPN[iAnalysisType]->ProjectionY(Form("%s_Cent_%.0f_%.0f_copy",fHistPN[iAnalysisType]->GetName(),centrMin,centrMax),binMinX,binMaxX));
842 TH1D *hTemp3 = dynamic_cast<TH1D *>(fHistNN[iAnalysisType]->ProjectionY(Form("%s_Cent_%.0f_%.0f",fHistNN[iAnalysisType]->GetName(),centrMin,centrMax),binMinX,binMaxX));
843 TH1D *hTemp4 = dynamic_cast<TH1D *>(fHistPP[iAnalysisType]->ProjectionY(Form("%s_Cent_%.0f_%.0f",fHistPP[iAnalysisType]->GetName(),centrMin,centrMax),binMinX,binMaxX));
844 TH1D *hTemp5 = dynamic_cast<TH1D *>(fHistN[iAnalysisType]->ProjectionY(Form("%s_Cent_%.0f_%.0f",fHistN[iAnalysisType]->GetName(),centrMin,centrMax),binMinX,binMaxX));
845 TH1D *hTemp6 = dynamic_cast<TH1D *>(fHistP[iAnalysisType]->ProjectionY(Form("%s_Cent_%.0f_%.0f",fHistP[iAnalysisType]->GetName(),centrMin,centrMax),binMinX,binMaxX));
847 if((hTemp1)&&(hTemp2)&&(hTemp3)&&(hTemp4)) {
852 hTemp1->Add(hTemp3,-2.);
853 hTemp1->Scale(1./hTemp5->GetEntries());
854 hTemp2->Add(hTemp4,-2.);
855 hTemp2->Scale(1./hTemp6->GetEntries());
856 gHistBalanceFunctionHistogram->Add(hTemp1,hTemp2,1.,1.);
857 gHistBalanceFunctionHistogram->Scale(0.5/fP2Step[iAnalysisType]);
860 // do the acceptance correction (only for Eta and etaWindow > 0)
861 if(iAnalysisType == kEta && etaWindow > 0){
862 for(Int_t iBin = 0; iBin < gHistBalanceFunctionHistogram->GetNbinsX(); iBin++){
864 Double_t notCorrected = gHistBalanceFunctionHistogram->GetBinContent(iBin+1);
865 Double_t corrected = notCorrected / (1 - (gHistBalanceFunctionHistogram->GetBinCenter(iBin+1))/ etaWindow );
866 gHistBalanceFunctionHistogram->SetBinContent(iBin+1, corrected);
871 PrintResults(iAnalysisType,gHistBalanceFunctionHistogram);
873 return gHistBalanceFunctionHistogram;