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