]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSHits2FastRecPoints.C
Methos SetType added
[u/mrichter/AliRoot.git] / ITS / AliITSHits2FastRecPoints.C
1 TFile* AccessFile(TString inFile="galice.root", TString acctype="R");
2
3 void AliITSHits2FastRecPoints (Int_t evNumber1=0,Int_t evNumber2=0, TString inFile = "galice.root", TString outFile="galice.root", Int_t nsignal=25, Int_t size=-1) 
4 {
5   /////////////////////////////////////////////////////////////////////////
6   //   
7   //   This macro creates fast recpoints, optionally on a separate file
8   //   
9   /////////////////////////////////////////////////////////////////////////
10
11   // Dynamically link some shared libs
12
13   if (gClassTable->GetID("AliRun") < 0) {
14     gROOT->LoadMacro("loadlibs.C");
15     loadlibs();
16   } else {
17     delete gAlice;
18     gAlice=0;
19   }
20
21   // Connect the Root Galice file containing Geometry, Kine and Hits
22   TFile *file;
23   if(outFile.Data() == inFile.Data()){
24     file = AccessFile(inFile,"U");
25   }
26   else {
27     file = AccessFile(inFile);
28   }
29
30   TFile *file2 = 0;   // possible output file for TreeR
31   if(!(outFile.Data() == inFile.Data())){
32       // open output file and create TreeR on it
33     file2 = gAlice->InitTreeFile("R",outFile);
34   }
35
36   AliITS *ITS  = (AliITS*) gAlice->GetModule("ITS");
37   if (!ITS) return;
38
39   // Set the simulation model
40
41   for (Int_t i=0;i<3;i++) {
42     ITS->SetSimulationModel(i,new AliITSsimulationFastPoints());
43   }
44    
45
46   //
47   // Event Loop
48   //
49
50   Int_t nbgr_ev=0;
51   TStopwatch timer;
52
53   cout << "Creating fast reconstructed points from hits for the ITS..." << endl;
54
55   for (int ev=evNumber1; ev<= evNumber2; ev++) {
56     cout << "...working on event "<< ev << " ..." << endl;
57     Int_t nparticles = gAlice->GetEvent(ev);
58     cout << "event         " <<ev<<endl;
59     cout << "nparticles  " <<nparticles<<endl;
60     gAlice->SetEvent(ev);
61     if(!gAlice->TreeR() && file2 == 0) gAlice-> MakeTree("R");
62     if(!gAlice->TreeR() && file2 != 0) gAlice->MakeTree("R",file2);
63     ITS->MakeBranch("RF");
64     if (ev < evNumber1) continue;
65     if (nparticles <= 0) return;
66
67     Int_t bgr_ev=Int_t(ev/nsignal);
68     //printf("bgr_ev %d\n",bgr_ev);
69     timer.Start();
70     ITS->HitsToFastRecPoints(ev,bgr_ev,size," ","All"," ");
71     timer.Stop(); timer.Print();
72   } // event loop 
73
74   delete gAlice; gAlice=0;
75   file->Close();
76 }
77
78
79 //-------------------------------------------------------------------
80 TFile * AccessFile(TString FileName, TString acctype){
81
82   // Function used to open the input file and fetch the AliRun object
83
84   TFile *retfil = 0;
85   TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(FileName);
86   if (file) {file->Close(); delete file; file = 0;}
87   if(acctype.Contains("U")){
88     file = new TFile(FileName,"update");
89   }
90   if(acctype.Contains("N") && !file){
91     file = new TFile(FileName,"recreate");
92   }
93   if(!file) file = new TFile(FileName);   // default readonly
94   if (!file->IsOpen()) {
95         cerr<<"Can't open "<<FileName<<" !" << endl;
96         return retfil;
97   } 
98
99   // Get AliRun object from file or return if not on file
100   if (gAlice) {delete gAlice; gAlice = 0;}
101   gAlice = (AliRun*)file->Get("gAlice");
102   if (!gAlice) {
103         cerr << "AliRun object not found on file"<< endl;
104         return retfil;
105   } 
106   return file;
107 }