]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSFindClustersV2.C
Replace FXS fileID from DASSD_DB_results to CALIBRATION for compatibility with the...
[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 gAlice->GetRunLoader();
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");
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       rawreader->SelectEvents(7);
82     }
83     //    if (!fEquipIdMap.IsNull() && fRawReader)fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
84     Reconstruct(rl,rawreader,opt);
85   }
86   else {
87     cout<< "Starting from DIGITS \n";
88     Reconstruct(rl,opt);
89   }
90
91   return 0;
92 }
93
94 //________________________________________________________________________
95 void Reconstruct(AliRunLoader* runLoader,Option_t *opt){
96 // reconstruct clusters starting from DIGITS
97 // MC truth if available is used to label clusters according to the particles
98 // originating them
99
100   AliITSLoader* loader = 
101                 static_cast<AliITSLoader*>(runLoader->GetLoader("ITSLoader"));
102   if (!loader) {
103     Error("Reconstruct", "ITS loader not found");
104     return;
105   }
106   AliITSDetTypeRec* rec = new AliITSDetTypeRec();
107   rec->SetITSgeom(loader->GetITSgeom());
108   rec->SetDefaults();
109
110   runLoader->LoadKinematics();
111   TTree *trK=(TTree*)runLoader->TreeK();
112   if(trK){
113     cout<<"kine tree found - MC labels will be used in RP \n";
114     if(runLoader->LoadgAlice())gAlice = runLoader->GetAliRun();
115   }
116   else{
117     cout<<"kine tree not found - MC labels will not b used\n";
118   }
119
120   Int_t nEvents = runLoader->GetNumberOfEvents();
121   // loop on the events
122
123   for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
124     runLoader->GetEvent(iEvent);
125     cout<<">>>>>>>   Processing event number "<<iEvent+1<<endl;
126     loader->LoadRecPoints("update");
127     loader->CleanRecPoints();
128     loader->MakeRecPointsContainer();
129     loader->LoadDigits("read");
130     TTree *tR = loader->TreeR();
131     TTree *tD = loader->TreeD();
132     if(!tR){
133       cout<<"Tree R pointer not found - Abort \n";
134       break;
135     }
136     rec->SetTreeAddressD(tD);
137     rec->MakeBranch(tR,"R");
138     rec->SetTreeAddressR(tR);
139     rec->DigitsToRecPoints(tD,tR,0,opt,kTRUE);    
140     loader->WriteRecPoints("OVERWRITE");
141     loader->UnloadRecPoints();
142     loader->UnloadDigits();
143     runLoader->UnloadKinematics();
144   }
145
146 }
147
148 void Reconstruct(AliRunLoader* runLoader, AliRawReader *rawreader,Option_t *opt){
149 // reconstruct clusters starting from raw data (root or DATE format)
150 // MC truth if available is used to label clusters according to the particles
151 // originating them
152
153   AliITSLoader* loader = static_cast<AliITSLoader*>(runLoader->GetLoader("ITSLoader"));
154   if (!loader) {
155     Error("Reconstruct", "ITS loader not found");
156     return;
157   }
158   AliITSDetTypeRec* rec = new AliITSDetTypeRec();
159   rec->SetITSgeom(loader->GetITSgeom());
160   rec->SetDefaults();
161   // direct clusterfinding starting from raw data is implemented only 
162   // in AliITSClusterFinderV2
163   rec->SetDefaultClusterFindersV2(kTRUE);
164
165   runLoader->LoadKinematics();
166   TTree *trK=(TTree*)runLoader->TreeK();
167   if(trK){
168     cout<<"kine tree found - MC labels will be used in RP \n";
169     if(runLoader->LoadgAlice())gAlice = runLoader->GetAliRun();
170   }
171   Int_t nEvents = runLoader->GetNumberOfEvents();
172   rawreader->RewindEvents();
173   for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
174     rawreader->NextEvent();
175     runLoader->GetEvent(iEvent);
176     cout<<">>>>>>>   Processing event number: "<<iEvent<<endl;
177     loader->LoadRecPoints("update");
178     loader->CleanRecPoints();
179     loader->MakeRecPointsContainer();
180     TTree *tR = loader->TreeR();
181     if(!tR){
182       cout<<"Tree R pointer not found - Abort \n";
183       break;
184     }
185     rec->DigitsToRecPoints(rawreader,tR,opt);
186     rec->ResetRecPoints();
187     rec->ResetClusters();    
188     loader->WriteRecPoints("OVERWRITE");
189     loader->UnloadRecPoints();
190     loader->UnloadDigits();
191     runLoader->UnloadKinematics();
192   }
193
194 }