]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/UNICOR/AliUnicorAnalCorrel.cxx
Fixing warnings
[u/mrichter/AliRoot.git] / PWG2 / UNICOR / AliUnicorAnalCorrel.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> 2005
17
18 //=============================================================================
19 // two-particle correlation 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 <TMath.h>
27 #include <TRandom2.h>
28 #include "AliUnicorEvent.h"
29 #include "AliUnicorHN.h"
30 #include "AliUnicorAnalCorrel.h"
31
32 ClassImp(AliUnicorAnalCorrel)
33  
34 //=============================================================================
35 AliUnicorAnalCorrel::AliUnicorAnalCorrel(Char_t *nam, Double_t emi, Double_t ema, 
36                          Int_t pid0, Int_t pid1): 
37   AliUnicorAnal(nam), fPid0(pid0), fPid1(pid1), fMass0(0), fMass1(0), fPa()
38 {
39   // constructor
40   // emi and ema define the rapidity range for histogram
41
42   TParticlePDG *part0 = AliUnicorAnal::fgPDG.GetParticle(fPid0);
43   TParticlePDG *part1 = AliUnicorAnal::fgPDG.GetParticle(fPid1);
44   fMass0 = part0? part0->Mass() : 0;
45   fMass1 = part1? part1->Mass() : 0;
46   double pi = TMath::Pi();
47
48   // correlation function
49
50   //int npt = 7; double ptbins[]={0,0.1,0.2,0.3,0.4,0.5,0.7,1.0};
51   //int npt = 6; double ptbins[]={0,0.1,0.25,0.35,0.55,1.0,2.0}; // like Adam, except last bin split
52   int npt = 7; double ptbins[]={0,0.1,0.25,0.40,0.55,0.7,1.0,2.0}; // like Adam in Mar-2010, + first+last
53
54   double qbins[100];
55   for (int i=0;i<20;i++) qbins[i]=i*0.005;
56   for (int i=0;i<45;i++) qbins[20+i]=0.1+i*0.02;
57
58   TAxis *ax[8];
59   ax[0] = new TAxis(3,-0.5,2.5);ax[0]->SetTitle("trumixrot");
60   ax[1] = new TAxis(5,0,1.0);   ax[1]->SetTitle("centrality");
61   ax[2] = new TAxis(3,emi,ema); ax[2]->SetTitle("y");          // pair y
62   //  ax[3] = new TAxis(8,-pi,pi);  ax[3]->SetTitle("phi");      // wrt event plane
63   ax[3] = new TAxis(1,-pi,pi);  ax[3]->SetTitle("phi");        // wrt event plane
64   ax[4] = new TAxis(npt,ptbins);ax[4]->SetTitle("kt (GeV/c)"); // pair pt/2
65   ax[5] = new TAxis(8,0,pi);    ax[5]->SetTitle("q-theta");
66   ax[6] = new TAxis(16,-pi,pi); ax[6]->SetTitle("q-phi");
67   //ax[7] = new TAxis(64,qbins);  ax[7]->SetTitle("q (GeV/c)");
68   ax[7] = new TAxis(150,0,3.0); ax[7]->SetTitle("q (GeV/c)");
69   AliUnicorHN *pair = new AliUnicorHN("pair",8,ax);
70   for (int i=0; i<8; i++) delete ax[i];
71   fHistos.Add(pair);
72
73   // correlation function bin monitor (needed to get <kt> etc.)
74
75   ax[0] = new TAxis(5,0,1);      ax[0]->SetTitle("centrality");
76   ax[1] = new TAxis(30,emi,ema); ax[1]->SetTitle("y");           // pair y
77   ax[2] = new TAxis(100,0,2);    ax[2]->SetTitle("kt (GeV/c)");  // pair pt/2
78   AliUnicorHN *bimo = new AliUnicorHN("bimo",3,ax);
79   for (int i=0; i<3; i++) delete ax[i];
80   fHistos.Add(bimo);
81
82   // two-track resolution monitoring histogram
83
84   ax[0] = new TAxis(3,-0.5,2.5);    ax[0]->SetTitle("trumixrot");
85   ax[1] = new TAxis(2,-0.5,1.5);    ax[1]->SetTitle("cut applied");
86   ax[2] = new TAxis(npt,ptbins);    ax[2]->SetTitle("(pair pt)/2 (GeV)");
87   ax[3] = new TAxis(80,-0.02,0.02); ax[3]->SetTitle("dtheta");
88   ax[4] = new TAxis(80,-0.04,0.04); ax[4]->SetTitle("dphi");
89   AliUnicorHN *twot = new AliUnicorHN("twot",5,ax);
90   for (int i=0; i<5; i++) delete ax[i];
91   fHistos.Add(twot);
92
93   gROOT->cd();
94   printf("%s object named %s created\n",ClassName(),GetName());
95 }
96 //=============================================================================
97 void AliUnicorAnalCorrel::Process(Int_t tmr, const AliUnicorEvent * const ev0, const AliUnicorEvent * const ev1, Double_t phirot) 
98 {
99   // process pairs from one or two (if mixing) events
100   // tmr tells which histogram (bins) to fill: tru,mix,rot
101
102   // Could be possibly accelerated by checking the "good particle" only once 
103   // and caching the result. (Maybe the optimizer does it already.)
104
105   static TRandom2 ran;
106   AliUnicorHN *pair = (AliUnicorHN*) fHistos.At(0);
107   AliUnicorHN *bimo = (AliUnicorHN*) fHistos.At(1);
108   AliUnicorHN *twot = (AliUnicorHN*) fHistos.At(2);
109
110   // mixing-and-rotating-proof centrality and reaction plane angle
111   // (but not rotation-proof for rotation angles much different from 0 and 180)
112   // true and rotated pairs are within the triangle (j<i), mixed - all
113   // thus, proper rotation is either by 180, or by 170 AND 190, etc. 
114
115   double cent = (ev0->Centrality()+ev1->Centrality())/2.0;
116   double q0x,q0y,q1x,q1y;
117   ev0->RP(q0x,q0y);
118   ev1->RP(q1x,q1y); 
119   double rpphi = atan2(q0y+q1y,q0x+q1x);
120
121   // loop over pairs 
122
123   for (int i=0; i<ev0->NParticles(); i++) {
124     if (!ev0->ParticleGood(i,fPid0)) continue;
125     for (int j=0; j<ev1->NParticles(); j++) {
126       if (ev0==ev1 && j<i && fPid0==fPid1 ) continue; 
127       if (ev0==ev1 && j==i) continue; // beware, not even when rotated or non-identical
128       if (!ev1->ParticleGood(j,fPid1)) continue;
129       fPa.Set0(fMass0,ev0->ParticleP(i),ev0->ParticleTheta(i),ev0->ParticlePhi(i));
130       fPa.Set1(fMass1,ev1->ParticleP(j),ev1->ParticleTheta(j),ev1->ParticlePhi(j)+phirot);
131       if (ev0==ev1 && fPid0==fPid1 && ran.Rndm()>=0.5) fPa.Swap();
132       twot->Fill((double) tmr, 0.0, fPa.Pt()/2.0, fPa.DTheta(), fPa.DPhi(),1.0);
133       if (!ev0->PairGood(ev0->ParticleP(i),ev0->ParticleTheta(i),ev0->ParticlePhi(i),
134                          ev1->ParticleP(j),ev1->ParticleTheta(j),ev1->ParticlePhi(j)+phirot)) continue;
135       twot->Fill((double) tmr, 1.0, fPa.Pt()/2.0, fPa.DTheta(), fPa.DPhi(),1.0);
136       fPa.CalcLAB();
137       fPa.CalcPairCM();
138       //fPa.CalcLcmsCM();
139       if (fPa.QCM()==0) return; // should not be too frequent
140       double phi = TVector2::Phi_mpi_pi(fPa.Phi()-rpphi);
141       double weigth = 1.0;
142       //      if (tmr==0) weigth = 1.0+0.5*exp(-fPa.QCM()*fPa.QCM()*1*1/0.197/0.197); 
143       pair->Fill((double) tmr,           // 0 for tru, 1 for mix, 2 for rot
144                  cent,                   // centrality
145                  fPa.Rapidity(),         // pair rapidity
146                  phi,                    // pair phi wrt reaction plane
147                  fPa.Pt()/2.0,           // half of pair pt
148                  fPa.QCMTheta(),         // polar angle of Q
149                  fPa.QCMPhiOut(),        // azimuthal angle of Q w.r.t. out
150                  fPa.QCM(),              // |p2-p1| in c.m.s.
151                  weigth);                // weigth
152       if (tmr) continue;
153       if (fPa.QCM()>0.2) continue;
154       bimo->Fill(cent, fPa.Rapidity(), fPa.Pt()/2.0, 1.0);
155     }
156   }
157 }
158 //=============================================================================