]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/UNICOR/AliUnicorAnalPtfluc.cxx
Update to macros for wider multiplicity axis and fixed rapidity interval
[u/mrichter/AliRoot.git] / PWG2 / UNICOR / AliUnicorAnalPtfluc.cxx
1 /************************************************************************* 
2 * Copyright(c) 1998-2048, ALICE Experiment at CERN, All rights reserved. * 
3 *                                                                        * 
4 * Author: The ALICE Off-line Project.                                    * 
5 * Contributors are mentioned in the code where appropriate.              * 
6 *                                                                        * 
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 **************************************************************************/
15
16 // Author: Dariusz Miskowiec <mailto:d.miskowiec@gsi.de> 2008
17
18 //=============================================================================
19 // pt-fluctuations analyzer
20 // Loop over pairs and fill pair histograms. The first particle is always 
21 // taken from ev0 and the second from ev1 so for ev0=ev1 one is getting true 
22 // pairs, otherwise mixed ones. 
23 //=============================================================================
24
25 #include <TROOT.h>
26 #include <TRandom2.h>
27 #include <TVector2.h>
28 #include <TMath.h>
29 #include "AliUnicorEvent.h"
30 #include "AliUnicorHN.h"
31 #include "AliUnicorAnalPtfluc.h"
32
33 ClassImp(AliUnicorAnalPtfluc)
34
35 //=============================================================================
36 AliUnicorAnalPtfluc::AliUnicorAnalPtfluc(const char *nam, Int_t pid0, Int_t pid1) : 
37   AliUnicorAnal(nam), fPid0(pid0), fPid1(pid1) 
38 {
39   // constructor
40
41   TAxis *ax[5];
42   ax[0] = new TAxis(2,-0.5,1.5);   ax[0]->SetTitle("trumix");
43   ax[1] = new TAxis(9,0,0.9);      ax[1]->SetTitle("centrality");
44   ax[2] = new TAxis(6,-0.5,5.5);   ax[2]->SetTitle("n-pt0-pt1-pt00-pt11-pt01");
45   ax[3] = new TAxis(48,-180,180);  ax[3]->SetTitle("dphi (deg)");
46   ax[4] = new TAxis(40,-2,2);      ax[4]->SetTitle("deta");
47   AliUnicorHN *pair = new AliUnicorHN("pair",5,ax);
48   for (int i=0; i<5; i++) delete ax[i];
49   fHistos.Add(pair);
50   gROOT->cd();
51   printf("%s object named %s created\n",ClassName(),GetName());
52 }
53 //=============================================================================
54 void AliUnicorAnalPtfluc::Process(Int_t tmr, AliUnicorEvent * const ev0, AliUnicorEvent * const ev1) 
55 {
56   // process pairs from one or two (if mixing) events
57
58   double ptmin=0.1;  // GeV
59   double ptmax=1.5;  // GeV
60   double etamin=-9;  
61   double etamax=9;  
62
63   // mixing-and-rotating-proof centrality
64
65   double cent = (ev0->Centrality()+ev1->Centrality())/2.0;
66
67   // loop over pairs 
68
69   AliUnicorHN *pair = (AliUnicorHN*) fHistos.At(0);
70   static TRandom2 ran;
71   for (int i=0; i<ev0->NParticles(); i++) {
72     if (!ev0->ParticleGood(i,fPid0)) continue;
73     double eta0 = ev0->ParticleEta(i);
74     double phi0 = ev0->ParticlePhi(i);
75     double pt0 = ev0->ParticlePt(i);
76     if (eta0 < etamin) continue;
77     if (eta0 > etamax) continue;
78     if (pt0 < ptmin) continue;
79     if (pt0 > ptmax) continue;
80     for (int j=0; j<ev1->NParticles(); j++) {
81       if (ev0==ev1 && j==i) continue;
82       if (ev0==ev1 && j<i && fPid0==fPid1) continue;
83       if (!ev1->ParticleGood(j,fPid1)) continue;
84       double eta1 = ev1->ParticleEta(j);
85       double phi1 = ev1->ParticlePhi(j);
86       double pt1 = ev1->ParticlePt(j);
87       if (eta1 < etamin) continue;
88       if (eta1 > etamax) continue;
89       if (pt1 < ptmin) continue;
90       if (pt1 > ptmax) continue;
91       double deta = eta1-eta0;
92       double dphi = phi1-phi0;
93       // randomize order
94       if (ran.Rndm()<0.5) { 
95         double buf = pt0;
96         pt0 = pt1;
97         pt1 = buf;
98         deta = -deta;
99         dphi = -dphi;
100       }
101       dphi = TVector2::Phi_mpi_pi(dphi);
102       dphi*=TMath::RadToDeg();
103       pair->Fill((double) tmr, cent, 0.0, dphi, deta, 1.0);   // number of pairs
104       pair->Fill((double) tmr, cent, 1.0, dphi, deta, pt0);
105       pair->Fill((double) tmr, cent, 2.0, dphi, deta, pt1);
106       pair->Fill((double) tmr, cent, 3.0, dphi, deta, pt0*pt0);
107       pair->Fill((double) tmr, cent, 4.0, dphi, deta, pt1*pt1);
108       pair->Fill((double) tmr, cent, 5.0, dphi, deta, pt0*pt1);
109     }
110   }
111 }
112 //=============================================================================