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