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