]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSDetTypeSim.cxx
Renamed method in AliCDBManager (A.Colla)
[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   Bool_t deleteManager = kFALSE;
405   if(!AliCDBManager::Instance()->IsDefaultStorageSet()) {
406     AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT");
407     deleteManager = kTRUE;
408   }
409   AliCDBStorage *storage = AliCDBManager::Instance()->GetDefaultStorage();
410
411   AliCDBEntry *entrySPD = storage->Get("ITS/Calib/RespSPD", fRunNumber);
412   TObjArray *respSPD = (TObjArray *)entrySPD->GetObject();
413   entrySPD->SetObject(NULL);
414   entrySPD->SetOwner(kTRUE);
415   AliCDBEntry *entrySDD = storage->Get("ITS/Calib/RespSDD", fRunNumber);
416   TObjArray *respSDD = (TObjArray *)entrySDD->GetObject();
417   entrySDD->SetObject(NULL);
418   entrySDD->SetOwner(kTRUE);
419   AliCDBEntry *entrySSD = storage->Get("ITS/Calib/RespSSD", fRunNumber);
420   TObjArray *respSSD = (TObjArray *)entrySSD->GetObject();
421   entrySSD->SetObject(NULL);
422   entrySSD->SetOwner(kTRUE);
423   // DB entries are dleted. In this waymetadeta objects are deleted as well
424   delete entrySPD;
425   delete entrySDD;
426   delete entrySSD;
427   if(deleteManager){
428     AliCDBManager::Instance()->Destroy();
429     AliCDBManager::Instance()->UnsetDefaultStorage();
430     storage = 0;   // the storage is killed by AliCDBManager::Instance()->Destroy()
431   }
432
433   if ((! respSPD)||(! respSDD)||(! respSSD)) {
434     AliWarning("Can not get calibration from calibration database !");
435     return kFALSE;
436   }
437   fNMod[0] = respSPD->GetEntries();
438   fNMod[1] = respSDD->GetEntries();
439   fNMod[2] = respSSD->GetEntries();
440   AliInfo(Form("%i SPD, %i SDD and %i SSD in calibration database",
441                fNMod[0], fNMod[1], fNMod[2]));
442   
443   for (Int_t i=0; i<fNMod[0]; i++) {
444   AliITSresponse*  res = (AliITSresponse*) respSPD->At(i);
445   SetResponseModel(i, res);
446  }
447   for (Int_t i=0; i<fNMod[1]; i++) {
448     AliITSresponse* res = (AliITSresponse*) respSDD->At(i);
449     Int_t iMod = i + fNMod[0];
450     SetResponseModel(iMod, res);
451  }
452   for (Int_t i=0; i<fNMod[2]; i++) {
453      AliITSresponse* res = (AliITSresponse*) respSSD->At(i);
454     Int_t iMod = i + fNMod[0] + fNMod[1];
455     SetResponseModel(iMod, res);
456  }
457
458   return kTRUE;
459 }
460
461
462
463 //_______________________________________________________________________
464 void AliITSDetTypeSim::SetDefaultSimulation(){
465
466   //Set default simulation for detector type
467
468
469   if(fGeom==0){
470     Warning("SetDefaultSimulation","fGeom is 0!\n");
471     return;
472   }
473   if(fResponse==0){
474     Warning("SetDefaultSimulation","fResponse is 0!\n");
475     return;
476   }
477   if(fSegmentation ==0){
478     Warning("SetDefaultSimulation","fSegmentation is 0!\n");
479     for(Int_t i=0;i<fgkNdettypes;i++)SetDefaultSegmentation(i);
480   }
481   else {
482     for(Int_t i=0;i<fgkNdettypes;i++){
483       if(!GetSegmentationModel(i)){
484         Warning("SetDefaultSimulation","Segmentation not defined for det %d - Default taken\n");
485         SetDefaultSegmentation(i);
486       }
487     }
488   }
489   AliITSsimulation* sim;
490
491   for(Int_t idet=0;idet<fgkNdettypes;idet++){
492    //SPD
493     if(idet==0){
494       sim = GetSimulationModel(idet);
495       if(!sim){
496         sim = new AliITSsimulationSPD(this);
497         SetSimulationModel(idet,sim);
498       } else{
499         sim->SetDetType(this);
500         sim->SetSegmentationModel(0,(AliITSsegmentationSPD*)GetSegmentationModel(idet));
501         sim->Init();
502       }
503     }
504     //SDD
505     if(idet==1){
506       sim = GetSimulationModel(idet);
507       if(!sim){
508         sim = new AliITSsimulationSDD(this);
509         SetSimulationModel(idet,sim);
510       } else {
511         sim->SetDetType(this);
512         sim->SetSegmentationModel(1,(AliITSsegmentationSDD*)GetSegmentationModel(idet));
513         sim->Init();
514
515       }
516       
517     }
518     //SSD
519     if(idet==2){
520       sim = GetSimulationModel(idet);
521       if(!sim){
522         sim = new AliITSsimulationSSD(this);
523         SetSimulationModel(idet,sim);
524
525       } else{
526         sim->SetDetType(this);
527         sim->SetSegmentationModel(2,(AliITSsegmentationSSD*)GetSegmentationModel(idet));    
528         sim->Init();
529       }
530
531     }
532
533   }
534 }
535
536
537
538
539 //___________________________________________________________________
540 void AliITSDetTypeSim::SetTreeAddressS(TTree* treeS, Char_t* name){
541   // Set branch address for the ITS summable digits Trees.
542   
543   char branchname[30];
544
545   if(!treeS){
546     return;
547   }
548   if (fSDigits == 0x0){
549     fSDigits = new TClonesArray("AliITSpListItem",1000);
550   }
551   TBranch *branch;
552   sprintf(branchname,"%s",name);
553   branch = treeS->GetBranch(branchname);
554   if (branch) branch->SetAddress(&fSDigits);
555
556 }
557 //___________________________________________________________________
558 void AliITSDetTypeSim::SetTreeAddressD(TTree* treeD, Char_t* name){
559   // Set branch address for the digit Trees.
560   
561   const char *det[3] = {"SPD","SDD","SSD"};
562   TBranch *branch;
563   
564   char branchname[30];
565   
566   if(!treeD){
567     return;
568   }
569   if(!fDigits){
570     fDigits = new TObjArray(fgkNdettypes); 
571   }
572   for(Int_t i=0;i<fgkNdettypes;i++){
573     Char_t* digclass = GetDigitClassName(i);
574     if(digclass==0x0){
575       if(i==0) SetDigitClassName(i,"AliITSdigitSPD");
576       if(i==1) SetDigitClassName(i,"AliITSdigitSDD");
577       if(i==2) SetDigitClassName(i,"AliITSdigitSSD");
578       digclass = GetDigitClassName(i);
579     }
580     TString classn = digclass;
581     if(!(fDigits->At(i))){
582       fDigits->AddAt(new TClonesArray(classn.Data(),1000),i);
583     }else{
584       ResetDigits(i);
585     }
586     
587     if(fgkNdettypes==3) sprintf(branchname,"%sDigits%s",name,det[i]);
588     else sprintf(branchname,"%sDigits%d",name,i+1);
589     if(fDigits){
590       branch = treeD->GetBranch(branchname);
591       if(branch) branch->SetAddress(&((*fDigits)[i]));
592     }
593   }
594
595 }
596 //___________________________________________________________________
597 void AliITSDetTypeSim::ResetDigits(){
598   // Reset number of digits and the digits array for the ITS detector.
599   
600
601   if(!fDigits){
602     Error("ResetDigits","fDigits is null!");
603     return;
604   }
605   for(Int_t i=0;i<fgkNdettypes;i++){
606     ResetDigits(i);
607   }
608 }
609 //___________________________________________________________________
610 void AliITSDetTypeSim::ResetDigits(Int_t branch){
611   // Reset number of digits and the digits array for this branch.
612
613   if(fDigits->At(branch)){
614     ((TClonesArray*)fDigits->At(branch))->Clear();
615   }
616   if(fNDigits) fNDigits[branch]=0;
617
618 }
619
620
621
622 //_______________________________________________________________________
623 void AliITSDetTypeSim::SDigitsToDigits(Option_t* opt, Char_t* name){
624   // Standard Summable digits to Digits function.
625   if(!fGeom){
626     Warning("SDigitsToDigits","fGeom is null!!");
627     return;
628   }
629   
630   const char *all = strstr(opt,"All");
631   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
632                         strstr(opt,"SSD")};
633   if( !det[0] && !det[1] && !det[2] ) all = "All";
634   else all = 0;
635   static Bool_t setDef = kTRUE;
636   if(setDef) SetDefaultSimulation();
637   setDef = kFALSE;
638   
639   AliITSsimulation *sim =0;
640   TTree* trees = fLoader->TreeS();
641   if( !(trees && GetSDigits()) ){
642     Error("SDigits2Digits","Error: No trees or SDigits. Returning.");
643     return;
644   } 
645   sprintf(name,"%s",name);
646   TBranch* brchSDigits = trees->GetBranch(name);
647   
648   Int_t id;
649   for(Int_t module=0;module<fGeom->GetIndexMax();module++){
650      id = fGeom->GetModuleType(module);
651     if (!all && !det[id]) continue;
652     sim = (AliITSsimulation*)GetSimulationModel(id);
653     printf("module=%d name=%s\n",module,sim->ClassName());
654     if(!sim){
655       Error("SDigit2Digits","The simulation class was not "
656             "instanciated for module %d type %s!",module,
657             fGeom->GetModuleTypeName(module));
658       exit(1);
659     }
660     sim->InitSimulationModule(module,gAlice->GetEvNumber());
661     
662     fSDigits->Clear();
663     brchSDigits->GetEvent(module);
664     sim->AddSDigitsToModule(fSDigits,0);
665     sim->FinishSDigitiseModule();
666     fLoader->TreeD()->Fill();
667     ResetDigits();
668   }
669   fLoader->TreeD()->GetEntries();
670   fLoader->TreeD()->AutoSave();
671   fLoader->TreeD()->Reset();
672 }
673
674
675
676 //_________________________________________________________
677 void AliITSDetTypeSim::AddSumDigit(AliITSpListItem &sdig){
678   
679   //Adds the module full of summable digits to the summable digits tree.
680   TClonesArray &lsdig = *fSDigits;
681   new(lsdig[fNSDigits++]) AliITSpListItem(sdig);
682 }
683 //__________________________________________________________
684 void AliITSDetTypeSim::AddRealDigit(Int_t branch, Int_t *digits){
685   //   Add a real digit - as coming from data.
686   TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
687   new(ldigits[fNDigits[branch]++]) AliITSdigit(digits); 
688 }
689 //__________________________________________________________
690 void AliITSDetTypeSim::AddSimDigit(Int_t branch, AliITSdigit* d){
691   
692   //    Add a simulated digit.
693   TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
694   switch(branch){
695   case 0:
696     new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
697     break;
698   case 1:
699     new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
700     break;
701   case 2:
702     new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
703     break;
704   } 
705   
706 }
707
708 //______________________________________________________________________
709 void AliITSDetTypeSim::AddSimDigit(Int_t branch,Float_t phys,Int_t *digits,
710                                    Int_t *tracks,Int_t *hits,Float_t *charges){
711   //   Add a simulated digit to the list.
712
713   TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
714   AliITSresponseSDD *resp = 0;
715   switch(branch){
716   case 0:
717     new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(digits,tracks,hits);
718     break;
719   case 1:
720     resp = (AliITSresponseSDD*)GetResponseModel(fGeom->GetStartSDD());
721     new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(phys,digits,tracks,
722                                                    hits,charges,resp);
723     break;
724   case 2:
725     new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(digits,tracks,hits);
726     break;
727   } 
728 }
729
730
731
732 //______________________________________________________________________
733 void AliITSDetTypeSim::StoreCalibration(Int_t firstRun, Int_t lastRun, AliCDBMetaData &md) {
734
735   // Store calibration in the calibration database
736
737   // The database must be created in an external piece of code (i.e. 
738   // a configuration macro )
739
740   if(!AliCDBManager::Instance()->IsDefaultStorageSet()) {
741     //AliError("No storage set!");
742     AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT");
743     //return;
744   }
745
746   if (!fResponse) {
747     AliError("AliITSresponse classes are not defined - nothing done");
748     return;
749   }
750   AliCDBId idRespSPD("ITS/Calib/RespSPD",firstRun, lastRun);
751   AliCDBId idRespSDD("ITS/Calib/RespSDD",firstRun, lastRun);
752   AliCDBId idRespSSD("ITS/Calib/RespSSD",firstRun, lastRun);
753
754   TObjArray respSPD(fNMod[0]);
755   TObjArray respSDD(fNMod[1]-fNMod[0]);
756   TObjArray respSSD(fNMod[2]-fNMod[1]);
757   respSPD.SetOwner(kFALSE);
758   respSSD.SetOwner(kFALSE);
759   respSSD.SetOwner(kFALSE);
760
761   Int_t index[fgkNdettypes];
762   for (Int_t i = 0; i<fgkNdettypes; i++ ) {
763     index[i] = 0;
764     for (Int_t j = 0; j<=i; j++ )
765       index[i]+=fNMod[j];
766   }
767
768   for (Int_t i = 0; i<index[0]; i++ )
769     respSPD.Add(fResponse->At(i));
770
771   for (Int_t i = index[0]; i<index[1]; i++ )
772     respSDD.Add(fResponse->At(i));
773
774   for (Int_t i = index[1]; i<index[2]; i++ )
775     respSSD.Add(fResponse->At(i));
776
777   AliCDBManager::Instance()->GetDefaultStorage()->Put(&respSPD, idRespSPD, &md);
778
779   AliCDBManager::Instance()->GetDefaultStorage()->Put(&respSDD, idRespSDD, &md);
780   
781   AliCDBManager::Instance()->GetDefaultStorage()->Put(&respSSD, idRespSSD, &md);
782 }
783
784