]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSDetTypeSim.cxx
Bug fix
[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 "AliITSresponseSDD.h"
46 #include "AliITSCalibrationSDD.h"
47 #include "AliITSMapSDD.h"
48 #include "AliITSDriftSpeedArraySDD.h"
49 #include "AliITSCalibrationSSD.h"
50 #include "AliITSNoiseSSD.h"
51 #include "AliITSGainSSD.h"
52 #include "AliITSBadChannelsSSD.h"
53 #include "AliITSCalibrationSSD.h"
54 #include "AliITSsegmentationSPD.h"
55 #include "AliITSsegmentationSDD.h"
56 #include "AliITSsegmentationSSD.h"
57 #include "AliITSsimulation.h"
58 #include "AliITSsimulationSPD.h"
59 #include "AliITSsimulationSDD.h"
60 #include "AliITSsimulationSSD.h"
61
62
63 const Int_t AliITSDetTypeSim::fgkNdettypes = 3;
64 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSPD =  240;
65 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSDD =  260;
66 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSSD = 1698;
67
68 ClassImp(AliITSDetTypeSim)
69
70 //----------------------------------------------------------------------
71 AliITSDetTypeSim::AliITSDetTypeSim():
72 TObject(),
73 fSimulation(),   // [NDet]
74 fSegmentation(), // [NDet]
75 fCalibration(),     // [NMod]
76 fPreProcess(),   // [] e.g. Fill fHitModule with hits
77 fPostProcess(),  // [] e.g. Wright Raw data
78 fNSDigits(0),    //! number of SDigits
79 fSDigits(),      //! [NMod][NSDigits]
80 fNDigits(0),     //! number of Digits
81 fRunNumber(0),   //! Run number (to access DB)
82 fDigits(),       //! [NMod][NDigits]
83 fDDLMapSDD(0),
84 fHitClassName(), // String with Hit class name.
85 fSDigClassName(),// String with SDigit class name.
86 fDigClassName(), // String with digit class name.
87 fLoader(0),      // local pointer to loader
88 fFirstcall(kTRUE){ // flag
89     // Default Constructor
90     // Inputs:
91     //    none.
92     // Outputs:
93     //    none.
94     // Return:
95     //    A properly zero-ed AliITSDetTypeSim class.
96
97   fSimulation = new TObjArray(fgkNdettypes);
98   fSegmentation = new TObjArray(fgkNdettypes);
99   fSegmentation->SetOwner(kTRUE);
100   fSDigits = new TClonesArray("AliITSpListItem",1000);
101   fDigits = new TObjArray(fgkNdettypes);
102   fNDigits = new Int_t[fgkNdettypes];
103   fDDLMapSDD=new AliITSDDLModuleMapSDD();
104   fNMod[0] = fgkDefaultNModulesSPD;
105   fNMod[1] = fgkDefaultNModulesSDD;
106   fNMod[2] = fgkDefaultNModulesSSD;
107   SetRunNumber();
108 }
109 //----------------------------------------------------------------------
110 AliITSDetTypeSim::~AliITSDetTypeSim(){
111     // Destructor
112     // Inputs:
113     //    none.
114     // Outputs:
115     //    none.
116     // Return:
117     //    Nothing.
118
119     if(fSimulation){
120         fSimulation->Delete();
121         delete fSimulation;
122     }
123     fSimulation = 0;
124     if(fSegmentation){
125         fSegmentation->Delete();
126         delete fSegmentation;
127     }
128     fSegmentation = 0;
129     if(fCalibration && fRunNumber<0){
130         AliITSresponse* rspd = ((AliITSCalibration*)fCalibration->At(
131                                GetITSgeom()->GetStartSPD()))->GetResponse();
132         AliITSresponse* rsdd = ((AliITSCalibration*)fCalibration->At(
133                                GetITSgeom()->GetStartSDD()))->GetResponse();
134         AliITSresponse* rssd = ((AliITSCalibration*)fCalibration->At(
135                                GetITSgeom()->GetStartSSD()))->GetResponse();
136         if(rspd) delete rspd;
137         if(rsdd) delete rsdd;
138         if(rssd) delete rssd;
139         fCalibration->Delete();
140         delete fCalibration;
141     }
142     fCalibration = 0;
143     if(fPreProcess){
144         fPreProcess->Delete();
145         delete fPreProcess;
146     }
147     fPreProcess = 0;
148     if(fPostProcess){
149         fPostProcess->Delete();
150         delete fPostProcess;
151     }
152     fPostProcess = 0;
153     if(fDDLMapSDD) delete fDDLMapSDD;
154     if(fNDigits) delete [] fNDigits;
155     fNDigits = 0;
156     if (fLoader)fLoader->GetModulesFolder()->Remove(this);
157     fLoader = 0; // Not deleting it.
158     if (fSDigits) {
159         fSDigits->Delete();
160         delete fSDigits;
161     }
162     fSDigits=0;
163     if (fDigits) {
164       fDigits->Delete();
165       delete fDigits;
166     }
167     fDigits=0;
168 }
169 //----------------------------------------------------------------------
170 AliITSDetTypeSim::AliITSDetTypeSim(const AliITSDetTypeSim &source) : TObject(source),
171 fSimulation(source.fSimulation),   // [NDet]
172 fSegmentation(source.fSegmentation), // [NDet]
173 fCalibration(source.fCalibration),     // [NMod]
174 fPreProcess(source.fPreProcess),   // [] e.g. Fill fHitModule with hits
175 fPostProcess(source.fPostProcess),  // [] e.g. Wright Raw data
176 fNSDigits(source.fNSDigits),    //! number of SDigits
177 fSDigits(source.fSDigits),      //! [NMod][NSDigits]
178 fNDigits(source.fNDigits),     //! number of Digits
179 fRunNumber(source.fRunNumber),   //! Run number (to access DB)
180 fDigits(source.fDigits),       //! [NMod][NDigits]
181 fDDLMapSDD(source.fDDLMapSDD),
182 fHitClassName(source.fHitClassName), // String with Hit class name.
183 fSDigClassName(source.fSDigClassName),// String with SDigit class name.
184 fDigClassName(), // String with digit class name.
185 fLoader(source.fLoader),      // local pointer to loader
186 fFirstcall(source.fFirstcall)
187 {
188     // Copy Constructor for object AliITSDetTypeSim not allowed
189   for(Int_t i=0;i<fgkNdettypes;i++){
190     fDigClassName[i] = source.fDigClassName[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){
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){
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(GetITSgeom());
282     }else if(idet==1){
283         seg = new AliITSsegmentationSDD(GetITSgeom());
284     }else {
285         seg = new AliITSsegmentationSSD(GetITSgeom());
286     }
287     SetSegmentationModel(idet,seg);
288 }
289 //______________________________________________________________________
290 void AliITSDetTypeSim::SetSegmentationModel(Int_t dettype,
291                                             AliITSsegmentation *seg){
292    
293   //Set segmentation model for detector type
294   if(fSegmentation==0x0){
295     fSegmentation = new TObjArray(fgkNdettypes);
296     fSegmentation->SetOwner(kTRUE);
297   }
298   fSegmentation->AddAt(seg,dettype);
299 }
300 //______________________________________________________________________
301 AliITSsegmentation* AliITSDetTypeSim::GetSegmentationModel(Int_t dettype){
302   //Get segmentation model for detector type
303    
304    if(fSegmentation==0) {
305        Warning("GetSegmentationModel","fSegmentation is 0!");
306        return 0; 
307    } 
308    return (AliITSsegmentation*)(fSegmentation->At(dettype));
309 }
310 //_______________________________________________________________________
311 AliITSsegmentation* AliITSDetTypeSim::GetSegmentationModelByModule(Int_t module){
312     //Get segmentation model by module number
313     if(GetITSgeom()==0){
314         Warning("GetSegmentationModelByModule","GetITSgeom() is 0!");
315         return 0;
316     }     
317     return GetSegmentationModel(GetITSgeom()->GetModuleType(module));
318 }
319 //_______________________________________________________________________
320 void AliITSDetTypeSim::CreateCalibrationArray() {
321     //Create the container of calibration functions with correct size
322     if (fCalibration) {
323         Warning("CreateCalibration","pointer to calibration object exists\n");
324         fCalibration->Delete();
325         delete fCalibration;
326     }
327
328     Int_t nModTot = GetITSgeom()->GetIndexMax();
329     fCalibration = new TObjArray(nModTot);
330     fCalibration->SetOwner(kTRUE);
331     fCalibration->Clear();
332 }
333 //_______________________________________________________________________
334 void AliITSDetTypeSim::SetCalibrationModel(Int_t iMod, AliITSCalibration *resp){
335     //Set response model for modules
336
337     if (fCalibration==0) CreateCalibrationArray();
338  
339     if (fCalibration->At(iMod)!=0)
340         delete (AliITSCalibration*) fCalibration->At(iMod);
341     fCalibration->AddAt(resp, iMod);
342 }
343 //______________________________________________________________________
344 void AliITSDetTypeSim::ResetCalibrationArray(){
345     //resets response array
346     if(fCalibration && fRunNumber<0){  // if fRunNumber<0 fCalibration is owner
347         AliITSresponse* rspd = ((AliITSCalibration*)fCalibration->At(
348                                 GetITSgeom()->GetStartSPD()))->GetResponse();
349         AliITSresponse* rsdd = ((AliITSCalibration*)fCalibration->At(
350                                 GetITSgeom()->GetStartSDD()))->GetResponse();
351         AliITSresponse* rssd = ((AliITSCalibration*)fCalibration->At(
352                                 GetITSgeom()->GetStartSSD()))->GetResponse();
353         if(rspd) delete rspd;
354         if(rsdd) delete rsdd;
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   return (AliITSCalibration*)(fCalibration->At(iMod));
375 }
376 //_______________________________________________________________________
377 void AliITSDetTypeSim::SetDefaults(){
378     //Set defaults for segmentation and response
379
380     if(GetITSgeom()==0){
381         Warning("SetDefaults","GetITSgeom() is 0!");
382         return;
383     } // end if
384     if (fCalibration==0) {
385         CreateCalibrationArray();
386     } // end if
387
388     ResetSegmentation();
389     if(!GetCalibration()){AliFatal("Exit"); exit(0);}
390
391     for(Int_t idet=0;idet<fgkNdettypes;idet++){
392         //SPD
393         if(idet==0){
394             if(!GetSegmentationModel(idet)) SetDefaultSegmentation(idet);
395             const char *kData0=(GetCalibrationModel(GetITSgeom()->GetStartSPD()))->DataType();
396             if (strstr(kData0,"real")) {
397                 SetDigitClassName(idet,"AliITSdigit");
398             }else {
399                 SetDigitClassName(idet,"AliITSdigitSPD");
400             } // end if
401         } // end if idet==0
402         //SDD
403         if(idet==1){
404             if(!GetSegmentationModel(idet)) SetDefaultSegmentation(idet);
405             AliITSCalibrationSDD* rsp = 
406                 (AliITSCalibrationSDD*)GetCalibrationModel(
407                     GetITSgeom()->GetStartSDD());
408             const char *kopt = ((AliITSresponseSDD*)rsp->GetResponse())->
409                 ZeroSuppOption();
410             if((!strstr(kopt,"2D"))&&(!strstr(kopt,"1D"))) {
411                 SetDigitClassName(idet,"AliITSdigit");
412             }else {
413                 SetDigitClassName(idet,"AliITSdigitSDD");
414             } // end if
415         } // end if idet==1
416         //SSD
417         if(idet==2){
418             if(!GetSegmentationModel(idet))SetDefaultSegmentation(idet);
419             const char *kData2 = (GetCalibrationModel(
420                                   GetITSgeom()->GetStartSSD())->DataType());
421             if (strstr(kData2,"real")) {
422                 SetDigitClassName(idet,"AliITSdigit");
423             }else {
424                 SetDigitClassName(idet,"AliITSdigitSSD");
425             } // end if
426         } // end if idet==2
427     }// end for idet
428 }
429 //______________________________________________________________________
430 Bool_t AliITSDetTypeSim::GetCalibration() {
431   // Get Default calibration if a storage is not defined.
432
433   if(!fFirstcall){
434     AliITSCalibration* cal = GetCalibrationModel(0);
435     if(cal)return kTRUE;
436   }
437   else {
438     fFirstcall = kFALSE;
439   }
440
441   SetRunNumber((Int_t)AliCDBManager::Instance()->GetRun());
442   Int_t run=GetRunNumber();
443
444   Bool_t origCacheStatus = AliCDBManager::Instance()->GetCacheFlag();
445   Bool_t isCacheActive = kTRUE;
446   if(GetRunNumber()<0){
447     isCacheActive=kFALSE;
448     fCalibration->SetOwner(kTRUE);
449   }
450   else{
451     fCalibration->SetOwner(kFALSE);
452   }
453
454   AliCDBManager::Instance()->SetCacheFlag(isCacheActive);
455
456   AliCDBEntry *entrySPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDDead", run);
457   AliCDBEntry *entrySDD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSDD", run);
458   AliCDBEntry *drSpSDD = AliCDBManager::Instance()->Get("ITS/Calib/DriftSpeedSDD",run);
459   AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD",run);
460   AliCDBEntry *mapASDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsAnodeSDD",run);
461   AliCDBEntry *mapTSDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsTimeSDD",run);
462   // AliCDBEntry *entrySSD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSSD", run);
463   AliCDBEntry *entryNoiseSSD = AliCDBManager::Instance()->Get("ITS/Calib/NoiseSSD");
464   AliCDBEntry *entryGainSSD = AliCDBManager::Instance()->Get("ITS/Calib/GainSSD");
465   AliCDBEntry *entryBadChannelsSSD = AliCDBManager::Instance()->Get("ITS/Calib/BadChannelsSSD");
466
467   AliCDBEntry *entry2SPD = AliCDBManager::Instance()->Get("ITS/Calib/RespSPD", run);
468   AliCDBEntry *entry2SDD = AliCDBManager::Instance()->Get("ITS/Calib/RespSDD", run);
469   AliCDBEntry *entry2SSD = AliCDBManager::Instance()->Get("ITS/Calib/RespSSD", run);
470
471   if(!entrySPD || !entrySDD || !entryNoiseSSD || !entryGainSSD || !entryBadChannelsSSD || 
472      !entry2SPD || !entry2SDD || !entry2SSD || !drSpSDD || !ddlMapSDD || !mapASDD ||!mapTSDD){
473     AliFatal("Calibration object retrieval failed! ");
474     return kFALSE;
475   }     
476
477 //   if(!entrySPD || !entrySDD || !entryNoiseSSD || !entryGainSSD || !entryBadChannelsSSD || 
478 //      !entry2SPD || !entry2SDD || !entry2SSD){
479 //     AliFatal("Calibration object retrieval failed! ");
480 //     return kFALSE;
481 //   }          
482
483   TObjArray *calSPD = (TObjArray *)entrySPD->GetObject();
484   if(!isCacheActive)entrySPD->SetObject(NULL);
485   entrySPD->SetOwner(kTRUE);
486
487   AliITSresponseSPD *pSPD = (AliITSresponseSPD*)entry2SPD->GetObject();
488   if(!isCacheActive)entry2SPD->SetObject(NULL);
489   entry2SPD->SetOwner(kTRUE);
490    
491   TObjArray *calSDD = (TObjArray *)entrySDD->GetObject();
492   if(!isCacheActive)entrySDD->SetObject(NULL);
493   entrySDD->SetOwner(kTRUE);
494
495   AliITSresponseSDD *pSDD = (AliITSresponseSDD*)entry2SDD->GetObject();
496   if(!isCacheActive)entry2SDD->SetObject(NULL);
497   entry2SDD->SetOwner(kTRUE);
498
499   TObjArray *drSp = (TObjArray *)drSpSDD->GetObject();
500   if(!isCacheActive)drSpSDD->SetObject(NULL);
501   drSpSDD->SetOwner(kTRUE);
502
503   AliITSDDLModuleMapSDD *ddlsdd=(AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject();
504   if(!isCacheActive)ddlMapSDD->SetObject(NULL);
505   ddlMapSDD->SetOwner(kTRUE);
506
507   TObjArray *mapAn = (TObjArray *)mapASDD->GetObject();
508   if(!isCacheActive)mapASDD->SetObject(NULL);
509   mapASDD->SetOwner(kTRUE);
510
511   TObjArray *mapT = (TObjArray *)mapTSDD->GetObject();
512   if(!isCacheActive)mapTSDD->SetObject(NULL);
513   mapTSDD->SetOwner(kTRUE);
514
515   /*
516   TObjArray *calSSD = (TObjArray *)entrySSD->GetObject();
517   if(!isCacheActive)entrySSD->SetObject(NULL);
518   entrySSD->SetOwner(kTRUE);
519   */
520
521   TObjArray *noiseSSD = (TObjArray *)entryNoiseSSD->GetObject();
522   if(!isCacheActive)entryNoiseSSD->SetObject(NULL);
523   entryNoiseSSD->SetOwner(kTRUE);
524
525   TObjArray *gainSSD = (TObjArray *)entryGainSSD->GetObject();
526   if(!isCacheActive)entryGainSSD->SetObject(NULL);
527   entryGainSSD->SetOwner(kTRUE);
528
529   TObjArray *badchannelsSSD = (TObjArray *)entryBadChannelsSSD->GetObject();
530   if(!isCacheActive)entryBadChannelsSSD->SetObject(NULL);
531   entryBadChannelsSSD->SetOwner(kTRUE);
532
533   AliITSresponseSSD *pSSD = (AliITSresponseSSD*)entry2SSD->GetObject();
534   if(!isCacheActive)entry2SSD->SetObject(NULL);
535   entry2SSD->SetOwner(kTRUE);
536   
537   // DB entries are deleted. In this way metadeta objects are deleted as well
538   if(!isCacheActive){
539     delete entrySPD;
540     delete entrySDD;
541     delete entryNoiseSSD;
542     delete entryGainSSD;
543     delete entryBadChannelsSSD;
544     delete entry2SPD;
545     delete entry2SDD;
546     delete entry2SSD;
547     delete mapASDD;   
548     delete mapTSDD;
549     delete drSpSDD;
550     delete ddlMapSDD;
551   }
552   
553   AliCDBManager::Instance()->SetCacheFlag(origCacheStatus);
554
555   if ((!pSPD)||(!pSDD)||(!pSSD) || (!calSPD) || (!calSDD) || (!drSp) || (!ddlsdd)
556       || (!mapAn) || (!mapT) || (!noiseSSD)|| (!gainSSD)|| (!badchannelsSSD)) {
557     AliWarning("Can not get calibration from calibration database !");
558     return kFALSE;
559   }
560
561   fNMod[0] = calSPD->GetEntries();
562   fNMod[1] = calSDD->GetEntries();
563   fNMod[2] = noiseSSD->GetEntries();
564   AliInfo(Form("%i SPD, %i SDD and %i SSD in calibration database",
565                fNMod[0], fNMod[1], fNMod[2]));
566   AliITSCalibration* cal;
567   for (Int_t i=0; i<fNMod[0]; i++) {
568     cal = (AliITSCalibration*) calSPD->At(i);
569     cal->SetResponse(pSPD);
570     SetCalibrationModel(i, cal);
571   }
572   for (Int_t i=0; i<fNMod[1]; i++) {
573     cal = (AliITSCalibration*) calSDD->At(i);
574     cal->SetResponse(pSDD);
575     Int_t i0=2*i;
576     Int_t i1=1+2*i;
577     AliITSDriftSpeedArraySDD* arr0 = (AliITSDriftSpeedArraySDD*) drSp->At(i0);
578     AliITSMapSDD* ma0 = (AliITSMapSDD*)mapAn->At(i0);
579     AliITSMapSDD* mt0 = (AliITSMapSDD*)mapT->At(i0);
580     AliITSDriftSpeedArraySDD* arr1 = (AliITSDriftSpeedArraySDD*) drSp->At(i1);
581     AliITSMapSDD* ma1 = (AliITSMapSDD*)mapAn->At(i1);
582     AliITSMapSDD* mt1 = (AliITSMapSDD*)mapT->At(i1);
583     cal->SetDriftSpeed(0,arr0);
584     cal->SetDriftSpeed(1,arr1);
585     cal->SetMapA(0,ma0);
586     cal->SetMapA(1,ma1);
587     cal->SetMapT(0,mt0);
588     cal->SetMapT(1,mt1);
589     fDDLMapSDD->SetDDLMap(ddlsdd);
590     Int_t iMod = i + fNMod[0];
591     SetCalibrationModel(iMod, cal);
592  }
593   for (Int_t i=0; i<fNMod[2]; i++) {
594     AliITSCalibrationSSD *calibSSD = new AliITSCalibrationSSD();
595     calibSSD->SetResponse((AliITSresponse*)pSSD);
596     
597     AliITSNoiseSSD *noise = (AliITSNoiseSSD*) (noiseSSD->At(i));
598     calibSSD->SetNoise(noise);
599     AliITSGainSSD *gain = (AliITSGainSSD*) (gainSSD->At(i));
600     calibSSD->SetGain(gain);
601     AliITSBadChannelsSSD *bad = (AliITSBadChannelsSSD*) (badchannelsSSD->At(i));
602     calibSSD->SetBadChannels(bad);
603
604     Int_t iMod = i + fNMod[0] + fNMod[1];
605     SetCalibrationModel(iMod, calibSSD);
606
607  }
608   return kTRUE;
609 }
610 //_______________________________________________________________________
611 void AliITSDetTypeSim::SetDefaultSimulation(){
612   //Set default simulation for detector type
613
614   if(GetITSgeom()==0){
615     Warning("SetDefaultSimulation","GetITSgeom() is 0!");
616     return;
617   }
618   if(fCalibration==0){
619     Warning("SetDefaultSimulation","fCalibration is 0!");
620     return;
621   }
622   if(fSegmentation==0){
623     Warning("SetDefaultSimulation","fSegmentation is 0!");
624     for(Int_t i=0;i<fgkNdettypes;i++) SetDefaultSegmentation(i);
625   }else for(Int_t i=0;i<fgkNdettypes;i++) if(!GetSegmentationModel(i)){
626       Warning("SetDefaultSimulation",
627               "Segmentation not defined for det %d - Default taken\n!",i);
628       SetDefaultSegmentation(i);
629   }
630   AliITSsimulation* sim;
631
632   for(Int_t idet=0;idet<fgkNdettypes;idet++){
633    //SPD
634     if(idet==0){
635       sim = GetSimulationModel(idet); 
636       if(!sim){
637         sim = new AliITSsimulationSPD(this);
638         SetSimulationModel(idet,sim);
639       }
640     }
641     //SDD
642     if(idet==1){
643       sim = GetSimulationModel(idet);
644       if(!sim){
645         sim = new AliITSsimulationSDD(this);
646         SetSimulationModel(idet,sim);
647       }      
648     }
649     //SSD
650     if(idet==2){
651       sim = GetSimulationModel(idet);
652       if(!sim){
653         sim = new AliITSsimulationSSD(this);
654         SetSimulationModel(idet,sim);
655       }
656     }
657
658   }
659 }
660 //___________________________________________________________________
661 void AliITSDetTypeSim::SetTreeAddressS(TTree* treeS, Char_t* name){
662   // Set branch address for the ITS summable digits Trees.  
663   char branchname[30];
664
665   if(!treeS){
666     return;
667   }
668   if (fSDigits ==  0x0){
669     fSDigits = new TClonesArray("AliITSpListItem",1000);
670   }
671   TBranch *branch;
672   sprintf(branchname,"%s",name);
673   branch = treeS->GetBranch(branchname);
674   if (branch) branch->SetAddress(&fSDigits);
675
676 }
677 //___________________________________________________________________
678 void AliITSDetTypeSim::SetTreeAddressD(TTree* treeD, Char_t* name){
679   // Set branch address for the digit Trees.
680   
681   const char *det[3] = {"SPD","SDD","SSD"};
682   TBranch *branch;
683   
684   char branchname[30];
685   
686   if(!treeD){
687     return;
688   }
689   if(!fDigits){
690     fDigits = new TObjArray(fgkNdettypes); 
691   }
692   for(Int_t i=0;i<fgkNdettypes;i++){
693     Char_t* digclass = GetDigitClassName(i);
694     if(digclass==0x0){
695       if(i==0) SetDigitClassName(i,"AliITSdigitSPD");
696       if(i==1) SetDigitClassName(i,"AliITSdigitSDD");
697       if(i==2) SetDigitClassName(i,"AliITSdigitSSD");
698       digclass = GetDigitClassName(i);
699     }
700     TString classn = digclass;
701     if(!(fDigits->At(i))){
702       fDigits->AddAt(new TClonesArray(classn.Data(),1000),i);
703     }else{
704       ResetDigits(i);
705     }
706     
707     if(fgkNdettypes==3) sprintf(branchname,"%sDigits%s",name,det[i]);
708     else sprintf(branchname,"%sDigits%d",name,i+1);
709     if(fDigits){
710       branch = treeD->GetBranch(branchname);
711       if(branch) branch->SetAddress(&((*fDigits)[i]));
712     }
713   }
714
715 }
716 //___________________________________________________________________
717 void AliITSDetTypeSim::ResetDigits(){
718   // Reset number of digits and the digits array for the ITS detector.  
719
720   if(!fDigits){
721     Error("ResetDigits","fDigits is null!");
722     return;
723   }
724   for(Int_t i=0;i<fgkNdettypes;i++){
725     ResetDigits(i);
726   }
727 }
728 //___________________________________________________________________
729 void AliITSDetTypeSim::ResetDigits(Int_t branch){
730   // Reset number of digits and the digits array for this branch.
731
732   if(fDigits->At(branch)){
733     ((TClonesArray*)fDigits->At(branch))->Clear();
734   }
735   if(fNDigits) fNDigits[branch]=0;
736
737 }
738
739
740
741 //_______________________________________________________________________
742 void AliITSDetTypeSim::SDigitsToDigits(Option_t* opt, Char_t* name){
743   // Standard Summable digits to Digits function.
744   if(!GetITSgeom()){
745     Warning("SDigitsToDigits","GetITSgeom() is null!!");
746     return;
747   }
748   
749   const char *all = strstr(opt,"All");
750   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
751                         strstr(opt,"SSD")};
752   if( !det[0] && !det[1] && !det[2] ) all = "All";
753   else all = 0;
754   static Bool_t setDef = kTRUE;
755   if(setDef) SetDefaultSimulation();
756   setDef = kFALSE;
757   
758   AliITSsimulation *sim =0;
759   TTree* trees = fLoader->TreeS();
760   if( !(trees && GetSDigits()) ){
761     Error("SDigits2Digits","Error: No trees or SDigits. Returning.");
762     return;
763   } 
764   sprintf(name,"%s",name);
765   TBranch* brchSDigits = trees->GetBranch(name);
766   
767   Int_t id;
768   for(Int_t module=0;module<GetITSgeom()->GetIndexMax();module++){
769      id = GetITSgeom()->GetModuleType(module);
770     if (!all && !det[id]) continue;
771     sim = (AliITSsimulation*)GetSimulationModel(id);
772     if(!sim){
773       Error("SDigit2Digits","The simulation class was not "
774             "instanciated for module %d type %s!",module,
775             GetITSgeom()->GetModuleTypeName(module));
776       exit(1);
777     }
778     sim->InitSimulationModule(module,gAlice->GetEvNumber());
779     
780     fSDigits->Clear();
781     brchSDigits->GetEvent(module);
782     sim->AddSDigitsToModule(fSDigits,0);
783     sim->FinishSDigitiseModule();
784     fLoader->TreeD()->Fill();
785     ResetDigits();
786   }
787   fLoader->TreeD()->GetEntries();
788   fLoader->TreeD()->AutoSave();
789   fLoader->TreeD()->Reset();
790 }
791 //_________________________________________________________
792 void AliITSDetTypeSim::AddSumDigit(AliITSpListItem &sdig){  
793   //Adds the module full of summable digits to the summable digits tree.
794
795   TClonesArray &lsdig = *fSDigits;
796   new(lsdig[fNSDigits++]) AliITSpListItem(sdig);
797 }
798 //__________________________________________________________
799 void AliITSDetTypeSim::AddRealDigit(Int_t branch, Int_t *digits){
800   //   Add a real digit - as coming from data.
801
802   TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
803   new(ldigits[fNDigits[branch]++]) AliITSdigit(digits); 
804 }
805 //__________________________________________________________
806 void AliITSDetTypeSim::AddSimDigit(Int_t branch, AliITSdigit* d){  
807   //    Add a simulated digit.
808
809   TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
810   switch(branch){
811   case 0:
812     new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
813     break;
814   case 1:
815     new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
816     break;
817   case 2:
818     new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
819     break;
820   }  
821 }
822 //______________________________________________________________________
823 void AliITSDetTypeSim::AddSimDigit(Int_t branch,Float_t phys,Int_t *digits,
824                                    Int_t *tracks,Int_t *hits,Float_t *charges){
825   //   Add a simulated digit to the list.
826
827   TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
828   AliITSCalibrationSDD *resp = 0;
829   switch(branch){
830   case 0:
831     new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(digits,tracks,hits);
832     break;
833   case 1:
834     resp = (AliITSCalibrationSDD*)GetCalibrationModel(
835         GetITSgeom()->GetStartSDD());
836     new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(phys,digits,tracks,
837                                                    hits,charges,resp);
838     break;
839   case 2:
840     new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(digits,tracks,hits);
841     break;
842   } 
843 }
844 //______________________________________________________________________
845 void AliITSDetTypeSim::StoreCalibration(Int_t firstRun, Int_t lastRun,
846                                         AliCDBMetaData &md) {
847   // Store calibration in the calibration database
848   // The database must be created in an external piece of code (i.e. 
849   // a configuration macro )
850
851   if(!AliCDBManager::Instance()->IsDefaultStorageSet()) {
852     AliWarning("No storage set! Will use dummy one");
853     AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT");
854   }
855
856   if (!fCalibration) {
857     AliError("AliITSCalibration classes are not defined - nothing done");
858     return;
859   }
860   AliCDBId idRespSPD("ITS/Calib/SPDDead",firstRun, lastRun);
861   AliCDBId idRespSDD("ITS/Calib/CalibSDD",firstRun, lastRun);
862   AliCDBId idRespSSD("ITS/Calib/CalibSSD",firstRun, lastRun);
863
864   TObjArray respSPD(fNMod[0]);
865   TObjArray respSDD(fNMod[1]-fNMod[0]);
866   TObjArray respSSD(fNMod[2]-fNMod[1]);
867   respSPD.SetOwner(kFALSE);
868   respSSD.SetOwner(kFALSE);
869   respSSD.SetOwner(kFALSE);
870
871   Int_t index[fgkNdettypes];
872   for (Int_t i = 0; i<fgkNdettypes; i++ ) {
873     index[i] = 0;
874     for (Int_t j = 0; j<=i; j++ )
875       index[i]+=fNMod[j];
876   }
877
878   for (Int_t i = 0; i<index[0]; i++ )
879     respSPD.Add(fCalibration->At(i));
880
881   for (Int_t i = index[0]; i<index[1]; i++ )
882     respSDD.Add(fCalibration->At(i));
883
884   for (Int_t i = index[1]; i<index[2]; i++ )
885     respSSD.Add(fCalibration->At(i));
886
887   AliCDBManager::Instance()->Put(&respSPD, idRespSPD, &md);
888   AliCDBManager::Instance()->Put(&respSDD, idRespSDD, &md);
889   AliCDBManager::Instance()->Put(&respSSD, idRespSSD, &md);
890 }
891
892