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