]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/DoVerticesSPD.C
Fix for pre_VTH choice (A. Mastroserio)
[u/mrichter/AliRoot.git] / ITS / DoVerticesSPD.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include "TFile.h"
3 #include "TGeoManager.h"
4 #include "TInterpreter.h"
5 #include "TClassTable.h"
6 #include "TNtuple.h"
7 #include "TTree.h"
8 #include "TArrayF.h"
9 #include "TParticle.h"
10 #include "AliESD.h"
11 #include "AliGenPythiaEventHeader.h"
12 #include "AliTracker.h"
13 #include "AliHeader.h"
14 #include "AliITSLoader.h"
15 #include "AliVertexerTracks.h"
16 #include "AliCDBManager.h"
17 #include "AliGeomManager.h"
18 #include "AliITSDetTypeRec.h"
19 #include "AliITSVertexer3D.h"
20 #include "AliITSVertexerZ.h"
21 #include "AliESDVertex.h"
22 #include "AliESDEvent.h"
23 #include "AliRun.h"
24 #include "AliRunLoader.h"
25 #include "AliGenEventHeader.h"
26 #include "AliStack.h"
27 #endif
28
29 /*  $Id$    */
30
31 Bool_t DoVerticesSPD(Int_t pileupalgo=1, Int_t optdebug=0){
32
33   TFile *fint = new TFile("VertexSPD.root","recreate");
34   TNtuple *nt = new TNtuple("ntvert","vertices","xtrue:ytrue:ztrue:zZ:zdiffZ:zerrZ:ntrksZ:x3D:xdiff3D:xerr3D:y3D:ydiff3D:yerr3D:z3D:zdiff3D:zerr3D:ntrks3D:dndy:ntrklets:nrp1:ptyp:is3D:isTriggered");
35
36   if (gClassTable->GetID("AliRun") < 0) {
37     gInterpreter->ExecuteMacro("loadlibs.C");
38   }
39   // Set OCDB if needed
40   AliCDBManager* man = AliCDBManager::Instance();
41   if (!man->IsDefaultStorageSet()) {
42     printf("Setting a local default storage and run number 0\n");
43     man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
44     man->SetRun(0);
45   }
46   else {
47     printf("Using deafult storage \n");
48   }
49  
50   // retrives geometry 
51   if(!gGeoManager){
52     AliGeomManager::LoadGeometry("geometry.root");
53   }
54
55   AliRunLoader* runLoader = AliRunLoader::Open("galice.root");
56   if (!runLoader) {
57     Error("DoVertices", "getting run loader from file %s failed", 
58           "galice.root");
59     return kFALSE;
60   }
61   runLoader->LoadgAlice();
62   gAlice = runLoader->GetAliRun();
63   if (!gAlice) {
64     Error("DoVertices", "no galice object found");
65     return kFALSE;
66   }
67   runLoader->LoadKinematics();
68   runLoader->LoadHeader();
69   AliITSLoader* ITSloader =  (AliITSLoader*) runLoader->GetLoader("ITSLoader");
70   ITSloader->LoadRecPoints("read");
71
72  Int_t totev=runLoader->GetNumberOfEvents();
73   if(optdebug)  printf("Number of events= %d\n",totev);
74
75   TFile* esdFile = TFile::Open("AliESDs.root");
76   if (!esdFile || !esdFile->IsOpen()) {
77     printf("Error in opening ESD file");
78     return kFALSE;
79   }
80
81   AliESDEvent * esd = new AliESDEvent;
82   TTree* tree = (TTree*) esdFile->Get("esdTree");
83   if (!tree) {
84     printf("Error: no ESD tree found");
85     return kFALSE;
86   }
87   esd->ReadFromTree(tree);
88   tree->GetEvent(0);
89   Int_t fieldkG=(Int_t)(esd->GetMagneticField()+0.001);
90   
91   if(fieldkG==5){
92     TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", 2, 1., 1., 10., AliMagF::k5kG));
93   }else if(fieldkG==2){
94     TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", 2, 1., 1., 10., AliMagF::k2kG));
95   }else{
96     TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", 2, 0., 1., 10., AliMagF::k5kG));
97   }
98   AliMagF* fld = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
99   printf("Magnetic field set to %f\n",-fld->SolenoidField());
100
101   AliITSDetTypeRec* detTypeRec = new AliITSDetTypeRec();
102   Double_t xnom=0.,ynom=0.;
103   AliITSVertexerZ *vertz = new AliITSVertexerZ(xnom,ynom);
104   vertz->Init("default");
105   AliITSVertexer3D *vert3d = new AliITSVertexer3D();
106   vert3d->Init("default");
107   vert3d->SetWideFiducialRegion(40.,1.);
108   vert3d->SetPileupAlgo(pileupalgo);
109   vert3d->PrintStatus();
110   vertz->SetDetTypeRec(detTypeRec);
111   vert3d->SetDetTypeRec(detTypeRec);
112
113   /* uncomment these lines to use diamond constrain */
114 //   Double_t posdiam[3]={0.03,0.1,0.};
115 //   Double_t sigdiam[3]={0.01,0.01,10.0};
116 //   AliESDVertex* diam=new AliESDVertex(posdiam,sigdiam);
117 //   vertz->SetVtxStart(diam);
118 //   vert3d->SetVtxStart(diam);
119   /* end lines to be uncommented to use diamond constrain */
120
121   Int_t goodz=0,good3d=0;
122   // Trigger mask
123   ULong64_t spdFO = (1 << 14);
124   ULong64_t v0left = (1 << 11);
125   ULong64_t v0right = (1 << 12);
126
127   for (Int_t iEvent = 0; iEvent < totev; iEvent++) {
128     TArrayF mcVertex(3); 
129     runLoader->GetEvent(iEvent);
130     runLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
131     AliGenPythiaEventHeader *evh=(AliGenPythiaEventHeader*)runLoader->GetHeader()->GenEventHeader();
132     Int_t ptype = evh->ProcessType();
133     if(optdebug){
134       printf("==============================================================\n");
135       printf("\nEvent: %d ---- Process Type = %d \n",iEvent,ptype);
136     }
137
138     AliStack* stack = runLoader->Stack();
139     TTree *treek=(TTree*)runLoader->TreeK();
140     Int_t npart = (Int_t)treek->GetEntries();
141     if(optdebug) printf("particles  %d\n",npart);
142
143     Double_t dNchdy = 0.;
144
145    // loop on particles to get generated dN/dy
146     for(Int_t iPart=0; iPart<npart; iPart++) {
147       if(!stack->IsPhysicalPrimary(iPart)) continue;
148       TParticle* part = (TParticle*)stack->Particle(iPart);
149       if(part->GetPDG()->Charge() == 0) continue;
150       Double_t eta=part->Eta();
151
152       if(TMath::Abs(eta)<1.5) dNchdy+=1.; 
153     }
154     if(optdebug) printf(" dNch/dy = %f\n",dNchdy);
155  
156     tree->GetEvent(iEvent);
157
158     TTree* cltree = ITSloader->TreeR();
159     ULong64_t triggerMask=esd->GetTriggerMask();
160     // MB1: SPDFO || V0L || V0R
161     Bool_t eventTriggered = (triggerMask & spdFO || ((triggerMask & v0left) || (triggerMask & v0right)));
162     AliESDVertex* vtxz = vertz->FindVertexForCurrentEvent(cltree);
163     AliESDVertex* vtx3d = vert3d->FindVertexForCurrentEvent(cltree);
164     AliMultiplicity *alimult = vert3d->GetMultiplicity();
165     Int_t  ntrklets=0;
166     Int_t nrecp1=0;
167     if(alimult){
168       ntrklets=alimult->GetNumberOfTracklets() ;
169       nrecp1=ntrklets+alimult->GetNumberOfSingleClusters();
170     }
171        
172
173
174     TDirectory *current = gDirectory;
175     fint->cd();
176     Float_t xnt[23];
177
178     Double_t zz = 0.;
179     Double_t zresz = 0.;
180     Double_t zdiffz = 0.;
181
182     Int_t ntrkz = 0;
183     if(vtxz){
184  
185       ntrkz = vtxz->GetNContributors();
186       if(ntrkz>0)goodz++;
187       zz=vtxz->GetZv();
188       zresz =vtxz->GetZRes(); // microns
189       zdiffz = 10000.*(zz-mcVertex[2]); // microns
190     }
191
192     Double_t x3d = 0.;
193     Double_t xerr3d = 0.;
194     Double_t xdiff3d = 0.;
195
196     Double_t y3d = 0.;
197     Double_t yerr3d = 0.;
198     Double_t ydiff3d = 0.;
199
200     Double_t z3d = 0.;
201     Double_t zerr3d = 0.;
202     Double_t zdiff3d = 0.;
203
204     Int_t ntrk3d = 0;
205     Bool_t is3d=kFALSE;
206     if(vtx3d){
207
208       if(vtx3d->IsFromVertexer3D()) is3d=kTRUE;
209       ntrk3d = vtx3d->GetNContributors();
210       if(is3d && ntrk3d>0)good3d++;
211       x3d = vtx3d->GetXv();
212       xerr3d=vtx3d->GetXRes();
213       xdiff3d = 10000.*(x3d-mcVertex[0]);  // microns
214
215       y3d = vtx3d->GetYv();
216       yerr3d=vtx3d->GetYRes();
217       ydiff3d = 10000.*(y3d-mcVertex[1]);  // microns
218
219       z3d = vtx3d->GetZv();
220       zerr3d=vtx3d->GetZRes();
221       zdiff3d = 10000.*(z3d-mcVertex[2]);  // microns
222
223     }
224
225     xnt[0]=mcVertex[0];//x
226     xnt[1]=mcVertex[1];//y
227     xnt[2]=mcVertex[2];//z
228
229     xnt[3]=zz;
230     xnt[4]=zdiffz;
231     xnt[5]=zresz;
232     xnt[6]=ntrkz;
233
234     xnt[7]=x3d;
235     xnt[8]=xdiff3d;
236     xnt[9]=xerr3d;
237     xnt[10]=y3d;
238     xnt[11]=ydiff3d;
239     xnt[12]=yerr3d;
240     xnt[13]=z3d;
241     xnt[14]=zdiff3d;
242     xnt[15]=zerr3d;
243     xnt[16]=ntrk3d;
244
245     xnt[17]=dNchdy;
246     xnt[18]=float(ntrklets);
247     xnt[19]=float(nrecp1);
248     xnt[20]=float(ptype);
249     xnt[21]=float(is3d);
250     xnt[22]=float(eventTriggered);
251     nt->Fill(xnt);
252     current->cd();
253     
254     if(optdebug){
255       printf("\nVertexerZ:  \tz(cm)=%9.5f \tTracklets=%d \tztrue(cm)=%9.5f \tzdiff(um)=%8.2f\n",zz,ntrkz,mcVertex[2],zdiffz);
256       printf("Vertexer3D: \tz(cm)=%9.5f \tTracklets=%d \tztrue(cm)=%9.5f \tzdiff(um)=%8.2f\n",z3d,ntrk3d,mcVertex[2],zdiff3d);
257     }
258     
259   }
260   fint->cd();
261   nt->Write();
262   fint->Close();
263   delete fint;
264   if(optdebug){
265    printf("***********************************************\n");
266    printf("Number of Z vertices found= %d\n",goodz);
267    printf("efficiency=%f\n",float(goodz)/float(totev));
268    printf("Number of 3D vertices found= %d\n",good3d);
269    printf("efficiency=%f\n",float(good3d)/float(totev));
270   }
271   return kTRUE;
272 }