]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSFindClustersV2.C
Increase speed of SDD cluster finder (F. Prino)
[u/mrichter/AliRoot.git] / ITS / AliITSFindClustersV2.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2
3 #include <Riostream.h>
4 #include <TROOT.h>
5 #include <TString.h>
6 #include "AliRawReader.h"
7 #include "AliRawReaderDate.h"
8 #include "AliRawReaderRoot.h"
9 #include "AliRun.h"
10 #include "AliRunLoader.h"
11 #include "AliITSInitGeometry.h"
12 #include "AliITSgeom.h"
13 #include "AliITSLoader.h"
14 #include "AliCDBManager.h"
15 #include "AliITSDetTypeRec.h"
16 #include "AliGeomManager.h"
17
18 #endif
19
20 /*
21 $Id$ 
22 */
23 /***************************************************************
24  *  This macro performs the ITS local reconstruction           *
25  *  It is intended for special purposes and tests.             *
26  *  The reccomended way to reconstruct ITs is via the          *
27  *  class STEER/AliReconstruction                              *
28  *  Present version: M.Masera - Previous version: J. Belikov   *
29  ***************************************************************/
30
31 void Reconstruct(AliRunLoader* runLoader,Option_t *opt);
32 void Reconstruct(AliRunLoader* runLoader, AliRawReader *rawreader,Option_t *opt);
33
34 Int_t AliITSFindClustersV2(char *inputRawData = NULL,TString filename="galice.root",Option_t *opt="All"){
35   // if kine tree is available MC labels are used for rec points
36   // set opt equal to "SPD" or to "SDD" or to "SSD" do build
37   // rec points for individual subdetectors 
38   if (gAlice) {
39     delete AliRunLoader::Instance();
40     delete gAlice;
41     gAlice = NULL;
42   }
43
44   // Get geometry
45   AliGeomManager::LoadGeometry("geometry.root");
46
47   //Get Run loader and ITS loader - set geometry
48   AliRunLoader* rl = AliRunLoader::Open(filename.Data());
49
50   AliITSInitGeometry initgeom;
51   AliITSgeom *geom = initgeom.CreateAliITSgeom();
52   printf("Geometry name: %s \n",(initgeom.GetGeometryName()).Data());
53   AliITSLoader* loader = static_cast<AliITSLoader*>(rl->GetLoader("ITSLoader"));
54   if (!loader) {
55     Error("Init", "ITS loader not found");
56     return -1;
57   }
58   loader->SetITSgeom(geom);
59
60   // Set OCDB if needed
61   AliCDBManager* man = AliCDBManager::Instance();
62   if (!man->IsDefaultStorageSet()) {
63     printf("Setting a local default storage and run number 0\n");
64     man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
65     man->SetRun(0);
66   }
67   else {
68     printf("Using deafult storage \n");
69   }
70
71   AliRawReader *rawreader = NULL;
72   TString fileRaw(inputRawData);
73   if(!fileRaw.IsNull()){
74     if (fileRaw.EndsWith(".root")) {
75       cout<<"Raw data format - ROOT file \n"; 
76       rawreader = new AliRawReaderRoot(fileRaw); // root format
77     }
78     else {
79       cout<<"Raw data format - DATE file \n";
80       rawreader = new AliRawReaderDate(fileRaw);  // DATE format
81     }
82     //    if (!fEquipIdMap.IsNull() && fRawReader)fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
83     Reconstruct(rl,rawreader,opt);
84   }
85   else {
86     cout<< "Starting from DIGITS \n";
87     Reconstruct(rl,opt);
88   }
89
90   return 0;
91 }
92
93 //________________________________________________________________________
94 void Reconstruct(AliRunLoader* runLoader,Option_t *opt){
95 // reconstruct clusters starting from DIGITS
96 // MC truth if available is used to label clusters according to the particles
97 // originating them
98
99   AliITSLoader* loader = 
100                 static_cast<AliITSLoader*>(runLoader->GetLoader("ITSLoader"));
101   if (!loader) {
102     Error("Reconstruct", "ITS loader not found");
103     return;
104   }
105   AliITSDetTypeRec* rec = new AliITSDetTypeRec();
106   rec->SetITSgeom(loader->GetITSgeom());
107   rec->SetDefaults();
108
109   runLoader->LoadKinematics();
110   TTree *trK=(TTree*)runLoader->TreeK();
111   if(trK){
112     cout<<"kine tree found - MC labels will be used in RP \n";
113     if(runLoader->LoadgAlice())gAlice = runLoader->GetAliRun();
114   }
115   else{
116     cout<<"kine tree not found - MC labels will not b used\n";
117   }
118
119   Int_t nEvents = runLoader->GetNumberOfEvents();
120   // loop on the events
121
122   for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
123     runLoader->GetEvent(iEvent);
124     cout<<">>>>>>>   Processing event number "<<iEvent+1<<endl;
125     loader->LoadRecPoints("update");
126     loader->CleanRecPoints();
127     loader->MakeRecPointsContainer();
128     loader->LoadDigits("read");
129     TTree *tR = loader->TreeR();
130     TTree *tD = loader->TreeD();
131     if(!tR){
132       cout<<"Tree R pointer not found - Abort \n";
133       break;
134     }
135     rec->SetTreeAddressD(tD);
136     rec->MakeBranch(tR,"R");
137     rec->SetTreeAddressR(tR);
138     rec->DigitsToRecPoints(tD,tR,0,opt,kTRUE);    
139     loader->WriteRecPoints("OVERWRITE");
140     loader->UnloadRecPoints();
141     loader->UnloadDigits();
142     runLoader->UnloadKinematics();
143   }
144
145 }
146
147 void Reconstruct(AliRunLoader* runLoader, AliRawReader *rawreader,Option_t *opt){
148 // reconstruct clusters starting from raw data (root or DATE format)
149 // MC truth if available is used to label clusters according to the particles
150 // originating them
151
152   AliITSLoader* loader = static_cast<AliITSLoader*>(runLoader->GetLoader("ITSLoader"));
153   if (!loader) {
154     Error("Reconstruct", "ITS loader not found");
155     return;
156   }
157   AliITSDetTypeRec* rec = new AliITSDetTypeRec();
158   rec->SetITSgeom(loader->GetITSgeom());
159   rec->SetDefaults();
160   // direct clusterfinding starting from raw data is implemented only 
161   // in AliITSClusterFinderV2
162   rec->SetDefaultClusterFindersV2(kTRUE);
163
164   runLoader->LoadKinematics();
165   TTree *trK=(TTree*)runLoader->TreeK();
166   if(trK){
167     cout<<"kine tree found - MC labels will be used in RP \n";
168     if(runLoader->LoadgAlice())gAlice = runLoader->GetAliRun();
169   }
170   Int_t nEvents = runLoader->GetNumberOfEvents();
171   rawreader->RewindEvents();
172   for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
173     rawreader->NextEvent();
174     runLoader->GetEvent(iEvent);
175     cout<<">>>>>>>   Processing event number: "<<iEvent<<endl;
176     loader->LoadRecPoints("update");
177     loader->CleanRecPoints();
178     loader->MakeRecPointsContainer();
179     TTree *tR = loader->TreeR();
180     if(!tR){
181       cout<<"Tree R pointer not found - Abort \n";
182       break;
183     }
184     rec->DigitsToRecPoints(rawreader,tR,opt);
185     rec->ResetRecPoints();
186     rec->ResetClusters();    
187     loader->WriteRecPoints("OVERWRITE");
188     loader->UnloadRecPoints();
189     loader->UnloadDigits();
190     runLoader->UnloadKinematics();
191   }
192
193 }