]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSRecPointContainer.cxx
track labels added
[u/mrichter/AliRoot.git] / ITS / AliITSRecPointContainer.cxx
1 #include <TTree.h>
2 #include "AliLog.h"
3 #include "AliITSRecoParam.h"
4 #include "AliITSReconstructor.h"
5 #include "AliITSRecPointContainer.h"
6 #include "AliITSRecPoint.h"
7 #include "AliRunLoader.h"
8
9 ClassImp(AliITSRecPointContainer)
10
11 //////////////////////////////////////////////////////////////////////
12 // Class to store ITS RecPoints for the duration of                 //
13 // one event processing                                             //
14 // The container is cleared at each event and new RP                //
15 // are loaded from TTree                                            //
16 // Origin masera@to.infn.it  Nov. 12 2009                           //
17 //////////////////////////////////////////////////////////////////////
18
19 /* $Id$ */
20
21 AliITSRecPointContainer* AliITSRecPointContainer::fgInstance = 0x0;
22
23 //______________________________________________________________________
24 AliITSRecPointContainer::AliITSRecPointContainer(const AliITSRecoParam* krp):TObject(),
25 fSPDNModules(0),
26 fSDDNModules(0),
27 fSSDNModules(0),
28 fArray(),
29 fCurrentEve(-1000),
30 fNextEvent(-1000),
31 fActualSize(0),
32 fDet(""),
33 fStatusOK(kTRUE){
34   // Default constructor
35
36   for(Int_t i=0;i<6;i++)fNClusters[i]=0;
37   if(fgkNModules != AliITSgeomTGeo::GetNModules())AliError(Form("The total number of modules is not %d, but %d",fgkNModules,AliITSgeomTGeo::GetNModules()));
38
39   Int_t modperlay[6];
40   for(Int_t i=0;i<6;i++)modperlay[i]=AliITSgeomTGeo::GetNDetectors(1+i)*AliITSgeomTGeo::GetNLadders(1+i);
41   fSPDNModules=modperlay[0]+modperlay[1];
42   fSDDNModules=modperlay[2]+modperlay[3];
43   fSSDNModules=modperlay[4]+modperlay[5];
44   //  AliInfo(Form("Total modules: %d \n SPD modules=%d , SDD modules=%d, SSD modules=%d ",fgkNModules,fSPDNModules,fSDDNModules,fSSDNModules));
45
46   // kLimits[0:5] --> low fluw; kLimits[6,11] --> High flux
47   const Int_t kLimits[12]={25,25,20,20,10,10,300,300,200,200,100,100};
48   Int_t offset=0;
49   if(!krp){
50     AliWarning("AliITSRecoParam is missing. Using defaults");
51   }
52   else {
53     if(krp->GetEventSpecie() == AliRecoParam::kHighMult)offset=6;
54   }
55   Int_t maxval[6];
56   TString values="";
57   for(Int_t i=0;i<6;i++){
58     maxval[i]=kLimits[i+offset];
59     values+=maxval[i];
60     values+=" ";
61     if(i>0)modperlay[i]+=modperlay[i-1];
62   }
63   AliInfo(Form("Container created with sizes/layer: %s",values.Data()));
64   Int_t layer=0;
65   for(Int_t i=0;i<fgkNModules;i++){
66     if(i>=modperlay[layer])++layer;
67     fArray[i]=new TClonesArray("AliITSRecPoint",maxval[layer]);
68   }
69 }
70
71
72 //______________________________________________________________________
73 AliITSRecPointContainer::~AliITSRecPointContainer(){
74   // Destructor
75   for(Int_t i=0;i<fgkNModules;i++){
76     if(fArray[i]){
77       fArray[i]->Delete();
78       delete fArray[i];
79     }
80   }
81 }
82
83 //______________________________________________________________________
84 void AliITSRecPointContainer::CookEntries(){
85   // From the number of entries in TTree R, the number of ITS subdetectors
86   // active for the present run is inferred
87   if(fActualSize == fgkNModules)fDet="ALL SPD SDD SSD ";
88   if(fActualSize == fSPDNModules) fDet = "SPD ";
89   if(fActualSize == fSDDNModules) fDet = "SDD ";
90   if(fActualSize == fSSDNModules)fDet = "SSD ";
91   if(fActualSize == (fSPDNModules+fSDDNModules)) fDet = "SPD SDD ";
92   if(fActualSize == (fSPDNModules+fSSDNModules))fDet = "SPD SSD ";
93   if(fActualSize == (fSDDNModules+fSSDNModules))fDet = "SDD SSD ";
94   if((!fDet.Contains("SPD")) && (!fDet.Contains("SDD")) &&
95      (!fDet.Contains("SSD"))){
96     AliError(Form("The number of active modules %d does not correspond to any standard configuration of the detector",fActualSize));
97     fStatusOK = kFALSE;
98   }
99 }
100 //______________________________________________________________________
101 TClonesArray* AliITSRecPointContainer::FetchClusters(Int_t mod, TTree* tR){
102   // retrieves Recpoints for module mod (offline mode: the event number is
103   // retrieved via the AliRunLoader object)
104   // The actual access to the RP TTree is done as follows:
105   // If the AliRunLoader object exists, the event number is taken from it
106   // If not, the data member fNextEvent is used. 
107   // To set fNextEvent it is necessary to call PrepareToRead in advance.
108   // if this is never done, fNextEvent will have its default negative value
109   // and an error message will be delivered.
110   AliRunLoader* rl = AliRunLoader::Instance();
111   Int_t cureve;
112   if(rl){
113     cureve = rl->GetEventNumber();
114   }
115   else if(fNextEvent>=0){
116     cureve = fNextEvent;
117   }
118   else {
119     AliError("The RunLoader is not defined, PrepareToRead was not invoked. Revise calling sequence. Nothing done");
120     return NULL;
121   }
122   return FetchClusters(mod,tR,cureve);
123 }
124 //______________________________________________________________________
125 TClonesArray* AliITSRecPointContainer::FetchClusters(Int_t mod, TTree* tR,Int_t cureve){
126   // retrieves Recpoints for module mod
127   // cureve is the current event number. If it is different w.r.t.
128   // the event number stored in fCurrentEve, the recpoints are read from
129   // the TTree. Otherwise, the RP stored in memory are used. 
130   if(cureve != fCurrentEve){
131     fCurrentEve = cureve;
132     Reset();
133     TBranch *branch = NULL;
134     branch = tR->GetBranch("ITSRecPoints");
135     if(!branch){
136       AliError("Branch ITSRecPoints not found on ITS recpoints TTree");
137       fStatusOK = kFALSE;
138       return NULL;
139     }
140
141     fActualSize = branch->GetEntries();
142     CookEntries();
143     if(fDet.IsNull())return NULL;
144     // it is assumed that the filling order of the tree is SPD, SDD, SSD
145     // even if one or two subdetector are missing
146     Int_t modL1=AliITSgeomTGeo::GetNDetectors(1)*AliITSgeomTGeo::GetNLadders(1);
147     if(IsSPDActive()){
148       for(Int_t i=0;i<fSPDNModules;i++){
149         branch->SetAddress(&fArray[i]);
150         branch->GetEvent(i);
151         if(i<modL1){
152           fNClusters[0]+=fArray[i]->GetEntries();
153         }
154         else {
155           fNClusters[1]+=fArray[i]->GetEntries();
156         }
157       }
158     }
159     if(IsSDDActive()){
160       Int_t start=0;
161       if(IsSPDActive())start+=fSPDNModules;
162       Int_t modL3=AliITSgeomTGeo::GetNDetectors(3)*AliITSgeomTGeo::GetNLadders(3);
163       Int_t counter = fSPDNModules;
164       for(Int_t i=start;i<start+fSDDNModules;i++){
165         branch->SetAddress(&fArray[counter]);
166         ++counter;
167         branch->GetEvent(i);
168         if((i-start)<modL3){
169           fNClusters[2]+=fArray[i]->GetEntries();
170         }
171         else {
172           fNClusters[3]+=fArray[i]->GetEntries();
173         }
174       }
175     }
176     if(IsSSDActive()){
177       Int_t start=0;
178       if(IsSPDActive())start+=fSPDNModules;
179       if(IsSDDActive())start+=fSDDNModules;
180       Int_t modL5=AliITSgeomTGeo::GetNDetectors(5)*AliITSgeomTGeo::GetNLadders(5);
181       Int_t counter = fSPDNModules+fSDDNModules;
182       for(Int_t i=start;i<start+fSSDNModules;i++){
183         branch->SetAddress(&fArray[counter]);
184         ++counter;
185         branch->GetEvent(i);
186         if((i-start)<modL5){
187           fNClusters[4]+=fArray[i]->GetEntries();
188         }
189         else {
190           fNClusters[5]+=fArray[i]->GetEntries();
191         }
192       }
193     }
194   }
195
196   if(CheckBoundaries(mod)){
197     return fArray[mod];
198   }
199   else {
200     AliError(Form("Module %d is out of boundaries",mod));
201     return NULL;
202   }
203   
204 }
205 //______________________________________________________________________
206 UInt_t AliITSRecPointContainer::GetNClustersInLayer(Int_t lay, TTree* tR, Int_t eventN){
207   // returns the number of clusters for laier lay
208   // layers are numbered from 1 to 6
209   if(lay<1 || lay >6){
210     AliError(Form("Layer %d is out of range",lay));
211     return 0;
212   }
213   if(eventN>=0){
214     FetchClusters(0,tR,eventN);
215   }
216   else {
217     FetchClusters(0,tR);
218   }
219   return fNClusters[lay-1];
220 }
221 //______________________________________________________________________
222 UInt_t AliITSRecPointContainer::GetNClustersInLayerFast(Int_t lay) const {
223   // returns the number of clusters for laier lay
224   // layers are numbered from 1 to 6
225   // No checks are done on the event number: the numer of clusters 
226   // for the event stored in memory is returned
227   if(lay<1 || lay >6){
228     AliError(Form("Layer %d is out of range",lay));
229     return 0;
230   }
231   return fNClusters[lay-1];
232 }
233 //______________________________________________________________________
234 AliITSRecPointContainer* AliITSRecPointContainer::Instance(const AliITSRecoParam* kptr){
235   // returns AliITSRecPointContainer instance (singleton)
236   if(!fgInstance){
237     if(!kptr){
238       fgInstance =  new AliITSRecPointContainer(AliITSReconstructor::GetRecoParam());
239     }
240     else {
241     fgInstance = new AliITSRecPointContainer(kptr);
242     }
243   }
244   return fgInstance;
245 }
246
247 //______________________________________________________________________
248 void AliITSRecPointContainer::Reset(){
249   // Resets the status of the object
250   ClearClus(0,fgkNModules);
251   fDet="";
252   for(Int_t i=0;i<6;i++)fNClusters[i]=0;
253 }
254 //______________________________________________________________________
255 void AliITSRecPointContainer::ResetSPD(){
256   // Resets only the entries in fArray concerning SPD
257   // This method should be used with care only when the filling
258   // of the container is not done from the RP TTree. 
259   fCurrentEve = -1000;  // protection: if FetchClusters method will be used
260                           // after this call, an ccess to the RP TTree will
261                           // be forced
262   ClearClus(0,fSPDNModules);
263 }
264
265 //______________________________________________________________________
266 void AliITSRecPointContainer::ResetSDD(){
267   // Resets only the entries in fArray concerning SDD
268   // This method should be used with care only when the filling
269   // of the container is not done from the RP TTree. 
270   fCurrentEve = -1000;  // protection: if FetchClusters method will be used
271                           // after this call, an ccess to the RP TTree will
272                           // be forced
273   Int_t first = fSPDNModules;
274   Int_t last = first + fSDDNModules; 
275   ClearClus(first,last);
276 }
277
278 //______________________________________________________________________
279 void AliITSRecPointContainer::ResetSSD(){
280   // Resets only the entries in fArray concerning SSD
281   // This method should be used with care only when the filling
282   // of the container is not done from the RP TTree. 
283   fCurrentEve = -1000;  // protection: if FetchClusters method will be used
284                           // after this call, an ccess to the RP TTree will
285                           // be forced
286   Int_t first = fSPDNModules + fSDDNModules;
287   Int_t last = first + fSSDNModules; 
288   ClearClus(first,last);
289 }
290