1 /***************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
21 /////////////////////////////////////////////////////////////////////
22 // Base simulation functions for ITS //
25 /////////////////////////////////////////////////////////////////////
27 #include "TClonesArray.h"
28 #include "TObjArray.h"
33 #include "AliCDBManager.h"
35 #include "AliCDBStorage.h"
36 #include "AliCDBEntry.h"
37 #include "AliCDBMetaData.h"
38 #include "AliITSdigit.h"
39 #include "AliITSdigitSPD.h"
40 #include "AliITSdigitSDD.h"
41 #include "AliITSdigitSSD.h"
42 #include "AliITSgeom.h"
43 #include "AliITSDetTypeSim.h"
44 #include "AliITSpListItem.h"
45 #include "AliITSCalibration.h"
46 #include "AliITSCalibrationSDD.h"
47 #include "AliITSMapSDD.h"
48 #include "AliITSDriftSpeedArraySDD.h"
49 #include "AliITSDriftSpeedSDD.h"
50 #include "AliITSHLTforSDD.h"
51 #include "AliITSCalibrationSSD.h"
52 #include "AliITSNoiseSSDv2.h"
53 #include "AliITSGainSSDv2.h"
54 #include "AliITSBadChannelsSSDv2.h"
55 #include "AliITSNoiseSSD.h"
56 #include "AliITSGainSSD.h"
57 #include "AliITSBadChannelsSSD.h"
58 #include "AliITSCalibrationSSD.h"
59 #include "AliITSsegmentationSPD.h"
60 #include "AliITSsegmentationSDD.h"
61 #include "AliITSsegmentationSSD.h"
62 #include "AliITSsimulation.h"
63 #include "AliITSsimulationSPD.h"
64 #include "AliITSsimulationSDD.h"
65 #include "AliITSsimulationSSD.h"
66 #include "AliITSDDLModuleMapSDD.h"
68 const Int_t AliITSDetTypeSim::fgkNdettypes = 3;
69 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSPD = 240;
70 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSDD = 260;
71 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSSD = 1698;
73 ClassImp(AliITSDetTypeSim)
75 //----------------------------------------------------------------------
76 AliITSDetTypeSim::AliITSDetTypeSim():
78 fSimulation(), // [NDet]
79 fSegmentation(), // [NDet]
80 fCalibration(), // [NMod]
82 fNSDigits(0), //! number of SDigits
83 fSDigits("AliITSpListItem",1000),
84 fNDigits(0), //! number of Digits
85 fRunNumber(0), //! Run number (to access DB)
86 fDigits(), //! [NMod][NDigits]
89 fkDigClassName(), // String with digit class name.
90 fLoader(0), // local pointer to loader
92 fIsHLTmodeC(0){ // flag
93 // Default Constructor
99 // A properly zero-ed AliITSDetTypeSim class.
101 fSimulation = new TObjArray(fgkNdettypes);
102 fSegmentation = new TObjArray(fgkNdettypes);
103 fSegmentation->SetOwner(kTRUE);
104 fDigits = new TObjArray(fgkNdettypes);
105 fNDigits = new Int_t[fgkNdettypes];
106 fDDLMapSDD=new AliITSDDLModuleMapSDD();
107 fSimuPar= new AliITSSimuParam();
108 fSSDCalibration=new AliITSCalibrationSSD();
109 fNMod[0] = fgkDefaultNModulesSPD;
110 fNMod[1] = fgkDefaultNModulesSDD;
111 fNMod[2] = fgkDefaultNModulesSSD;
114 //----------------------------------------------------------------------
115 AliITSDetTypeSim::~AliITSDetTypeSim(){
125 fSimulation->Delete();
130 fSegmentation->Delete();
131 delete fSegmentation;
134 if(fCalibration && fRunNumber<0){
136 fCalibration->Delete();
140 if(fSSDCalibration) delete fSSDCalibration;
141 if(fSimuPar) delete fSimuPar;
142 if(fDDLMapSDD) delete fDDLMapSDD;
143 if(fNDigits) delete [] fNDigits;
145 if (fLoader)fLoader->GetModulesFolder()->Remove(this);
146 fLoader = 0; // Not deleting it.
154 //----------------------------------------------------------------------
155 AliITSDetTypeSim::AliITSDetTypeSim(const AliITSDetTypeSim &source) : TObject(source),
156 fSimulation(source.fSimulation), // [NDet]
157 fSegmentation(source.fSegmentation), // [NDet]
158 fCalibration(source.fCalibration), // [NMod]
159 fSSDCalibration(source.fSSDCalibration),
160 fNSDigits(source.fNSDigits), //! number of SDigits
161 fSDigits(*((TClonesArray*)source.fSDigits.Clone())),
162 fNDigits(source.fNDigits), //! number of Digits
163 fRunNumber(source.fRunNumber), //! Run number (to access DB)
164 fDigits(source.fDigits), //! [NMod][NDigits]
165 fSimuPar(source.fSimuPar),
166 fDDLMapSDD(source.fDDLMapSDD),
167 fkDigClassName(), // String with digit class name.
168 fLoader(source.fLoader), // local pointer to loader
169 fFirstcall(source.fFirstcall),
170 fIsHLTmodeC(source.fIsHLTmodeC)
172 // Copy Constructor for object AliITSDetTypeSim not allowed
173 for(Int_t i=0;i<fgkNdettypes;i++){
174 fkDigClassName[i] = source.fkDigClassName[i];
177 //----------------------------------------------------------------------
178 AliITSDetTypeSim& AliITSDetTypeSim::operator=(const AliITSDetTypeSim &source){
179 // The = operator for object AliITSDetTypeSim
181 this->~AliITSDetTypeSim();
182 new(this) AliITSDetTypeSim(source);
186 //______________________________________________________________________
187 void AliITSDetTypeSim::SetITSgeom(AliITSgeom *geom){
188 // Sets/replaces the existing AliITSgeom object kept in AliITSLoader
191 // AliITSgoem *geom The AliITSgeom object to be used.
197 Error("SetITSgeom","No pointer to loader - nothing done");
201 fLoader->SetITSgeom(geom); // protections in AliITSLoader::SetITSgeom
205 //______________________________________________________________________
206 void AliITSDetTypeSim::SetLoader(AliITSLoader *loader){
207 // Sets the local copy of the AliITSLoader, and passes on the
208 // AliITSgeom object as needed.
210 // AliITSLoader *loader pointer to AliITSLoader for local use
216 if(fLoader==loader) return; // Same do nothing
217 if(fLoader){ // alread have an existing loader
219 "Already have an exisiting loader ptr=%p Nothing done",
224 //______________________________________________________________________
225 void AliITSDetTypeSim::SetSimulationModel(Int_t dettype,AliITSsimulation *sim){
227 //Set simulation model for detector type
229 if(fSimulation==0) fSimulation = new TObjArray(fgkNdettypes);
230 fSimulation->AddAt(sim,dettype);
232 //______________________________________________________________________
233 AliITSsimulation* AliITSDetTypeSim::GetSimulationModel(Int_t dettype) const {
235 //Get simulation model for detector type
237 Warning("GetSimulationModel","fSimulation is 0!");
240 return (AliITSsimulation*)(fSimulation->At(dettype));
242 //______________________________________________________________________
243 AliITSsimulation* AliITSDetTypeSim::GetSimulationModelByModule(Int_t module) const {
245 //Get simulation model by module number
246 if(GetITSgeom()==0) {
247 Warning("GetSimulationModelByModule","GetITSgeom() is 0!");
251 return GetSimulationModel(GetITSgeom()->GetModuleType(module));
253 //_______________________________________________________________________
254 void AliITSDetTypeSim::SetDefaultSegmentation(Int_t idet){
255 // Set default segmentation model objects
256 AliITSsegmentation *seg;
258 if(fSegmentation==0x0){
259 fSegmentation = new TObjArray(fgkNdettypes);
260 fSegmentation->SetOwner(kTRUE);
262 if(GetSegmentationModel(idet))
263 delete (AliITSsegmentation*)fSegmentation->At(idet);
265 seg = new AliITSsegmentationSPD();
267 seg = new AliITSsegmentationSDD();
268 AliITSCalibrationSDD* cal=(AliITSCalibrationSDD*)GetCalibrationModel(fgkDefaultNModulesSPD+1);
269 if(cal->IsAMAt20MHz()){
270 seg->SetPadSize(seg->Dpz(0),20.);
271 seg->SetNPads(seg->Npz()/2,128);
274 seg = new AliITSsegmentationSSD();
276 SetSegmentationModel(idet,seg);
278 //______________________________________________________________________
279 void AliITSDetTypeSim::SetSegmentationModel(Int_t dettype,
280 AliITSsegmentation *seg){
282 //Set segmentation model for detector type
283 if(fSegmentation==0x0){
284 fSegmentation = new TObjArray(fgkNdettypes);
285 fSegmentation->SetOwner(kTRUE);
287 fSegmentation->AddAt(seg,dettype);
289 //______________________________________________________________________
290 AliITSsegmentation* AliITSDetTypeSim::GetSegmentationModel(Int_t dettype) const{
291 //Get segmentation model for detector type
293 if(fSegmentation==0) {
294 Warning("GetSegmentationModel","fSegmentation is 0!");
297 return (AliITSsegmentation*)(fSegmentation->At(dettype));
299 //_______________________________________________________________________
300 AliITSsegmentation* AliITSDetTypeSim::GetSegmentationModelByModule(Int_t module) const{
301 //Get segmentation model by module number
303 Warning("GetSegmentationModelByModule","GetITSgeom() is 0!");
306 return GetSegmentationModel(GetITSgeom()->GetModuleType(module));
308 //_______________________________________________________________________
309 void AliITSDetTypeSim::CreateCalibrationArray() {
310 //Create the container of calibration functions with correct size
312 Warning("CreateCalibration","pointer to calibration object exists\n");
313 fCalibration->Delete();
317 Int_t nModTot = GetITSgeom()->GetIndexMax();
318 fCalibration = new TObjArray(nModTot);
319 fCalibration->SetOwner(kTRUE);
320 fCalibration->Clear();
322 //_______________________________________________________________________
323 void AliITSDetTypeSim::SetCalibrationModel(Int_t iMod, AliITSCalibration *resp){
324 //Set response model for modules
326 if (fCalibration==0) CreateCalibrationArray();
328 if (fCalibration->At(iMod)!=0)
329 delete (AliITSCalibration*) fCalibration->At(iMod);
330 fCalibration->AddAt(resp, iMod);
332 //______________________________________________________________________
333 void AliITSDetTypeSim::ResetCalibrationArray(){
334 //resets response array
335 if(fCalibration && fRunNumber<0){ // if fRunNumber<0 fCalibration is owner
337 AliITSresponse* rspd = ((AliITSCalibration*)fCalibration->At(
338 GetITSgeom()->GetStartSPD()))->GetResponse();
339 AliITSresponse* rssd = ((AliITSCalibration*)fCalibration->At(
340 GetITSgeom()->GetStartSSD()))->GetResponse();
341 if(rspd) delete rspd;
342 if(rssd) delete rssd;
344 fCalibration->Clear();
345 }else if (fCalibration && fRunNumber>=0){
346 fCalibration->Clear();
349 //______________________________________________________________________
350 void AliITSDetTypeSim::ResetSegmentation(){
351 //Resets segmentation array
352 if(fSegmentation) fSegmentation->Clear();
354 //_______________________________________________________________________
355 AliITSCalibration* AliITSDetTypeSim::GetCalibrationModel(Int_t iMod) const {
356 //Get response model for module number iMod
358 if(fCalibration==0) {
359 AliError("fCalibration is 0!");
362 if(iMod<fgkDefaultNModulesSPD+fgkDefaultNModulesSDD){
363 return (AliITSCalibration*)fCalibration->At(iMod);
365 Int_t i=iMod-(fgkDefaultNModulesSPD+fgkDefaultNModulesSDD);
366 fSSDCalibration->SetModule(i);
367 return (AliITSCalibration*)fSSDCalibration;
371 //_______________________________________________________________________
372 void AliITSDetTypeSim::SetDefaults(){
373 //Set defaults for segmentation and response
376 Warning("SetDefaults","GetITSgeom() is 0!");
379 if (fCalibration==0) {
380 CreateCalibrationArray();
384 if(!GetCalibration()){AliFatal("Exit"); exit(0);}
386 SetDigitClassName(0,"AliITSdigitSPD");
387 SetDigitClassName(1,"AliITSdigitSDD");
388 SetDigitClassName(2,"AliITSdigitSSD");
390 for(Int_t idet=0;idet<fgkNdettypes;idet++){
391 if(!GetSegmentationModel(idet)) SetDefaultSegmentation(idet);
394 //______________________________________________________________________
395 Bool_t AliITSDetTypeSim::GetCalibration() {
396 // Get Default calibration if a storage is not defined.
399 AliITSCalibration* cal = GetCalibrationModel(0);
406 SetRunNumber((Int_t)AliCDBManager::Instance()->GetRun());
407 Int_t run=GetRunNumber();
409 Bool_t origCacheStatus = AliCDBManager::Instance()->GetCacheFlag();
410 Bool_t isCacheActive = kTRUE;
411 if(GetRunNumber()<0){
412 isCacheActive=kFALSE;
413 fCalibration->SetOwner(kTRUE);
416 fCalibration->SetOwner(kFALSE);
419 AliCDBManager::Instance()->SetCacheFlag(isCacheActive);
421 AliCDBEntry *entrySPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDDead", run);
422 AliCDBEntry *entrySDD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSDD", run);
423 AliCDBEntry *drSpSDD = AliCDBManager::Instance()->Get("ITS/Calib/DriftSpeedSDD",run);
424 AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD",run);
425 AliCDBEntry *hltforSDD = AliCDBManager::Instance()->Get("ITS/Calib/HLTforSDD");
426 //AliCDBEntry *mapASDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsAnodeSDD",run);
427 AliCDBEntry *mapTSDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsTimeSDD",run);
428 // AliCDBEntry *entrySSD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSSD", run);
429 AliCDBEntry *entryNoiseSSD = AliCDBManager::Instance()->Get("ITS/Calib/NoiseSSD");
430 AliCDBEntry *entryGainSSD = AliCDBManager::Instance()->Get("ITS/Calib/GainSSD");
431 AliCDBEntry *entryBadChannelsSSD = AliCDBManager::Instance()->Get("ITS/Calib/BadChannelsSSD");
433 if(!entrySPD || !entrySDD || !entryNoiseSSD || !entryGainSSD || !entryBadChannelsSSD ||
434 !drSpSDD || !ddlMapSDD || !hltforSDD || !mapTSDD){
435 AliFatal("Calibration object retrieval failed! ");
440 TObjArray *calSPD = (TObjArray *)entrySPD->GetObject();
441 if(!isCacheActive)entrySPD->SetObject(NULL);
442 entrySPD->SetOwner(kTRUE);
445 TObjArray *calSDD = (TObjArray *)entrySDD->GetObject();
446 if(!isCacheActive)entrySDD->SetObject(NULL);
447 entrySDD->SetOwner(kTRUE);
449 TObjArray *drSp = (TObjArray *)drSpSDD->GetObject();
450 if(!isCacheActive)drSpSDD->SetObject(NULL);
451 drSpSDD->SetOwner(kTRUE);
453 AliITSDDLModuleMapSDD *ddlsdd=(AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject();
454 if(!isCacheActive)ddlMapSDD->SetObject(NULL);
455 ddlMapSDD->SetOwner(kTRUE);
457 AliITSHLTforSDD* hltsdd=(AliITSHLTforSDD*)hltforSDD->GetObject();
458 if(!isCacheActive)hltforSDD->SetObject(NULL);
459 hltforSDD->SetOwner(kTRUE);
461 // TObjArray *mapAn = (TObjArray *)mapASDD->GetObject();
462 // if(!isCacheActive)mapASDD->SetObject(NULL);
463 // mapASDD->SetOwner(kTRUE);
465 TObjArray *mapT = (TObjArray *)mapTSDD->GetObject();
466 if(!isCacheActive)mapTSDD->SetObject(NULL);
467 mapTSDD->SetOwner(kTRUE);
470 TObjArray *calSSD = (TObjArray *)entrySSD->GetObject();
471 if(!isCacheActive)entrySSD->SetObject(NULL);
472 entrySSD->SetOwner(kTRUE);
475 TObject *emptyssd = 0; TString ssdobjectname = 0;
476 AliITSNoiseSSDv2 *noiseSSD = new AliITSNoiseSSDv2();
477 emptyssd = (TObject *)entryNoiseSSD->GetObject();
478 ssdobjectname = emptyssd->GetName();
479 if(ssdobjectname=="TObjArray") {
480 TObjArray *noiseSSDOld = (TObjArray *)entryNoiseSSD->GetObject();
481 ReadOldSSDNoise(noiseSSDOld, noiseSSD);
483 else if(ssdobjectname=="AliITSNoiseSSDv2")
484 noiseSSD = (AliITSNoiseSSDv2 *)entryNoiseSSD->GetObject();
485 if(!isCacheActive)entryNoiseSSD->SetObject(NULL);
486 entryNoiseSSD->SetOwner(kTRUE);
488 AliITSGainSSDv2 *gainSSD = new AliITSGainSSDv2();
489 emptyssd = (TObject *)entryGainSSD->GetObject();
490 ssdobjectname = emptyssd->GetName();
491 if(ssdobjectname=="Gain") {
492 TObjArray *gainSSDOld = (TObjArray *)entryGainSSD->GetObject();
493 ReadOldSSDGain(gainSSDOld, gainSSD);
495 else if(ssdobjectname=="AliITSGainSSDv2")
496 gainSSD = (AliITSGainSSDv2 *)entryGainSSD->GetObject();
497 if(!isCacheActive)entryGainSSD->SetObject(NULL);
498 entryGainSSD->SetOwner(kTRUE);
500 AliITSBadChannelsSSDv2 *badChannelsSSD = new AliITSBadChannelsSSDv2();
501 emptyssd = (TObject *)entryBadChannelsSSD->GetObject();
502 ssdobjectname = emptyssd->GetName();
503 if(ssdobjectname=="TObjArray") {
504 TObjArray *badChannelsSSDOld = (TObjArray *)entryBadChannelsSSD->GetObject();
505 ReadOldSSDBadChannels(badChannelsSSDOld, badChannelsSSD);
507 else if(ssdobjectname=="AliITSBadChannelsSSDv2")
508 badChannelsSSD = (AliITSBadChannelsSSDv2*)entryBadChannelsSSD->GetObject();
509 if(!isCacheActive)entryBadChannelsSSD->SetObject(NULL);
510 entryBadChannelsSSD->SetOwner(kTRUE);
512 /*AliITSNoiseSSDv2 *noiseSSD = (AliITSNoiseSSDv2 *)entryNoiseSSD->GetObject();
513 if(!isCacheActive)entryNoiseSSD->SetObject(NULL);
514 entryNoiseSSD->SetOwner(kTRUE);
516 AliITSGainSSDv2 *gainSSD = (AliITSGainSSDv2 *)entryGainSSD->GetObject();
517 if(!isCacheActive)entryGainSSD->SetObject(NULL);
518 entryGainSSD->SetOwner(kTRUE);
520 AliITSBadChannelsSSDv2 *badchannelsSSD =
521 (AliITSBadChannelsSSDv2 *)entryBadChannelsSSD->GetObject();
522 if(!isCacheActive)entryBadChannelsSSD->SetObject(NULL);
523 entryBadChannelsSSD->SetOwner(kTRUE);*/
525 // DB entries are deleted. In this way metadeta objects are deleted as well
529 delete entryNoiseSSD;
531 delete entryBadChannelsSSD;
539 AliCDBManager::Instance()->SetCacheFlag(origCacheStatus);
541 if ((!calSPD) || (!calSDD) || (!drSp) || (!ddlsdd)
542 || (!hltsdd) || (!mapT) || (!noiseSSD)|| (!gainSSD)|| (!badChannelsSSD)) {
543 AliWarning("Can not get calibration from calibration database !");
547 fNMod[0] = calSPD->GetEntries();
548 fNMod[1] = calSDD->GetEntries();
549 // fNMod[2] = noiseSSD->GetEntries();
550 AliInfo(Form("%i SPD, %i SDD and %i SSD in calibration database",
551 fNMod[0], fNMod[1], fNMod[2]));
552 AliITSCalibration* cal;
553 for (Int_t i=0; i<fNMod[0]; i++) {
554 cal = (AliITSCalibration*) calSPD->At(i);
555 SetCalibrationModel(i, cal);
558 fDDLMapSDD->SetDDLMap(ddlsdd);
559 fIsHLTmodeC=hltsdd->IsHLTmodeC();
561 for (Int_t i=0; i<fgkDefaultNModulesSDD; i++) {
563 Int_t iMod = i + fgkDefaultNModulesSPD;
564 fDDLMapSDD->FindInDDLMap(iMod,iddl,icarlos);
566 AliITSCalibrationSDD* calsdddead=new AliITSCalibrationSDD();
567 calsdddead->SetBad();
568 AliITSDriftSpeedSDD* driftspdef = new AliITSDriftSpeedSDD();
569 AliITSDriftSpeedArraySDD* arrdrsp=new AliITSDriftSpeedArraySDD(1);
570 arrdrsp->AddDriftSpeed(driftspdef);
571 calsdddead->SetDriftSpeed(0,arrdrsp);
572 calsdddead->SetDriftSpeed(1,arrdrsp);
573 SetCalibrationModel(iMod, calsdddead);
574 AliWarning(Form("SDD module %d not present in DDL map: set it as dead",iMod));
576 cal = (AliITSCalibration*) calSDD->At(i);
579 AliITSDriftSpeedArraySDD* arr0 = (AliITSDriftSpeedArraySDD*) drSp->At(i0);
580 // AliITSMapSDD* ma0 = (AliITSMapSDD*)mapAn->At(i0);
581 AliITSMapSDD* mt0 = (AliITSMapSDD*)mapT->At(i0);
582 AliITSDriftSpeedArraySDD* arr1 = (AliITSDriftSpeedArraySDD*) drSp->At(i1);
583 // AliITSMapSDD* ma1 = (AliITSMapSDD*)mapAn->At(i1);
584 AliITSMapSDD* mt1 = (AliITSMapSDD*)mapT->At(i1);
585 cal->SetDriftSpeed(0,arr0);
586 cal->SetDriftSpeed(1,arr1);
587 // cal->SetMapA(0,ma0);
588 // cal->SetMapA(1,ma1);
591 SetCalibrationModel(iMod, cal);
595 fSSDCalibration->SetNoise(noiseSSD);
596 fSSDCalibration->SetGain(gainSSD);
597 fSSDCalibration->SetBadChannels(badChannelsSSD);
598 //fSSDCalibration->FillBadChipMap();
604 //_______________________________________________________________________
605 void AliITSDetTypeSim::SetDefaultSimulation(){
606 //Set default simulation for detector type
609 Warning("SetDefaultSimulation","GetITSgeom() is 0!");
613 Warning("SetDefaultSimulation","fCalibration is 0!");
616 if(fSegmentation==0){
617 Warning("SetDefaultSimulation","fSegmentation is 0!");
618 for(Int_t i=0;i<fgkNdettypes;i++) SetDefaultSegmentation(i);
619 }else for(Int_t i=0;i<fgkNdettypes;i++) if(!GetSegmentationModel(i)){
620 Warning("SetDefaultSimulation",
621 "Segmentation not defined for det %d - Default taken\n!",i);
622 SetDefaultSegmentation(i);
624 AliITSsimulation* sim;
626 for(Int_t idet=0;idet<fgkNdettypes;idet++){
629 sim = GetSimulationModel(idet);
631 sim = new AliITSsimulationSPD(this);
632 SetSimulationModel(idet,sim);
637 sim = GetSimulationModel(idet);
639 sim = new AliITSsimulationSDD(this);
640 SetSimulationModel(idet,sim);
645 sim = GetSimulationModel(idet);
647 sim = new AliITSsimulationSSD(this);
648 SetSimulationModel(idet,sim);
654 //___________________________________________________________________
655 void AliITSDetTypeSim::SetTreeAddressS(TTree* treeS, const Char_t* name){
656 // Set branch address for the ITS summable digits Trees.
663 sprintf(branchname,"%s",name);
664 branch = treeS->GetBranch(branchname);
665 TClonesArray *sdigi = &fSDigits;
666 if (branch) branch->SetAddress(&sdigi);
669 //___________________________________________________________________
670 void AliITSDetTypeSim::SetTreeAddressD(TTree* treeD, const Char_t* name){
671 // Set branch address for the digit Trees.
673 const char *det[3] = {"SPD","SDD","SSD"};
682 fDigits = new TObjArray(fgkNdettypes);
684 for(Int_t i=0;i<fgkNdettypes;i++){
685 const Char_t* digclass = GetDigitClassName(i);
687 if(i==0) SetDigitClassName(i,"AliITSdigitSPD");
688 if(i==1) SetDigitClassName(i,"AliITSdigitSDD");
689 if(i==2) SetDigitClassName(i,"AliITSdigitSSD");
690 digclass = GetDigitClassName(i);
692 TString classn = digclass;
693 if(!(fDigits->At(i))){
694 fDigits->AddAt(new TClonesArray(classn.Data(),1000),i);
699 if(fgkNdettypes==3) sprintf(branchname,"%sDigits%s",name,det[i]);
700 else sprintf(branchname,"%sDigits%d",name,i+1);
702 branch = treeD->GetBranch(branchname);
703 if(branch) branch->SetAddress(&((*fDigits)[i]));
708 //___________________________________________________________________
709 void AliITSDetTypeSim::ResetDigits(){
710 // Reset number of digits and the digits array for the ITS detector.
713 Error("ResetDigits","fDigits is null!");
716 for(Int_t i=0;i<fgkNdettypes;i++){
720 //___________________________________________________________________
721 void AliITSDetTypeSim::ResetDigits(Int_t branch){
722 // Reset number of digits and the digits array for this branch.
724 if(fDigits->At(branch)){
725 ((TClonesArray*)fDigits->At(branch))->Clear();
727 if(fNDigits) fNDigits[branch]=0;
733 //_______________________________________________________________________
734 void AliITSDetTypeSim::SDigitsToDigits(Option_t* opt, Char_t* name){
735 // Standard Summable digits to Digits function.
737 Warning("SDigitsToDigits","GetITSgeom() is null!!");
741 const char *all = strstr(opt,"All");
742 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
744 if( !det[0] && !det[1] && !det[2] ) all = "All";
746 static Bool_t setDef = kTRUE;
747 if(setDef) SetDefaultSimulation();
750 AliITSsimulation *sim =0;
751 TTree* trees = fLoader->TreeS();
752 if( !(trees && GetSDigits()) ){
753 Error("SDigits2Digits","Error: No trees or SDigits. Returning.");
756 sprintf(name,"%s",name);
757 TBranch* brchSDigits = trees->GetBranch(name);
760 for(Int_t module=0;module<GetITSgeom()->GetIndexMax();module++){
761 id = GetITSgeom()->GetModuleType(module);
762 if (!all && !det[id]) continue;
763 sim = (AliITSsimulation*)GetSimulationModel(id);
765 Error("SDigit2Digits","The simulation class was not "
766 "instanciated for module %d type %s!",module,
767 GetITSgeom()->GetModuleTypeName(module));
770 sim->InitSimulationModule(module,gAlice->GetEvNumber());
773 brchSDigits->GetEvent(module);
774 sim->AddSDigitsToModule(&fSDigits,0);
775 sim->FinishSDigitiseModule();
776 fLoader->TreeD()->Fill();
779 fLoader->TreeD()->GetEntries();
780 fLoader->TreeD()->AutoSave();
781 fLoader->TreeD()->Reset();
783 //_________________________________________________________
784 void AliITSDetTypeSim::AddSumDigit(AliITSpListItem &sdig){
785 //Adds the module full of summable digits to the summable digits tree.
787 new(fSDigits[fNSDigits++]) AliITSpListItem(sdig);
789 //__________________________________________________________
790 void AliITSDetTypeSim::AddSimDigit(Int_t branch, const AliITSdigit* d){
791 // Add a simulated digit.
793 TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
796 new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
799 new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
802 new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
806 //______________________________________________________________________
807 void AliITSDetTypeSim::AddSimDigit(Int_t branch,Float_t phys,Int_t *digits,
808 Int_t *tracks,Int_t *hits,Float_t *charges,
810 // Add a simulated digit to the list.
812 TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
815 new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(digits,tracks,hits);
818 new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(phys,digits,tracks,
819 hits,charges,sigexpanded);
822 new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(digits,tracks,hits);
826 //______________________________________________________________________
827 void AliITSDetTypeSim::ReadOldSSDNoise(const TObjArray *array,
828 AliITSNoiseSSDv2 *noiseSSD) {
829 //Reads the old SSD calibration object and converts it to the new format
830 const Int_t fgkSSDSTRIPSPERMODULE = 1536;
831 const Int_t fgkSSDPSIDESTRIPSPERMODULE = 768;
833 Int_t gNMod = array->GetEntries();
834 cout<<"Converting old calibration object for noise..."<<endl;
837 Double_t noise = 0.0;
838 for (Int_t iModule = 0; iModule < gNMod; iModule++) {
839 AliITSNoiseSSD *noiseModule = (AliITSNoiseSSD*) (array->At(iModule));
840 for(Int_t iStrip = 0; iStrip < fgkSSDSTRIPSPERMODULE; iStrip++) {
841 noise = (iStrip < fgkSSDPSIDESTRIPSPERMODULE) ? noiseModule->GetNoiseP(iStrip) : noiseModule->GetNoiseN(1535 - iStrip);
842 if(iStrip < fgkSSDPSIDESTRIPSPERMODULE)
843 noiseSSD->AddNoiseP(iModule,iStrip,noise);
844 if(iStrip >= fgkSSDPSIDESTRIPSPERMODULE)
845 noiseSSD->AddNoiseN(iModule,1535 - iStrip,noise);
850 //______________________________________________________________________
851 void AliITSDetTypeSim::ReadOldSSDBadChannels(const TObjArray *array,
852 AliITSBadChannelsSSDv2 *badChannelsSSD) {
853 //Reads the old SSD calibration object and converts it to the new format
854 Int_t nMod = array->GetEntries();
855 cout<<"Converting old calibration object for bad channels..."<<endl;
856 for (Int_t iModule = 0; iModule < nMod; iModule++) {
857 //for (Int_t iModule = 0; iModule < 1; iModule++) {
858 AliITSBadChannelsSSD *bad = (AliITSBadChannelsSSD*) (array->At(iModule));
859 TArrayI arrayPSide = bad->GetBadPChannelsList();
860 for(Int_t iPCounter = 0; iPCounter < arrayPSide.GetSize(); iPCounter++)
861 badChannelsSSD->AddBadChannelP(iModule,
863 (Char_t)arrayPSide.At(iPCounter));
865 TArrayI arrayNSide = bad->GetBadNChannelsList();
866 for(Int_t iNCounter = 0; iNCounter < arrayNSide.GetSize(); iNCounter++)
867 badChannelsSSD->AddBadChannelN(iModule,
869 (Char_t)arrayNSide.At(iNCounter));
874 //______________________________________________________________________
875 void AliITSDetTypeSim::ReadOldSSDGain(const TObjArray *array,
876 AliITSGainSSDv2 *gainSSD) {
877 //Reads the old SSD calibration object and converts it to the new format
879 Int_t nMod = array->GetEntries();
880 cout<<"Converting old calibration object for gain..."<<endl;
883 for (Int_t iModule = 0; iModule < nMod; iModule++) {
884 AliITSGainSSD *gainModule = (AliITSGainSSD*) (array->At(iModule));
885 TArrayF arrayPSide = gainModule->GetGainP();
886 for(Int_t iPCounter = 0; iPCounter < arrayPSide.GetSize(); iPCounter++)
887 gainSSD->AddGainP(iModule,
889 arrayPSide.At(iPCounter));
890 TArrayF arrayNSide = gainModule->GetGainN();
891 for(Int_t iNCounter = 0; iNCounter < arrayNSide.GetSize(); iNCounter++)
892 gainSSD->AddGainN(iModule,
894 arrayNSide.At(iNCounter));