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