]> git.uio.no Git - u/mrichter/AliRoot.git/blob - UNICOR/AliDAnalCorrel.cxx
ALICE version of the Universal Correlation analysis package UNICOR.
[u/mrichter/AliRoot.git] / UNICOR / AliDAnalCorrel.cxx
1 // Author: Dariusz Miskowiec <mailto:d.miskowiec@gsi.de> 2005
2
3 //=============================================================================
4 // two-particle correlation analyzer
5 //=============================================================================
6
7 #include <TROOT.h>
8 #include <TRandom2.h>
9 #include "AliDEvent.h"
10 #include "AliDHN.h"
11 #include "AliDAnalCorrel.h"
12
13 ClassImp(AliDAnalCorrel)
14  
15 //=============================================================================
16 AliDAnalCorrel::AliDAnalCorrel(Char_t *nam, Double_t emi, Double_t ema, 
17                          Int_t pid0, Int_t pid1): 
18   AliDAnal(nam), fPid0(pid0), fPid1(pid1), fMass0(0), fMass1(0), fPa()
19 {
20   // constructor
21   // emi and ema define the rapidity range for histogram
22
23   TParticlePDG *part0 = AliDAnal::fgPDG.GetParticle(fPid0);
24   TParticlePDG *part1 = AliDAnal::fgPDG.GetParticle(fPid1);
25   fMass0 = part0? part0->Mass() : 0;
26   fMass1 = part1? part1->Mass() : 0;
27
28   double pi = TMath::Pi();
29   TAxis *ax[8];
30   ax[0] = new TAxis(3,-0.5,2.5);ax[0]->SetTitle("trumixrot");
31   ax[1] = new TAxis(5,0,0.5);   ax[1]->SetTitle("centrality");
32   ax[2] = new TAxis(3,emi,ema); ax[2]->SetTitle("pair y");
33   ax[3] = new TAxis(8,-pi,pi);  ax[3]->SetTitle("pair phi"); // wrt event plane
34   double a0[]={0,0.1,0.2,0.3,0.4,0.5,0.7,1.0};
35   ax[4] = new TAxis(7,a0);      ax[4]->SetTitle("(pair pt)/2 (GeV)");
36   ax[5] = new TAxis(8,0,pi);    ax[5]->SetTitle("q-theta");
37   ax[6] = new TAxis(16,-pi,pi); ax[6]->SetTitle("q-phi");
38   double a1[100];
39   for (int i=0;i<20;i++) a1[i]=i*0.005;
40   for (int i=0;i<45;i++) a1[20+i]=0.1+i*0.02;
41   ax[7] = new TAxis(64,a1);     ax[7]->SetTitle("q (GeV/c)");
42   AliDHN *pair = new AliDHN("pair",8,ax);
43   for (int i=0; i<8; i++) delete ax[i];
44   fHistos.Add(pair);
45   gROOT->cd();
46   printf("%s object named %s created\n",ClassName(),GetName());
47 }
48 //=============================================================================
49 void AliDAnalCorrel::Process(Int_t tmr, AliDEvent *ev0, AliDEvent *ev1, Double_t phirot) 
50 {
51   // process pairs from one or two (if mixing) events
52   // tmr tells which histogram (bins) to fill: tru,mix,rot
53
54   static TRandom2 ran;
55   AliDHN *pair = (AliDHN*) fHistos.At(0);
56
57   // mixing-and-rotating-proof centrality and reaction plane angle
58   // (but not rotation-proof for rotation angles much different from 0 and 180)
59   // true and rotated pairs are within the triangle (j<i), mixed - all
60   // thus, proper rotation is either by 180, or by 170 AND 190, etc. 
61
62   double cent = (ev0->Centrality()+ev1->Centrality())/2.0;
63   double q0x,q0y,q1x,q1y;
64   ev0->RP(q0x,q0y);
65   ev1->RP(q1x,q1y); 
66   double rpphi = atan2(q0y+q1y,q0x+q1x);
67
68   // loop over pairs 
69
70   for (int i=0; i<ev0->NParticles(); i++) {
71     if (!ev0->ParticleGood(i,fPid0)) continue;
72     for (int j=0; j<ev1->NParticles(); j++) {
73       if (ev0==ev1 && j<i && fPid0==fPid1 ) continue; 
74       if (ev0==ev1 && j==i) continue; // beware, not even when rotated or non-identical
75       if (!ev1->ParticleGood(j,fPid1)) continue;
76       if (!ev0->PairGood(ev0->ParticleP(i),ev0->ParticleTheta(i),ev0->ParticlePhi(i),
77                          ev1->ParticleP(j),ev1->ParticleTheta(j),ev1->ParticlePhi(j)+phirot)) continue;
78       fPa.Set0(fMass0,ev0->ParticleP(i),ev0->ParticleTheta(i),ev0->ParticlePhi(i));
79       fPa.Set1(fMass1,ev1->ParticleP(j),ev1->ParticleTheta(j),ev1->ParticlePhi(j)+phirot);
80       if (ev0==ev1 && fPid0==fPid1 && ran.Rndm()>=0.5) fPa.Swap();
81       fPa.CalcLAB();
82       fPa.CalcPairCM();
83       if (fPa.QCM()==0) return; // should not be too frequent
84       double phi = TVector2::Phi_mpi_pi(fPa.Phi()-rpphi);
85       pair->Fill(tmr,                    // 0 for tru, 1 for mix, 2 for rot
86                  cent,                   // centrality
87                  fPa.Rapidity(),         // pair rapidity
88                  phi,                    // pair phi wrt reaction plane
89                  fPa.Pt()/2.0,           // half of pair pt
90                  fPa.QCMTheta(),         // polar angle of Q
91                  fPa.QCMPhiOut(),        // azimuthal angle of Q w.r.t. out
92                  fPa.QCM(),              // |p2-p1| in c.m.s.
93                  1.0);                   // weigth
94     }
95   }
96 }
97 //=============================================================================