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