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