]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSDetTypeSim.cxx
enhanced HOMER reader to work on normal buffer
[u/mrichter/AliRoot.git] / ITS / AliITSDetTypeSim.cxx
1 /***************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors 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 /////////////////////////////////////////////////////////////////////
22 // Base simulation functions for ITS                               //
23 //                                                                 //
24 //                                                                 //
25 /////////////////////////////////////////////////////////////////////          
26 #include "TBranch.h"
27 #include "TClonesArray.h"
28 #include "TObjArray.h"
29 #include "TTree.h"
30
31 #include "AliRun.h"
32
33 #include "AliCDBManager.h"
34 #include "AliCDBId.h"
35 #include "AliCDBStorage.h"
36 #include "AliCDBEntry.h"
37 #include "AliCDBMetaData.h"
38
39 #include "AliITSdigit.h"
40 #include "AliITSdigitSPD.h"
41 #include "AliITSdigitSDD.h"
42 #include "AliITSdigitSSD.h"
43 #include "AliITSDetTypeSim.h"
44 #include "AliITSpListItem.h"
45 #include "AliITSresponseSDD.h"
46 #include "AliITSCalibrationSDD.h"
47 #include "AliITSMapSDD.h"
48 #include "AliITSCalibrationSSD.h"
49 #include "AliITSNoiseSSD.h"
50 #include "AliITSGainSSD.h"
51 #include "AliITSBadChannelsSSD.h"
52 #include "AliITSCalibrationSSD.h"
53 #include "AliITSsegmentationSPD.h"
54 #include "AliITSsegmentationSDD.h"
55 #include "AliITSsegmentationSSD.h"
56 #include "AliITSsimulation.h"
57 #include "AliITSsimulationSPD.h"
58 #include "AliITSsimulationSDD.h"
59 #include "AliITSsimulationSSD.h"
60
61
62 const Int_t AliITSDetTypeSim::fgkNdettypes = 3;
63 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSPD =  240;
64 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSDD =  260;
65 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSSD = 1698;
66
67 ClassImp(AliITSDetTypeSim)
68
69 //----------------------------------------------------------------------
70 AliITSDetTypeSim::AliITSDetTypeSim():
71 TObject(),
72 fSimulation(),   // [NDet]
73 fSegmentation(), // [NDet]
74 fCalibration(),     // [NMod]
75 fPreProcess(),   // [] e.g. Fill fHitModule with hits
76 fPostProcess(),  // [] e.g. Wright Raw data
77 fNSDigits(0),    //! number of SDigits
78 fSDigits(),      //! [NMod][NSDigits]
79 fNDigits(0),     //! number of Digits
80 fRunNumber(0),   //! Run number (to access DB)
81 fDigits(),       //! [NMod][NDigits]
82 fHitClassName(), // String with Hit class name.
83 fSDigClassName(),// String with SDigit class name.
84 fDigClassName(), // String with digit class name.
85 fLoader(0),      // local pointer to loader
86 fFirstcall(kTRUE){ // flag
87     // Default Constructor
88     // Inputs:
89     //    none.
90     // Outputs:
91     //    none.
92     // Return:
93     //    A properly zero-ed AliITSDetTypeSim class.
94
95   fSimulation = new TObjArray(fgkNdettypes);
96   fSegmentation = new TObjArray(fgkNdettypes);
97   fSegmentation->SetOwner(kTRUE);
98   fSDigits = new TClonesArray("AliITSpListItem",1000);
99   fDigits = new TObjArray(fgkNdettypes);
100   fNDigits = new Int_t[fgkNdettypes];
101   fNMod[0] = fgkDefaultNModulesSPD;
102   fNMod[1] = fgkDefaultNModulesSDD;
103   fNMod[2] = fgkDefaultNModulesSSD;
104   SetRunNumber();
105 }
106 //----------------------------------------------------------------------
107 AliITSDetTypeSim::~AliITSDetTypeSim(){
108     // Destructor
109     // Inputs:
110     //    none.
111     // Outputs:
112     //    none.
113     // Return:
114     //    Nothing.
115
116     if(fSimulation){
117         fSimulation->Delete();
118         delete fSimulation;
119     }
120     fSimulation = 0;
121     if(fSegmentation){
122         fSegmentation->Delete();
123         delete fSegmentation;
124     }
125     fSegmentation = 0;
126     if(fCalibration && fRunNumber<0){
127         AliITSresponse* rspd = ((AliITSCalibration*)fCalibration->At(
128                                GetITSgeom()->GetStartSPD()))->GetResponse();
129         AliITSresponse* rsdd = ((AliITSCalibration*)fCalibration->At(
130                                GetITSgeom()->GetStartSDD()))->GetResponse();
131         AliITSresponse* rssd = ((AliITSCalibration*)fCalibration->At(
132                                GetITSgeom()->GetStartSSD()))->GetResponse();
133         if(rspd) delete rspd;
134         if(rsdd) delete rsdd;
135         if(rssd) delete rssd;
136         fCalibration->Delete();
137         delete fCalibration;
138     }
139     fCalibration = 0;
140     if(fPreProcess){
141         fPreProcess->Delete();
142         delete fPreProcess;
143     }
144     fPreProcess = 0;
145     if(fPostProcess){
146         fPostProcess->Delete();
147         delete fPostProcess;
148     }
149     fPostProcess = 0;
150     if(fNDigits) delete [] fNDigits;
151     fNDigits = 0;
152     if (fLoader)fLoader->GetModulesFolder()->Remove(this);
153     fLoader = 0; // Not deleting it.
154     if (fSDigits) {
155         fSDigits->Delete();
156         delete fSDigits;
157     }
158     fSDigits=0;
159     if (fDigits) {
160       fDigits->Delete();
161       delete fDigits;
162     }
163     fDigits=0;
164 }
165 //----------------------------------------------------------------------
166 AliITSDetTypeSim::AliITSDetTypeSim(const AliITSDetTypeSim &source) : TObject(source),
167 fSimulation(source.fSimulation),   // [NDet]
168 fSegmentation(source.fSegmentation), // [NDet]
169 fCalibration(source.fCalibration),     // [NMod]
170 fPreProcess(source.fPreProcess),   // [] e.g. Fill fHitModule with hits
171 fPostProcess(source.fPostProcess),  // [] e.g. Wright Raw data
172 fNSDigits(source.fNSDigits),    //! number of SDigits
173 fSDigits(source.fSDigits),      //! [NMod][NSDigits]
174 fNDigits(source.fNDigits),     //! number of Digits
175 fRunNumber(source.fRunNumber),   //! Run number (to access DB)
176 fDigits(source.fDigits),       //! [NMod][NDigits]
177 fHitClassName(source.fHitClassName), // String with Hit class name.
178 fSDigClassName(source.fSDigClassName),// String with SDigit class name.
179 fDigClassName(), // String with digit class name.
180 fLoader(source.fLoader),      // local pointer to loader
181 fFirstcall(source.fFirstcall)
182 {
183     // Copy Constructor for object AliITSDetTypeSim not allowed
184   for(Int_t i=0;i<fgkNdettypes;i++){
185     fDigClassName[i] = source.fDigClassName[i];
186   }
187 }
188 //----------------------------------------------------------------------
189 AliITSDetTypeSim& AliITSDetTypeSim::operator=(const AliITSDetTypeSim &source){
190     // The = operator for object AliITSDetTypeSim
191  
192   this->~AliITSDetTypeSim();
193   new(this) AliITSDetTypeSim(source);
194   return *this;
195 }
196
197 //______________________________________________________________________
198 void AliITSDetTypeSim::SetITSgeom(AliITSgeom *geom){
199     // Sets/replaces the existing AliITSgeom object kept in AliITSLoader
200     // 
201     // Inputs:
202     //   AliITSgoem   *geom  The AliITSgeom object to be used.
203     // Output:
204     //   none.
205     // Return:
206     //   none.
207   if(!fLoader){
208     Error("SetITSgeom","No pointer to loader - nothing done");
209     return;
210   }
211   else {
212     fLoader->SetITSgeom(geom);  // protections in AliITSLoader::SetITSgeom
213   }
214  
215 }
216 //______________________________________________________________________
217 void AliITSDetTypeSim::SetLoader(AliITSLoader *loader){
218     // Sets the local copy of the AliITSLoader, and passes on the
219     // AliITSgeom object as needed.
220     // Inputs
221     //   AliITSLoader  *loader pointer to AliITSLoader for local use
222     // Outputs:
223     //   none.
224     // Return:
225     //  none.
226
227     if(fLoader==loader) return; // Same do nothing
228     if(fLoader){ // alread have an existing loader
229         Error("SetLoader",
230                 "Already have an exisiting loader ptr=%p Nothing done",
231                 fLoader);
232     } // end if
233     fLoader = loader;
234 }
235 //______________________________________________________________________
236 void AliITSDetTypeSim::SetSimulationModel(Int_t dettype,AliITSsimulation *sim){
237
238   //Set simulation model for detector type
239
240   if(fSimulation==0) fSimulation = new TObjArray(fgkNdettypes);
241   fSimulation->AddAt(sim,dettype);
242 }
243 //______________________________________________________________________
244 AliITSsimulation* AliITSDetTypeSim::GetSimulationModel(Int_t dettype){
245
246   //Get simulation model for detector type
247   if(fSimulation==0)  {
248     Warning("GetSimulationModel","fSimulation is 0!");
249     return 0;     
250   }
251   return (AliITSsimulation*)(fSimulation->At(dettype));
252 }
253 //______________________________________________________________________
254 AliITSsimulation* AliITSDetTypeSim::GetSimulationModelByModule(Int_t module){
255
256   //Get simulation model by module number
257   if(GetITSgeom()==0) {
258     Warning("GetSimulationModelByModule","GetITSgeom() is 0!");
259     return 0;
260   }
261   
262   return GetSimulationModel(GetITSgeom()->GetModuleType(module));
263 }
264 //_______________________________________________________________________
265 void AliITSDetTypeSim::SetDefaultSegmentation(Int_t idet){
266     // Set default segmentation model objects
267     AliITSsegmentation *seg;
268
269     if(fSegmentation==0x0){
270         fSegmentation = new TObjArray(fgkNdettypes);
271         fSegmentation->SetOwner(kTRUE);
272     }
273     if(GetSegmentationModel(idet))
274         delete (AliITSsegmentation*)fSegmentation->At(idet);
275     if(idet==0){
276         seg = new AliITSsegmentationSPD(GetITSgeom());
277     }else if(idet==1){
278         seg = new AliITSsegmentationSDD(GetITSgeom());
279     }else {
280         seg = new AliITSsegmentationSSD(GetITSgeom());
281     }
282     SetSegmentationModel(idet,seg);
283 }
284 //______________________________________________________________________
285 void AliITSDetTypeSim::SetSegmentationModel(Int_t dettype,
286                                             AliITSsegmentation *seg){
287    
288   //Set segmentation model for detector type
289   if(fSegmentation==0x0){
290     fSegmentation = new TObjArray(fgkNdettypes);
291     fSegmentation->SetOwner(kTRUE);
292   }
293   fSegmentation->AddAt(seg,dettype);
294 }
295 //______________________________________________________________________
296 AliITSsegmentation* AliITSDetTypeSim::GetSegmentationModel(Int_t dettype){
297   //Get segmentation model for detector type
298    
299    if(fSegmentation==0) {
300        Warning("GetSegmentationModel","fSegmentation is 0!");
301        return 0; 
302    } 
303    return (AliITSsegmentation*)(fSegmentation->At(dettype));
304 }
305 //_______________________________________________________________________
306 AliITSsegmentation* AliITSDetTypeSim::GetSegmentationModelByModule(Int_t module){
307     //Get segmentation model by module number
308     if(GetITSgeom()==0){
309         Warning("GetSegmentationModelByModule","GetITSgeom() is 0!");
310         return 0;
311     }     
312     return GetSegmentationModel(GetITSgeom()->GetModuleType(module));
313 }
314 //_______________________________________________________________________
315 void AliITSDetTypeSim::CreateCalibrationArray() {
316     //Create the container of calibration functions with correct size
317     if (fCalibration) {
318         Warning("CreateCalibration","pointer to calibration object exists\n");
319         fCalibration->Delete();
320         delete fCalibration;
321     }
322
323     Int_t nModTot = GetITSgeom()->GetIndexMax();
324     fCalibration = new TObjArray(nModTot);
325     fCalibration->SetOwner(kTRUE);
326     fCalibration->Clear();
327 }
328 //_______________________________________________________________________
329 void AliITSDetTypeSim::SetCalibrationModel(Int_t iMod, AliITSCalibration *resp){
330     //Set response model for modules
331
332     if (fCalibration==0) CreateCalibrationArray();
333  
334     if (fCalibration->At(iMod)!=0)
335         delete (AliITSCalibration*) fCalibration->At(iMod);
336     fCalibration->AddAt(resp, iMod);
337 }
338 //______________________________________________________________________
339 void AliITSDetTypeSim::ResetCalibrationArray(){
340     //resets response array
341     if(fCalibration && fRunNumber<0){  // if fRunNumber<0 fCalibration is owner
342         AliITSresponse* rspd = ((AliITSCalibration*)fCalibration->At(
343                                 GetITSgeom()->GetStartSPD()))->GetResponse();
344         AliITSresponse* rsdd = ((AliITSCalibration*)fCalibration->At(
345                                 GetITSgeom()->GetStartSDD()))->GetResponse();
346         AliITSresponse* rssd = ((AliITSCalibration*)fCalibration->At(
347                                 GetITSgeom()->GetStartSSD()))->GetResponse();
348         if(rspd) delete rspd;
349         if(rsdd) delete rsdd;
350         if(rssd) delete rssd;
351         fCalibration->Clear();
352     }else if (fCalibration && fRunNumber>=0){
353         fCalibration->Clear();
354     }
355 }
356 //______________________________________________________________________
357 void AliITSDetTypeSim::ResetSegmentation(){
358     //Resets segmentation array
359     if(fSegmentation) fSegmentation->Clear();
360 }
361 //_______________________________________________________________________
362 AliITSCalibration* AliITSDetTypeSim::GetCalibrationModel(Int_t iMod){
363     //Get response model for module number iMod 
364  
365     if(fCalibration==0) {
366         AliError("fCalibration is 0!");
367         return 0; 
368     }
369   return (AliITSCalibration*)(fCalibration->At(iMod));
370 }
371 //_______________________________________________________________________
372 void AliITSDetTypeSim::SetDefaults(){
373     //Set defaults for segmentation and response
374
375     if(GetITSgeom()==0){
376         Warning("SetDefaults","GetITSgeom() is 0!");
377         return;
378     } // end if
379     if (fCalibration==0) {
380         CreateCalibrationArray();
381     } // end if
382
383     ResetSegmentation();
384     if(!GetCalibration()){AliFatal("Exit"); exit(0);}
385
386     for(Int_t idet=0;idet<fgkNdettypes;idet++){
387         //SPD
388         if(idet==0){
389             if(!GetSegmentationModel(idet)) SetDefaultSegmentation(idet);
390             const char *kData0=(GetCalibrationModel(GetITSgeom()->GetStartSPD()))->DataType();
391             if (strstr(kData0,"real")) {
392                 SetDigitClassName(idet,"AliITSdigit");
393             }else {
394                 SetDigitClassName(idet,"AliITSdigitSPD");
395             } // end if
396         } // end if idet==0
397         //SDD
398         if(idet==1){
399             if(!GetSegmentationModel(idet)) SetDefaultSegmentation(idet);
400             AliITSCalibrationSDD* rsp = 
401                 (AliITSCalibrationSDD*)GetCalibrationModel(
402                     GetITSgeom()->GetStartSDD());
403             const char *kopt = ((AliITSresponseSDD*)rsp->GetResponse())->
404                 ZeroSuppOption();
405             if((!strstr(kopt,"2D"))&&(!strstr(kopt,"1D"))) {
406                 SetDigitClassName(idet,"AliITSdigit");
407             }else {
408                 SetDigitClassName(idet,"AliITSdigitSDD");
409             } // end if
410         } // end if idet==1
411         //SSD
412         if(idet==2){
413             if(!GetSegmentationModel(idet))SetDefaultSegmentation(idet);
414             const char *kData2 = (GetCalibrationModel(
415                                   GetITSgeom()->GetStartSSD())->DataType());
416             if (strstr(kData2,"real")) {
417                 SetDigitClassName(idet,"AliITSdigit");
418             }else {
419                 SetDigitClassName(idet,"AliITSdigitSSD");
420             } // end if
421         } // end if idet==2
422     }// end for idet
423 }
424 //______________________________________________________________________
425 Bool_t AliITSDetTypeSim::GetCalibration() {
426   // Get Default calibration if a storage is not defined.
427
428   if(!fFirstcall){
429     AliITSCalibration* cal = GetCalibrationModel(0);
430     if(cal)return kTRUE;
431   }
432   else {
433     fFirstcall = kFALSE;
434   }
435
436   SetRunNumber((Int_t)AliCDBManager::Instance()->GetRun());
437   Int_t run=GetRunNumber();
438
439   Bool_t origCacheStatus = AliCDBManager::Instance()->GetCacheFlag();
440   Bool_t isCacheActive = kTRUE;
441   if(GetRunNumber()<0){
442     isCacheActive=kFALSE;
443     fCalibration->SetOwner(kTRUE);
444   }
445   else{
446     fCalibration->SetOwner(kFALSE);
447   }
448
449   AliCDBManager::Instance()->SetCacheFlag(isCacheActive);
450
451   AliCDBEntry *entrySPD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSPD", run);
452   AliCDBEntry *entrySDD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSDD", run);
453   AliCDBEntry *mapASDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsAnodeSDD",run);
454   AliCDBEntry *mapTSDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsTimeSDD",run);
455   // AliCDBEntry *entrySSD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSSD", run);
456   AliCDBEntry *entryNoiseSSD = AliCDBManager::Instance()->Get("ITS/Calib/NoiseSSD");
457   AliCDBEntry *entryGainSSD = AliCDBManager::Instance()->Get("ITS/Calib/GainSSD");
458   AliCDBEntry *entryBadChannelsSSD = AliCDBManager::Instance()->Get("ITS/Calib/BadChannelsSSD");
459
460   AliCDBEntry *entry2SPD = AliCDBManager::Instance()->Get("ITS/Calib/RespSPD", run);
461   AliCDBEntry *entry2SDD = AliCDBManager::Instance()->Get("ITS/Calib/RespSDD", run);
462   AliCDBEntry *entry2SSD = AliCDBManager::Instance()->Get("ITS/Calib/RespSSD", run);
463
464   if(!entrySPD || !entrySDD || !entryNoiseSSD || !entryGainSSD || !entryBadChannelsSSD || 
465      !entry2SPD || !entry2SDD || !entry2SSD || !mapASDD ||!mapTSDD){
466     AliFatal("Calibration object retrieval failed! ");
467     return kFALSE;
468   }     
469
470 //   if(!entrySPD || !entrySDD || !entryNoiseSSD || !entryGainSSD || !entryBadChannelsSSD || 
471 //      !entry2SPD || !entry2SDD || !entry2SSD){
472 //     AliFatal("Calibration object retrieval failed! ");
473 //     return kFALSE;
474 //   }          
475
476   TObjArray *calSPD = (TObjArray *)entrySPD->GetObject();
477   if(!isCacheActive)entrySPD->SetObject(NULL);
478   entrySPD->SetOwner(kTRUE);
479
480   AliITSresponseSPD *pSPD = (AliITSresponseSPD*)entry2SPD->GetObject();
481   if(!isCacheActive)entry2SPD->SetObject(NULL);
482   entry2SPD->SetOwner(kTRUE);
483    
484   TObjArray *calSDD = (TObjArray *)entrySDD->GetObject();
485   if(!isCacheActive)entrySDD->SetObject(NULL);
486   entrySDD->SetOwner(kTRUE);
487
488   AliITSresponseSDD *pSDD = (AliITSresponseSDD*)entry2SDD->GetObject();
489   if(!isCacheActive)entry2SDD->SetObject(NULL);
490   entry2SDD->SetOwner(kTRUE);
491
492   TObjArray *mapAn = (TObjArray *)mapASDD->GetObject();
493   if(!isCacheActive)mapASDD->SetObject(NULL);
494   mapASDD->SetOwner(kTRUE);
495
496   TObjArray *mapT = (TObjArray *)mapTSDD->GetObject();
497   if(!isCacheActive)mapTSDD->SetObject(NULL);
498   mapTSDD->SetOwner(kTRUE);
499
500   /*
501   TObjArray *calSSD = (TObjArray *)entrySSD->GetObject();
502   if(!isCacheActive)entrySSD->SetObject(NULL);
503   entrySSD->SetOwner(kTRUE);
504   */
505
506   TObjArray *noiseSSD = (TObjArray *)entryNoiseSSD->GetObject();
507   if(!isCacheActive)entryNoiseSSD->SetObject(NULL);
508   entryNoiseSSD->SetOwner(kTRUE);
509
510   TObjArray *gainSSD = (TObjArray *)entryGainSSD->GetObject();
511   if(!isCacheActive)entryGainSSD->SetObject(NULL);
512   entryGainSSD->SetOwner(kTRUE);
513
514   TObjArray *badchannelsSSD = (TObjArray *)entryBadChannelsSSD->GetObject();
515   if(!isCacheActive)entryBadChannelsSSD->SetObject(NULL);
516   entryBadChannelsSSD->SetOwner(kTRUE);
517
518   AliITSresponseSSD *pSSD = (AliITSresponseSSD*)entry2SSD->GetObject();
519   if(!isCacheActive)entry2SSD->SetObject(NULL);
520   entry2SSD->SetOwner(kTRUE);
521   
522   // DB entries are deleted. In this way metadeta objects are deleted as well
523   if(!isCacheActive){
524     delete entrySPD;
525     delete entrySDD;
526     delete entryNoiseSSD;
527     delete entryGainSSD;
528     delete entryBadChannelsSSD;
529     delete entry2SPD;
530     delete entry2SDD;
531     delete entry2SSD;
532     delete mapASDD;   
533     delete mapTSDD;
534   }
535   
536   AliCDBManager::Instance()->SetCacheFlag(origCacheStatus);
537
538   if ((!pSPD)||(!pSDD)||(!pSSD) || (!calSPD) || (!calSDD) 
539       || (!mapAn) || (!mapT) || (!noiseSSD)|| (!gainSSD)|| (!badchannelsSSD)) {
540     AliWarning("Can not get calibration from calibration database !");
541     return kFALSE;
542   }
543
544   fNMod[0] = calSPD->GetEntries();
545   fNMod[1] = calSDD->GetEntries();
546   fNMod[2] = noiseSSD->GetEntries();
547   AliInfo(Form("%i SPD, %i SDD and %i SSD in calibration database",
548                fNMod[0], fNMod[1], fNMod[2]));
549   AliITSCalibration* cal;
550   for (Int_t i=0; i<fNMod[0]; i++) {
551     cal = (AliITSCalibration*) calSPD->At(i);
552     cal->SetResponse(pSPD);
553     SetCalibrationModel(i, cal);
554  }
555   Int_t nn1=0;
556   Int_t nn2=0;
557   for (Int_t i=0; i<fNMod[1]; i++) {
558     cal = (AliITSCalibration*) calSDD->At(i);
559     cal->SetResponse(pSDD);
560     AliITSMapSDD* m1 = (AliITSMapSDD*)mapAn->At(nn1);nn1++;
561     AliITSMapSDD* m2 = (AliITSMapSDD*)mapAn->At(nn1);nn1++;
562     AliITSMapSDD* m3 = (AliITSMapSDD*)mapT->At(nn2);nn2++;
563     AliITSMapSDD* m4 = (AliITSMapSDD*)mapT->At(nn2);nn2++; 
564     cal->SetMapA(0,m1);
565     cal->SetMapA(1,m2);
566     cal->SetMapT(0,m3);
567     cal->SetMapT(1,m4);
568     Int_t iMod = i + fNMod[0];
569     SetCalibrationModel(iMod, cal);
570  }
571   for (Int_t i=0; i<fNMod[2]; i++) {
572     AliITSCalibrationSSD *calibSSD = new AliITSCalibrationSSD();
573     calibSSD->SetResponse((AliITSresponse*)pSSD);
574     
575     AliITSNoiseSSD *noise = (AliITSNoiseSSD*) (noiseSSD->At(i));
576     calibSSD->SetNoise(noise);
577     AliITSGainSSD *gain = (AliITSGainSSD*) (gainSSD->At(i));
578     calibSSD->SetGain(gain);
579     AliITSBadChannelsSSD *bad = (AliITSBadChannelsSSD*) (badchannelsSSD->At(i));
580     calibSSD->SetBadChannels(bad);
581
582     Int_t iMod = i + fNMod[0] + fNMod[1];
583     SetCalibrationModel(iMod, calibSSD);
584
585  }
586   return kTRUE;
587 }
588 //_______________________________________________________________________
589 void AliITSDetTypeSim::SetDefaultSimulation(){
590   //Set default simulation for detector type
591
592   if(GetITSgeom()==0){
593     Warning("SetDefaultSimulation","GetITSgeom() is 0!");
594     return;
595   }
596   if(fCalibration==0){
597     Warning("SetDefaultSimulation","fCalibration is 0!");
598     return;
599   }
600   if(fSegmentation==0){
601     Warning("SetDefaultSimulation","fSegmentation is 0!");
602     for(Int_t i=0;i<fgkNdettypes;i++) SetDefaultSegmentation(i);
603   }else for(Int_t i=0;i<fgkNdettypes;i++) if(!GetSegmentationModel(i)){
604       Warning("SetDefaultSimulation",
605               "Segmentation not defined for det %d - Default taken\n!",i);
606       SetDefaultSegmentation(i);
607   }
608   AliITSsimulation* sim;
609
610   for(Int_t idet=0;idet<fgkNdettypes;idet++){
611    //SPD
612     if(idet==0){
613       sim = GetSimulationModel(idet); 
614       if(!sim){
615         sim = new AliITSsimulationSPD(this);
616         SetSimulationModel(idet,sim);
617       }
618     }
619     //SDD
620     if(idet==1){
621       sim = GetSimulationModel(idet);
622       if(!sim){
623         sim = new AliITSsimulationSDD(this);
624         SetSimulationModel(idet,sim);
625       }      
626     }
627     //SSD
628     if(idet==2){
629       sim = GetSimulationModel(idet);
630       if(!sim){
631         sim = new AliITSsimulationSSD(this);
632         SetSimulationModel(idet,sim);
633       }
634     }
635
636   }
637 }
638 //___________________________________________________________________
639 void AliITSDetTypeSim::SetTreeAddressS(TTree* treeS, Char_t* name){
640   // Set branch address for the ITS summable digits Trees.  
641   char branchname[30];
642
643   if(!treeS){
644     return;
645   }
646   if (fSDigits ==  0x0){
647     fSDigits = new TClonesArray("AliITSpListItem",1000);
648   }
649   TBranch *branch;
650   sprintf(branchname,"%s",name);
651   branch = treeS->GetBranch(branchname);
652   if (branch) branch->SetAddress(&fSDigits);
653
654 }
655 //___________________________________________________________________
656 void AliITSDetTypeSim::SetTreeAddressD(TTree* treeD, Char_t* name){
657   // Set branch address for the digit Trees.
658   
659   const char *det[3] = {"SPD","SDD","SSD"};
660   TBranch *branch;
661   
662   char branchname[30];
663   
664   if(!treeD){
665     return;
666   }
667   if(!fDigits){
668     fDigits = new TObjArray(fgkNdettypes); 
669   }
670   for(Int_t i=0;i<fgkNdettypes;i++){
671     Char_t* digclass = GetDigitClassName(i);
672     if(digclass==0x0){
673       if(i==0) SetDigitClassName(i,"AliITSdigitSPD");
674       if(i==1) SetDigitClassName(i,"AliITSdigitSDD");
675       if(i==2) SetDigitClassName(i,"AliITSdigitSSD");
676       digclass = GetDigitClassName(i);
677     }
678     TString classn = digclass;
679     if(!(fDigits->At(i))){
680       fDigits->AddAt(new TClonesArray(classn.Data(),1000),i);
681     }else{
682       ResetDigits(i);
683     }
684     
685     if(fgkNdettypes==3) sprintf(branchname,"%sDigits%s",name,det[i]);
686     else sprintf(branchname,"%sDigits%d",name,i+1);
687     if(fDigits){
688       branch = treeD->GetBranch(branchname);
689       if(branch) branch->SetAddress(&((*fDigits)[i]));
690     }
691   }
692
693 }
694 //___________________________________________________________________
695 void AliITSDetTypeSim::ResetDigits(){
696   // Reset number of digits and the digits array for the ITS detector.  
697
698   if(!fDigits){
699     Error("ResetDigits","fDigits is null!");
700     return;
701   }
702   for(Int_t i=0;i<fgkNdettypes;i++){
703     ResetDigits(i);
704   }
705 }
706 //___________________________________________________________________
707 void AliITSDetTypeSim::ResetDigits(Int_t branch){
708   // Reset number of digits and the digits array for this branch.
709
710   if(fDigits->At(branch)){
711     ((TClonesArray*)fDigits->At(branch))->Clear();
712   }
713   if(fNDigits) fNDigits[branch]=0;
714
715 }
716
717
718
719 //_______________________________________________________________________
720 void AliITSDetTypeSim::SDigitsToDigits(Option_t* opt, Char_t* name){
721   // Standard Summable digits to Digits function.
722   if(!GetITSgeom()){
723     Warning("SDigitsToDigits","GetITSgeom() is null!!");
724     return;
725   }
726   
727   const char *all = strstr(opt,"All");
728   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
729                         strstr(opt,"SSD")};
730   if( !det[0] && !det[1] && !det[2] ) all = "All";
731   else all = 0;
732   static Bool_t setDef = kTRUE;
733   if(setDef) SetDefaultSimulation();
734   setDef = kFALSE;
735   
736   AliITSsimulation *sim =0;
737   TTree* trees = fLoader->TreeS();
738   if( !(trees && GetSDigits()) ){
739     Error("SDigits2Digits","Error: No trees or SDigits. Returning.");
740     return;
741   } 
742   sprintf(name,"%s",name);
743   TBranch* brchSDigits = trees->GetBranch(name);
744   
745   Int_t id;
746   for(Int_t module=0;module<GetITSgeom()->GetIndexMax();module++){
747      id = GetITSgeom()->GetModuleType(module);
748     if (!all && !det[id]) continue;
749     sim = (AliITSsimulation*)GetSimulationModel(id);
750     if(!sim){
751       Error("SDigit2Digits","The simulation class was not "
752             "instanciated for module %d type %s!",module,
753             GetITSgeom()->GetModuleTypeName(module));
754       exit(1);
755     }
756     sim->InitSimulationModule(module,gAlice->GetEvNumber());
757     
758     fSDigits->Clear();
759     brchSDigits->GetEvent(module);
760     sim->AddSDigitsToModule(fSDigits,0);
761     sim->FinishSDigitiseModule();
762     fLoader->TreeD()->Fill();
763     ResetDigits();
764   }
765   fLoader->TreeD()->GetEntries();
766   fLoader->TreeD()->AutoSave();
767   fLoader->TreeD()->Reset();
768 }
769 //_________________________________________________________
770 void AliITSDetTypeSim::AddSumDigit(AliITSpListItem &sdig){  
771   //Adds the module full of summable digits to the summable digits tree.
772
773   TClonesArray &lsdig = *fSDigits;
774   new(lsdig[fNSDigits++]) AliITSpListItem(sdig);
775 }
776 //__________________________________________________________
777 void AliITSDetTypeSim::AddRealDigit(Int_t branch, Int_t *digits){
778   //   Add a real digit - as coming from data.
779
780   TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
781   new(ldigits[fNDigits[branch]++]) AliITSdigit(digits); 
782 }
783 //__________________________________________________________
784 void AliITSDetTypeSim::AddSimDigit(Int_t branch, AliITSdigit* d){  
785   //    Add a simulated digit.
786
787   TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
788   switch(branch){
789   case 0:
790     new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
791     break;
792   case 1:
793     new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
794     break;
795   case 2:
796     new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
797     break;
798   }  
799 }
800 //______________________________________________________________________
801 void AliITSDetTypeSim::AddSimDigit(Int_t branch,Float_t phys,Int_t *digits,
802                                    Int_t *tracks,Int_t *hits,Float_t *charges){
803   //   Add a simulated digit to the list.
804
805   TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
806   AliITSCalibrationSDD *resp = 0;
807   switch(branch){
808   case 0:
809     new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(digits,tracks,hits);
810     break;
811   case 1:
812     resp = (AliITSCalibrationSDD*)GetCalibrationModel(
813         GetITSgeom()->GetStartSDD());
814     new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(phys,digits,tracks,
815                                                    hits,charges,resp);
816     break;
817   case 2:
818     new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(digits,tracks,hits);
819     break;
820   } 
821 }
822 //______________________________________________________________________
823 void AliITSDetTypeSim::StoreCalibration(Int_t firstRun, Int_t lastRun,
824                                         AliCDBMetaData &md) {
825   // Store calibration in the calibration database
826   // The database must be created in an external piece of code (i.e. 
827   // a configuration macro )
828
829   if(!AliCDBManager::Instance()->IsDefaultStorageSet()) {
830     AliWarning("No storage set! Will use dummy one");
831     AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT");
832   }
833
834   if (!fCalibration) {
835     AliError("AliITSCalibration classes are not defined - nothing done");
836     return;
837   }
838   AliCDBId idRespSPD("ITS/Calib/CalibSPD",firstRun, lastRun);
839   AliCDBId idRespSDD("ITS/Calib/CalibSDD",firstRun, lastRun);
840   AliCDBId idRespSSD("ITS/Calib/CalibSSD",firstRun, lastRun);
841
842   TObjArray respSPD(fNMod[0]);
843   TObjArray respSDD(fNMod[1]-fNMod[0]);
844   TObjArray respSSD(fNMod[2]-fNMod[1]);
845   respSPD.SetOwner(kFALSE);
846   respSSD.SetOwner(kFALSE);
847   respSSD.SetOwner(kFALSE);
848
849   Int_t index[fgkNdettypes];
850   for (Int_t i = 0; i<fgkNdettypes; i++ ) {
851     index[i] = 0;
852     for (Int_t j = 0; j<=i; j++ )
853       index[i]+=fNMod[j];
854   }
855
856   for (Int_t i = 0; i<index[0]; i++ )
857     respSPD.Add(fCalibration->At(i));
858
859   for (Int_t i = index[0]; i<index[1]; i++ )
860     respSDD.Add(fCalibration->At(i));
861
862   for (Int_t i = index[1]; i<index[2]; i++ )
863     respSSD.Add(fCalibration->At(i));
864
865   AliCDBManager::Instance()->Put(&respSPD, idRespSPD, &md);
866   AliCDBManager::Instance()->Put(&respSDD, idRespSDD, &md);
867   AliCDBManager::Instance()->Put(&respSSD, idRespSSD, &md);
868 }
869
870