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