]> git.uio.no Git - u/mrichter/AliRoot.git/blob - UNICOR/AliDAnalSingle.cxx
Removing -Wno-long-double and adding -single_module to the compilation options (macosx)
[u/mrichter/AliRoot.git] / UNICOR / AliDAnalSingle.cxx
1 // Author: Dariusz Miskowiec <mailto:d.miskowiec@gsi.de> 2007
2
3 //=============================================================================
4 // single particle analyzer
5 //=============================================================================
6
7 #include <cmath>
8 #include <TROOT.h>
9 #include <TMath.h>
10 #include <TAxis.h>
11 #include <TParticlePDG.h>
12 #include "AliDHN.h"
13 #include "AliDEvent.h"
14 #include "AliDAnalSingle.h"
15
16 ClassImp(AliDAnalSingle)
17
18 //=============================================================================
19 AliDAnalSingle::AliDAnalSingle(Char_t *nam, Double_t emi, Double_t ema, Int_t pid) : 
20   AliDAnal(nam), fPid(pid), fMass(0.0) 
21 {
22   // constructor
23   // emi and ema define the rapidity range for histograms
24
25   fPid = pid;
26   TParticlePDG *part = AliDAnal::fgPDG.GetParticle(fPid);
27   fMass = part? part->Mass() : 0;
28
29   double pi = TMath::Pi();
30   TAxis *ax[10];
31   ax[0] = new TAxis(30,-1,1);    ax[0]->SetTitle("vertex z");
32   ax[1] = new TAxis(80,emi,ema); ax[1]->SetTitle("eta");
33   ax[2] = new TAxis(90,-pi,pi);  ax[2]->SetTitle("phi");
34   AliDHN *zep = new AliDHN("zep",3,ax);
35   for (int i=0; i<3; i++) delete ax[i];
36
37   ax[0] = new TAxis(20,0,1);     ax[0]->SetTitle("centrality");
38   ax[1] = new TAxis(80,emi,ema); ax[1]->SetTitle("y");
39   ax[2] = new TAxis(80,0,2);     ax[2]->SetTitle("pt (GeV)");
40   AliDHN *cyp = new AliDHN("cyp",3,ax);
41   for (int i=0; i<3; i++) delete ax[i];
42
43   ax[0] = new TAxis(10,emi,ema); ax[0]->SetTitle("eta");
44   ax[1] = new TAxis(150,0,3);    ax[1]->SetTitle("p (GeV)");
45   ax[2] = new TAxis(150,0.5,3.5);ax[2]->SetTitle("sqrt(dedx (mips))");
46   AliDHN *epd = new AliDHN("epd",3,ax);
47   for (int i=0; i<3; i++) delete ax[i];
48
49   fHistos.Add(zep);
50   fHistos.Add(cyp);
51   fHistos.Add(epd);
52   gROOT->cd();
53   printf("%s object named %s created\n",ClassName(),GetName());
54 }
55 //=============================================================================
56 void AliDAnalSingle::Process(AliDEvent *ev) 
57 {
58   // fill single particle histograms
59
60   AliDHN *zep = (AliDHN*) fHistos.At(0);
61   AliDHN *cyp = (AliDHN*) fHistos.At(1);
62   AliDHN *epd = (AliDHN*) fHistos.At(2);
63   for (int i=0; i<ev->NParticles(); i++) {
64     if (!ev->ParticleGood(i,fPid)) continue;
65     zep->Fill(ev->Zver(),ev->ParticleEta(i),ev->ParticlePhi(i),1.0);
66     double y = fMass>0? ev->ParticleY(i,fMass) : ev->ParticleEta(i);
67     cyp->Fill(ev->Centrality(),y,ev->ParticlePt(i),1.0);
68     epd->Fill(ev->ParticleEta(i),ev->ParticleP(i),sqrt(ev->ParticleDedx(i)),1.0);
69   }
70 }
71 //=============================================================================