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