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