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