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