]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliV0FindVertices.C
Example macro for the creation of tags (P.Christakoglou)
[u/mrichter/AliRoot.git] / ITS / AliV0FindVertices.C
CommitLineData
a62f4fcc 1/****************************************************************************
2 * Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch *
3 ****************************************************************************/
4
566abf75 5#if !defined(__CINT__) || defined(__MAKECINT__)
548b1270 6 #include "Riostream.h"
ca28c5f5 7 #include "AliV0vertexer.h"
8 #include "TFile.h"
a62f4fcc 9 #include "TKey.h"
ca28c5f5 10 #include "TStopwatch.h"
566abf75 11
12 #include "AliRun.h"
c84a5e9e 13 #include "AliTracker.h"
a62f4fcc 14 #include "AliMagF.h"
15 #include "AliESD.h"
566abf75 16 #include "AliRunLoader.h"
ca28c5f5 17#endif
18
566abf75 19extern AliRun *gAlice;
ca28c5f5 20
566abf75 21Int_t AliV0FindVertices(Int_t nev=5) {
22 cerr<<"Looking for V0 vertices...\n";
ca28c5f5 23
566abf75 24 if (gAlice) {
25 delete gAlice->GetRunLoader();
26 delete gAlice;
27 gAlice=0;
28 }
a62f4fcc 29
566abf75 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 }
a62f4fcc 35
36 if (rl->LoadgAlice()) {
37 cerr<<"AliV0FindVertices.C : LoadgAlice returned error"<<endl;
38 delete rl;
39 return 3;
566abf75 40 }
ca28c5f5 41
c84a5e9e 42 // Magnetic field
43 AliTracker::SetFieldMap(gAlice->Field(),1); // 1 means uniform magnetic field
a62f4fcc 44
ca28c5f5 45 Double_t cuts[]={33, // max. allowed chi2
04b2a5f1 46 0.16,// min. allowed negative daughter's impact parameter
47 0.05,// min. allowed positive daughter's impact parameter
ca28c5f5 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;
5d390eda 54 AliV0vertexer vtxer(cuts);
55 Int_t rc=0;
566abf75 56 if (nev>rl->GetNumberOfEvents()) nev=rl->GetNumberOfEvents();
a62f4fcc 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());
5d390eda 68 for (Int_t i=0; i<nev; i++) {
a62f4fcc 69 itsf->cd();
70 if ((key=(TKey*)next())==0) break;
71 cerr<<"Processing event number: "<<i<<endl;
72 AliESD *event=(AliESD*)key->ReadObj();
566abf75 73
a62f4fcc 74 //Double_t vtx[3]={0.,0.,0.}; vtxer.SetVertex(vtx); // primary vertex (cm)
566abf75 75
a62f4fcc 76 rc=vtxer.Tracks2V0vertices(event);
566abf75 77
a62f4fcc 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;
5d390eda 88 }
ca28c5f5 89 timer.Stop(); timer.Print();
90
a62f4fcc 91 itsf->Close();
92 v0f->Close();
93
566abf75 94 delete rl;
ca28c5f5 95
96 return rc;
97}