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