]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSDetTypeSim.cxx
update (alberto):
[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 "AliITSCalibrationSSD.h"
48 #include "AliITSsegmentationSPD.h"
49 #include "AliITSsegmentationSDD.h"
50 #include "AliITSsegmentationSSD.h"
51 #include "AliITSsimulation.h"
52 #include "AliITSsimulationSPD.h"
53 #include "AliITSsimulationSDD.h"
54 #include "AliITSsimulationSSD.h"
55
56
57 const Int_t AliITSDetTypeSim::fgkNdettypes = 3;
58 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSPD =  240;
59 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSDD =  260;
60 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSSD = 1698;
61
62 ClassImp(AliITSDetTypeSim)
63
64 //----------------------------------------------------------------------
65 AliITSDetTypeSim::AliITSDetTypeSim():
66 TObject(),
67 fSimulation(),   // [NDet]
68 fSegmentation(), // [NDet]
69 fCalibration(),     // [NMod]
70 fPreProcess(),   // [] e.g. Fill fHitModule with hits
71 fPostProcess(),  // [] e.g. Wright Raw data
72 fNSDigits(0),    //! number of SDigits
73 fSDigits(),      //! [NMod][NSDigits]
74 fNDigits(0),     //! number of Digits
75 fRunNumber(0),   //! Run number (to access DB)
76 fDigits(),       //! [NMod][NDigits]
77 fHitClassName(), // String with Hit class name.
78 fSDigClassName(),// String with SDigit class name.
79 fDigClassName(), // String with digit class name.
80 fLoader(0),      // local pointer to loader
81 fFirstcall(kTRUE){ // flag
82     // Default Constructor
83     // Inputs:
84     //    none.
85     // Outputs:
86     //    none.
87     // Return:
88     //    A properly zero-ed AliITSDetTypeSim class.
89
90   fSimulation = new TObjArray(fgkNdettypes);
91   fSegmentation = new TObjArray(fgkNdettypes);
92   fSegmentation->SetOwner(kTRUE);
93   fSDigits = new TClonesArray("AliITSpListItem",1000);
94   fDigits = new TObjArray(fgkNdettypes);
95   fNDigits = new Int_t[fgkNdettypes];
96   fNMod[0] = fgkDefaultNModulesSPD;
97   fNMod[1] = fgkDefaultNModulesSDD;
98   fNMod[2] = fgkDefaultNModulesSSD;
99   SetRunNumber();
100 }
101 //----------------------------------------------------------------------
102 AliITSDetTypeSim::~AliITSDetTypeSim(){
103     // Destructor
104     // Inputs:
105     //    none.
106     // Outputs:
107     //    none.
108     // Return:
109     //    Nothing.
110
111     if(fSimulation){
112         fSimulation->Delete();
113         delete fSimulation;
114     }
115     fSimulation = 0;
116     if(fSegmentation){
117         fSegmentation->Delete();
118         delete fSegmentation;
119     }
120     fSegmentation = 0;
121     if(fCalibration && fRunNumber<0){
122         AliITSresponse* rspd = ((AliITSCalibration*)fCalibration->At(
123                                GetITSgeom()->GetStartSPD()))->GetResponse();
124         AliITSresponse* rsdd = ((AliITSCalibration*)fCalibration->At(
125                                GetITSgeom()->GetStartSDD()))->GetResponse();
126         AliITSresponse* rssd = ((AliITSCalibration*)fCalibration->At(
127                                GetITSgeom()->GetStartSSD()))->GetResponse();
128         if(rspd) delete rspd;
129         if(rsdd) delete rsdd;
130         if(rssd) delete rssd;
131         fCalibration->Delete();
132         delete fCalibration;
133     }
134     fCalibration = 0;
135     if(fPreProcess){
136         fPreProcess->Delete();
137         delete fPreProcess;
138     }
139     fPreProcess = 0;
140     if(fPostProcess){
141         fPostProcess->Delete();
142         delete fPostProcess;
143     }
144     fPostProcess = 0;
145     if(fNDigits) delete [] fNDigits;
146     fNDigits = 0;
147     if (fLoader)fLoader->GetModulesFolder()->Remove(this);
148     fLoader = 0; // Not deleting it.
149     if (fSDigits) {
150         fSDigits->Delete();
151         delete fSDigits;
152     }
153     fSDigits=0;
154     if (fDigits) {
155       fDigits->Delete();
156       delete fDigits;
157     }
158     fDigits=0;
159 }
160 //----------------------------------------------------------------------
161 AliITSDetTypeSim::AliITSDetTypeSim(const AliITSDetTypeSim &source) : TObject(source){
162     // Copy Constructor for object AliITSDetTypeSim not allowed
163
164     if(this==&source) return;
165     Error("Copy constructor",
166           "You are not allowed to make a copy of the AliITSDetTypeSim");
167     exit(1);
168 }
169 //----------------------------------------------------------------------
170 AliITSDetTypeSim& AliITSDetTypeSim::operator=(const AliITSDetTypeSim &source){
171     // The = operator for object AliITSDetTypeSim
172  
173     if(&source==this) return *this;
174     Error("operator=",
175           "You are not allowed to make a copy of the AliITSDetTypeSIm");
176     exit(1);
177     return *this;
178 }
179
180 //______________________________________________________________________
181 void AliITSDetTypeSim::SetITSgeom(AliITSgeom *geom){
182     // Sets/replaces the existing AliITSgeom object kept in AliITSLoader
183     // 
184     // Inputs:
185     //   AliITSgoem   *geom  The AliITSgeom object to be used.
186     // Output:
187     //   none.
188     // Return:
189     //   none.
190   if(!fLoader){
191     Error("SetITSgeom","No pointer to loader - nothing done");
192     return;
193   }
194   else {
195     fLoader->SetITSgeom(geom);  // protections in AliITSLoader::SetITSgeom
196   }
197  
198 }
199 //______________________________________________________________________
200 void AliITSDetTypeSim::SetLoader(AliITSLoader *loader){
201     // Sets the local copy of the AliITSLoader, and passes on the
202     // AliITSgeom object as needed.
203     // Inputs
204     //   AliITSLoader  *loader pointer to AliITSLoader for local use
205     // Outputs:
206     //   none.
207     // Return:
208     //  none.
209
210     if(fLoader==loader) return; // Same do nothing
211     if(fLoader){ // alread have an existing loader
212         Error("SetLoader",
213                 "Already have an exisiting loader ptr=%p Nothing done",
214                 fLoader);
215     } // end if
216     fLoader = loader;
217 }
218 //______________________________________________________________________
219 void AliITSDetTypeSim::SetSimulationModel(Int_t dettype,AliITSsimulation *sim){
220
221   //Set simulation model for detector type
222
223   if(fSimulation==0) fSimulation = new TObjArray(fgkNdettypes);
224   fSimulation->AddAt(sim,dettype);
225 }
226 //______________________________________________________________________
227 AliITSsimulation* AliITSDetTypeSim::GetSimulationModel(Int_t dettype){
228
229   //Get simulation model for detector type
230   if(fSimulation==0)  {
231     Warning("GetSimulationModel","fSimulation is 0!");
232     return 0;     
233   }
234   return (AliITSsimulation*)(fSimulation->At(dettype));
235 }
236 //______________________________________________________________________
237 AliITSsimulation* AliITSDetTypeSim::GetSimulationModelByModule(Int_t module){
238
239   //Get simulation model by module number
240   if(GetITSgeom()==0) {
241     Warning("GetSimulationModelByModule","GetITSgeom() is 0!");
242     return 0;
243   }
244   
245   return GetSimulationModel(GetITSgeom()->GetModuleType(module));
246 }
247 //_______________________________________________________________________
248 void AliITSDetTypeSim::SetDefaultSegmentation(Int_t idet){
249     // Set default segmentation model objects
250     AliITSsegmentation *seg;
251
252     if(fSegmentation==0x0){
253         fSegmentation = new TObjArray(fgkNdettypes);
254         fSegmentation->SetOwner(kTRUE);
255     }
256     if(GetSegmentationModel(idet))
257         delete (AliITSsegmentation*)fSegmentation->At(idet);
258     if(idet==0){
259         seg = new AliITSsegmentationSPD(GetITSgeom());
260     }else if(idet==1){
261         AliITSCalibration* res=GetCalibrationModel(
262             GetITSgeom()->GetStartSDD());
263         seg = new AliITSsegmentationSDD(GetITSgeom(),res);
264     }else {
265         seg = new AliITSsegmentationSSD(GetITSgeom());
266     }
267     SetSegmentationModel(idet,seg);
268 }
269 //______________________________________________________________________
270 void AliITSDetTypeSim::SetSegmentationModel(Int_t dettype,
271                                             AliITSsegmentation *seg){
272    
273   //Set segmentation model for detector type
274   if(fSegmentation==0x0){
275     fSegmentation = new TObjArray(fgkNdettypes);
276     fSegmentation->SetOwner(kTRUE);
277   }
278   fSegmentation->AddAt(seg,dettype);
279 }
280 //______________________________________________________________________
281 AliITSsegmentation* AliITSDetTypeSim::GetSegmentationModel(Int_t dettype){
282   //Get segmentation model for detector type
283    
284    if(fSegmentation==0) {
285        Warning("GetSegmentationModel","fSegmentation is 0!");
286        return 0; 
287    } 
288    return (AliITSsegmentation*)(fSegmentation->At(dettype));
289 }
290 //_______________________________________________________________________
291 AliITSsegmentation* AliITSDetTypeSim::GetSegmentationModelByModule(Int_t module){
292     //Get segmentation model by module number
293     if(GetITSgeom()==0){
294         Warning("GetSegmentationModelByModule","GetITSgeom() is 0!");
295         return 0;
296     }     
297     return GetSegmentationModel(GetITSgeom()->GetModuleType(module));
298 }
299 //_______________________________________________________________________
300 void AliITSDetTypeSim::CreateCalibrationArray() {
301     //Create the container of calibration functions with correct size
302     if (fCalibration) {
303         Warning("CreateCalibration","pointer to calibration object exists\n");
304         fCalibration->Delete();
305         delete fCalibration;
306     }
307
308     Int_t nModTot = GetITSgeom()->GetIndexMax();
309     fCalibration = new TObjArray(nModTot);
310     fCalibration->SetOwner(kTRUE);
311     fCalibration->Clear();
312 }
313 //_______________________________________________________________________
314 void AliITSDetTypeSim::SetCalibrationModel(Int_t iMod, AliITSCalibration *resp){
315     //Set response model for modules
316
317     if (fCalibration==0) CreateCalibrationArray();
318  
319     if (fCalibration->At(iMod)!=0)
320         delete (AliITSCalibration*) fCalibration->At(iMod);
321     fCalibration->AddAt(resp, iMod);
322 }
323 //______________________________________________________________________
324 void AliITSDetTypeSim::ResetCalibrationArray(){
325     //resets response array
326     if(fCalibration && fRunNumber<0){  // if fRunNumber<0 fCalibration is owner
327         AliITSresponse* rspd = ((AliITSCalibration*)fCalibration->At(
328                                 GetITSgeom()->GetStartSPD()))->GetResponse();
329         AliITSresponse* rsdd = ((AliITSCalibration*)fCalibration->At(
330                                 GetITSgeom()->GetStartSDD()))->GetResponse();
331         AliITSresponse* rssd = ((AliITSCalibration*)fCalibration->At(
332                                 GetITSgeom()->GetStartSSD()))->GetResponse();
333         if(rspd) delete rspd;
334         if(rsdd) delete rsdd;
335         if(rssd) delete rssd;
336         fCalibration->Clear();
337     }else if (fCalibration && fRunNumber>=0){
338         fCalibration->Clear();
339     }
340 }
341 //______________________________________________________________________
342 void AliITSDetTypeSim::ResetSegmentation(){
343     //Resets segmentation array
344     if(fSegmentation) fSegmentation->Clear();
345 }
346 //_______________________________________________________________________
347 AliITSCalibration* AliITSDetTypeSim::GetCalibrationModel(Int_t iMod){
348     //Get response model for module number iMod 
349  
350     if(fCalibration==0) {
351         AliError("fCalibration is 0!");
352         return 0; 
353     }
354   return (AliITSCalibration*)(fCalibration->At(iMod));
355 }
356 //_______________________________________________________________________
357 void AliITSDetTypeSim::SetDefaults(){
358     //Set defaults for segmentation and response
359
360     if(GetITSgeom()==0){
361         Warning("SetDefaults","GetITSgeom() is 0!");
362         return;
363     } // end if
364     if (fCalibration==0) {
365         CreateCalibrationArray();
366     } // end if
367
368     ResetSegmentation();
369     if(!GetCalibration()){AliFatal("Exit"); exit(0);}
370
371     for(Int_t idet=0;idet<fgkNdettypes;idet++){
372         //SPD
373         if(idet==0){
374             if(!GetSegmentationModel(idet)) SetDefaultSegmentation(idet);
375             const char *kData0=(GetCalibrationModel(GetITSgeom()->GetStartSPD()))->DataType();
376             if (strstr(kData0,"real")) {
377                 SetDigitClassName(idet,"AliITSdigit");
378             }else {
379                 SetDigitClassName(idet,"AliITSdigitSPD");
380             } // end if
381         } // end if idet==0
382         //SDD
383         if(idet==1){
384             if(!GetSegmentationModel(idet)) SetDefaultSegmentation(idet);
385             AliITSCalibrationSDD* rsp = 
386                 (AliITSCalibrationSDD*)GetCalibrationModel(
387                     GetITSgeom()->GetStartSDD());
388             const char *kopt = ((AliITSresponseSDD*)rsp->GetResponse())->
389                 ZeroSuppOption();
390             if((!strstr(kopt,"2D"))&&(!strstr(kopt,"1D"))) {
391                 SetDigitClassName(idet,"AliITSdigit");
392             }else {
393                 SetDigitClassName(idet,"AliITSdigitSDD");
394             } // end if
395         } // end if idet==1
396         //SSD
397         if(idet==2){
398             if(!GetSegmentationModel(idet))SetDefaultSegmentation(idet);
399             const char *kData2 = (GetCalibrationModel(
400                                   GetITSgeom()->GetStartSSD())->DataType());
401             if (strstr(kData2,"real")) {
402                 SetDigitClassName(idet,"AliITSdigit");
403             }else {
404                 SetDigitClassName(idet,"AliITSdigitSSD");
405             } // end if
406         } // end if idet==2
407     }// end for idet
408 }
409 //______________________________________________________________________
410 Bool_t AliITSDetTypeSim::GetCalibration() {
411   // Get Default calibration if a storage is not defined.
412
413   if(!fFirstcall){
414     AliITSCalibration* cal = GetCalibrationModel(0);
415     if(cal)return kTRUE;
416   }
417   else {
418     fFirstcall = kFALSE;
419   }
420
421   SetRunNumber((Int_t)AliCDBManager::Instance()->GetRun());
422   Int_t run=GetRunNumber();
423   if(run<0)run=0;   // if the run number is not yet set, use fake run # 0
424
425   Bool_t origCacheStatus = AliCDBManager::Instance()->GetCacheFlag();
426   Bool_t isCacheActive = kTRUE;
427   if(GetRunNumber()<0){
428     isCacheActive=kFALSE;
429     fCalibration->SetOwner(kTRUE);
430   }
431   else{
432     fCalibration->SetOwner(kFALSE);
433   }
434
435   AliCDBManager::Instance()->SetCacheFlag(isCacheActive);
436
437   AliCDBEntry *entrySPD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSPD", run);
438   AliCDBEntry *entrySDD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSDD", run);
439   AliCDBEntry *entrySSD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSSD", run);
440   AliCDBEntry *entry2SPD = AliCDBManager::Instance()->Get("ITS/Calib/RespSPD", run);
441   AliCDBEntry *entry2SDD = AliCDBManager::Instance()->Get("ITS/Calib/RespSDD", run);
442   AliCDBEntry *entry2SSD = AliCDBManager::Instance()->Get("ITS/Calib/RespSSD", run);
443
444   if(!entrySPD || !entrySDD || !entrySSD || !entry2SPD || !entry2SDD || !entry2SSD){
445         AliWarning("Calibration object retrieval failed! Dummy calibration will be used.");
446         AliCDBStorage *localStor = 
447                 AliCDBManager::Instance()->GetStorage("local://$ALICE_ROOT");
448         
449         entrySPD = localStor->Get("ITS/Calib/CalibSPD", run);
450         entrySDD = localStor->Get("ITS/Calib/CalibSDD", run);
451         entrySSD = localStor->Get("ITS/Calib/CalibSSD", run);
452         entry2SPD = localStor->Get("ITS/Calib/RespSPD", run);
453         entry2SDD = localStor->Get("ITS/Calib/RespSDD", run);
454         entry2SSD = localStor->Get("ITS/Calib/RespSSD", run);
455   }
456
457   if(!entrySPD || !entrySDD || !entrySSD || !entry2SPD || !entry2SDD || !entry2SSD){
458     AliError("Calibration data was not found in $ALICE_ROOT!");
459     return kFALSE;
460   }
461
462   TObjArray *calSPD = (TObjArray *)entrySPD->GetObject();
463   if(!isCacheActive)entrySPD->SetObject(NULL);
464   entrySPD->SetOwner(kTRUE);
465
466   AliITSresponseSPD *pSPD = (AliITSresponseSPD*)entry2SPD->GetObject();
467   if(!isCacheActive)entry2SPD->SetObject(NULL);
468   entry2SPD->SetOwner(kTRUE);
469    
470   TObjArray *calSDD = (TObjArray *)entrySDD->GetObject();
471   if(!isCacheActive)entrySDD->SetObject(NULL);
472   entrySDD->SetOwner(kTRUE);
473
474   AliITSresponseSDD *pSDD = (AliITSresponseSDD*)entry2SDD->GetObject();
475   if(!isCacheActive)entry2SDD->SetObject(NULL);
476   entry2SDD->SetOwner(kTRUE);
477
478   TObjArray *calSSD = (TObjArray *)entrySSD->GetObject();
479   if(!isCacheActive)entrySSD->SetObject(NULL);
480   entrySSD->SetOwner(kTRUE);
481
482   AliITSresponseSSD *pSSD = (AliITSresponseSSD*)entry2SSD->GetObject();
483   if(!isCacheActive)entry2SSD->SetObject(NULL);
484   entry2SSD->SetOwner(kTRUE);
485   
486   // DB entries are deleted. In this way metadeta objects are deleted as well
487   if(!isCacheActive){
488     delete entrySPD;
489     delete entrySDD;
490     delete entrySSD;
491     delete entry2SPD;
492     delete entry2SDD;
493     delete entry2SSD;
494   }
495   
496   AliCDBManager::Instance()->SetCacheFlag(origCacheStatus);
497
498   if ((!pSPD)||(!pSDD)||(!pSSD)) {
499     AliWarning("Can not get calibration from calibration database !");
500     return kFALSE;
501   }
502   if ((! calSPD)||(! calSDD)||(! calSSD)) {
503     AliWarning("Can not get calibration from calibration database !");
504     return kFALSE;
505   }
506   fNMod[0] = calSPD->GetEntries();
507   fNMod[1] = calSDD->GetEntries();
508   fNMod[2] = calSSD->GetEntries();
509   AliInfo(Form("%i SPD, %i SDD and %i SSD in calibration database",
510                fNMod[0], fNMod[1], fNMod[2]));
511   AliITSCalibration* cal;
512   for (Int_t i=0; i<fNMod[0]; i++) {
513     cal = (AliITSCalibration*) calSPD->At(i);
514     cal->SetResponse(pSPD);
515     SetCalibrationModel(i, cal);
516  }
517   for (Int_t i=0; i<fNMod[1]; i++) {
518     cal = (AliITSCalibration*) calSDD->At(i);
519     cal->SetResponse(pSDD);
520     Int_t iMod = i + fNMod[0];
521     SetCalibrationModel(iMod, cal);
522  }
523   for (Int_t i=0; i<fNMod[2]; i++) {
524     cal = (AliITSCalibration*) calSSD->At(i);
525     cal->SetResponse(pSSD);
526     Int_t iMod = i + fNMod[0] + fNMod[1];
527     SetCalibrationModel(iMod, cal);
528  }
529   return kTRUE;
530 }
531 //_______________________________________________________________________
532 void AliITSDetTypeSim::SetDefaultSimulation(){
533   //Set default simulation for detector type
534
535   if(GetITSgeom()==0){
536     Warning("SetDefaultSimulation","GetITSgeom() is 0!");
537     return;
538   }
539   if(fCalibration==0){
540     Warning("SetDefaultSimulation","fCalibration is 0!");
541     return;
542   }
543   if(fSegmentation==0){
544     Warning("SetDefaultSimulation","fSegmentation is 0!");
545     for(Int_t i=0;i<fgkNdettypes;i++) SetDefaultSegmentation(i);
546   }else for(Int_t i=0;i<fgkNdettypes;i++) if(!GetSegmentationModel(i)){
547       Warning("SetDefaultSimulation",
548               "Segmentation not defined for det %d - Default taken\n!",i);
549       SetDefaultSegmentation(i);
550   }
551   AliITSsimulation* sim;
552
553   for(Int_t idet=0;idet<fgkNdettypes;idet++){
554    //SPD
555     if(idet==0){
556       sim = GetSimulationModel(idet); 
557       if(!sim){
558         sim = new AliITSsimulationSPD(this);
559         SetSimulationModel(idet,sim);
560       }
561     }
562     //SDD
563     if(idet==1){
564       sim = GetSimulationModel(idet);
565       if(!sim){
566         sim = new AliITSsimulationSDD(this);
567         SetSimulationModel(idet,sim);
568       }      
569     }
570     //SSD
571     if(idet==2){
572       sim = GetSimulationModel(idet);
573       if(!sim){
574         sim = new AliITSsimulationSSD(this);
575         SetSimulationModel(idet,sim);
576       }
577     }
578
579   }
580 }
581 //___________________________________________________________________
582 void AliITSDetTypeSim::SetTreeAddressS(TTree* treeS, Char_t* name){
583   // Set branch address for the ITS summable digits Trees.  
584   char branchname[30];
585
586   if(!treeS){
587     return;
588   }
589   if (fSDigits ==  0x0){
590     fSDigits = new TClonesArray("AliITSpListItem",1000);
591   }
592   TBranch *branch;
593   sprintf(branchname,"%s",name);
594   branch = treeS->GetBranch(branchname);
595   if (branch) branch->SetAddress(&fSDigits);
596
597 }
598 //___________________________________________________________________
599 void AliITSDetTypeSim::SetTreeAddressD(TTree* treeD, Char_t* name){
600   // Set branch address for the digit Trees.
601   
602   const char *det[3] = {"SPD","SDD","SSD"};
603   TBranch *branch;
604   
605   char branchname[30];
606   
607   if(!treeD){
608     return;
609   }
610   if(!fDigits){
611     fDigits = new TObjArray(fgkNdettypes); 
612   }
613   for(Int_t i=0;i<fgkNdettypes;i++){
614     Char_t* digclass = GetDigitClassName(i);
615     if(digclass==0x0){
616       if(i==0) SetDigitClassName(i,"AliITSdigitSPD");
617       if(i==1) SetDigitClassName(i,"AliITSdigitSDD");
618       if(i==2) SetDigitClassName(i,"AliITSdigitSSD");
619       digclass = GetDigitClassName(i);
620     }
621     TString classn = digclass;
622     if(!(fDigits->At(i))){
623       fDigits->AddAt(new TClonesArray(classn.Data(),1000),i);
624     }else{
625       ResetDigits(i);
626     }
627     
628     if(fgkNdettypes==3) sprintf(branchname,"%sDigits%s",name,det[i]);
629     else sprintf(branchname,"%sDigits%d",name,i+1);
630     if(fDigits){
631       branch = treeD->GetBranch(branchname);
632       if(branch) branch->SetAddress(&((*fDigits)[i]));
633     }
634   }
635
636 }
637 //___________________________________________________________________
638 void AliITSDetTypeSim::ResetDigits(){
639   // Reset number of digits and the digits array for the ITS detector.  
640
641   if(!fDigits){
642     Error("ResetDigits","fDigits is null!");
643     return;
644   }
645   for(Int_t i=0;i<fgkNdettypes;i++){
646     ResetDigits(i);
647   }
648 }
649 //___________________________________________________________________
650 void AliITSDetTypeSim::ResetDigits(Int_t branch){
651   // Reset number of digits and the digits array for this branch.
652
653   if(fDigits->At(branch)){
654     ((TClonesArray*)fDigits->At(branch))->Clear();
655   }
656   if(fNDigits) fNDigits[branch]=0;
657
658 }
659
660
661
662 //_______________________________________________________________________
663 void AliITSDetTypeSim::SDigitsToDigits(Option_t* opt, Char_t* name){
664   // Standard Summable digits to Digits function.
665   if(!GetITSgeom()){
666     Warning("SDigitsToDigits","GetITSgeom() is null!!");
667     return;
668   }
669   
670   const char *all = strstr(opt,"All");
671   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
672                         strstr(opt,"SSD")};
673   if( !det[0] && !det[1] && !det[2] ) all = "All";
674   else all = 0;
675   static Bool_t setDef = kTRUE;
676   if(setDef) SetDefaultSimulation();
677   setDef = kFALSE;
678   
679   AliITSsimulation *sim =0;
680   TTree* trees = fLoader->TreeS();
681   if( !(trees && GetSDigits()) ){
682     Error("SDigits2Digits","Error: No trees or SDigits. Returning.");
683     return;
684   } 
685   sprintf(name,"%s",name);
686   TBranch* brchSDigits = trees->GetBranch(name);
687   
688   Int_t id;
689   for(Int_t module=0;module<GetITSgeom()->GetIndexMax();module++){
690      id = GetITSgeom()->GetModuleType(module);
691     if (!all && !det[id]) continue;
692     sim = (AliITSsimulation*)GetSimulationModel(id);
693     if(!sim){
694       Error("SDigit2Digits","The simulation class was not "
695             "instanciated for module %d type %s!",module,
696             GetITSgeom()->GetModuleTypeName(module));
697       exit(1);
698     }
699     sim->InitSimulationModule(module,gAlice->GetEvNumber());
700     
701     fSDigits->Clear();
702     brchSDigits->GetEvent(module);
703     sim->AddSDigitsToModule(fSDigits,0);
704     sim->FinishSDigitiseModule();
705     fLoader->TreeD()->Fill();
706     ResetDigits();
707   }
708   fLoader->TreeD()->GetEntries();
709   fLoader->TreeD()->AutoSave();
710   fLoader->TreeD()->Reset();
711 }
712 //_________________________________________________________
713 void AliITSDetTypeSim::AddSumDigit(AliITSpListItem &sdig){  
714   //Adds the module full of summable digits to the summable digits tree.
715
716   TClonesArray &lsdig = *fSDigits;
717   new(lsdig[fNSDigits++]) AliITSpListItem(sdig);
718 }
719 //__________________________________________________________
720 void AliITSDetTypeSim::AddRealDigit(Int_t branch, Int_t *digits){
721   //   Add a real digit - as coming from data.
722
723   TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
724   new(ldigits[fNDigits[branch]++]) AliITSdigit(digits); 
725 }
726 //__________________________________________________________
727 void AliITSDetTypeSim::AddSimDigit(Int_t branch, AliITSdigit* d){  
728   //    Add a simulated digit.
729
730   TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
731   switch(branch){
732   case 0:
733     new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
734     break;
735   case 1:
736     new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
737     break;
738   case 2:
739     new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
740     break;
741   }  
742 }
743 //______________________________________________________________________
744 void AliITSDetTypeSim::AddSimDigit(Int_t branch,Float_t phys,Int_t *digits,
745                                    Int_t *tracks,Int_t *hits,Float_t *charges){
746   //   Add a simulated digit to the list.
747
748   TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
749   AliITSCalibrationSDD *resp = 0;
750   switch(branch){
751   case 0:
752     new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(digits,tracks,hits);
753     break;
754   case 1:
755     resp = (AliITSCalibrationSDD*)GetCalibrationModel(
756         GetITSgeom()->GetStartSDD());
757     new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(phys,digits,tracks,
758                                                    hits,charges,resp);
759     break;
760   case 2:
761     new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(digits,tracks,hits);
762     break;
763   } 
764 }
765 //______________________________________________________________________
766 void AliITSDetTypeSim::StoreCalibration(Int_t firstRun, Int_t lastRun,
767                                         AliCDBMetaData &md) {
768   // Store calibration in the calibration database
769   // The database must be created in an external piece of code (i.e. 
770   // a configuration macro )
771
772   if(!AliCDBManager::Instance()->IsDefaultStorageSet()) {
773     AliWarning("No storage set! Will use dummy one");
774     AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT");
775   }
776
777   if (!fCalibration) {
778     AliError("AliITSCalibration classes are not defined - nothing done");
779     return;
780   }
781   AliCDBId idRespSPD("ITS/Calib/CalibSPD",firstRun, lastRun);
782   AliCDBId idRespSDD("ITS/Calib/CalibSDD",firstRun, lastRun);
783   AliCDBId idRespSSD("ITS/Calib/CalibSSD",firstRun, lastRun);
784
785   TObjArray respSPD(fNMod[0]);
786   TObjArray respSDD(fNMod[1]-fNMod[0]);
787   TObjArray respSSD(fNMod[2]-fNMod[1]);
788   respSPD.SetOwner(kFALSE);
789   respSSD.SetOwner(kFALSE);
790   respSSD.SetOwner(kFALSE);
791
792   Int_t index[fgkNdettypes];
793   for (Int_t i = 0; i<fgkNdettypes; i++ ) {
794     index[i] = 0;
795     for (Int_t j = 0; j<=i; j++ )
796       index[i]+=fNMod[j];
797   }
798
799   for (Int_t i = 0; i<index[0]; i++ )
800     respSPD.Add(fCalibration->At(i));
801
802   for (Int_t i = index[0]; i<index[1]; i++ )
803     respSDD.Add(fCalibration->At(i));
804
805   for (Int_t i = index[1]; i<index[2]; i++ )
806     respSSD.Add(fCalibration->At(i));
807
808   AliCDBManager::Instance()->Put(&respSPD, idRespSPD, &md);
809   AliCDBManager::Instance()->Put(&respSDD, idRespSDD, &md);
810   AliCDBManager::Instance()->Put(&respSSD, idRespSSD, &md);
811 }
812
813