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