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