]>
Commit | Line | Data |
---|---|---|
621688e4 | 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 | ||
47 | double pi = TMath::Pi(); | |
48 | TAxis *ax[8]; | |
49 | ax[0] = new TAxis(3,-0.5,2.5);ax[0]->SetTitle("trumixrot"); | |
50 | ax[1] = new TAxis(5,0,0.5); ax[1]->SetTitle("centrality"); | |
51 | ax[2] = new TAxis(3,emi,ema); ax[2]->SetTitle("pair y"); | |
52 | //ax[3] = new TAxis(8,-pi,pi); ax[3]->SetTitle("pair phi"); // wrt event plane | |
53 | ax[3] = new TAxis(1,-pi,pi); ax[3]->SetTitle("pair phi"); // wrt event plane | |
54 | double a0[]={0,0.1,0.2,0.3,0.4,0.5,0.7,1.0}; | |
55 | ax[4] = new TAxis(7,a0); ax[4]->SetTitle("(pair pt)/2 (GeV)"); | |
56 | ax[5] = new TAxis(8,0,pi); ax[5]->SetTitle("q-theta"); | |
57 | ax[6] = new TAxis(16,-pi,pi); ax[6]->SetTitle("q-phi"); | |
58 | double a1[100]; | |
59 | for (int i=0;i<20;i++) a1[i]=i*0.005; | |
60 | for (int i=0;i<45;i++) a1[20+i]=0.1+i*0.02; | |
61 | ax[7] = new TAxis(64,a1); ax[7]->SetTitle("q (GeV/c)"); | |
62 | AliUnicorHN *pair = new AliUnicorHN("pair",8,ax); | |
63 | for (int i=0; i<8; i++) delete ax[i]; | |
64 | fHistos.Add(pair); | |
65 | gROOT->cd(); | |
66 | printf("%s object named %s created\n",ClassName(),GetName()); | |
67 | } | |
68 | //============================================================================= | |
69 | void AliUnicorAnalCorrel::Process(Int_t tmr, AliUnicorEvent *ev0, AliUnicorEvent *ev1, Double_t phirot) | |
70 | { | |
71 | // process pairs from one or two (if mixing) events | |
72 | // tmr tells which histogram (bins) to fill: tru,mix,rot | |
73 | ||
74 | // Could be possibly accelerated by checking the "good particle" only once | |
75 | // and caching the result. (Maybe the optimizer does it already.) | |
76 | ||
77 | static TRandom2 ran; | |
78 | AliUnicorHN *pair = (AliUnicorHN*) fHistos.At(0); | |
79 | ||
80 | // mixing-and-rotating-proof centrality and reaction plane angle | |
81 | // (but not rotation-proof for rotation angles much different from 0 and 180) | |
82 | // true and rotated pairs are within the triangle (j<i), mixed - all | |
83 | // thus, proper rotation is either by 180, or by 170 AND 190, etc. | |
84 | ||
85 | double cent = (ev0->Centrality()+ev1->Centrality())/2.0; | |
86 | double q0x,q0y,q1x,q1y; | |
87 | ev0->RP(q0x,q0y); | |
88 | ev1->RP(q1x,q1y); | |
89 | double rpphi = atan2(q0y+q1y,q0x+q1x); | |
90 | ||
91 | // loop over pairs | |
92 | ||
93 | for (int i=0; i<ev0->NParticles(); i++) { | |
94 | if (!ev0->ParticleGood(i,fPid0)) continue; | |
95 | for (int j=0; j<ev1->NParticles(); j++) { | |
96 | if (ev0==ev1 && j<i && fPid0==fPid1 ) continue; | |
97 | if (ev0==ev1 && j==i) continue; // beware, not even when rotated or non-identical | |
98 | if (!ev1->ParticleGood(j,fPid1)) continue; | |
99 | if (!ev0->PairGood(ev0->ParticleP(i),ev0->ParticleTheta(i),ev0->ParticlePhi(i), | |
100 | ev1->ParticleP(j),ev1->ParticleTheta(j),ev1->ParticlePhi(j)+phirot)) continue; | |
101 | fPa.Set0(fMass0,ev0->ParticleP(i),ev0->ParticleTheta(i),ev0->ParticlePhi(i)); | |
102 | fPa.Set1(fMass1,ev1->ParticleP(j),ev1->ParticleTheta(j),ev1->ParticlePhi(j)+phirot); | |
103 | if (ev0==ev1 && fPid0==fPid1 && ran.Rndm()>=0.5) fPa.Swap(); | |
104 | fPa.CalcLAB(); | |
105 | fPa.CalcPairCM(); | |
106 | if (fPa.QCM()==0) return; // should not be too frequent | |
107 | double phi = TVector2::Phi_mpi_pi(fPa.Phi()-rpphi); | |
108 | pair->Fill(tmr, // 0 for tru, 1 for mix, 2 for rot | |
109 | cent, // centrality | |
110 | fPa.Rapidity(), // pair rapidity | |
111 | phi, // pair phi wrt reaction plane | |
112 | fPa.Pt()/2.0, // half of pair pt | |
113 | fPa.QCMTheta(), // polar angle of Q | |
114 | fPa.QCMPhiOut(), // azimuthal angle of Q w.r.t. out | |
115 | fPa.QCM(), // |p2-p1| in c.m.s. | |
116 | 1.0); // weigth | |
117 | } | |
118 | } | |
119 | } | |
120 | //============================================================================= |