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"
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 "AliITSCalibrationSDD.h"
46 #include "AliITSMapSDD.h"
47 #include "AliITSDriftSpeedArraySDD.h"
48 #include "AliITSDriftSpeedSDD.h"
49 #include "AliITSCalibrationSSD.h"
50 #include "AliITSNoiseSSDv2.h"
51 #include "AliITSGainSSDv2.h"
52 #include "AliITSBadChannelsSSDv2.h"
53 #include "AliITSNoiseSSD.h"
54 #include "AliITSGainSSD.h"
55 #include "AliITSBadChannelsSSD.h"
56 #include "AliITSCalibrationSSD.h"
57 #include "AliITSsegmentationSPD.h"
58 #include "AliITSsegmentationSDD.h"
59 #include "AliITSsegmentationSSD.h"
60 #include "AliITSsimulation.h"
61 #include "AliITSsimulationSPD.h"
62 #include "AliITSsimulationSDD.h"
63 #include "AliITSsimulationSSD.h"
66 const Int_t AliITSDetTypeSim::fgkNdettypes = 3;
67 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSPD = 240;
68 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSDD = 260;
69 const Int_t AliITSDetTypeSim::fgkDefaultNModulesSSD = 1698;
71 ClassImp(AliITSDetTypeSim)
73 //----------------------------------------------------------------------
74 AliITSDetTypeSim::AliITSDetTypeSim():
76 fSimulation(), // [NDet]
77 fSegmentation(), // [NDet]
78 fCalibration(), // [NMod]
80 fPreProcess(), // [] e.g. Fill fHitModule with hits
81 fPostProcess(), // [] e.g. Wright Raw data
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 fHitClassName(), // String with Hit class name.
90 fSDigClassName(),// String with SDigit class name.
91 fDigClassName(), // String with digit class name.
92 fLoader(0), // local pointer to loader
93 fFirstcall(kTRUE){ // flag
94 // Default Constructor
100 // A properly zero-ed AliITSDetTypeSim class.
102 fSimulation = new TObjArray(fgkNdettypes);
103 fSegmentation = new TObjArray(fgkNdettypes);
104 fSegmentation->SetOwner(kTRUE);
105 fDigits = new TObjArray(fgkNdettypes);
106 fNDigits = new Int_t[fgkNdettypes];
107 fDDLMapSDD=new AliITSDDLModuleMapSDD();
108 fSimuPar= new AliITSSimuParam();
109 fSSDCalibration=new AliITSCalibrationSSD();
110 fNMod[0] = fgkDefaultNModulesSPD;
111 fNMod[1] = fgkDefaultNModulesSDD;
112 fNMod[2] = fgkDefaultNModulesSSD;
115 //----------------------------------------------------------------------
116 AliITSDetTypeSim::~AliITSDetTypeSim(){
126 fSimulation->Delete();
131 fSegmentation->Delete();
132 delete fSegmentation;
135 if(fCalibration && fRunNumber<0){
136 AliITSresponse* rspd = ((AliITSCalibration*)fCalibration->At(
137 GetITSgeom()->GetStartSPD()))->GetResponse();
138 AliITSresponse* rssd = ((AliITSCalibration*)fCalibration->At(
139 GetITSgeom()->GetStartSSD()))->GetResponse();
140 if(rspd) delete rspd;
141 if(rssd) delete rssd;
142 fCalibration->Delete();
146 if(fSSDCalibration) delete fSSDCalibration;
148 fPreProcess->Delete();
153 fPostProcess->Delete();
157 if(fSimuPar) delete fSimuPar;
158 if(fDDLMapSDD) delete fDDLMapSDD;
159 if(fNDigits) delete [] fNDigits;
161 if (fLoader)fLoader->GetModulesFolder()->Remove(this);
162 fLoader = 0; // Not deleting it.
170 //----------------------------------------------------------------------
171 AliITSDetTypeSim::AliITSDetTypeSim(const AliITSDetTypeSim &source) : TObject(source),
172 fSimulation(source.fSimulation), // [NDet]
173 fSegmentation(source.fSegmentation), // [NDet]
174 fCalibration(source.fCalibration), // [NMod]
175 fSSDCalibration(source.fSSDCalibration),
176 fPreProcess(source.fPreProcess), // [] e.g. Fill fHitModule with hits
177 fPostProcess(source.fPostProcess), // [] e.g. Wright Raw data
178 fNSDigits(source.fNSDigits), //! number of SDigits
179 fSDigits(*((TClonesArray*)source.fSDigits.Clone())),
180 fNDigits(source.fNDigits), //! number of Digits
181 fRunNumber(source.fRunNumber), //! Run number (to access DB)
182 fDigits(source.fDigits), //! [NMod][NDigits]
183 fSimuPar(source.fSimuPar),
184 fDDLMapSDD(source.fDDLMapSDD),
185 fHitClassName(source.fHitClassName), // String with Hit class name.
186 fSDigClassName(source.fSDigClassName),// String with SDigit class name.
187 fDigClassName(), // String with digit class name.
188 fLoader(source.fLoader), // local pointer to loader
189 fFirstcall(source.fFirstcall)
191 // Copy Constructor for object AliITSDetTypeSim not allowed
192 for(Int_t i=0;i<fgkNdettypes;i++){
193 fDigClassName[i] = source.fDigClassName[i];
196 //----------------------------------------------------------------------
197 AliITSDetTypeSim& AliITSDetTypeSim::operator=(const AliITSDetTypeSim &source){
198 // The = operator for object AliITSDetTypeSim
200 this->~AliITSDetTypeSim();
201 new(this) AliITSDetTypeSim(source);
205 //______________________________________________________________________
206 void AliITSDetTypeSim::SetITSgeom(AliITSgeom *geom){
207 // Sets/replaces the existing AliITSgeom object kept in AliITSLoader
210 // AliITSgoem *geom The AliITSgeom object to be used.
216 Error("SetITSgeom","No pointer to loader - nothing done");
220 fLoader->SetITSgeom(geom); // protections in AliITSLoader::SetITSgeom
224 //______________________________________________________________________
225 void AliITSDetTypeSim::SetLoader(AliITSLoader *loader){
226 // Sets the local copy of the AliITSLoader, and passes on the
227 // AliITSgeom object as needed.
229 // AliITSLoader *loader pointer to AliITSLoader for local use
235 if(fLoader==loader) return; // Same do nothing
236 if(fLoader){ // alread have an existing loader
238 "Already have an exisiting loader ptr=%p Nothing done",
243 //______________________________________________________________________
244 void AliITSDetTypeSim::SetSimulationModel(Int_t dettype,AliITSsimulation *sim){
246 //Set simulation model for detector type
248 if(fSimulation==0) fSimulation = new TObjArray(fgkNdettypes);
249 fSimulation->AddAt(sim,dettype);
251 //______________________________________________________________________
252 AliITSsimulation* AliITSDetTypeSim::GetSimulationModel(Int_t dettype){
254 //Get simulation model for detector type
256 Warning("GetSimulationModel","fSimulation is 0!");
259 return (AliITSsimulation*)(fSimulation->At(dettype));
261 //______________________________________________________________________
262 AliITSsimulation* AliITSDetTypeSim::GetSimulationModelByModule(Int_t module){
264 //Get simulation model by module number
265 if(GetITSgeom()==0) {
266 Warning("GetSimulationModelByModule","GetITSgeom() is 0!");
270 return GetSimulationModel(GetITSgeom()->GetModuleType(module));
272 //_______________________________________________________________________
273 void AliITSDetTypeSim::SetDefaultSegmentation(Int_t idet){
274 // Set default segmentation model objects
275 AliITSsegmentation *seg;
277 if(fSegmentation==0x0){
278 fSegmentation = new TObjArray(fgkNdettypes);
279 fSegmentation->SetOwner(kTRUE);
281 if(GetSegmentationModel(idet))
282 delete (AliITSsegmentation*)fSegmentation->At(idet);
284 seg = new AliITSsegmentationSPD();
286 seg = new AliITSsegmentationSDD();
287 AliITSCalibrationSDD* cal=(AliITSCalibrationSDD*)GetCalibrationModel(fgkDefaultNModulesSPD+1);
288 if(cal->IsAMAt20MHz()){
289 seg->SetPadSize(seg->Dpz(0),20.);
290 seg->SetNPads(seg->Npz()/2,128);
293 seg = new AliITSsegmentationSSD();
295 SetSegmentationModel(idet,seg);
297 //______________________________________________________________________
298 void AliITSDetTypeSim::SetSegmentationModel(Int_t dettype,
299 AliITSsegmentation *seg){
301 //Set segmentation model for detector type
302 if(fSegmentation==0x0){
303 fSegmentation = new TObjArray(fgkNdettypes);
304 fSegmentation->SetOwner(kTRUE);
306 fSegmentation->AddAt(seg,dettype);
308 //______________________________________________________________________
309 AliITSsegmentation* AliITSDetTypeSim::GetSegmentationModel(Int_t dettype){
310 //Get segmentation model for detector type
312 if(fSegmentation==0) {
313 Warning("GetSegmentationModel","fSegmentation is 0!");
316 return (AliITSsegmentation*)(fSegmentation->At(dettype));
318 //_______________________________________________________________________
319 AliITSsegmentation* AliITSDetTypeSim::GetSegmentationModelByModule(Int_t module){
320 //Get segmentation model by module number
322 Warning("GetSegmentationModelByModule","GetITSgeom() is 0!");
325 return GetSegmentationModel(GetITSgeom()->GetModuleType(module));
327 //_______________________________________________________________________
328 void AliITSDetTypeSim::CreateCalibrationArray() {
329 //Create the container of calibration functions with correct size
331 Warning("CreateCalibration","pointer to calibration object exists\n");
332 fCalibration->Delete();
336 Int_t nModTot = GetITSgeom()->GetIndexMax();
337 fCalibration = new TObjArray(nModTot);
338 fCalibration->SetOwner(kTRUE);
339 fCalibration->Clear();
341 //_______________________________________________________________________
342 void AliITSDetTypeSim::SetCalibrationModel(Int_t iMod, AliITSCalibration *resp){
343 //Set response model for modules
345 if (fCalibration==0) CreateCalibrationArray();
347 if (fCalibration->At(iMod)!=0)
348 delete (AliITSCalibration*) fCalibration->At(iMod);
349 fCalibration->AddAt(resp, iMod);
351 //______________________________________________________________________
352 void AliITSDetTypeSim::ResetCalibrationArray(){
353 //resets response array
354 if(fCalibration && fRunNumber<0){ // if fRunNumber<0 fCalibration is owner
355 AliITSresponse* rspd = ((AliITSCalibration*)fCalibration->At(
356 GetITSgeom()->GetStartSPD()))->GetResponse();
357 AliITSresponse* rssd = ((AliITSCalibration*)fCalibration->At(
358 GetITSgeom()->GetStartSSD()))->GetResponse();
359 if(rspd) delete rspd;
360 if(rssd) delete rssd;
361 fCalibration->Clear();
362 }else if (fCalibration && fRunNumber>=0){
363 fCalibration->Clear();
366 //______________________________________________________________________
367 void AliITSDetTypeSim::ResetSegmentation(){
368 //Resets segmentation array
369 if(fSegmentation) fSegmentation->Clear();
371 //_______________________________________________________________________
372 AliITSCalibration* AliITSDetTypeSim::GetCalibrationModel(Int_t iMod){
373 //Get response model for module number iMod
375 if(fCalibration==0) {
376 AliError("fCalibration is 0!");
379 if(iMod<fgkDefaultNModulesSPD+fgkDefaultNModulesSDD){
380 return (AliITSCalibration*)fCalibration->At(iMod);
382 Int_t i=iMod-(fgkDefaultNModulesSPD+fgkDefaultNModulesSDD);
383 fSSDCalibration->SetModule(i);
384 return (AliITSCalibration*)fSSDCalibration;
388 //_______________________________________________________________________
389 void AliITSDetTypeSim::SetDefaults(){
390 //Set defaults for segmentation and response
393 Warning("SetDefaults","GetITSgeom() is 0!");
396 if (fCalibration==0) {
397 CreateCalibrationArray();
401 if(!GetCalibration()){AliFatal("Exit"); exit(0);}
403 SetDigitClassName(0,"AliITSdigitSPD");
404 SetDigitClassName(1,"AliITSdigitSDD");
405 SetDigitClassName(2,"AliITSdigitSSD");
407 for(Int_t idet=0;idet<fgkNdettypes;idet++){
408 if(!GetSegmentationModel(idet)) SetDefaultSegmentation(idet);
411 //______________________________________________________________________
412 Bool_t AliITSDetTypeSim::GetCalibration() {
413 // Get Default calibration if a storage is not defined.
416 AliITSCalibration* cal = GetCalibrationModel(0);
423 SetRunNumber((Int_t)AliCDBManager::Instance()->GetRun());
424 Int_t run=GetRunNumber();
426 Bool_t origCacheStatus = AliCDBManager::Instance()->GetCacheFlag();
427 Bool_t isCacheActive = kTRUE;
428 if(GetRunNumber()<0){
429 isCacheActive=kFALSE;
430 fCalibration->SetOwner(kTRUE);
433 fCalibration->SetOwner(kFALSE);
436 AliCDBManager::Instance()->SetCacheFlag(isCacheActive);
438 AliCDBEntry *entrySPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDDead", run);
439 AliCDBEntry *entrySDD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSDD", run);
440 AliCDBEntry *drSpSDD = AliCDBManager::Instance()->Get("ITS/Calib/DriftSpeedSDD",run);
441 AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD",run);
442 AliCDBEntry *mapASDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsAnodeSDD",run);
443 AliCDBEntry *mapTSDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsTimeSDD",run);
444 // AliCDBEntry *entrySSD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSSD", run);
445 AliCDBEntry *entryNoiseSSD = AliCDBManager::Instance()->Get("ITS/Calib/NoiseSSD");
446 AliCDBEntry *entryGainSSD = AliCDBManager::Instance()->Get("ITS/Calib/GainSSD");
447 AliCDBEntry *entryBadChannelsSSD = AliCDBManager::Instance()->Get("ITS/Calib/BadChannelsSSD");
449 AliCDBEntry *entry2SPD = AliCDBManager::Instance()->Get("ITS/Calib/RespSPD", run);
450 AliCDBEntry *entry2SSD = AliCDBManager::Instance()->Get("ITS/Calib/RespSSD", run);
452 if(!entrySPD || !entrySDD || !entryNoiseSSD || !entryGainSSD || !entryBadChannelsSSD ||
453 !entry2SPD || !entry2SSD || !drSpSDD || !ddlMapSDD || !mapASDD ||!mapTSDD){
454 AliFatal("Calibration object retrieval failed! ");
458 // if(!entrySPD || !entrySDD || !entryNoiseSSD || !entryGainSSD || !entryBadChannelsSSD ||
459 // !entry2SPD || !entry2SSD){
460 // AliFatal("Calibration object retrieval failed! ");
464 TObjArray *calSPD = (TObjArray *)entrySPD->GetObject();
465 if(!isCacheActive)entrySPD->SetObject(NULL);
466 entrySPD->SetOwner(kTRUE);
468 AliITSresponseSPD *pSPD = (AliITSresponseSPD*)entry2SPD->GetObject();
469 if(!isCacheActive)entry2SPD->SetObject(NULL);
470 entry2SPD->SetOwner(kTRUE);
472 TObjArray *calSDD = (TObjArray *)entrySDD->GetObject();
473 if(!isCacheActive)entrySDD->SetObject(NULL);
474 entrySDD->SetOwner(kTRUE);
476 TObjArray *drSp = (TObjArray *)drSpSDD->GetObject();
477 if(!isCacheActive)drSpSDD->SetObject(NULL);
478 drSpSDD->SetOwner(kTRUE);
480 AliITSDDLModuleMapSDD *ddlsdd=(AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject();
481 if(!isCacheActive)ddlMapSDD->SetObject(NULL);
482 ddlMapSDD->SetOwner(kTRUE);
484 TObjArray *mapAn = (TObjArray *)mapASDD->GetObject();
485 if(!isCacheActive)mapASDD->SetObject(NULL);
486 mapASDD->SetOwner(kTRUE);
488 TObjArray *mapT = (TObjArray *)mapTSDD->GetObject();
489 if(!isCacheActive)mapTSDD->SetObject(NULL);
490 mapTSDD->SetOwner(kTRUE);
493 TObjArray *calSSD = (TObjArray *)entrySSD->GetObject();
494 if(!isCacheActive)entrySSD->SetObject(NULL);
495 entrySSD->SetOwner(kTRUE);
498 TObject *emptyssd = 0; TString ssdobjectname = 0;
499 AliITSNoiseSSDv2 *noiseSSD = new AliITSNoiseSSDv2();
500 emptyssd = (TObject *)entryNoiseSSD->GetObject();
501 ssdobjectname = emptyssd->GetName();
502 if(ssdobjectname=="TObjArray") {
503 TObjArray *noiseSSDOld = (TObjArray *)entryNoiseSSD->GetObject();
504 ReadOldSSDNoise(noiseSSDOld, noiseSSD);
506 else if(ssdobjectname=="AliITSNoiseSSDv2")
507 noiseSSD = (AliITSNoiseSSDv2 *)entryNoiseSSD->GetObject();
508 if(!isCacheActive)entryNoiseSSD->SetObject(NULL);
509 entryNoiseSSD->SetOwner(kTRUE);
511 AliITSGainSSDv2 *gainSSD = new AliITSGainSSDv2();
512 emptyssd = (TObject *)entryGainSSD->GetObject();
513 ssdobjectname = emptyssd->GetName();
514 if(ssdobjectname=="Gain") {
515 TObjArray *gainSSDOld = (TObjArray *)entryGainSSD->GetObject();
516 ReadOldSSDGain(gainSSDOld, gainSSD);
518 else if(ssdobjectname=="AliITSGainSSDv2")
519 gainSSD = (AliITSGainSSDv2 *)entryGainSSD->GetObject();
520 if(!isCacheActive)entryGainSSD->SetObject(NULL);
521 entryGainSSD->SetOwner(kTRUE);
523 AliITSBadChannelsSSDv2 *badChannelsSSD = new AliITSBadChannelsSSDv2();
524 emptyssd = (TObject *)entryBadChannelsSSD->GetObject();
525 ssdobjectname = emptyssd->GetName();
526 if(ssdobjectname=="TObjArray") {
527 TObjArray *badChannelsSSDOld = (TObjArray *)entryBadChannelsSSD->GetObject();
528 ReadOldSSDBadChannels(badChannelsSSDOld, badChannelsSSD);
530 else if(ssdobjectname=="AliITSBadChannelsSSDv2")
531 badChannelsSSD = (AliITSBadChannelsSSDv2*)entryBadChannelsSSD->GetObject();
532 if(!isCacheActive)entryBadChannelsSSD->SetObject(NULL);
533 entryBadChannelsSSD->SetOwner(kTRUE);
535 /*AliITSNoiseSSDv2 *noiseSSD = (AliITSNoiseSSDv2 *)entryNoiseSSD->GetObject();
536 if(!isCacheActive)entryNoiseSSD->SetObject(NULL);
537 entryNoiseSSD->SetOwner(kTRUE);
539 AliITSGainSSDv2 *gainSSD = (AliITSGainSSDv2 *)entryGainSSD->GetObject();
540 if(!isCacheActive)entryGainSSD->SetObject(NULL);
541 entryGainSSD->SetOwner(kTRUE);
543 AliITSBadChannelsSSDv2 *badchannelsSSD =
544 (AliITSBadChannelsSSDv2 *)entryBadChannelsSSD->GetObject();
545 if(!isCacheActive)entryBadChannelsSSD->SetObject(NULL);
546 entryBadChannelsSSD->SetOwner(kTRUE);*/
548 AliITSresponseSSD *pSSD = (AliITSresponseSSD*)entry2SSD->GetObject();
549 if(!isCacheActive)entry2SSD->SetObject(NULL);
550 entry2SSD->SetOwner(kTRUE);
552 // DB entries are deleted. In this way metadeta objects are deleted as well
556 delete entryNoiseSSD;
558 delete entryBadChannelsSSD;
567 AliCDBManager::Instance()->SetCacheFlag(origCacheStatus);
569 if ((!pSPD)||(!pSSD) || (!calSPD) || (!calSDD) || (!drSp) || (!ddlsdd)
570 || (!mapAn) || (!mapT) || (!noiseSSD)|| (!gainSSD)|| (!badChannelsSSD)) {
571 AliWarning("Can not get calibration from calibration database !");
575 fNMod[0] = calSPD->GetEntries();
576 fNMod[1] = calSDD->GetEntries();
577 // fNMod[2] = noiseSSD->GetEntries();
578 AliInfo(Form("%i SPD, %i SDD and %i SSD in calibration database",
579 fNMod[0], fNMod[1], fNMod[2]));
580 AliITSCalibration* cal;
581 for (Int_t i=0; i<fNMod[0]; i++) {
582 cal = (AliITSCalibration*) calSPD->At(i);
583 cal->SetResponse(pSPD);
584 SetCalibrationModel(i, cal);
587 fDDLMapSDD->SetDDLMap(ddlsdd);
589 for (Int_t i=0; i<fgkDefaultNModulesSDD; i++) {
591 Int_t iMod = i + fgkDefaultNModulesSPD;
592 fDDLMapSDD->FindInDDLMap(iMod,iddl,icarlos);
594 AliITSCalibrationSDD* calsdddead=new AliITSCalibrationSDD();
595 calsdddead->SetBad();
596 AliITSDriftSpeedSDD* driftspdef = new AliITSDriftSpeedSDD();
597 AliITSDriftSpeedArraySDD* arrdrsp=new AliITSDriftSpeedArraySDD(1);
598 arrdrsp->AddDriftSpeed(driftspdef);
599 calsdddead->SetDriftSpeed(0,arrdrsp);
600 calsdddead->SetDriftSpeed(1,arrdrsp);
601 SetCalibrationModel(iMod, calsdddead);
602 AliWarning(Form("SDD module %d not present in DDL map: set it as dead",iMod));
604 cal = (AliITSCalibration*) calSDD->At(i);
607 AliITSDriftSpeedArraySDD* arr0 = (AliITSDriftSpeedArraySDD*) drSp->At(i0);
608 AliITSMapSDD* ma0 = (AliITSMapSDD*)mapAn->At(i0);
609 AliITSMapSDD* mt0 = (AliITSMapSDD*)mapT->At(i0);
610 AliITSDriftSpeedArraySDD* arr1 = (AliITSDriftSpeedArraySDD*) drSp->At(i1);
611 AliITSMapSDD* ma1 = (AliITSMapSDD*)mapAn->At(i1);
612 AliITSMapSDD* mt1 = (AliITSMapSDD*)mapT->At(i1);
613 cal->SetDriftSpeed(0,arr0);
614 cal->SetDriftSpeed(1,arr1);
619 SetCalibrationModel(iMod, cal);
623 fSSDCalibration->SetResponse((AliITSresponse*)pSSD);
624 fSSDCalibration->SetNoise(noiseSSD);
625 fSSDCalibration->SetGain(gainSSD);
626 fSSDCalibration->SetBadChannels(badChannelsSSD);
627 //fSSDCalibration->FillBadChipMap();
631 for (Int_t i=0; i<fNMod[2]; i++) {
632 AliITSCalibrationSSD *calibSSD = new AliITSCalibrationSSD();
633 calibSSD->SetResponse((AliITSresponse*)pSSD);
635 AliITSNoiseSSD *noise = (AliITSNoiseSSD*) (noiseSSD->At(i));
636 calibSSD->SetNoise(noise);
637 AliITSGainSSD *gain = (AliITSGainSSD*) (gainSSD->At(i));
638 calibSSD->SetGain(gain);
639 AliITSBadChannelsSSD *bad = (AliITSBadChannelsSSD*) (badchannelsSSD->At(i));
640 calibSSD->SetBadChannels(bad);
642 Int_t iMod = i + fgkDefaultNModulesSPD + fgkDefaultNModulesSDD;
643 SetCalibrationModel(iMod, calibSSD);
651 //_______________________________________________________________________
652 void AliITSDetTypeSim::SetDefaultSimulation(){
653 //Set default simulation for detector type
656 Warning("SetDefaultSimulation","GetITSgeom() is 0!");
660 Warning("SetDefaultSimulation","fCalibration is 0!");
663 if(fSegmentation==0){
664 Warning("SetDefaultSimulation","fSegmentation is 0!");
665 for(Int_t i=0;i<fgkNdettypes;i++) SetDefaultSegmentation(i);
666 }else for(Int_t i=0;i<fgkNdettypes;i++) if(!GetSegmentationModel(i)){
667 Warning("SetDefaultSimulation",
668 "Segmentation not defined for det %d - Default taken\n!",i);
669 SetDefaultSegmentation(i);
671 AliITSsimulation* sim;
673 for(Int_t idet=0;idet<fgkNdettypes;idet++){
676 sim = GetSimulationModel(idet);
678 sim = new AliITSsimulationSPD(this);
679 SetSimulationModel(idet,sim);
684 sim = GetSimulationModel(idet);
686 sim = new AliITSsimulationSDD(this);
687 SetSimulationModel(idet,sim);
692 sim = GetSimulationModel(idet);
694 sim = new AliITSsimulationSSD(this);
695 SetSimulationModel(idet,sim);
701 //___________________________________________________________________
702 void AliITSDetTypeSim::SetTreeAddressS(TTree* treeS, Char_t* name){
703 // Set branch address for the ITS summable digits Trees.
710 sprintf(branchname,"%s",name);
711 branch = treeS->GetBranch(branchname);
712 TClonesArray *sdigi = &fSDigits;
713 if (branch) branch->SetAddress(&sdigi);
716 //___________________________________________________________________
717 void AliITSDetTypeSim::SetTreeAddressD(TTree* treeD, Char_t* name){
718 // Set branch address for the digit Trees.
720 const char *det[3] = {"SPD","SDD","SSD"};
729 fDigits = new TObjArray(fgkNdettypes);
731 for(Int_t i=0;i<fgkNdettypes;i++){
732 Char_t* digclass = GetDigitClassName(i);
734 if(i==0) SetDigitClassName(i,"AliITSdigitSPD");
735 if(i==1) SetDigitClassName(i,"AliITSdigitSDD");
736 if(i==2) SetDigitClassName(i,"AliITSdigitSSD");
737 digclass = GetDigitClassName(i);
739 TString classn = digclass;
740 if(!(fDigits->At(i))){
741 fDigits->AddAt(new TClonesArray(classn.Data(),1000),i);
746 if(fgkNdettypes==3) sprintf(branchname,"%sDigits%s",name,det[i]);
747 else sprintf(branchname,"%sDigits%d",name,i+1);
749 branch = treeD->GetBranch(branchname);
750 if(branch) branch->SetAddress(&((*fDigits)[i]));
755 //___________________________________________________________________
756 void AliITSDetTypeSim::ResetDigits(){
757 // Reset number of digits and the digits array for the ITS detector.
760 Error("ResetDigits","fDigits is null!");
763 for(Int_t i=0;i<fgkNdettypes;i++){
767 //___________________________________________________________________
768 void AliITSDetTypeSim::ResetDigits(Int_t branch){
769 // Reset number of digits and the digits array for this branch.
771 if(fDigits->At(branch)){
772 ((TClonesArray*)fDigits->At(branch))->Clear();
774 if(fNDigits) fNDigits[branch]=0;
780 //_______________________________________________________________________
781 void AliITSDetTypeSim::SDigitsToDigits(Option_t* opt, Char_t* name){
782 // Standard Summable digits to Digits function.
784 Warning("SDigitsToDigits","GetITSgeom() is null!!");
788 const char *all = strstr(opt,"All");
789 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
791 if( !det[0] && !det[1] && !det[2] ) all = "All";
793 static Bool_t setDef = kTRUE;
794 if(setDef) SetDefaultSimulation();
797 AliITSsimulation *sim =0;
798 TTree* trees = fLoader->TreeS();
799 if( !(trees && GetSDigits()) ){
800 Error("SDigits2Digits","Error: No trees or SDigits. Returning.");
803 sprintf(name,"%s",name);
804 TBranch* brchSDigits = trees->GetBranch(name);
807 for(Int_t module=0;module<GetITSgeom()->GetIndexMax();module++){
808 id = GetITSgeom()->GetModuleType(module);
809 if (!all && !det[id]) continue;
810 sim = (AliITSsimulation*)GetSimulationModel(id);
812 Error("SDigit2Digits","The simulation class was not "
813 "instanciated for module %d type %s!",module,
814 GetITSgeom()->GetModuleTypeName(module));
817 sim->InitSimulationModule(module,gAlice->GetEvNumber());
820 brchSDigits->GetEvent(module);
821 sim->AddSDigitsToModule(&fSDigits,0);
822 sim->FinishSDigitiseModule();
823 fLoader->TreeD()->Fill();
826 fLoader->TreeD()->GetEntries();
827 fLoader->TreeD()->AutoSave();
828 fLoader->TreeD()->Reset();
830 //_________________________________________________________
831 void AliITSDetTypeSim::AddSumDigit(AliITSpListItem &sdig){
832 //Adds the module full of summable digits to the summable digits tree.
834 new(fSDigits[fNSDigits++]) AliITSpListItem(sdig);
836 //__________________________________________________________
837 void AliITSDetTypeSim::AddSimDigit(Int_t branch, AliITSdigit* d){
838 // Add a simulated digit.
840 TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
843 new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
846 new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
849 new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
853 //______________________________________________________________________
854 void AliITSDetTypeSim::AddSimDigit(Int_t branch,Float_t phys,Int_t *digits,
855 Int_t *tracks,Int_t *hits,Float_t *charges,
857 // Add a simulated digit to the list.
859 TClonesArray &ldigits = *((TClonesArray*)fDigits->At(branch));
862 new(ldigits[fNDigits[branch]++]) AliITSdigitSPD(digits,tracks,hits);
865 new(ldigits[fNDigits[branch]++]) AliITSdigitSDD(phys,digits,tracks,
866 hits,charges,sigexpanded);
869 new(ldigits[fNDigits[branch]++]) AliITSdigitSSD(digits,tracks,hits);
873 //______________________________________________________________________
874 void AliITSDetTypeSim::StoreCalibration(Int_t firstRun, Int_t lastRun,
875 AliCDBMetaData &md) {
876 // Store calibration in the calibration database
877 // The database must be created in an external piece of code (i.e.
878 // a configuration macro )
880 if(!AliCDBManager::Instance()->IsDefaultStorageSet()) {
881 AliWarning("No storage set! Will use dummy one");
882 AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT");
886 AliError("AliITSCalibration classes are not defined - nothing done");
889 AliCDBId idRespSPD("ITS/Calib/SPDDead",firstRun, lastRun);
890 AliCDBId idRespSDD("ITS/Calib/CalibSDD",firstRun, lastRun);
891 AliCDBId idRespSSD("ITS/Calib/CalibSSD",firstRun, lastRun);
893 TObjArray respSPD(fNMod[0]);
894 TObjArray respSDD(fNMod[1]-fNMod[0]);
895 TObjArray respSSD(fNMod[2]-fNMod[1]);
896 respSPD.SetOwner(kFALSE);
897 respSSD.SetOwner(kFALSE);
898 respSSD.SetOwner(kFALSE);
900 Int_t index[fgkNdettypes];
901 for (Int_t i = 0; i<fgkNdettypes; i++ ) {
903 for (Int_t j = 0; j<=i; j++ )
907 for (Int_t i = 0; i<index[0]; i++ )
908 respSPD.Add(fCalibration->At(i));
910 for (Int_t i = index[0]; i<index[1]; i++ )
911 respSDD.Add(fCalibration->At(i));
913 for (Int_t i = index[1]; i<index[2]; i++ )
914 respSSD.Add(fCalibration->At(i));
916 AliCDBManager::Instance()->Put(&respSPD, idRespSPD, &md);
917 AliCDBManager::Instance()->Put(&respSDD, idRespSDD, &md);
918 AliCDBManager::Instance()->Put(&respSSD, idRespSSD, &md);
921 //______________________________________________________________________
922 void AliITSDetTypeSim::ReadOldSSDNoise(TObjArray *array,
923 AliITSNoiseSSDv2 *noiseSSD) {
924 //Reads the old SSD calibration object and converts it to the new format
925 const Int_t fgkSSDSTRIPSPERMODULE = 1536;
926 const Int_t fgkSSDPSIDESTRIPSPERMODULE = 768;
928 Int_t fNMod = array->GetEntries();
929 cout<<"Converting old calibration object for noise..."<<endl;
932 Double_t noise = 0.0;
933 for (Int_t iModule = 0; iModule < fNMod; iModule++) {
934 AliITSNoiseSSD *noiseModule = (AliITSNoiseSSD*) (array->At(iModule));
935 for(Int_t iStrip = 0; iStrip < fgkSSDSTRIPSPERMODULE; iStrip++) {
936 noise = (iStrip < fgkSSDPSIDESTRIPSPERMODULE) ? noiseModule->GetNoiseP(iStrip) : noiseModule->GetNoiseN(1535 - iStrip);
937 if(iStrip < fgkSSDPSIDESTRIPSPERMODULE)
938 noiseSSD->AddNoiseP(iModule,iStrip,noise);
939 if(iStrip >= fgkSSDPSIDESTRIPSPERMODULE)
940 noiseSSD->AddNoiseN(iModule,1535 - iStrip,noise);
945 //______________________________________________________________________
946 void AliITSDetTypeSim::ReadOldSSDBadChannels(TObjArray *array,
947 AliITSBadChannelsSSDv2 *badChannelsSSD) {
948 //Reads the old SSD calibration object and converts it to the new format
949 Int_t fNMod = array->GetEntries();
950 cout<<"Converting old calibration object for bad channels..."<<endl;
951 for (Int_t iModule = 0; iModule < fNMod; iModule++) {
952 //for (Int_t iModule = 0; iModule < 1; iModule++) {
953 AliITSBadChannelsSSD *bad = (AliITSBadChannelsSSD*) (array->At(iModule));
954 TArrayI arrayPSide = bad->GetBadPChannelsList();
955 for(Int_t iPCounter = 0; iPCounter < arrayPSide.GetSize(); iPCounter++)
956 badChannelsSSD->AddBadChannelP(iModule,
958 (Char_t)arrayPSide.At(iPCounter));
960 TArrayI arrayNSide = bad->GetBadNChannelsList();
961 for(Int_t iNCounter = 0; iNCounter < arrayNSide.GetSize(); iNCounter++)
962 badChannelsSSD->AddBadChannelN(iModule,
964 (Char_t)arrayNSide.At(iNCounter));
969 //______________________________________________________________________
970 void AliITSDetTypeSim::ReadOldSSDGain(TObjArray *array,
971 AliITSGainSSDv2 *gainSSD) {
972 //Reads the old SSD calibration object and converts it to the new format
974 Int_t fNMod = array->GetEntries();
975 cout<<"Converting old calibration object for gain..."<<endl;
978 for (Int_t iModule = 0; iModule < fNMod; iModule++) {
979 AliITSGainSSD *gainModule = (AliITSGainSSD*) (array->At(iModule));
980 TArrayF arrayPSide = gainModule->GetGainP();
981 for(Int_t iPCounter = 0; iPCounter < arrayPSide.GetSize(); iPCounter++)
982 gainSSD->AddGainP(iModule,
984 arrayPSide.At(iPCounter));
985 TArrayF arrayNSide = gainModule->GetGainN();
986 for(Int_t iNCounter = 0; iNCounter < arrayNSide.GetSize(); iNCounter++)
987 gainSSD->AddGainN(iModule,
989 arrayNSide.At(iNCounter));