]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSDetTypeRec.cxx
Additional protection in the destructor: when you have a chain of calls returning...
[u/mrichter/AliRoot.git] / ITS / AliITSDetTypeRec.cxx
1 /***************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Conributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /*
17  $Id$
18 */
19
20 ////////////////////////////////////////////////////////////////////////
21 // This class defines the "Standard" reconstruction for the ITS       // 
22 // detector.                                                          //
23 //                                                                    //
24 ////////////////////////////////////////////////////////////////////////
25 #include "TObjArray.h"
26 #include "TTree.h"
27
28 #include "AliCDBManager.h"
29 #include "AliCDBStorage.h"
30 #include "AliCDBEntry.h"
31 #include "AliITSClusterFinder.h"
32 #include "AliITSClusterFinderV2.h"
33 #include "AliITSClusterFinderV2SPD.h"
34 #include "AliITSClusterFinderV2SDD.h"
35 #include "AliITSClusterFinderV2SSD.h"
36 #include "AliITSClusterFinderSPD.h"
37 #include "AliITSClusterFinderSDD.h"
38 #include "AliITSClusterFinderSSD.h"
39 #include "AliITSDetTypeRec.h"
40 #include "AliITSgeom.h"
41 #include "AliITSRawCluster.h"
42 #include "AliITSRawClusterSPD.h"
43 #include "AliITSRawClusterSDD.h"
44 #include "AliITSRawClusterSSD.h"
45 #include "AliITSRecPoint.h"
46 #include "AliITSReconstructor.h"
47 #include "AliITSRecoParam.h"
48 #include "AliITSCalibrationSDD.h"
49 #include "AliITSMapSDD.h"
50 #include "AliITSDriftSpeedArraySDD.h"
51 #include "AliITSDriftSpeedSDD.h"
52 #include "AliITSCalibrationSSD.h"
53 #include "AliITSNoiseSSD.h"
54 #include "AliITSGainSSD.h"
55 #include "AliITSBadChannelsSSD.h"
56 #include "AliITSPedestalSSD.h"
57 #include "AliITSsegmentationSPD.h"
58 #include "AliITSsegmentationSDD.h"
59 #include "AliITSsegmentationSSD.h"
60 #include "AliLog.h"
61
62
63 const Int_t AliITSDetTypeRec::fgkNdettypes = 3;
64 const Int_t AliITSDetTypeRec::fgkDefaultNModulesSPD =  240;
65 const Int_t AliITSDetTypeRec::fgkDefaultNModulesSDD =  260;
66 const Int_t AliITSDetTypeRec::fgkDefaultNModulesSSD = 1698;
67
68 ClassImp(AliITSDetTypeRec)
69
70 //________________________________________________________________
71 AliITSDetTypeRec::AliITSDetTypeRec(): TObject(),
72 fNMod(0),
73 fITSgeom(0),
74 fReconstruction(0),
75 fSegmentation(0),
76 fCalibration(0),
77 fSSDCalibration(0),
78 fSPDDead(0),
79 fPreProcess(0),
80 fPostProcess(0),
81 fDigits(0),
82 fDDLMapSDD(0),
83 fRespSDD(0),
84 fNdtype(0),
85 fCtype(0),
86 fNctype(0),
87 fRecPoints(0),
88 fNRecPoints(0),
89 fSelectedVertexer(),
90 fFirstcall(kTRUE){
91     // Standard Constructor
92     // Inputs:
93     //    none.
94     // Outputs:
95     //    none.
96     // Return:
97     //   
98
99   fReconstruction = new TObjArray(fgkNdettypes);
100   fDigits = new TObjArray(fgkNdettypes);
101   for(Int_t i=0; i<3; i++){
102     fClusterClassName[i]=0;
103     fDigClassName[i]=0;
104     fRecPointClassName[i]=0;
105   }
106   fSSDCalibration=new AliITSCalibrationSSD();
107   fNdtype = new Int_t[fgkNdettypes];
108   fCtype = new TObjArray(fgkNdettypes);
109   fNctype = new Int_t[fgkNdettypes];
110   fNMod = new Int_t [fgkNdettypes];
111   fNMod[0] = fgkDefaultNModulesSPD;
112   fNMod[1] = fgkDefaultNModulesSDD;
113   fNMod[2] = fgkDefaultNModulesSSD;
114   fRecPoints = new TClonesArray("AliITSRecPoint",3000);
115   fNRecPoints = 0;
116   
117   for(Int_t i=0;i<fgkNdettypes;i++){
118     fNdtype[i]=0;
119     fNctype[i]=0;
120   }
121   
122   SelectVertexer(" "); 
123 }
124
125 //______________________________________________________________________
126 AliITSDetTypeRec::AliITSDetTypeRec(const AliITSDetTypeRec & rec):TObject(rec),
127 fNMod(rec.fNMod),
128 fITSgeom(rec.fITSgeom),
129 fReconstruction(rec.fReconstruction),
130 fSegmentation(rec.fSegmentation),
131 fCalibration(rec.fCalibration),
132 fSSDCalibration(rec.fSSDCalibration),
133 fSPDDead(rec.fSPDDead),
134 fPreProcess(rec.fPreProcess),
135 fPostProcess(rec.fPostProcess),
136 fDigits(rec.fDigits),
137 fDDLMapSDD(rec.fDDLMapSDD),
138 fRespSDD(rec.fRespSDD),
139 fNdtype(rec.fNdtype),
140 fCtype(rec.fCtype),
141 fNctype(rec.fNctype),
142 fRecPoints(rec.fRecPoints),
143 fNRecPoints(rec.fNRecPoints),
144 fSelectedVertexer(rec.fSelectedVertexer),
145 fFirstcall(rec.fFirstcall)
146 {
147
148   // Copy constructor. 
149
150 }
151 //______________________________________________________________________
152 AliITSDetTypeRec& AliITSDetTypeRec::operator=(const AliITSDetTypeRec& source){
153     // Assignment operator. 
154     this->~AliITSDetTypeRec();
155     new(this) AliITSDetTypeRec(source);
156     return *this;
157
158 }
159
160 //_____________________________________________________________________
161 AliITSDetTypeRec::~AliITSDetTypeRec(){
162   //Destructor
163
164   if(fReconstruction){
165     fReconstruction->Delete();
166     delete fReconstruction;
167     fReconstruction = 0;
168   }
169   if(fSegmentation){
170     fSegmentation->Delete();
171     delete fSegmentation;
172     fSegmentation = 0;
173   }
174   if(fCalibration){
175     if(!(AliCDBManager::Instance()->GetCacheFlag())) {
176       AliITSCalibration * calSPD = (AliITSCalibration*)fCalibration->At(GetITSgeom()->GetStartSPD()); 
177       if (calSPD) {
178         AliITSresponse* rspd = calSPD->GetResponse();    
179         if(rspd) delete rspd;
180       }
181
182       AliITSCalibration * calSSD = (AliITSCalibration*)fCalibration->At(GetITSgeom()->GetStartSSD());
183       if (calSSD) {
184         AliITSresponse* rssd = calSSD->GetResponse();
185         if(rssd) delete rssd;
186       }
187       fCalibration->Delete();
188       delete fCalibration;
189       fCalibration = 0;
190       if(fRespSDD) delete fRespSDD;
191       if(fDDLMapSDD) delete fDDLMapSDD;
192    }
193   }
194   if(fSSDCalibration) delete fSSDCalibration;
195    if(fSPDDead){
196     if(!(AliCDBManager::Instance()->GetCacheFlag())) {
197       fSPDDead->Delete();
198       delete fSPDDead;
199       fSPDDead = 0;
200     }
201   }  
202   if(fPreProcess) delete fPreProcess;
203   if(fPostProcess) delete fPostProcess;
204   if(fDigits){
205     fDigits->Delete();
206     delete fDigits;
207     fDigits=0;
208   }
209   if(fRecPoints){
210     fRecPoints->Delete();
211     delete fRecPoints;
212     fRecPoints=0;
213   }
214   if(fCtype) {
215     fCtype->Delete();
216     delete fCtype;
217     fCtype = 0;
218   }
219   delete [] fNctype;
220   delete [] fNdtype;
221   delete [] fNMod;
222   
223   if (fITSgeom) delete fITSgeom;
224  
225 }
226
227 //___________________________________________________________________
228 void AliITSDetTypeRec::SetReconstructionModel(Int_t dettype,AliITSClusterFinder *clf){
229
230   //Set reconstruction model for detector type
231
232   if(fReconstruction==0) fReconstruction = new TObjArray(fgkNdettypes);
233   if(fReconstruction->At(dettype)!=0) delete fReconstruction->At(dettype);
234   fReconstruction->AddAt(clf,dettype);
235 }
236 //______________________________________________________________________
237 AliITSClusterFinder* AliITSDetTypeRec::GetReconstructionModel(Int_t dettype){
238
239   //Get reconstruction model for detector type
240   if(fReconstruction==0)  {
241     Warning("GetReconstructionModel","fReconstruction is 0!");
242     return 0;     
243   }
244   return (AliITSClusterFinder*)fReconstruction->At(dettype);
245 }
246
247 //______________________________________________________________________
248 void AliITSDetTypeRec::SetSegmentationModel(Int_t dettype,AliITSsegmentation *seg){
249    
250   //Set segmentation model for detector type
251   
252   if(fSegmentation==0) fSegmentation = new TObjArray(fgkNdettypes);
253   if(fSegmentation->At(dettype)!=0) delete fSegmentation->At(dettype);
254   fSegmentation->AddAt(seg,dettype);
255
256 }
257 //______________________________________________________________________
258 AliITSsegmentation* AliITSDetTypeRec::GetSegmentationModel(Int_t dettype){
259
260   //Get segmentation model for detector type
261    
262    if(fSegmentation==0) {
263      Warning("GetSegmentationModel","fSegmentation is 0!");
264      return 0; 
265    } 
266    return (AliITSsegmentation*)fSegmentation->At(dettype);
267
268 }
269 //_______________________________________________________________________
270 void AliITSDetTypeRec::SetCalibrationModel(Int_t iMod, AliITSCalibration *cal){
271
272   //Set calibration (response) for the module iMod of type dettype
273   if (fCalibration==0) {
274     fCalibration = new TObjArray(GetITSgeom()->GetIndexMax());
275     fCalibration->SetOwner(kTRUE);
276     fCalibration->Clear();
277   }
278
279   if (fCalibration->At(iMod) != 0)
280     delete (AliITSCalibration*) fCalibration->At(iMod);
281   fCalibration->AddAt(cal,iMod);
282
283 }
284 //_______________________________________________________________________
285 void AliITSDetTypeRec::SetSPDDeadModel(Int_t iMod, AliITSCalibration *cal){
286
287   //Set dead pixel info for the SPD module iMod
288   if (fSPDDead==0) {
289     fSPDDead = new TObjArray(fgkDefaultNModulesSPD);
290     fSPDDead->SetOwner(kTRUE);
291     fSPDDead->Clear();
292   }
293
294   if (fSPDDead->At(iMod) != 0)
295     delete (AliITSCalibration*) fSPDDead->At(iMod);
296   fSPDDead->AddAt(cal,iMod);
297 }
298 //_______________________________________________________________________
299 AliITSCalibration* AliITSDetTypeRec::GetCalibrationModel(Int_t iMod){
300   
301   //Get calibration model for module type
302   
303   if(fCalibration==0) {
304     Warning("GetalibrationModel","fCalibration is 0!");
305     return 0; 
306   }  
307
308   if(iMod<fgkDefaultNModulesSPD+fgkDefaultNModulesSDD){
309     return (AliITSCalibration*)fCalibration->At(iMod);
310   }else{
311     Int_t i=iMod-(fgkDefaultNModulesSPD+fgkDefaultNModulesSDD);
312     fSSDCalibration->SetModule(i);
313     return (AliITSCalibration*)fSSDCalibration;
314   }
315
316 }
317 //_______________________________________________________________________
318 AliITSCalibration* AliITSDetTypeRec::GetSPDDeadModel(Int_t iMod){
319   
320   //Get SPD dead for module iMod
321   
322   if(fSPDDead==0) {
323     AliWarning("fSPDDead is 0!");
324     return 0; 
325   }  
326
327   return (AliITSCalibration*)fSPDDead->At(iMod);
328 }
329
330 //______________________________________________________________________
331 void AliITSDetTypeRec::SetTreeAddressD(TTree *treeD){
332     // Set branch address for the tree of digits.
333
334     const char *det[4] = {"SPD","SDD","SSD","ITS"};
335     TBranch *branch;
336     Char_t* digclass;
337     Int_t i;
338     char branchname[30];
339
340     if(!treeD) return;
341     if (fDigits == 0x0) fDigits = new TObjArray(fgkNdettypes);
342     for (i=0; i<fgkNdettypes; i++) {
343         digclass = GetDigitClassName(i);
344         if(!(fDigits->At(i))) {
345             fDigits->AddAt(new TClonesArray(digclass,1000),i);
346         }else{
347             ResetDigits(i);
348         } 
349         if (fgkNdettypes==3) sprintf(branchname,"%sDigits%s",det[3],det[i]);
350         else  sprintf(branchname,"%sDigits%d",det[3],i+1);
351         if (fDigits) {
352             branch = treeD->GetBranch(branchname);
353             if (branch) branch->SetAddress(&((*fDigits)[i]));
354         } 
355     } 
356 }
357
358 //_______________________________________________________________________
359 TBranch* AliITSDetTypeRec::MakeBranchInTree(TTree *tree, const char* name, 
360                                        const char *classname, 
361                                        void* address,Int_t size, 
362                                        Int_t splitlevel)
363
364 //
365 // Makes branch in given tree and diverts them to a separate file
366 // 
367 //
368 //
369     
370   if (tree == 0x0) {
371     Error("MakeBranchInTree","Making Branch %s Tree is NULL",name);
372     return 0x0;
373   }
374   TBranch *branch = tree->GetBranch(name);
375   if (branch) {  
376     return branch;
377   }
378   if (classname){
379     branch = tree->Branch(name,classname,address,size,splitlevel);
380   }
381   else {
382     branch = tree->Bronch(name, "TClonesArray", address, size, splitlevel);
383   }
384   
385   return branch;
386 }
387
388 //____________________________________________________________________
389 void AliITSDetTypeRec::SetDefaults(){
390   
391   //Set defaults for segmentation and response
392
393   if(!GetITSgeom()){
394     Warning("SetDefaults","null pointer to AliITSgeomGeom !");
395     return;
396   }
397
398   AliITSsegmentation* seg;
399   if(!GetCalibration()) {AliFatal("Exit");exit(0);}  
400
401   for(Int_t dettype=0;dettype<fgkNdettypes;dettype++){
402     if(dettype==0){
403       seg = new AliITSsegmentationSPD();
404       SetSegmentationModel(dettype,seg);
405       SetDigitClassName(dettype,"AliITSdigitSPD");
406       SetClusterClassName(dettype,"AliITSRawClusterSPD");
407
408     }
409     if(dettype==1){
410       seg = new AliITSsegmentationSDD();
411       SetSegmentationModel(dettype,seg);
412       SetDigitClassName(dettype,"AliITSdigitSDD");
413       SetClusterClassName(dettype,"AliITSRawClusterSDD");
414     }
415     if(dettype==2){
416       AliITSsegmentationSSD* seg2 = new AliITSsegmentationSSD();
417       SetSegmentationModel(dettype,seg2);
418       SetDigitClassName(dettype,"AliITSdigitSSD");
419       SetClusterClassName(dettype,"AliITSRawClusterSSD");
420     }
421   }
422   
423 }
424 //______________________________________________________________________
425 Bool_t AliITSDetTypeRec::GetCalibration() {
426   // Get Default calibration if a storage is not defined.
427
428   if(!fFirstcall){
429     AliITSCalibration* cal = GetCalibrationModel(0);
430     if(cal)return kTRUE;
431   }
432   else {
433     fFirstcall = kFALSE;
434   }
435
436   //  SetRunNumber((Int_t)AliCDBManager::Instance()->GetRun());
437   //  Int_t run=GetRunNumber();
438
439   Bool_t cacheStatus = AliCDBManager::Instance()->GetCacheFlag();
440   if (fCalibration==0) {
441     fCalibration = new TObjArray(GetITSgeom()->GetIndexMax());
442     fCalibration->SetOwner(!cacheStatus);
443     fCalibration->Clear();
444   }
445
446   // dead pixel are not used for local reconstruction
447   AliCDBEntry *entrySPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDNoisy");
448   AliCDBEntry *deadSPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDDead");
449   AliCDBEntry *entrySDD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSDD");
450  
451  //  AliCDBEntry *entrySSD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSSD");
452   AliCDBEntry *entryNoiseSSD = AliCDBManager::Instance()->Get("ITS/Calib/NoiseSSD");
453   AliCDBEntry *entryPedestalSSD = AliCDBManager::Instance()->Get("ITS/Calib/PedestalSSD");
454   AliCDBEntry *entryGainSSD = AliCDBManager::Instance()->Get("ITS/Calib/GainSSD");
455   AliCDBEntry *entryBadChannelsSSD = AliCDBManager::Instance()->Get("ITS/Calib/BadChannelsSSD");
456   // Entry for the AliITSRecoParam object
457   AliCDBEntry *entryRP = AliCDBManager::Instance()->Get("ITS/Calib/RecoParam/");
458   AliCDBEntry *entry2SPD = AliCDBManager::Instance()->Get("ITS/Calib/RespSPD");
459   AliCDBEntry *entry2SDD = AliCDBManager::Instance()->Get("ITS/Calib/RespSDD");
460   AliCDBEntry *entry2SSD = AliCDBManager::Instance()->Get("ITS/Calib/RespSSD");
461   AliCDBEntry *drSpSDD = AliCDBManager::Instance()->Get("ITS/Calib/DriftSpeedSDD");
462   AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD");
463   AliCDBEntry *mapASDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsAnodeSDD");
464   AliCDBEntry *mapTSDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsTimeSDD");
465
466   if(!entrySPD || !deadSPD || !entrySDD || !entryNoiseSSD || !entryGainSSD || 
467      !entryPedestalSSD || !entryBadChannelsSSD || 
468      !entry2SPD || !entry2SDD || !entry2SSD || !drSpSDD || !ddlMapSDD || !mapASDD || !mapTSDD ){
469     AliFatal("Calibration object retrieval failed! ");
470     return kFALSE;
471   }     
472
473   TObjArray *calSPD = (TObjArray *)entrySPD->GetObject();
474   if(!cacheStatus)entrySPD->SetObject(NULL);
475   entrySPD->SetOwner(kTRUE);
476  
477   TObjArray *caldeadSPD = (TObjArray *)deadSPD->GetObject();
478   if(!cacheStatus)deadSPD->SetObject(NULL);
479   deadSPD->SetOwner(kTRUE);
480
481   AliITSresponseSPD *pSPD = (AliITSresponseSPD*)entry2SPD->GetObject();
482   if(!cacheStatus)entry2SPD->SetObject(NULL);
483   entry2SPD->SetOwner(kTRUE);
484     
485   TObjArray *calSDD = (TObjArray *)entrySDD->GetObject();
486   if(!cacheStatus)entrySDD->SetObject(NULL);
487   entrySDD->SetOwner(kTRUE);
488  
489   AliITSresponseSDD *pSDD = (AliITSresponseSDD*)entry2SDD->GetObject();
490   if(!cacheStatus)entry2SDD->SetObject(NULL);
491   entry2SDD->SetOwner(kTRUE);
492
493   TObjArray *drSp = (TObjArray *)drSpSDD->GetObject();
494   if(!cacheStatus)drSpSDD->SetObject(NULL);
495   drSpSDD->SetOwner(kTRUE);
496
497   AliITSDDLModuleMapSDD *ddlsdd=(AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject();
498   if(!cacheStatus)ddlMapSDD->SetObject(NULL);
499   ddlMapSDD->SetOwner(kTRUE);
500
501   TObjArray *mapAn = (TObjArray *)mapASDD->GetObject();
502   if(!cacheStatus)mapASDD->SetObject(NULL);
503   mapASDD->SetOwner(kTRUE);
504
505   TObjArray *mapT = (TObjArray *)mapTSDD->GetObject();
506   if(!cacheStatus)mapTSDD->SetObject(NULL);
507   mapTSDD->SetOwner(kTRUE);
508
509   AliITSNoiseSSD *noiseSSD = (AliITSNoiseSSD *)entryNoiseSSD->GetObject();
510   if(!cacheStatus)entryNoiseSSD->SetObject(NULL);
511   entryNoiseSSD->SetOwner(kTRUE);
512
513   AliITSPedestalSSD *pedestalSSD = 
514     (AliITSPedestalSSD*)entryPedestalSSD->GetObject();
515   if(!cacheStatus)entryPedestalSSD->SetObject(NULL);
516   entryPedestalSSD->SetOwner(kTRUE);
517
518   AliITSGainSSD *gainSSD = (AliITSGainSSD *)entryGainSSD->GetObject();
519   if(!cacheStatus)entryGainSSD->SetObject(NULL);
520   entryGainSSD->SetOwner(kTRUE);
521
522   AliITSBadChannelsSSD *badchannelsSSD = 
523     (AliITSBadChannelsSSD*)entryBadChannelsSSD->GetObject();
524   if(!cacheStatus)entryBadChannelsSSD->SetObject(NULL);
525   entryBadChannelsSSD->SetOwner(kTRUE);
526
527   AliITSresponseSSD *pSSD = (AliITSresponseSSD*)entry2SSD->GetObject();
528   if(!cacheStatus)entry2SSD->SetObject(NULL);
529   entry2SSD->SetOwner(kTRUE);
530
531   if(!AliITSReconstructor::GetRecoParam()){
532     if(!entryRP) AliFatal("Calibration object (RecoParam) retrieval from OCDB failed! Hint: as an alternative you can set it in your reconstruction macro ");
533     AliITSRecoParam *rp = (AliITSRecoParam*)entryRP->GetObject();
534     if(!cacheStatus)entryRP->SetObject(NULL);
535     entryRP->SetOwner(kTRUE);
536     AliITSReconstructor::SetRecoParam(rp);
537   }
538   else {
539     AliWarning("AliITSRecoParam object has been already set in AliITSReconstructor. The OCDB instance will not be used\n");
540   }
541
542
543   // DB entries are deleted. In this way metadeta objects are deleted as well
544   if(!cacheStatus){
545     delete entrySPD;
546     delete deadSPD;
547     delete entrySDD;
548     delete entryNoiseSSD;
549     delete entryPedestalSSD;
550     delete entryGainSSD;
551     delete entryBadChannelsSSD;
552     delete entry2SPD;
553     delete entry2SDD;
554     delete entry2SSD;
555     delete mapASDD;
556     delete mapTSDD;
557     delete drSpSDD;
558     delete ddlMapSDD;
559   }
560
561   if ((!pSPD)||(!pSDD)||(!pSSD) || (!calSPD) || (!caldeadSPD) ||(!calSDD) || (!drSp) || (!ddlsdd)
562       || (!mapAn) || (!mapT) || (!noiseSSD)|| (!gainSSD)|| (!badchannelsSSD)) {
563     AliWarning("Can not get calibration from calibration database !");
564     return kFALSE;
565   }
566
567   fNMod[0] = calSPD->GetEntries();
568   fNMod[1] = calSDD->GetEntries();
569   //fNMod[2] = noiseSSD->GetEntries();
570   AliInfo(Form("%i SPD, %i SDD and %i SSD in calibration database",
571                fNMod[0], fNMod[1], fNMod[2]));
572   AliITSCalibration* cal;
573   for (Int_t i=0; i<fNMod[0]; i++) {
574     cal = (AliITSCalibration*) calSPD->At(i);
575     cal->SetResponse((AliITSresponse*)pSPD);
576     SetCalibrationModel(i, cal);
577     cal = (AliITSCalibration*) caldeadSPD->At(i);
578     SetSPDDeadModel(i, cal);
579   }
580
581   fDDLMapSDD=ddlsdd;
582   fRespSDD=pSDD;
583   for(Int_t iddl=0; iddl<AliITSDDLModuleMapSDD::GetNDDLs(); iddl++){
584     for(Int_t icar=0; icar<AliITSDDLModuleMapSDD::GetNModPerDDL();icar++){
585       Int_t iMod=fDDLMapSDD->GetModuleNumber(iddl,icar);
586       if(iMod==-1) continue;
587       Int_t i=iMod - fgkDefaultNModulesSPD;
588       cal = (AliITSCalibration*) calSDD->At(i);
589       Int_t i0=2*i;
590       Int_t i1=1+2*i;
591       AliITSDriftSpeedArraySDD* arr0 = (AliITSDriftSpeedArraySDD*) drSp->At(i0);
592       AliITSMapSDD* ma0 = (AliITSMapSDD*)mapAn->At(i0);
593       AliITSMapSDD* mt0 = (AliITSMapSDD*)mapT->At(i0);
594       AliITSDriftSpeedArraySDD* arr1 = (AliITSDriftSpeedArraySDD*) drSp->At(i1);
595       AliITSMapSDD* ma1 = (AliITSMapSDD*)mapAn->At(i1);
596       AliITSMapSDD* mt1 = (AliITSMapSDD*)mapT->At(i1);
597       cal->SetDriftSpeed(0,arr0);
598       cal->SetDriftSpeed(1,arr1);
599       cal->SetMapA(0,ma0);
600       cal->SetMapA(1,ma1);
601       cal->SetMapT(0,mt0);
602       cal->SetMapT(1,mt1);
603       SetCalibrationModel(iMod, cal);
604     }
605   }
606
607   fSSDCalibration->SetResponse((AliITSresponse*)pSSD);
608   fSSDCalibration->SetNoise(noiseSSD);
609   fSSDCalibration->SetPedestal(pedestalSSD);
610   fSSDCalibration->SetGain(gainSSD);
611   fSSDCalibration->SetBadChannels(badchannelsSSD);
612   //fSSDCalibration->FillBadChipMap();
613
614
615   /*
616   for (Int_t i=0; i<fNMod[2]; i++) {
617
618     AliITSCalibrationSSD *calibSSD = new AliITSCalibrationSSD();
619     calibSSD->SetResponse((AliITSresponse*)pSSD);
620     
621     AliITSNoiseSSD *noise = (AliITSNoiseSSD*) (noiseSSD->At(i));
622     calibSSD->SetNoise(noise);
623     AliITSPedestalSSD *pedestal = (AliITSPedestalSSD*) (pedestalSSD->At(i));
624     calibSSD->SetPedestal(pedestal);
625     AliITSGainSSD *gain = (AliITSGainSSD*) (gainSSD->At(i));
626     calibSSD->SetGain(gain);
627     AliITSBadChannelsSSD *bad = (AliITSBadChannelsSSD*) (badchannelsSSD->At(i));
628     calibSSD->SetBadChannels(bad);
629     calibSSD->FillBadChipMap();
630
631     Int_t iMod = i + fgkDefaultNModulesSPD + fgkDefaultNModulesSDD;
632     SetCalibrationModel(iMod, calibSSD);
633  }
634   */
635
636   return kTRUE;
637 }
638
639
640 //________________________________________________________________
641 void AliITSDetTypeRec::SetDefaultClusterFinders(){
642   
643   //set defaults for standard cluster finder
644
645   if(!GetITSgeom()){
646     Warning("SetDefaults","null pointer to AliITSgeom!");
647     return;
648   }
649
650   AliITSClusterFinder *clf; 
651
652   for(Int_t dettype=0;dettype<fgkNdettypes;dettype++){
653     //SPD
654     if(dettype==0){
655       if(!GetReconstructionModel(dettype)){
656         TClonesArray *dig0 = DigitsAddress(0);
657         TClonesArray *rec0 = ClustersAddress(0);
658         clf = new AliITSClusterFinderSPD(this,dig0,rec0);
659         SetReconstructionModel(dettype,clf);
660
661       }
662     }
663    
664     //SDD
665     if(dettype==1){
666       if(!GetReconstructionModel(dettype)){
667         TClonesArray *dig1 = DigitsAddress(1);
668         TClonesArray *rec1 = ClustersAddress(1);
669         clf = new AliITSClusterFinderSDD(this,dig1,rec1);
670         SetReconstructionModel(dettype,clf);
671       }
672
673     }
674     //SSD
675     if(dettype==2){
676       if(!GetReconstructionModel(dettype)){
677         TClonesArray* dig2 = DigitsAddress(2);
678         clf = new AliITSClusterFinderSSD(this,dig2);
679         SetReconstructionModel(dettype,clf);
680       }
681     }
682
683  }
684  
685   
686 }
687
688 //________________________________________________________________
689 void AliITSDetTypeRec::SetDefaultClusterFindersV2(Bool_t rawdata){
690
691   //Set defaults for cluster finder V2
692
693   if(!GetITSgeom()){
694     Warning("SetDefaults","Null pointer to AliITSgeom !");
695     return;
696   }
697
698   AliITSClusterFinder *clf; 
699
700   for(Int_t dettype=0;dettype<fgkNdettypes;dettype++){
701     //SPD
702     if(dettype==0){
703       if(!GetReconstructionModel(dettype)){
704         clf = new AliITSClusterFinderV2SPD(this);
705         clf->InitGeometry();
706         if(!rawdata) clf->SetDigits(DigitsAddress(0));
707         SetReconstructionModel(dettype,clf);
708
709       }
710     }
711     //SDD
712     if(dettype==1){
713       if(!GetReconstructionModel(dettype)){
714         clf = new AliITSClusterFinderV2SDD(this);
715         clf->InitGeometry();
716         if(!rawdata) clf->SetDigits(DigitsAddress(1));
717         SetReconstructionModel(dettype,clf);
718       }
719
720     }
721
722     //SSD
723     if(dettype==2){
724       if(!GetReconstructionModel(dettype)){
725         clf = new AliITSClusterFinderV2SSD(this);
726         clf->InitGeometry();
727         if(!rawdata) clf->SetDigits(DigitsAddress(2));
728         SetReconstructionModel(dettype,clf);
729       }
730     }
731
732  }
733    
734 }
735 //______________________________________________________________________
736 void AliITSDetTypeRec::MakeBranch(TTree* tree, Option_t* option){
737
738   //Creates branches for clusters and recpoints
739   Bool_t cR = (strstr(option,"R")!=0);
740   Bool_t cRF = (strstr(option,"RF")!=0);
741   
742   if(cRF)cR = kFALSE;
743
744   if(cR) MakeBranchR(tree);
745   if(cRF) MakeBranchRF(tree);
746
747 }
748
749 //___________________________________________________________________
750 void AliITSDetTypeRec::AddCluster(Int_t id, AliITSRawCluster *c){
751
752   // Adds a raw cluster to the list
753   TClonesArray &lc = *((TClonesArray*)fCtype->At(id));  
754   switch(id){
755   case 0:
756     new(lc[fNctype[id]++]) AliITSRawClusterSPD(*((AliITSRawClusterSPD*)c));
757     break;
758   case 1:
759     new(lc[fNctype[id]++]) AliITSRawClusterSDD(*((AliITSRawClusterSDD*)c));
760     break;
761   case 2:
762     new(lc[fNctype[id]++]) AliITSRawClusterSSD(*((AliITSRawClusterSSD*)c));
763     break;
764   } 
765 }
766 //___________________________________________________________________
767 void AliITSDetTypeRec::ResetDigits(){
768   // Reset number of digits and the digits array for the ITS detector.
769   
770   if(!fDigits) return;
771   for(Int_t i=0;i<fgkNdettypes;i++){
772     ResetDigits(i);
773   }
774 }
775 //___________________________________________________________________
776 void AliITSDetTypeRec::ResetDigits(Int_t branch){
777   // Reset number of digits and the digits array for this branch.
778   
779   if(fDigits->At(branch)) ((TClonesArray*)fDigits->At(branch))->Clear();
780   if(fNdtype) fNdtype[branch]=0;
781
782 }
783
784 //__________________________________________________________________
785 void AliITSDetTypeRec::ResetClusters(){
786
787   //Resets number of clusters and the cluster array 
788   for(Int_t i=0;i<fgkNdettypes;i++){
789     ResetClusters(i);
790   }
791 }
792
793 //__________________________________________________________________
794 void AliITSDetTypeRec::ResetClusters(Int_t i){
795
796   //Resets number of clusters and the cluster array for this branch
797
798   if (fCtype->At(i))    ((TClonesArray*)fCtype->At(i))->Clear();
799   if (fNctype)  fNctype[i]=0;
800 }
801 //__________________________________________________________________
802 void AliITSDetTypeRec::MakeBranchR(TTree *treeR, Option_t *opt){
803
804   //Creates tree branches for recpoints
805   // Inputs:
806   //      cont char *file  File name where RecPoints branch is to be written
807   //                       to. If blank it write the SDigits to the same
808   //                       file in which the Hits were found.
809
810   Int_t buffsz = 4000;
811   char branchname[30];
812
813   // only one branch for rec points for all detector types
814   Bool_t oFast= (strstr(opt,"Fast")!=0);
815   
816   Char_t detname[10] = "ITS";
817  
818   
819   if(oFast){
820     sprintf(branchname,"%sRecPointsF",detname);
821   } else {
822     sprintf(branchname,"%sRecPoints",detname);
823   }
824   
825   if(!fRecPoints)fRecPoints = new TClonesArray("AliITSRecPoint",1000);
826   if (treeR)
827     MakeBranchInTree(treeR,branchname,0,&fRecPoints,buffsz,99);
828 }
829 //______________________________________________________________________
830 void AliITSDetTypeRec::SetTreeAddressR(TTree *treeR){
831     // Set branch address for the Reconstructed points Trees.
832     // Inputs:
833     //      TTree *treeR   Tree containing the RecPoints.
834     // Outputs:
835     //      none.
836     // Return:
837
838     char branchname[30];
839     Char_t namedet[10]="ITS";
840
841     if(!treeR) return;
842     if(fRecPoints==0x0) fRecPoints = new TClonesArray("AliITSRecPoint",1000);
843     TBranch *branch;
844     sprintf(branchname,"%sRecPoints",namedet);
845     branch = treeR->GetBranch(branchname);
846     if (branch) {
847       branch->SetAddress(&fRecPoints);
848     }else {
849       sprintf(branchname,"%sRecPointsF",namedet);
850       branch = treeR->GetBranch(branchname);
851       if (branch) {
852         branch->SetAddress(&fRecPoints);
853       }
854
855     }
856 }
857 //____________________________________________________________________
858 void AliITSDetTypeRec::AddRecPoint(const AliITSRecPoint &r){
859     // Add a reconstructed space point to the list
860     // Inputs:
861     //      const AliITSRecPoint &r RecPoint class to be added to the tree
862     //                              of reconstructed points TreeR.
863     // Outputs:
864     //      none.
865     // Return:
866     //      none.
867
868     TClonesArray &lrecp = *fRecPoints;
869     new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
870 }
871
872 //______________________________________________________________________
873 void AliITSDetTypeRec::DigitsToRecPoints(TTree *treeD,TTree *treeR,Int_t lastentry,Option_t *opt, Bool_t v2){
874   // cluster finding and reconstruction of space points
875   // the condition below will disappear when the geom class will be
876   // initialized for all versions - for the moment it is only for v5 !
877   // 7 is the SDD beam test version
878   // Inputs:
879   //      TTree *treeD     Digits tree
880   //      TTree *treeR     Clusters tree
881   //      Int_t lastentry  Offset for module when not all of the modules
882   //                       are processed.
883   //      Option_t *opt    String indicating which ITS sub-detectors should
884   //                       be processed. If ="All" then all of the ITS
885   //                       sub detectors are processed.
886
887   const char *all = strstr(opt,"All");
888   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
889                         strstr(opt,"SSD")};
890   if(!v2) {
891     SetDefaultClusterFinders();
892     AliInfo("Original cluster finder has been selected\n");
893   }
894   else   { 
895     SetDefaultClusterFindersV2();
896     AliInfo("V2 cluster finder has been selected \n");
897   }
898
899   AliITSClusterFinder *rec     = 0;
900   Int_t id,module,first=0;
901   for(module=0;module<GetITSgeom()->GetIndexMax();module++){
902       id       = GetITSgeom()->GetModuleType(module);
903       if (!all && !det[id]) continue;
904       if(det[id]) first = GetITSgeom()->GetStartDet(id);
905       rec = (AliITSClusterFinder*)GetReconstructionModel(id);
906       TClonesArray *itsDigits  = DigitsAddress(id);
907       if (!rec)
908           AliFatal("The reconstruction class was not instanciated!");
909       ResetDigits();  // MvL: Not sure we neeed this when rereading anyways
910       if (all) {
911           treeD->GetEvent(lastentry+module);
912       }else {
913           treeD->GetEvent(lastentry+(module-first));
914       }
915       Int_t ndigits = itsDigits->GetEntriesFast();
916       if(ndigits>0){
917         rec->SetDetTypeRec(this);
918         rec->SetDigits(DigitsAddress(id));
919         //      rec->SetClusters(ClustersAddress(id));
920         rec->FindRawClusters(module);
921       } // end if
922       treeR->Fill();
923       ResetRecPoints();
924       ResetClusters();
925   } 
926 }
927 //______________________________________________________________________
928 void AliITSDetTypeRec::DigitsToRecPoints(AliRawReader* rawReader,TTree *treeR,Option_t *opt){
929   // cluster finding and reconstruction of space points
930   // the condition below will disappear when the geom class will be
931   // initialized for all versions - for the moment it is only for v5 !
932   // 7 is the SDD beam test version
933   // Inputs:
934   //      AliRawReader *rawReader  Pointer to the raw-data reader
935   //      TTree *treeR             Clusters tree
936   // Outputs:
937   //      none.
938   // Return:
939   //      none.
940   const char *all = strstr(opt,"All");
941   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
942                         strstr(opt,"SSD")};
943   AliITSClusterFinderV2 *rec     = 0;
944   Int_t id=0;
945
946   TClonesArray *array=new TClonesArray("AliITSRecPoint",1000);
947   TBranch *branch = treeR->Branch("ITSRecPoints",&array);
948   delete array;
949  
950   TClonesArray** clusters = new TClonesArray*[GetITSgeom()->GetIndexMax()]; 
951   for (Int_t iModule = 0; iModule < GetITSgeom()->GetIndexMax(); iModule++) {
952     clusters[iModule] = NULL;
953   }
954   for(id=0;id<3;id++){
955     if (!all && !det[id]) continue;
956     rec = (AliITSClusterFinderV2*)GetReconstructionModel(id);
957     if (!rec)
958       AliFatal("The reconstruction class was not instanciated");
959     rec->SetDetTypeRec(this);
960     rec->RawdataToClusters(rawReader,clusters);    
961   } 
962   Int_t nClusters =0;
963   TClonesArray *emptyArray=new TClonesArray("AliITSRecPoint");
964   for(Int_t iModule=0;iModule<GetITSgeom()->GetIndexMax();iModule++){
965     id = GetITSgeom()->GetModuleType(iModule);
966     if (!all && !det[id]) continue;
967     array = clusters[iModule];
968     if(!array){
969       AliDebug(1,Form("data for module %d missing!",iModule));
970       array = emptyArray;
971     }
972     branch->SetAddress(&array);
973     treeR->Fill();
974     nClusters+=array->GetEntriesFast();
975
976     if (array != emptyArray) {
977       array->Delete();
978       delete array;
979     }
980   }
981   delete emptyArray;
982
983   delete[] clusters;
984   Info("DigitsToRecPoints", "total number of found recpoints in ITS: %d\n", 
985        nClusters);
986   
987 }
988
989