]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/FEMTOSCOPY/UNICOR/AliUnicorAnalHighpt.cxx
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGCF / FEMTOSCOPY / UNICOR / AliUnicorAnalHighpt.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 // high-pt 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. Overlaps a bit with AliUnicorAnalCorrel. 
23 //=============================================================================
24
25 #include <TROOT.h>
26 #include <TMath.h>
27 #include <TVector2.h>
28 #include "AliUnicorEvent.h"
29 #include "AliUnicorHN.h"
30 #include "AliUnicorAnalHighpt.h"
31
32 ClassImp(AliUnicorAnalHighpt)
33  
34 //=============================================================================
35 AliUnicorAnalHighpt::AliUnicorAnalHighpt(const char *nam, Double_t emi, Double_t ema, Int_t pid0, 
36                          Int_t pid1): AliUnicorAnal(nam), fPid0(pid0), fPid1(pid1)
37 {
38   // constructor
39   // emi and ema define the rapidity range for histogram
40
41   Double_t ewi = ema-emi; // width of the pseudorapidity range
42   double pi = TMath::Pi();
43   TAxis *ax[9];
44   ax[0] = new TAxis(2,-0.5,1.5);   ax[0]->SetTitle("trumix");
45   ax[1] = new TAxis(2,-0.5,1.5);   ax[1]->SetTitle("weigth"); // 1 or ass pt
46   ax[2] = new TAxis(5,0,0.5);      ax[2]->SetTitle("centrality");
47   ax[3] = new TAxis(3,emi,ema);    ax[3]->SetTitle("trig eta");
48   ax[4] = new TAxis(1,-pi,pi);     ax[4]->SetTitle("trig phi"); // w.r.t. RP
49   //  ax[4] = new TAxis(8,-pi,pi);     ax[4]->SetTitle("trig phi"); // w.r.t. RP
50   double a5[]={0,2.5,3,4,6,8,10,15,20,30,50,100};
51   ax[5] = new TAxis(11,a5);        ax[5]->SetTitle("trig pt (GeV)");
52   ax[6] = new TAxis(20,-ewi,ewi);  ax[6]->SetTitle("ass eta - trig eta"); 
53   ax[7] = new TAxis(36,-1,2*pi-1); ax[7]->SetTitle("ass phi - trig phi"); 
54   ax[8] = new TAxis(10,0,1);       ax[8]->SetTitle("ass pt / trig pt");
55
56   AliUnicorHN *pair = new AliUnicorHN("pair",9,ax);
57   for (int i=0; i<9; i++) delete ax[i];
58   fHistos.Add(pair);
59   gROOT->cd();
60   printf("%s object named %s created\n",ClassName(),GetName());
61 }
62 //=============================================================================
63 void AliUnicorAnalHighpt::Process(const AliUnicorEvent * const ev0, const AliUnicorEvent * const ev1) 
64 {
65   // process pairs from one or two (if mixing) events
66   // true pairs are within the triangle (j<i), mixed - all
67
68   // mixing-proof centrality and reaction plane angle
69
70   double cent = (ev0->Centrality()+ev1->Centrality())/2.0;
71   double q0x=0,q0y=0,q1x=0,q1y=0;
72   ev0->RP(q0x,q0y);
73   ev1->RP(q1x,q1y); 
74   double rpphi = atan2(q0y+q1y,q0x+q1x);
75
76   // loop over pairs 
77
78   int tmr = (ev0!=ev1); // 0 for true, 1 for mixed
79   AliUnicorHN *pair = (AliUnicorHN*) fHistos.At(0);
80   for (int i=0; i<ev0->NParticles(); i++) {
81     if (!ev0->ParticleGood(i,fPid0)) continue;
82     double eta0 = ev0->ParticleEta(i);
83     double phi0 = ev0->ParticlePhi(i);
84     double pt0 = ev0->ParticlePt(i);
85     for (int j=0; j<ev1->NParticles(); j++) {
86       if (ev0==ev1 && j==i) continue; // same particle
87       if (ev0==ev1 && j<i && fPid0==fPid1) continue;
88       if (!ev1->ParticleGood(j,fPid1)) continue;
89       double eta1 = ev1->ParticleEta(j);
90       double phi1 = ev1->ParticlePhi(j);
91       double pt1 = ev1->ParticlePt(j);
92       int order = (pt0>pt1); 
93       double deta = order*(eta1-eta0); // ass-trig
94       double dphi = order*(phi1-phi0);
95       dphi = TVector2::Phi_mpi_pi(dphi);
96       if (dphi<-1) dphi+=2*TMath::Pi(); 
97       double relphi = TVector2::Phi_mpi_pi(phi0-rpphi);  // wrt. reaction plane
98       pair->Fill((double)tmr,0.0,cent,eta0,relphi,pt0,deta,dphi,pt1/pt0,1.0); // number of pairs
99       pair->Fill((double)tmr,1.0,cent,eta0,relphi,pt0,deta,dphi,pt1/pt0,pt1); // weigthed with ass pt
100     }
101   }
102 }
103 //=============================================================================