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