]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliV0FindVertices.C
Now the full chain includes raw data.
[u/mrichter/AliRoot.git] / ITS / AliV0FindVertices.C
1 /****************************************************************************
2  *           Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch                 *
3  ****************************************************************************/
4
5 #if !defined(__CINT__) || defined(__MAKECINT__)
6   #include "Riostream.h"
7   #include "AliV0vertexer.h"
8   #include "TFile.h"
9   #include "TKey.h"
10   #include "TStopwatch.h"
11
12   #include "AliRun.h"
13   #include "AliTracker.h"
14   #include "AliMagF.h"
15   #include "AliESD.h"
16   #include "AliRunLoader.h"
17 #endif
18
19 extern AliRun *gAlice;
20
21 Int_t AliV0FindVertices(Int_t nev=5) {
22    cerr<<"Looking for V0 vertices...\n";
23
24    if (gAlice) {
25       delete gAlice->GetRunLoader();
26       delete gAlice; 
27       gAlice=0;
28    } 
29
30    AliRunLoader* rl = AliRunLoader::Open("galice.root");
31    if (rl == 0x0) {
32       cerr<<"AliV0FindVertices.C : Can not open session RL=NULL"<< endl;
33       return 1;
34    }
35
36    if (rl->LoadgAlice()) {
37       cerr<<"AliV0FindVertices.C : LoadgAlice returned error"<<endl;
38       delete rl;
39       return 3;
40    }
41
42    // Magnetic field
43    AliTracker::SetFieldMap(gAlice->Field(),1); // 1 means uniform magnetic field
44        
45    Double_t cuts[]={33,  // max. allowed chi2
46                     0.16,// min. allowed negative daughter's impact parameter 
47                     0.05,// min. allowed positive daughter's impact parameter 
48                     0.080,// max. allowed DCA between the daughter tracks
49                     0.998,// max. allowed cosine of V0's pointing angle
50                     0.9,  // min. radius of the fiducial volume
51                     2.9   // max. radius of the fiducial volume
52                    };
53    TStopwatch timer;
54    AliV0vertexer vtxer(cuts);
55    Int_t rc=0;
56    if (nev>rl->GetNumberOfEvents()) nev=rl->GetNumberOfEvents();
57
58    TFile *v0f=TFile::Open("AliESDv0.root","RECREATE");
59    if ((!v0f)||(!v0f->IsOpen())) {
60       cerr<<"Can't AliESDv0.root !\n"; return 1;
61    }
62    TFile *itsf=TFile::Open("AliESDits.root");
63    if ((!itsf)||(!itsf->IsOpen())) {
64       cerr<<"Can't AliESDits.root !\n"; return 1;
65    }
66    TKey *key=0;
67    TIter next(itsf->GetListOfKeys());
68    for (Int_t i=0; i<nev; i++) {
69      itsf->cd();
70      if ((key=(TKey*)next())==0) break;
71      cerr<<"Processing event number: "<<i<<endl;
72      AliESD *event=(AliESD*)key->ReadObj();
73
74      //Double_t vtx[3]={0.,0.,0.}; vtxer.SetVertex(vtx); // primary vertex (cm)
75
76      rc=vtxer.Tracks2V0vertices(event);
77
78      if (rc==0) {
79         Char_t ename[100]; 
80         sprintf(ename,"%d",i);
81         v0f->cd();
82         if (!event->Write(ename)) rc++;
83      } 
84      if (rc) {
85         cerr<<"Something bad happened...\n";
86      }
87      delete event;
88    }
89    timer.Stop(); timer.Print();
90     
91    itsf->Close();
92    v0f->Close();
93
94    delete rl;
95
96    return rc;
97 }