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