]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSDetTypeSim.cxx
Fixes for bug #71984 ENUM_AS_BOOLEAN problem reported by Coverity
[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 #include "AliITSdigit.h"
39 #include "AliITSdigitSPD.h"
40 #include "AliITSdigitSDD.h"
41 #include "AliITSdigitSSD.h"
42 #include "AliITSgeom.h"
43 #include "AliITSDetTypeSim.h"
44 #include "AliITSpListItem.h"
45 #include "AliITSCalibration.h"
46 #include "AliITSCalibrationSDD.h"
47 #include "AliITSMapSDD.h"
48 #include "AliITSCorrMapSDD.h"
49 #include "AliITSDriftSpeedArraySDD.h"
50 #include "AliITSDriftSpeedSDD.h"
51 #include "AliITSCalibrationSSD.h"
52 #include "AliITSNoiseSSDv2.h"
53 #include "AliITSGainSSDv2.h"
54 #include "AliITSBadChannelsSSDv2.h"
55 #include "AliITSNoiseSSD.h"
56 #include "AliITSGainSSD.h"
57 #include "AliITSBadChannelsSSD.h"
58 #include "AliITSCalibrationSSD.h"
59 #include "AliITSsegmentationSPD.h"
60 #include "AliITSsegmentationSDD.h"
61 #include "AliITSsegmentationSSD.h"
62 #include "AliITSsimulation.h"
63 #include "AliITSsimulationSPD.h"
64 #include "AliITSsimulationSDD.h"
65 #include "AliITSsimulationSSD.h"
66 #include "AliITSDDLModuleMapSDD.h"
67 #include "AliITSTriggerConditions.h"
68 #include "AliBaseLoader.h"
69
70 const Int_t AliITSDetTypeSim::fgkNdettypes = 3;
71 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSPD =  240;
72 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSDD =  260;
73 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSSD = 1698;
74
75 ClassImp(AliITSDetTypeSim)
76
77 //----------------------------------------------------------------------
78 AliITSDetTypeSim::AliITSDetTypeSim():
79 TObject(),
80 fSimulation(),   // [NDet]
81 fSegmentation(), // [NDet]
82 fCalibration(),     // [NMod]
83 fSSDCalibration(0),
84 fSPDNoisy(0),
85 fNSDigits(0),    //! number of SDigits
86 fSDigits("AliITSpListItem",1000),   
87 fNDigits(0),     //! number of Digits
88 fRunNumber(0),   //! Run number (to access DB)
89 fDigits(),       //! [NMod][NDigits]
90 fSimuPar(0),
91 fDDLMapSDD(0),
92 fAveGainSDD(0),
93 fkDigClassName(), // String with digit class name.
94 fLoader(0),      // local pointer to loader
95 fFirstcall(kTRUE),
96 fFOGenerator(),
97 fTriggerConditions(NULL)
98
99     // Default Constructor
100     // Inputs:
101     //    none.
102     // Outputs:
103     //    none.
104     // Return:
105     //    A properly zero-ed AliITSDetTypeSim class.
106
107   fSimulation = new TObjArray(fgkNdettypes);
108   fSegmentation = new TObjArray(fgkNdettypes);
109   fSegmentation->SetOwner(kTRUE);
110   fDigits = new TObjArray(fgkNdettypes);
111   fNDigits = new Int_t[fgkNdettypes];
112   fDDLMapSDD=new AliITSDDLModuleMapSDD();
113   fSimuPar= new AliITSSimuParam();
114   fSSDCalibration=new AliITSCalibrationSSD();
115   fNMod[0] = fgkDefaultNModulesSPD;
116   fNMod[1] = fgkDefaultNModulesSDD;
117   fNMod[2] = fgkDefaultNModulesSSD;
118   SetRunNumber();
119 }
120 //----------------------------------------------------------------------
121 AliITSDetTypeSim::~AliITSDetTypeSim(){
122     // Destructor
123     // Inputs:
124     //    none.
125     // Outputs:
126     //    none.
127     // Return:
128     //    Nothing.
129
130     if(fSimulation){
131         fSimulation->Delete();
132         delete fSimulation;
133     }
134     fSimulation = 0;
135     if(fSegmentation){
136         fSegmentation->Delete();
137         delete fSegmentation;
138     }
139     fSegmentation = 0;
140     if(fCalibration && fRunNumber<0){
141
142         fCalibration->Delete();
143         delete fCalibration;
144     }
145     fCalibration = 0;
146     if(fSSDCalibration) {
147       if(!(AliCDBManager::Instance()->GetCacheFlag())) {
148         delete fSSDCalibration;
149         fSSDCalibration = NULL;
150       }
151     }
152     if(fSPDNoisy){
153     if(!(AliCDBManager::Instance()->GetCacheFlag())) {
154        fSPDNoisy->Delete();
155        delete fSPDNoisy;
156        fSPDNoisy = 0;
157       }
158     }
159     if(fSimuPar) delete fSimuPar;
160     if(fDDLMapSDD) delete fDDLMapSDD;
161     if(fNDigits) delete [] fNDigits;
162     fNDigits = 0;
163     if (fLoader)fLoader->GetModulesFolder()->Remove(this);
164     fLoader = 0; // Not deleting it.
165     fSDigits.Delete();
166     if (fDigits) {
167       fDigits->Delete();
168       delete fDigits;
169     }
170     fDigits=0;
171 }
172 //----------------------------------------------------------------------
173 AliITSDetTypeSim::AliITSDetTypeSim(const AliITSDetTypeSim &source) : TObject(source),
174 fSimulation(source.fSimulation),   // [NDet]
175 fSegmentation(source.fSegmentation), // [NDet]
176 fCalibration(source.fCalibration),     // [NMod]
177 fSSDCalibration(source.fSSDCalibration),
178 fSPDNoisy(source.fSPDNoisy),
179 fNSDigits(source.fNSDigits),    //! number of SDigits
180 fSDigits(*((TClonesArray*)source.fSDigits.Clone())),
181 fNDigits(source.fNDigits),     //! number of Digits
182 fRunNumber(source.fRunNumber),   //! Run number (to access DB)
183 fDigits(source.fDigits),       //! [NMod][NDigits]
184 fSimuPar(source.fSimuPar),
185 fDDLMapSDD(source.fDDLMapSDD),
186 fAveGainSDD(source.fAveGainSDD),
187 fkDigClassName(), // String with digit class name.
188 fLoader(source.fLoader),      // local pointer to loader
189 fFirstcall(source.fFirstcall),
190 fFOGenerator(source.fFOGenerator),
191 fTriggerConditions(source.fTriggerConditions) 
192 {
193     // Copy Constructor for object AliITSDetTypeSim not allowed
194   for(Int_t i=0;i<fgkNdettypes;i++){
195     fkDigClassName[i] = source.fkDigClassName[i];
196   }
197 }
198 //----------------------------------------------------------------------
199 AliITSDetTypeSim& AliITSDetTypeSim::operator=(const AliITSDetTypeSim &source){
200     // The = operator for object AliITSDetTypeSim
201  
202   this->~AliITSDetTypeSim();
203   new(this) AliITSDetTypeSim(source);
204   return *this;
205 }
206
207 //______________________________________________________________________
208 void AliITSDetTypeSim::SetITSgeom(AliITSgeom *geom){
209     // Sets/replaces the existing AliITSgeom object kept in AliITSLoader
210     // 
211     // Inputs:
212     //   AliITSgoem   *geom  The AliITSgeom object to be used.
213     // Output:
214     //   none.
215     // Return:
216     //   none.
217   if(!fLoader){
218     Error("SetITSgeom","No pointer to loader - nothing done");
219     return;
220   }
221   else {
222     fLoader->SetITSgeom(geom);  // protections in AliITSLoader::SetITSgeom
223   }
224  
225 }
226 //______________________________________________________________________
227 void AliITSDetTypeSim::SetLoader(AliITSLoader *loader){
228     // Sets the local copy of the AliITSLoader, and passes on the
229     // AliITSgeom object as needed.
230     // Inputs
231     //   AliITSLoader  *loader pointer to AliITSLoader for local use
232     // Outputs:
233     //   none.
234     // Return:
235     //  none.
236
237     if(fLoader==loader) return; // Same do nothing
238     if(fLoader){ // alread have an existing loader
239         Error("SetLoader",
240                 "Already have an exisiting loader ptr=%p Nothing done",
241                 fLoader);
242     } // end if
243     fLoader = loader;
244 }
245 //______________________________________________________________________
246 void AliITSDetTypeSim::SetSimulationModel(Int_t dettype,AliITSsimulation *sim){
247
248   //Set simulation model for detector type
249
250   if(fSimulation==0) fSimulation = new TObjArray(fgkNdettypes);
251   fSimulation->AddAt(sim,dettype);
252 }
253 //______________________________________________________________________
254 AliITSsimulation* AliITSDetTypeSim::GetSimulationModel(Int_t dettype) const { 
255
256   //Get simulation model for detector type
257   if(fSimulation==0)  {
258     Warning("GetSimulationModel","fSimulation is 0!");
259     return 0;     
260   }
261   return (AliITSsimulation*)(fSimulation->At(dettype));
262 }
263 //______________________________________________________________________
264 AliITSsimulation* AliITSDetTypeSim::GetSimulationModelByModule(Int_t module) const {
265
266   //Get simulation model by module number
267   if(GetITSgeom()==0) {
268     Warning("GetSimulationModelByModule","GetITSgeom() is 0!");
269     return 0;
270   }
271   
272   return GetSimulationModel(GetITSgeom()->GetModuleType(module));
273 }
274 //_______________________________________________________________________
275 void AliITSDetTypeSim::SetDefaultSegmentation(Int_t idet){
276     // Set default segmentation model objects
277     AliITSsegmentation *seg;
278
279     if(fSegmentation==0x0){
280         fSegmentation = new TObjArray(fgkNdettypes);
281         fSegmentation->SetOwner(kTRUE);
282     }
283     if(GetSegmentationModel(idet))
284         delete (AliITSsegmentation*)fSegmentation->At(idet);
285     if(idet==0){
286         seg = new AliITSsegmentationSPD();
287     }else if(idet==1){
288       seg = new AliITSsegmentationSDD();
289       AliITSCalibrationSDD* cal=(AliITSCalibrationSDD*)GetCalibrationModel(fgkDefaultNModulesSPD+1);
290       if(cal->IsAMAt20MHz()){ 
291         seg->SetPadSize(seg->Dpz(0),20.);
292         seg->SetNPads(seg->Npz()/2,128);
293       }
294     }else {
295         seg = new AliITSsegmentationSSD();
296     }
297     SetSegmentationModel(idet,seg);
298 }
299 //______________________________________________________________________
300 void AliITSDetTypeSim::SetSegmentationModel(Int_t dettype,
301                                             AliITSsegmentation *seg){
302    
303   //Set segmentation model for detector type
304   if(fSegmentation==0x0){
305     fSegmentation = new TObjArray(fgkNdettypes);
306     fSegmentation->SetOwner(kTRUE);
307   }
308   fSegmentation->AddAt(seg,dettype);
309 }
310 //______________________________________________________________________
311 AliITSsegmentation* AliITSDetTypeSim::GetSegmentationModel(Int_t dettype) const{
312   //Get segmentation model for detector type
313    
314    if(fSegmentation==0) {
315        Warning("GetSegmentationModel","fSegmentation is 0!");
316        return 0; 
317    } 
318    return (AliITSsegmentation*)(fSegmentation->At(dettype));
319 }
320 //_______________________________________________________________________
321 AliITSsegmentation* AliITSDetTypeSim::GetSegmentationModelByModule(Int_t module) const{
322     //Get segmentation model by module number
323     if(GetITSgeom()==0){
324         Warning("GetSegmentationModelByModule","GetITSgeom() is 0!");
325         return 0;
326     }     
327     return GetSegmentationModel(GetITSgeom()->GetModuleType(module));
328 }
329 //_______________________________________________________________________
330 void AliITSDetTypeSim::CreateCalibrationArray() {
331     //Create the container of calibration functions with correct size
332     if (fCalibration) {
333         Warning("CreateCalibration","pointer to calibration object exists\n");
334         fCalibration->Delete();
335         delete fCalibration;
336     }
337
338     Int_t nModTot = GetITSgeom()->GetIndexMax();
339     fCalibration = new TObjArray(nModTot);
340     fCalibration->SetOwner(kTRUE);
341     fCalibration->Clear();
342 }
343 //_______________________________________________________________________
344 void AliITSDetTypeSim::SetCalibrationModel(Int_t iMod, AliITSCalibration *resp){
345     //Set response model for modules
346
347     if (fCalibration==0) CreateCalibrationArray();
348  
349     if (fCalibration->At(iMod)!=0)
350         delete (AliITSCalibration*) fCalibration->At(iMod);
351     fCalibration->AddAt(resp, iMod);
352 }
353 //_______________________________________________________________________
354 void AliITSDetTypeSim::SetSPDNoisyModel(Int_t iMod, AliITSCalibration *cal){
355   //Set noisy pixel info for the SPD module iMod
356   if (fSPDNoisy==0) {
357     fSPDNoisy = new TObjArray(fgkDefaultNModulesSPD);
358     fSPDNoisy->SetOwner(kTRUE);
359     fSPDNoisy->Clear();
360   }
361
362   if (fSPDNoisy->At(iMod) != 0)
363     delete (AliITSCalibration*) fSPDNoisy->At(iMod);
364   fSPDNoisy->AddAt(cal,iMod);
365 }
366 //______________________________________________________________________
367 void AliITSDetTypeSim::ResetCalibrationArray(){
368     //resets response array
369     if(fCalibration && fRunNumber<0){  // if fRunNumber<0 fCalibration is owner
370       /*
371         AliITSresponse* rspd = ((AliITSCalibration*)fCalibration->At(
372                                 GetITSgeom()->GetStartSPD()))->GetResponse();
373         AliITSresponse* rssd = ((AliITSCalibration*)fCalibration->At(
374                                 GetITSgeom()->GetStartSSD()))->GetResponse();
375         if(rspd) delete rspd;
376         if(rssd) delete rssd;
377       */
378         fCalibration->Clear();
379     }else if (fCalibration && fRunNumber>=0){
380         fCalibration->Clear();
381     }
382 }
383 //______________________________________________________________________
384 void AliITSDetTypeSim::ResetSegmentation(){
385     //Resets segmentation array
386     if(fSegmentation) fSegmentation->Clear();
387 }
388 //_______________________________________________________________________
389 AliITSCalibration* AliITSDetTypeSim::GetCalibrationModel(Int_t iMod) const {
390     //Get response model for module number iMod 
391  
392     if(fCalibration==0) {
393         AliError("fCalibration is 0!");
394         return 0; 
395     }
396   if(iMod<fgkDefaultNModulesSPD+fgkDefaultNModulesSDD){
397     return (AliITSCalibration*)fCalibration->At(iMod);
398   }else{
399     Int_t i=iMod-(fgkDefaultNModulesSPD+fgkDefaultNModulesSDD);
400     fSSDCalibration->SetModule(i);
401     return (AliITSCalibration*)fSSDCalibration;
402   }
403
404 }
405 //_______________________________________________________________________
406 AliITSCalibration* AliITSDetTypeSim::GetSPDNoisyModel(Int_t iMod) const {
407   //Get SPD noisy calib for module iMod 
408   if(fSPDNoisy==0) {
409     AliWarning("fSPDNoisy is 0!");
410     return 0; 
411   }
412   return (AliITSCalibration*)fSPDNoisy->At(iMod);
413 }
414 //_______________________________________________________________________
415 void AliITSDetTypeSim::SetDefaults(){
416     //Set defaults for segmentation and response
417
418     if(GetITSgeom()==0){
419         Warning("SetDefaults","GetITSgeom() is 0!");
420         return;
421     } // end if
422     if (fCalibration==0) {
423         CreateCalibrationArray();
424     } // end if
425
426     ResetSegmentation();
427     if(!GetCalibration()){AliFatal("Exit"); exit(0);}
428
429     SetDigitClassName(0,"AliITSdigitSPD");
430     SetDigitClassName(1,"AliITSdigitSDD");
431     SetDigitClassName(2,"AliITSdigitSSD");
432
433     for(Int_t idet=0;idet<fgkNdettypes;idet++){
434       if(!GetSegmentationModel(idet)) SetDefaultSegmentation(idet);
435     }
436 }
437 //______________________________________________________________________
438 Bool_t AliITSDetTypeSim::GetCalibration() {
439   // Get Default calibration if a storage is not defined.
440
441   if(!fFirstcall){
442     AliITSCalibration* cal = GetCalibrationModel(0);
443     if(cal)return kTRUE;
444   }
445   else {
446     fFirstcall = kFALSE;
447   }
448
449   SetRunNumber((Int_t)AliCDBManager::Instance()->GetRun());
450   Int_t run=GetRunNumber();
451
452   Bool_t origCacheStatus = AliCDBManager::Instance()->GetCacheFlag();
453   Bool_t isCacheActive = kTRUE;
454   if(GetRunNumber()<0){
455     isCacheActive=kFALSE;
456     fCalibration->SetOwner(kTRUE);
457   }
458   else{
459     fCalibration->SetOwner(kFALSE);
460   }
461
462   AliCDBManager::Instance()->SetCacheFlag(isCacheActive);
463
464   AliCDBEntry *deadSPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDDead", run);
465   AliCDBEntry *noisySPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDNoisy", run);
466   AliCDBEntry *foEffSPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDFOEfficiency", run);
467   AliCDBEntry *foNoiSPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDFONoise", run);
468   AliCDBEntry *entrySDD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSDD", run);
469   AliCDBEntry *drSpSDD = AliCDBManager::Instance()->Get("ITS/Calib/DriftSpeedSDD",run);
470   AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD",run);
471   //AliCDBEntry *mapASDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsAnodeSDD",run);
472   AliCDBEntry *mapTSDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsTimeSDD",run);
473   // AliCDBEntry *entrySSD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSSD", run);
474   AliCDBEntry *entryNoiseSSD = AliCDBManager::Instance()->Get("ITS/Calib/NoiseSSD");
475   AliCDBEntry *entryGainSSD = AliCDBManager::Instance()->Get("ITS/Calib/GainSSD");
476   AliCDBEntry *entryBadChannelsSSD = AliCDBManager::Instance()->Get("ITS/Calib/BadChannelsSSD");
477
478 if(!deadSPD || !noisySPD || !foEffSPD || !foNoiSPD 
479      || !entrySDD || !entryNoiseSSD || !entryGainSSD || !entryBadChannelsSSD 
480      || !drSpSDD || !ddlMapSDD || !mapTSDD){
481     AliFatal("Calibration object retrieval failed! ");
482     return kFALSE;
483   }             
484         
485
486   TObjArray *calDeadSPD = (TObjArray*) deadSPD->GetObject();
487   if (!isCacheActive) deadSPD->SetObject(NULL);
488   deadSPD->SetOwner(kTRUE);
489
490   TObjArray *calNoisySPD = (TObjArray*) noisySPD->GetObject();
491   if (!isCacheActive) noisySPD->SetObject(NULL);
492   noisySPD->SetOwner(kTRUE);
493
494   AliITSFOEfficiencySPD *calFoEffSPD = (AliITSFOEfficiencySPD*) foEffSPD->GetObject();
495   if (!isCacheActive) foEffSPD->SetObject(NULL);
496   foEffSPD->SetOwner(kTRUE);
497
498   AliITSFONoiseSPD *calFoNoiSPD = (AliITSFONoiseSPD*) foNoiSPD->GetObject();
499   if (!isCacheActive) foNoiSPD->SetObject(NULL);
500   foNoiSPD->SetOwner(kTRUE);
501    
502   TObjArray *calSDD = (TObjArray *)entrySDD->GetObject();
503   if(!isCacheActive)entrySDD->SetObject(NULL);
504   entrySDD->SetOwner(kTRUE);
505
506   TObjArray *drSp = (TObjArray *)drSpSDD->GetObject();
507   if(!isCacheActive)drSpSDD->SetObject(NULL);
508   drSpSDD->SetOwner(kTRUE);
509
510   AliITSDDLModuleMapSDD *ddlsdd=(AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject();
511   if(!isCacheActive)ddlMapSDD->SetObject(NULL);
512   ddlMapSDD->SetOwner(kTRUE);
513
514 //   TObjArray *mapAn = (TObjArray *)mapASDD->GetObject();
515 //   if(!isCacheActive)mapASDD->SetObject(NULL);
516 //   mapASDD->SetOwner(kTRUE);
517
518   TObjArray *mapT = (TObjArray *)mapTSDD->GetObject();
519   if(!isCacheActive)mapTSDD->SetObject(NULL);
520   mapTSDD->SetOwner(kTRUE);
521
522   /*
523   TObjArray *calSSD = (TObjArray *)entrySSD->GetObject();
524   if(!isCacheActive)entrySSD->SetObject(NULL);
525   entrySSD->SetOwner(kTRUE);
526   */
527
528   TObject *emptyssd = 0; TString ssdobjectname;
529   AliITSNoiseSSDv2 *noiseSSD = NULL; 
530   emptyssd = (TObject *)entryNoiseSSD->GetObject();
531   ssdobjectname = emptyssd->GetName();
532   if(ssdobjectname=="TObjArray") {
533     noiseSSD = new AliITSNoiseSSDv2(); 
534     TObjArray *noiseSSDOld = (TObjArray *)entryNoiseSSD->GetObject();
535     ReadOldSSDNoise(noiseSSDOld, noiseSSD);
536   }
537   else if(ssdobjectname=="AliITSNoiseSSDv2")
538     noiseSSD = (AliITSNoiseSSDv2 *)entryNoiseSSD->GetObject();
539   if(!isCacheActive)entryNoiseSSD->SetObject(NULL);
540   entryNoiseSSD->SetOwner(kTRUE);
541
542   AliITSGainSSDv2 *gainSSD = NULL;
543   emptyssd = (TObject *)entryGainSSD->GetObject();
544   ssdobjectname = emptyssd->GetName();
545   if(ssdobjectname=="Gain") {
546     TObjArray *gainSSDOld = (TObjArray *)entryGainSSD->GetObject();
547     gainSSD = new AliITSGainSSDv2();
548     ReadOldSSDGain(gainSSDOld, gainSSD);
549   }
550   else if(ssdobjectname=="AliITSGainSSDv2")
551     gainSSD = (AliITSGainSSDv2 *)entryGainSSD->GetObject();
552   if(!isCacheActive)entryGainSSD->SetObject(NULL);
553   entryGainSSD->SetOwner(kTRUE);
554
555   AliITSBadChannelsSSDv2 *badChannelsSSD = NULL;
556   emptyssd = (TObject *)entryBadChannelsSSD->GetObject();
557   ssdobjectname = emptyssd->GetName();
558   if(ssdobjectname=="TObjArray") {
559     TObjArray *badChannelsSSDOld = (TObjArray *)entryBadChannelsSSD->GetObject();
560     badChannelsSSD = new AliITSBadChannelsSSDv2();
561     ReadOldSSDBadChannels(badChannelsSSDOld, badChannelsSSD);
562   }
563   else if(ssdobjectname=="AliITSBadChannelsSSDv2")
564     badChannelsSSD = (AliITSBadChannelsSSDv2*)entryBadChannelsSSD->GetObject();
565   if(!isCacheActive)entryBadChannelsSSD->SetObject(NULL);
566   entryBadChannelsSSD->SetOwner(kTRUE);
567
568   /*AliITSNoiseSSDv2 *noiseSSD = (AliITSNoiseSSDv2 *)entryNoiseSSD->GetObject();
569   if(!isCacheActive)entryNoiseSSD->SetObject(NULL);
570   entryNoiseSSD->SetOwner(kTRUE);
571
572   AliITSGainSSDv2 *gainSSD = (AliITSGainSSDv2 *)entryGainSSD->GetObject();
573   if(!isCacheActive)entryGainSSD->SetObject(NULL);
574   entryGainSSD->SetOwner(kTRUE);
575
576   AliITSBadChannelsSSDv2 *badchannelsSSD = 
577     (AliITSBadChannelsSSDv2 *)entryBadChannelsSSD->GetObject();
578   if(!isCacheActive)entryBadChannelsSSD->SetObject(NULL);
579   entryBadChannelsSSD->SetOwner(kTRUE);*/
580
581   // DB entries are deleted. In this way metadeta objects are deleted as well
582   if(!isCacheActive){
583     delete deadSPD;
584     delete noisySPD;
585     delete foEffSPD;
586     delete foNoiSPD;
587     delete entrySDD;
588     delete entryNoiseSSD;
589     delete entryGainSSD;
590     delete entryBadChannelsSSD;
591 //    delete mapASDD;   
592     delete mapTSDD;
593     delete drSpSDD;
594     delete ddlMapSDD;
595   }
596   
597   AliCDBManager::Instance()->SetCacheFlag(origCacheStatus);
598
599  if ((!calDeadSPD) || (!calNoisySPD) || (!calFoEffSPD) || (!calFoNoiSPD) 
600       || (!calSDD) || (!drSp) || (!ddlsdd)  
601       || (!mapT) || (!noiseSSD)|| (!gainSSD)|| (!badChannelsSSD)) {
602     AliWarning("Can not get calibration from calibration database !");
603     return kFALSE;
604   }
605
606
607   fNMod[0] = calDeadSPD->GetEntries();
608   fNMod[1] = calSDD->GetEntries();
609   //  fNMod[2] = noiseSSD->GetEntries();
610   AliDebug(1,Form("%i SPD, %i SDD and %i SSD in calibration database",
611                fNMod[0], fNMod[1], fNMod[2]));
612   AliITSCalibration* cal;
613   for (Int_t i=0; i<fNMod[0]; i++) {
614     cal = (AliITSCalibration*) calDeadSPD->At(i);
615     SetCalibrationModel(i, cal);
616     cal = (AliITSCalibration*) calNoisySPD->At(i);
617     SetSPDNoisyModel(i, cal);
618   }
619
620   fFOGenerator.SetEfficiency(calFoEffSPD); // this cal object is used only by the generator
621   fFOGenerator.SetNoise(calFoNoiSPD); // this cal object is used only by the generator
622   
623   fDDLMapSDD->SetDDLMap(ddlsdd);
624   Float_t avegain=0.;
625   Float_t nGdAnodes=0;
626   Bool_t oldMapFormat=kFALSE;
627   TObject* objmap=(TObject*)mapT->At(0);
628   TString cname(objmap->ClassName());
629   if(cname.CompareTo("AliITSMapSDD")==0){ 
630     oldMapFormat=kTRUE;
631     AliInfo("SDD Maps converted to new format");
632   }
633
634   for (Int_t i=0; i<fgkDefaultNModulesSDD; i++) {
635     Int_t iddl,icarlos;
636     Int_t iMod = i + fgkDefaultNModulesSPD;
637     fDDLMapSDD->FindInDDLMap(iMod,iddl,icarlos);
638     if(iddl<0){ 
639       AliITSCalibrationSDD* calsdddead=new AliITSCalibrationSDD();
640       calsdddead->SetBad();      
641       AliITSDriftSpeedSDD* driftspdef = new AliITSDriftSpeedSDD();
642       AliITSDriftSpeedArraySDD* arrdrsp=new AliITSDriftSpeedArraySDD(1);
643       arrdrsp->AddDriftSpeed(driftspdef);
644       calsdddead->SetDriftSpeed(0,arrdrsp);
645       calsdddead->SetDriftSpeed(1,arrdrsp);
646       SetCalibrationModel(iMod, calsdddead);
647       AliWarning(Form("SDD module %d not present in DDL map: set it as dead",iMod));
648     }else{
649       cal = (AliITSCalibration*) calSDD->At(i);
650       for(Int_t iAnode=0;iAnode< ((AliITSCalibrationSDD*)cal)->NOfAnodes(); iAnode++){
651         if(((AliITSCalibrationSDD*)cal)->IsBadChannel(iAnode)) continue;
652         avegain+= ((AliITSCalibrationSDD*)cal)->GetChannelGain(iAnode);
653         nGdAnodes++;
654       }
655       Int_t i0=2*i;
656       Int_t i1=1+2*i;
657       AliITSDriftSpeedArraySDD* arr0 = (AliITSDriftSpeedArraySDD*) drSp->At(i0);
658       AliITSDriftSpeedArraySDD* arr1 = (AliITSDriftSpeedArraySDD*) drSp->At(i1);
659
660       AliITSCorrMapSDD* mt0 = 0;
661       AliITSCorrMapSDD* mt1 = 0;
662       if(oldMapFormat){ 
663         AliITSMapSDD* oldmap0=(AliITSMapSDD*)mapT->At(i0);
664         AliITSMapSDD* oldmap1=(AliITSMapSDD*)mapT->At(i1);
665         mt0=oldmap0->ConvertToNewFormat();
666         mt1=oldmap1->ConvertToNewFormat();
667       }else{
668         mt0=(AliITSCorrMapSDD*)mapT->At(i0);
669         mt1=(AliITSCorrMapSDD*)mapT->At(i1);
670       }
671       cal->SetDriftSpeed(0,arr0);
672       cal->SetDriftSpeed(1,arr1);
673       cal->SetMapT(0,mt0);
674       cal->SetMapT(1,mt1);
675       SetCalibrationModel(iMod, cal);
676     }
677   }
678   if(nGdAnodes) fAveGainSDD=avegain/nGdAnodes;
679   AliDebug(3,Form("SDD average gain=%f\n",fAveGainSDD));
680   fSSDCalibration->SetNoise(noiseSSD);
681   fSSDCalibration->SetGain(gainSSD);
682   fSSDCalibration->SetBadChannels(badChannelsSSD);
683   //fSSDCalibration->FillBadChipMap();
684
685
686
687   return kTRUE;
688 }
689 //_______________________________________________________________________
690 void AliITSDetTypeSim::SetDefaultSimulation(){
691   //Set default simulation for detector type
692
693   if(GetITSgeom()==0){
694     Warning("SetDefaultSimulation","GetITSgeom() is 0!");
695     return;
696   }
697   if(fCalibration==0){
698     Warning("SetDefaultSimulation","fCalibration is 0!");
699     return;
700   }
701   if(fSegmentation==0){
702     Warning("SetDefaultSimulation","fSegmentation is 0!");
703     for(Int_t i=0;i<fgkNdettypes;i++) SetDefaultSegmentation(i);
704   }else for(Int_t i=0;i<fgkNdettypes;i++) if(!GetSegmentationModel(i)){
705       Warning("SetDefaultSimulation",
706               "Segmentation not defined for det %d - Default taken\n!",i);
707       SetDefaultSegmentation(i);
708   }
709   AliITSsimulation* sim;
710
711   for(Int_t idet=0;idet<fgkNdettypes;idet++){
712    //SPD
713     if(idet==0){
714       sim = GetSimulationModel(idet); 
715       if(!sim){
716         sim = new AliITSsimulationSPD(this);
717         SetSimulationModel(idet,sim);
718       }
719     }
720     //SDD
721     if(idet==1){
722       sim = GetSimulationModel(idet);
723       if(!sim){
724         sim = new AliITSsimulationSDD(this);
725         SetSimulationModel(idet,sim);
726       }      
727     }
728     //SSD
729     if(idet==2){
730       sim = GetSimulationModel(idet);
731       if(!sim){
732         sim = new AliITSsimulationSSD(this);
733         SetSimulationModel(idet,sim);
734       }
735     }
736
737   }
738 }
739 //___________________________________________________________________
740 void AliITSDetTypeSim::SetTreeAddressS(TTree* treeS, const Char_t* name){
741   // Set branch address for the ITS summable digits Trees.  
742   char branchname[30];
743
744   if(!treeS){
745     return;
746   }
747   TBranch *branch;
748   sprintf(branchname,"%s",name);
749   branch = treeS->GetBranch(branchname);
750   TClonesArray *sdigi = &fSDigits;
751   if (branch) branch->SetAddress(&sdigi);
752
753 }
754 //___________________________________________________________________
755 void AliITSDetTypeSim::SetTreeAddressD(TTree* treeD, const Char_t* name){
756   // Set branch address for the digit Trees.
757   
758   const char *det[3] = {"SPD","SDD","SSD"};
759   TBranch *branch;
760   
761   char branchname[30];
762   
763   if(!treeD){
764     return;
765   }
766   if(!fDigits){
767     fDigits = new TObjArray(fgkNdettypes); 
768   }
769   for(Int_t i=0;i<fgkNdettypes;i++){
770     const Char_t* digclass = GetDigitClassName(i);
771     if(digclass==0x0){
772       if(i==0) SetDigitClassName(i,"AliITSdigitSPD");
773       if(i==1) SetDigitClassName(i,"AliITSdigitSDD");
774       if(i==2) SetDigitClassName(i,"AliITSdigitSSD");
775       digclass = GetDigitClassName(i);
776     }
777     TString classn = digclass;
778     if(!(fDigits->At(i))){
779       fDigits->AddAt(new TClonesArray(classn.Data(),1000),i);
780     }else{
781       ResetDigits(i);
782     }
783     
784     if(fgkNdettypes==3) sprintf(branchname,"%sDigits%s",name,det[i]);
785     else sprintf(branchname,"%sDigits%d",name,i+1);
786     if(fDigits){
787       branch = treeD->GetBranch(branchname);
788       if(branch) branch->SetAddress(&((*fDigits)[i]));
789     }
790   }
791
792 }
793 //___________________________________________________________________
794 void AliITSDetTypeSim::ResetDigits(){
795   // Reset number of digits and the digits array for the ITS detector.  
796
797   if(!fDigits){
798     Error("ResetDigits","fDigits is null!");
799     return;
800   }
801   for(Int_t i=0;i<fgkNdettypes;i++){
802     ResetDigits(i);
803   }
804 }
805 //___________________________________________________________________
806 void AliITSDetTypeSim::ResetDigits(Int_t branch){
807   // Reset number of digits and the digits array for this branch.
808
809   if(fDigits->At(branch)){
810     ((TClonesArray*)fDigits->At(branch))->Clear();
811   }
812   if(fNDigits) fNDigits[branch]=0;
813
814 }
815
816
817
818 //_______________________________________________________________________
819 void AliITSDetTypeSim::SDigitsToDigits(Option_t* opt, Char_t* name){
820   // Standard Summable digits to Digits function.
821   if(!GetITSgeom()){
822     Warning("SDigitsToDigits","GetITSgeom() is null!!");
823     return;
824   }
825   
826   const char *all = strstr(opt,"All");
827   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
828                         strstr(opt,"SSD")};
829   if( !det[0] && !det[1] && !det[2] ) all = "All";
830   else all = 0;
831   static Bool_t setDef = kTRUE;
832   if(setDef) SetDefaultSimulation();
833   setDef = kFALSE;
834   
835   AliITSsimulation *sim =0;
836   TTree* trees = fLoader->TreeS();
837   if( !(trees && GetSDigits()) ){
838     Error("SDigits2Digits","Error: No trees or SDigits. Returning.");
839     return;
840   } 
841   sprintf(name,"%s",name);
842   TBranch* brchSDigits = trees->GetBranch(name);
843   
844   Int_t id;
845   for(Int_t module=0;module<GetITSgeom()->GetIndexMax();module++){
846      id = GetITSgeom()->GetModuleType(module);
847     if (!all && !det[id]) continue;
848     sim = (AliITSsimulation*)GetSimulationModel(id);
849     if(!sim){
850       Error("SDigit2Digits","The simulation class was not "
851             "instanciated for module %d type %s!",module,
852             GetITSgeom()->GetModuleTypeName(module));
853       exit(1);
854     }
855     sim->InitSimulationModule(module,gAlice->GetEvNumber());
856     
857     fSDigits.Clear();
858     brchSDigits->GetEvent(module);
859     sim->AddSDigitsToModule(&fSDigits,0);
860     sim->FinishSDigitiseModule();
861     fLoader->TreeD()->Fill();
862     ResetDigits();
863   }
864   fLoader->TreeD()->GetEntries();
865   fLoader->TreeD()->AutoSave();
866   fLoader->TreeD()->Reset();
867 }
868 //_________________________________________________________
869 void AliITSDetTypeSim::AddSumDigit(AliITSpListItem &sdig){  
870   //Adds the module full of summable digits to the summable digits tree.
871
872   new(fSDigits[fNSDigits++]) AliITSpListItem(sdig);
873 }
874 //__________________________________________________________
875 void AliITSDetTypeSim::AddSimDigit(Int_t branch, const AliITSdigit* d){  
876   //    Add a simulated digit.
877
878   TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
879   switch(branch){
880   case 0:
881     new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
882     break;
883   case 1:
884     new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
885     break;
886   case 2:
887     new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
888     break;
889   }  
890 }
891 //______________________________________________________________________
892 void AliITSDetTypeSim::AddSimDigit(Int_t branch,Float_t phys,Int_t *digits,
893                                    Int_t *tracks,Int_t *hits,Float_t *charges, 
894                                    Int_t sigexpanded){
895   //   Add a simulated digit to the list.
896
897   TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
898   switch(branch){
899   case 0:
900     new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(digits,tracks,hits);
901     break;
902   case 1:
903     new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(phys,digits,tracks,
904                                                    hits,charges,sigexpanded);
905     break;
906   case 2:
907     new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(digits,tracks,hits);
908     break;
909   } 
910 }
911 //______________________________________________________________________
912 void AliITSDetTypeSim::ReadOldSSDNoise(const TObjArray *array, 
913                                        AliITSNoiseSSDv2 *noiseSSD) {
914   //Reads the old SSD calibration object and converts it to the new format
915   const Int_t fgkSSDSTRIPSPERMODULE = 1536;
916   const Int_t fgkSSDPSIDESTRIPSPERMODULE = 768;
917
918   Int_t gNMod = array->GetEntries();
919   cout<<"Converting old calibration object for noise..."<<endl;
920
921   //NOISE
922   Double_t noise = 0.0;
923   for (Int_t iModule = 0; iModule < gNMod; iModule++) {
924     AliITSNoiseSSD *noiseModule = (AliITSNoiseSSD*) (array->At(iModule));
925     for(Int_t iStrip = 0; iStrip < fgkSSDSTRIPSPERMODULE; iStrip++) {
926       noise = (iStrip < fgkSSDPSIDESTRIPSPERMODULE) ? noiseModule->GetNoiseP(iStrip) : noiseModule->GetNoiseN(1535 - iStrip);
927       if(iStrip < fgkSSDPSIDESTRIPSPERMODULE)
928         noiseSSD->AddNoiseP(iModule,iStrip,noise);
929       if(iStrip >= fgkSSDPSIDESTRIPSPERMODULE)
930         noiseSSD->AddNoiseN(iModule,1535 - iStrip,noise);
931     }//loop over strips
932   }//loop over modules      
933 }
934
935 //______________________________________________________________________
936 void AliITSDetTypeSim::ReadOldSSDBadChannels(const TObjArray *array, 
937                                              AliITSBadChannelsSSDv2 *badChannelsSSD) {
938   //Reads the old SSD calibration object and converts it to the new format
939   Int_t nMod = array->GetEntries();
940   cout<<"Converting old calibration object for bad channels..."<<endl;
941   for (Int_t iModule = 0; iModule < nMod; iModule++) {
942     //for (Int_t iModule = 0; iModule < 1; iModule++) {
943     AliITSBadChannelsSSD *bad = (AliITSBadChannelsSSD*) (array->At(iModule));
944     TArrayI arrayPSide = bad->GetBadPChannelsList();
945     for(Int_t iPCounter = 0; iPCounter < arrayPSide.GetSize(); iPCounter++) 
946       badChannelsSSD->AddBadChannelP(iModule,
947                                      iPCounter,
948                                      (Char_t)arrayPSide.At(iPCounter));
949         
950     TArrayI arrayNSide = bad->GetBadNChannelsList();
951     for(Int_t iNCounter = 0; iNCounter < arrayNSide.GetSize(); iNCounter++) 
952       badChannelsSSD->AddBadChannelN(iModule,
953                                      iNCounter,
954                                      (Char_t)arrayNSide.At(iNCounter));
955     
956   }//loop over modules      
957 }
958
959 //______________________________________________________________________
960 void AliITSDetTypeSim::ReadOldSSDGain(const TObjArray *array, 
961                                       AliITSGainSSDv2 *gainSSD) {
962   //Reads the old SSD calibration object and converts it to the new format
963
964   Int_t nMod = array->GetEntries();
965   cout<<"Converting old calibration object for gain..."<<endl;
966
967   //GAIN
968   for (Int_t iModule = 0; iModule < nMod; iModule++) {
969     AliITSGainSSD *gainModule = (AliITSGainSSD*) (array->At(iModule));
970     TArrayF arrayPSide = gainModule->GetGainP();
971     for(Int_t iPCounter = 0; iPCounter < arrayPSide.GetSize(); iPCounter++)
972       gainSSD->AddGainP(iModule,
973                         iPCounter,
974                         arrayPSide.At(iPCounter));
975     TArrayF arrayNSide = gainModule->GetGainN();
976     for(Int_t iNCounter = 0; iNCounter < arrayNSide.GetSize(); iNCounter++)
977       gainSSD->AddGainN(iModule,
978                         iNCounter,
979                         arrayNSide.At(iNCounter));
980   }//loop over modules 
981 }
982 //______________________________________________________________________
983 void AliITSDetTypeSim::ProcessSPDDigitForFastOr(UInt_t module, UInt_t colM, UInt_t rowM) {
984   // Processes wether a single fired pixel will give rise to a fast-or signal
985   fFOGenerator.ProcessPixelHitM(module,colM,rowM);
986 }
987 //_______________________________________________________________________
988 AliITSTriggerConditions* AliITSDetTypeSim::GetTriggerConditions() {
989   // Get Pixel Trigger Conditions (separate method since it is used only when simulating trigger)
990   if (fTriggerConditions==NULL) { // read from db
991     fRunNumber = ((Int_t)AliCDBManager::Instance()->GetRun());
992     Bool_t origCacheStatus = AliCDBManager::Instance()->GetCacheFlag();
993     Bool_t isCacheActive;
994     if (fRunNumber<0) isCacheActive=kFALSE;
995     else              isCacheActive=kTRUE;
996     AliCDBManager::Instance()->SetCacheFlag(isCacheActive);
997     AliCDBEntry *pitCond = AliCDBManager::Instance()->Get("TRIGGER/SPD/PITConditions", fRunNumber);
998     if (!pitCond) {
999       AliError("Trigger conditions retrieval failed! ");
1000       return NULL;
1001     }
1002     fTriggerConditions = (AliITSTriggerConditions*) pitCond->GetObject();
1003     if (!isCacheActive) pitCond->SetObject(NULL);
1004     pitCond->SetOwner(kTRUE);
1005     if (!isCacheActive) {
1006       delete pitCond;
1007     }
1008     AliCDBManager::Instance()->SetCacheFlag(origCacheStatus);
1009     if (fTriggerConditions==NULL) {
1010       AliWarning("fTriggerConditions is NULL!");
1011     }
1012   }
1013   return fTriggerConditions;
1014 }
1015 //_______________________________________________________________________
1016 void AliITSDetTypeSim::WriteFOSignals() {
1017   // write fo signals to event
1018
1019   if (!fLoader) {
1020     AliError("ITS loader is NULL.");
1021     return;
1022   }
1023
1024   if(!fLoader->TreeD()){
1025    AliError("No TreeD available");
1026    return;
1027   }
1028   TTree *tree = fLoader->TreeD();
1029   AliITSFOSignalsSPD *foSignals = new AliITSFOSignalsSPD(*GetFOSignals()); 
1030   tree->GetUserInfo()->Add(foSignals);
1031 }
1032