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