1 /***************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Conributors 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 **************************************************************************/
20 ////////////////////////////////////////////////////////////////////////
21 // This class defines the "Standard" reconstruction for the ITS //
24 ////////////////////////////////////////////////////////////////////////
25 #include "TObjArray.h"
28 #include "AliCDBManager.h"
29 #include "AliCDBEntry.h"
30 #include "AliITSClusterFinder.h"
31 #include "AliITSClusterFinderV2SPD.h"
32 #include "AliITSClusterFinderV2SDD.h"
33 #include "AliITSClusterFinderV2SSD.h"
34 #include "AliITSDetTypeRec.h"
35 #include "AliITSDDLModuleMapSDD.h"
36 #include "AliITSRecPoint.h"
37 #include "AliITSCalibrationSDD.h"
38 #include "AliITSMapSDD.h"
39 #include "AliITSHLTforSDD.h"
40 #include "AliITSCalibrationSSD.h"
41 #include "AliITSNoiseSSDv2.h"
42 #include "AliITSGainSSDv2.h"
43 #include "AliITSBadChannelsSSDv2.h"
44 #include "AliITSNoiseSSD.h"
45 #include "AliITSGainSSD.h"
46 #include "AliITSBadChannelsSSD.h"
47 #include "AliITSresponseSDD.h"
48 #include "AliITSsegmentationSPD.h"
49 #include "AliITSsegmentationSDD.h"
50 #include "AliITSsegmentationSSD.h"
52 #include "AliITSRawStreamSPD.h"
53 #include "AliITSTriggerConditions.h"
54 #include "AliITSFOSignalsSPD.h"
55 #include "AliRunLoader.h"
56 #include "AliDataLoader.h"
57 #include "AliITSLoader.h"
59 class AliITSDriftSpeedArraySDD;
60 class AliITSCorrMapSDD;
61 class AliITSRecoParam;
63 const Int_t AliITSDetTypeRec::fgkNdettypes = 3;
64 const Int_t AliITSDetTypeRec::fgkDefaultNModulesSPD = 240;
65 const Int_t AliITSDetTypeRec::fgkDefaultNModulesSDD = 260;
66 const Int_t AliITSDetTypeRec::fgkDefaultNModulesSSD = 1698;
68 ClassImp(AliITSDetTypeRec)
70 //________________________________________________________________
71 AliITSDetTypeRec::AliITSDetTypeRec(): TObject(),
79 fTriggerConditions(0),
90 fFastOrFiredMap(1200){
91 // Standard Constructor
99 fReconstruction = new TObjArray(fgkNdettypes);
100 fDigits = new TObjArray(fgkNdettypes);
101 for(Int_t i=0; i<3; i++){
104 fSSDCalibration=new AliITSCalibrationSSD();
105 fNMod = new Int_t [fgkNdettypes];
106 fNMod[0] = fgkDefaultNModulesSPD;
107 fNMod[1] = fgkDefaultNModulesSDD;
108 fNMod[2] = fgkDefaultNModulesSSD;
109 fRecPoints = new TClonesArray("AliITSRecPoint",3000);
115 //______________________________________________________________________
116 AliITSDetTypeRec::AliITSDetTypeRec(const AliITSDetTypeRec & rec):TObject(rec),
118 fITSgeom(rec.fITSgeom),
119 fReconstruction(rec.fReconstruction),
120 fSegmentation(rec.fSegmentation),
121 fCalibration(rec.fCalibration),
122 fSSDCalibration(rec.fSSDCalibration),
123 fSPDDead(rec.fSPDDead),
124 fTriggerConditions(rec.fTriggerConditions),
125 fDigits(rec.fDigits),
126 fFOSignals(rec.fFOSignals),
127 fDDLMapSDD(rec.fDDLMapSDD),
128 fRespSDD(rec.fRespSDD),
129 fAveGainSDD(rec.fAveGainSDD),
130 fIsHLTmodeC(rec.fIsHLTmodeC),
131 fRecPoints(rec.fRecPoints),
132 fNRecPoints(rec.fNRecPoints),
133 fFirstcall(rec.fFirstcall),
134 fLoadOnlySPDCalib(rec.fLoadOnlySPDCalib),
135 fFastOrFiredMap(rec.fFastOrFiredMap){
140 //______________________________________________________________________
141 AliITSDetTypeRec& AliITSDetTypeRec::operator=(const AliITSDetTypeRec& source){
142 // Assignment operator.
143 this->~AliITSDetTypeRec();
144 new(this) AliITSDetTypeRec(source);
149 //_____________________________________________________________________
150 AliITSDetTypeRec::~AliITSDetTypeRec(){
154 fReconstruction->Delete();
155 delete fReconstruction;
159 fSegmentation->Delete();
160 delete fSegmentation;
164 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
165 fCalibration->Delete();
168 if(fRespSDD) delete fRespSDD;
169 if(fDDLMapSDD) delete fDDLMapSDD;
172 if(fSSDCalibration) delete fSSDCalibration;
174 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
180 if(fTriggerConditions){
181 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
182 fTriggerConditions->Delete();
183 delete fTriggerConditions;
184 fTriggerConditions = 0;
193 fRecPoints->Delete();
199 if (fITSgeom) delete fITSgeom;
203 //___________________________________________________________________
204 void AliITSDetTypeRec::SetReconstructionModel(Int_t dettype,AliITSClusterFinder *clf){
206 //Set reconstruction model for detector type
208 if(fReconstruction==0) fReconstruction = new TObjArray(fgkNdettypes);
209 if(fReconstruction->At(dettype)!=0) delete fReconstruction->At(dettype);
210 fReconstruction->AddAt(clf,dettype);
212 //______________________________________________________________________
213 AliITSClusterFinder* AliITSDetTypeRec::GetReconstructionModel(Int_t dettype) const{
215 //Get reconstruction model for detector type
216 if(fReconstruction==0) {
217 Warning("GetReconstructionModel","fReconstruction is 0!");
220 return (AliITSClusterFinder*)fReconstruction->At(dettype);
223 //______________________________________________________________________
224 void AliITSDetTypeRec::SetSegmentationModel(Int_t dettype,AliITSsegmentation *seg){
226 //Set segmentation model for detector type
228 if(fSegmentation==0) fSegmentation = new TObjArray(fgkNdettypes);
229 if(fSegmentation->At(dettype)!=0) delete fSegmentation->At(dettype);
230 fSegmentation->AddAt(seg,dettype);
233 //______________________________________________________________________
234 AliITSsegmentation* AliITSDetTypeRec::GetSegmentationModel(Int_t dettype) const {
236 //Get segmentation model for detector type
238 if(fSegmentation==0) {
239 Warning("GetSegmentationModel","fSegmentation is 0!");
242 return (AliITSsegmentation*)fSegmentation->At(dettype);
245 //_______________________________________________________________________
246 void AliITSDetTypeRec::SetCalibrationModel(Int_t iMod, AliITSCalibration *cal){
248 //Set calibration (response) for the module iMod of type dettype
249 if (fCalibration==0) {
250 fCalibration = new TObjArray(GetITSgeom()->GetIndexMax());
251 fCalibration->SetOwner(kTRUE);
252 fCalibration->Clear();
255 if (fCalibration->At(iMod) != 0)
256 delete (AliITSCalibration*) fCalibration->At(iMod);
257 fCalibration->AddAt(cal,iMod);
260 //_______________________________________________________________________
261 void AliITSDetTypeRec::SetSPDDeadModel(Int_t iMod, AliITSCalibration *cal){
263 //Set dead pixel info for the SPD module iMod
265 fSPDDead = new TObjArray(fgkDefaultNModulesSPD);
266 fSPDDead->SetOwner(kTRUE);
270 if (fSPDDead->At(iMod) != 0)
271 delete (AliITSCalibration*) fSPDDead->At(iMod);
272 fSPDDead->AddAt(cal,iMod);
274 //_______________________________________________________________________
275 AliITSCalibration* AliITSDetTypeRec::GetCalibrationModel(Int_t iMod) const {
277 //Get calibration model for module type
279 if(fCalibration==0) {
280 Warning("GetalibrationModel","fCalibration is 0!");
284 if(iMod<fgkDefaultNModulesSPD+fgkDefaultNModulesSDD){
285 return (AliITSCalibration*)fCalibration->At(iMod);
287 Int_t i=iMod-(fgkDefaultNModulesSPD+fgkDefaultNModulesSDD);
288 fSSDCalibration->SetModule(i);
289 return (AliITSCalibration*)fSSDCalibration;
293 //_______________________________________________________________________
294 AliITSCalibration* AliITSDetTypeRec::GetSPDDeadModel(Int_t iMod) const {
296 //Get SPD dead for module iMod
299 AliWarning("fSPDDead is 0!");
302 return (AliITSCalibration*)fSPDDead->At(iMod);
304 //_______________________________________________________________________
305 AliITSTriggerConditions* AliITSDetTypeRec::GetTriggerConditions() const {
306 //Get Pixel Trigger Conditions
307 if (fTriggerConditions==0) {
308 AliWarning("fTriggerConditions is 0!");
310 return fTriggerConditions;
312 //______________________________________________________________________
313 void AliITSDetTypeRec::SetTreeAddressD(TTree* const treeD){
314 // Set branch address for the tree of digits.
316 const char *det[4] = {"SPD","SDD","SSD","ITS"};
318 const Char_t* digclass;
323 if (fDigits == 0x0) fDigits = new TObjArray(fgkNdettypes);
324 for (i=0; i<fgkNdettypes; i++) {
325 digclass = GetDigitClassName(i);
326 if(!(fDigits->At(i))) {
327 fDigits->AddAt(new TClonesArray(digclass,1000),i);
332 if (fgkNdettypes==3) sprintf(branchname,"%sDigits%s",det[3],det[i]);
333 else sprintf(branchname,"%sDigits%d",det[3],i+1);
335 branch = treeD->GetBranch(branchname);
336 if (branch) branch->SetAddress(&((*fDigits)[i]));
342 //_______________________________________________________________________
343 TBranch* AliITSDetTypeRec::MakeBranchInTree(TTree* const tree,
344 const char* name, const char *classname,
345 void* address,Int_t size,Int_t splitlevel)
348 // Makes branch in given tree and diverts them to a separate file
354 Error("MakeBranchInTree","Making Branch %s Tree is NULL",name);
357 TBranch *branch = tree->GetBranch(name);
362 branch = tree->Branch(name,classname,address,size,splitlevel);
365 branch = tree->Bronch(name, "TClonesArray", address, size, splitlevel);
371 //____________________________________________________________________
372 void AliITSDetTypeRec::SetDefaults(){
374 //Set defaults for segmentation and response
377 Warning("SetDefaults","null pointer to AliITSgeomGeom !");
381 AliITSsegmentation* seg;
382 if(!GetCalibration()) {AliFatal("Exit");exit(0);}
384 for(Int_t dettype=0;dettype<fgkNdettypes;dettype++){
386 seg = new AliITSsegmentationSPD();
387 SetSegmentationModel(dettype,seg);
388 SetDigitClassName(dettype,"AliITSdigitSPD");
390 if(fLoadOnlySPDCalib==kFALSE){
392 seg = new AliITSsegmentationSDD();
393 AliITSCalibrationSDD* cal=(AliITSCalibrationSDD*)GetCalibrationModel(fgkDefaultNModulesSPD+1);
394 if(cal->IsAMAt20MHz()){
395 seg->SetPadSize(seg->Dpz(0),20.);
396 seg->SetNPads(seg->Npz()/2,128);
398 SetSegmentationModel(dettype,seg);
399 SetDigitClassName(dettype,"AliITSdigitSDD");
403 AliITSsegmentationSSD* seg2 = new AliITSsegmentationSSD();
404 SetSegmentationModel(dettype,seg2);
405 SetDigitClassName(dettype,"AliITSdigitSSD");
409 //______________________________________________________________________
410 Bool_t AliITSDetTypeRec::GetCalibration() {
411 // Get Default calibration if a storage is not defined.
414 AliITSCalibration* cal = GetCalibrationModel(0);
420 // SetRunNumber((Int_t)AliCDBManager::Instance()->GetRun());
421 // Int_t run=GetRunNumber();
423 Bool_t cacheStatus = AliCDBManager::Instance()->GetCacheFlag();
424 if (fCalibration==0) {
425 fCalibration = new TObjArray(GetITSgeom()->GetIndexMax());
426 fCalibration->SetOwner(!cacheStatus);
427 fCalibration->Clear();
430 Bool_t retCode=GetCalibrationSPD(cacheStatus);
431 if(retCode==kFALSE) return kFALSE;
433 if(fLoadOnlySPDCalib==kFALSE){
434 retCode=GetCalibrationSDD(cacheStatus);
435 if(retCode==kFALSE) return kFALSE;
436 retCode=GetCalibrationSSD(cacheStatus);
437 if(retCode==kFALSE) return kFALSE;
440 AliInfo(Form("%i SPD, %i SDD and %i SSD in calibration database",
441 fNMod[0], fNMod[1], fNMod[2]));
444 //______________________________________________________________________
445 Bool_t AliITSDetTypeRec::GetCalibrationSPD(Bool_t cacheStatus) {
446 // Get SPD calibration objects from OCDB
447 // dead pixel are not used for local reconstruction
450 AliCDBEntry *noisySPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDNoisy");
451 AliCDBEntry *deadSPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDDead");
452 AliCDBEntry *pitCond = AliCDBManager::Instance()->Get("ITS/Calib/PITConditions");
453 if(!noisySPD || !deadSPD || !pitCond ){
454 AliFatal("SPD Calibration object retrieval failed! ");
458 TObjArray *calNoisySPD = (TObjArray*) noisySPD->GetObject();
459 if (!cacheStatus) noisySPD->SetObject(NULL);
460 noisySPD->SetOwner(kTRUE);
462 TObjArray *calDeadSPD = (TObjArray*) deadSPD->GetObject();
463 if (!cacheStatus) deadSPD->SetObject(NULL);
464 deadSPD->SetOwner(kTRUE);
466 AliITSTriggerConditions *calPitCond = (AliITSTriggerConditions*) pitCond->GetObject();
467 if (!cacheStatus) pitCond->SetObject(NULL);
468 pitCond->SetOwner(kTRUE);
475 if ((!calNoisySPD) || (!calDeadSPD) || (!calPitCond)){
476 AliWarning("Can not get SPD calibration from calibration database !");
480 fNMod[0] = calNoisySPD->GetEntries();
482 AliITSCalibration* cal;
483 for (Int_t i=0; i<fNMod[0]; i++) {
484 cal = (AliITSCalibration*) calNoisySPD->At(i);
485 SetCalibrationModel(i, cal);
486 cal = (AliITSCalibration*) calDeadSPD->At(i);
487 SetSPDDeadModel(i, cal);
489 fTriggerConditions = calPitCond;
494 //______________________________________________________________________
495 Bool_t AliITSDetTypeRec::GetCalibrationSDD(Bool_t cacheStatus) {
496 // Get SDD calibration objects from OCDB
498 AliCDBEntry *entrySDD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSDD");
499 AliCDBEntry *entry2SDD = AliCDBManager::Instance()->Get("ITS/Calib/RespSDD");
500 AliCDBEntry *drSpSDD = AliCDBManager::Instance()->Get("ITS/Calib/DriftSpeedSDD");
501 AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD");
502 AliCDBEntry *hltforSDD = AliCDBManager::Instance()->Get("ITS/Calib/HLTforSDD");
503 // AliCDBEntry *mapASDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsAnodeSDD");
504 AliCDBEntry *mapTSDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsTimeSDD");
506 if(!entrySDD || !entry2SDD || !drSpSDD || !ddlMapSDD || !hltforSDD || !mapTSDD ){
507 AliFatal("SDD Calibration object retrieval failed! ");
513 TObjArray *calSDD = (TObjArray *)entrySDD->GetObject();
514 if(!cacheStatus)entrySDD->SetObject(NULL);
515 entrySDD->SetOwner(kTRUE);
517 AliITSresponseSDD *pSDD = (AliITSresponseSDD*)entry2SDD->GetObject();
518 if(!cacheStatus)entry2SDD->SetObject(NULL);
519 entry2SDD->SetOwner(kTRUE);
521 TObjArray *drSp = (TObjArray *)drSpSDD->GetObject();
522 if(!cacheStatus)drSpSDD->SetObject(NULL);
523 drSpSDD->SetOwner(kTRUE);
525 AliITSDDLModuleMapSDD *ddlsdd=(AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject();
526 if(!cacheStatus)ddlMapSDD->SetObject(NULL);
527 ddlMapSDD->SetOwner(kTRUE);
529 AliITSHLTforSDD* hltsdd=(AliITSHLTforSDD*)hltforSDD->GetObject();
530 if(!cacheStatus)hltforSDD->SetObject(NULL);
531 hltforSDD->SetOwner(kTRUE);
533 // TObjArray *mapAn = (TObjArray *)mapASDD->GetObject();
534 // if(!cacheStatus)mapASDD->SetObject(NULL);
535 // mapASDD->SetOwner(kTRUE);
537 TObjArray *mapT = (TObjArray *)mapTSDD->GetObject();
538 if(!cacheStatus)mapTSDD->SetObject(NULL);
539 mapTSDD->SetOwner(kTRUE);
542 // DB entries are deleted. In this way metadeta objects are deleted as well
553 if ((!pSDD)||(!calSDD) || (!drSp) || (!ddlsdd) || (!hltsdd) || (!mapT) ){
554 AliWarning("Can not get SDD calibration from calibration database !");
558 fNMod[1] = calSDD->GetEntries();
562 fIsHLTmodeC=hltsdd->IsHLTmodeC();
563 AliITSCalibration* cal;
566 Bool_t oldMapFormat=kFALSE;
567 TObject* objmap=(TObject*)mapT->At(0);
568 TString cname(objmap->ClassName());
569 if(cname.CompareTo("AliITSMapSDD")==0){
571 AliInfo("SDD Maps converted to new format");
573 for(Int_t iddl=0; iddl<AliITSDDLModuleMapSDD::GetNDDLs(); iddl++){
574 for(Int_t icar=0; icar<AliITSDDLModuleMapSDD::GetNModPerDDL();icar++){
575 Int_t iMod=fDDLMapSDD->GetModuleNumber(iddl,icar);
576 if(iMod==-1) continue;
577 Int_t i=iMod - fgkDefaultNModulesSPD;
578 cal = (AliITSCalibration*) calSDD->At(i);
581 for(Int_t iAnode=0;iAnode< ((AliITSCalibrationSDD*)cal)->NOfAnodes(); iAnode++){
582 if(((AliITSCalibrationSDD*)cal)->IsBadChannel(iAnode)) continue;
583 avegain+= ((AliITSCalibrationSDD*)cal)->GetChannelGain(iAnode);
586 AliITSDriftSpeedArraySDD* arr0 = (AliITSDriftSpeedArraySDD*) drSp->At(i0);
587 AliITSDriftSpeedArraySDD* arr1 = (AliITSDriftSpeedArraySDD*) drSp->At(i1);
589 AliITSCorrMapSDD* mt0 = 0;
590 AliITSCorrMapSDD* mt1 = 0;
592 AliITSMapSDD* oldmap0=(AliITSMapSDD*)mapT->At(i0);
593 AliITSMapSDD* oldmap1=(AliITSMapSDD*)mapT->At(i1);
594 mt0=oldmap0->ConvertToNewFormat();
595 mt1=oldmap1->ConvertToNewFormat();
597 mt0=(AliITSCorrMapSDD*)mapT->At(i0);
598 mt1=(AliITSCorrMapSDD*)mapT->At(i1);
600 cal->SetDriftSpeed(0,arr0);
601 cal->SetDriftSpeed(1,arr1);
604 SetCalibrationModel(iMod, cal);
607 if(nGdAnodes) fAveGainSDD=avegain/nGdAnodes;
612 //______________________________________________________________________
613 Bool_t AliITSDetTypeRec::GetCalibrationSSD(Bool_t cacheStatus) {
614 // Get SSD calibration objects from OCDB
615 // AliCDBEntry *entrySSD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSSD");
617 AliCDBEntry *entryNoiseSSD = AliCDBManager::Instance()->Get("ITS/Calib/NoiseSSD");
618 AliCDBEntry *entryGainSSD = AliCDBManager::Instance()->Get("ITS/Calib/GainSSD");
619 AliCDBEntry *entryBadChannelsSSD = AliCDBManager::Instance()->Get("ITS/Calib/BadChannelsSSD");
621 if(!entryNoiseSSD || !entryGainSSD || !entryBadChannelsSSD){
622 AliFatal("SSD Calibration object retrieval failed! ");
626 TObject *emptyssd = 0; TString ssdobjectname;
627 AliITSNoiseSSDv2 *noiseSSD = new AliITSNoiseSSDv2();
628 emptyssd = (TObject *)entryNoiseSSD->GetObject();
629 ssdobjectname = emptyssd->GetName();
630 if(ssdobjectname=="TObjArray") {
631 TObjArray *noiseSSDOld = (TObjArray *)entryNoiseSSD->GetObject();
632 ReadOldSSDNoise(noiseSSDOld, noiseSSD);
634 else if(ssdobjectname=="AliITSNoiseSSDv2")
635 noiseSSD = (AliITSNoiseSSDv2 *)entryNoiseSSD->GetObject();
636 if(!cacheStatus)entryNoiseSSD->SetObject(NULL);
637 entryNoiseSSD->SetOwner(kTRUE);
639 AliITSGainSSDv2 *gainSSD = new AliITSGainSSDv2();
640 emptyssd = (TObject *)entryGainSSD->GetObject();
641 ssdobjectname = emptyssd->GetName();
642 if(ssdobjectname=="Gain") {
643 TObjArray *gainSSDOld = (TObjArray *)entryGainSSD->GetObject();
644 ReadOldSSDGain(gainSSDOld, gainSSD);
646 else if(ssdobjectname=="AliITSGainSSDv2")
647 gainSSD = (AliITSGainSSDv2 *)entryGainSSD->GetObject();
648 if(!cacheStatus)entryGainSSD->SetObject(NULL);
649 entryGainSSD->SetOwner(kTRUE);
651 AliITSBadChannelsSSDv2 *badChannelsSSD = new AliITSBadChannelsSSDv2();
652 emptyssd = (TObject *)entryBadChannelsSSD->GetObject();
653 ssdobjectname = emptyssd->GetName();
654 if(ssdobjectname=="TObjArray") {
655 TObjArray *badChannelsSSDOld = (TObjArray *)entryBadChannelsSSD->GetObject();
656 ReadOldSSDBadChannels(badChannelsSSDOld, badChannelsSSD);
658 else if(ssdobjectname=="AliITSBadChannelsSSDv2")
659 badChannelsSSD = (AliITSBadChannelsSSDv2*)entryBadChannelsSSD->GetObject();
660 if(!cacheStatus)entryBadChannelsSSD->SetObject(NULL);
661 entryBadChannelsSSD->SetOwner(kTRUE);
663 // DB entries are deleted. In this way metadeta objects are deleted as well
665 delete entryNoiseSSD;
667 delete entryBadChannelsSSD;
670 if ((!noiseSSD)|| (!gainSSD)|| (!badChannelsSSD)) {
671 AliWarning("Can not get SSD calibration from calibration database !");
675 fSSDCalibration->SetNoise(noiseSSD);
676 fSSDCalibration->SetGain(gainSSD);
677 fSSDCalibration->SetBadChannels(badChannelsSSD);
678 //fSSDCalibration->FillBadChipMap();
683 //________________________________________________________________
684 void AliITSDetTypeRec::SetDefaultClusterFindersV2(Bool_t rawdata){
686 //Set defaults for cluster finder V2
689 Warning("SetDefaults","Null pointer to AliITSgeom !");
693 AliITSClusterFinder *clf;
695 for(Int_t dettype=0;dettype<fgkNdettypes;dettype++){
698 if(!GetReconstructionModel(dettype)){
699 clf = new AliITSClusterFinderV2SPD(this);
701 if(!rawdata) clf->SetDigits(DigitsAddress(0));
702 SetReconstructionModel(dettype,clf);
708 if(!GetReconstructionModel(dettype)){
709 clf = new AliITSClusterFinderV2SDD(this);
711 if(!rawdata) clf->SetDigits(DigitsAddress(1));
712 SetReconstructionModel(dettype,clf);
719 if(!GetReconstructionModel(dettype)){
720 clf = new AliITSClusterFinderV2SSD(this);
722 if(!rawdata) clf->SetDigits(DigitsAddress(2));
723 SetReconstructionModel(dettype,clf);
730 //______________________________________________________________________
731 void AliITSDetTypeRec::MakeBranch(TTree* tree, Option_t* option){
733 //Creates branches for clusters and recpoints
734 Bool_t cR = (strstr(option,"R")!=0);
735 Bool_t cRF = (strstr(option,"RF")!=0);
739 if(cR) MakeBranchR(tree);
740 if(cRF) MakeBranchRF(tree);
744 //___________________________________________________________________
745 void AliITSDetTypeRec::ResetDigits(){
746 // Reset number of digits and the digits array for the ITS detector.
749 for(Int_t i=0;i<fgkNdettypes;i++){
753 //___________________________________________________________________
754 void AliITSDetTypeRec::ResetDigits(Int_t branch){
755 // Reset number of digits and the digits array for this branch.
757 if(fDigits->At(branch)) ((TClonesArray*)fDigits->At(branch))->Clear();
760 //__________________________________________________________________
761 void AliITSDetTypeRec::MakeBranchR(TTree *treeR, Option_t *opt){
763 //Creates tree branches for recpoints
765 // cont char *file File name where RecPoints branch is to be written
766 // to. If blank it write the SDigits to the same
767 // file in which the Hits were found.
772 // only one branch for rec points for all detector types
773 Bool_t oFast= (strstr(opt,"Fast")!=0);
775 Char_t detname[10] = "ITS";
779 sprintf(branchname,"%sRecPointsF",detname);
781 sprintf(branchname,"%sRecPoints",detname);
784 if(!fRecPoints)fRecPoints = new TClonesArray("AliITSRecPoint",1000);
786 MakeBranchInTree(treeR,branchname,0,&fRecPoints,buffsz,99);
788 //______________________________________________________________________
789 void AliITSDetTypeRec::SetTreeAddressR(TTree* const treeR){
790 // Set branch address for the Reconstructed points Trees.
792 // TTree *treeR Tree containing the RecPoints.
798 Char_t namedet[10]="ITS";
801 if(fRecPoints==0x0) fRecPoints = new TClonesArray("AliITSRecPoint",1000);
803 sprintf(branchname,"%sRecPoints",namedet);
804 branch = treeR->GetBranch(branchname);
806 branch->SetAddress(&fRecPoints);
809 sprintf(branchname,"%sRecPointsF",namedet);
810 branch = treeR->GetBranch(branchname);
812 branch->SetAddress(&fRecPoints);
816 //____________________________________________________________________
817 void AliITSDetTypeRec::AddRecPoint(const AliITSRecPoint &r){
818 // Add a reconstructed space point to the list
820 // const AliITSRecPoint &r RecPoint class to be added to the tree
821 // of reconstructed points TreeR.
827 TClonesArray &lrecp = *fRecPoints;
828 new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
831 //______________________________________________________________________
832 void AliITSDetTypeRec::DigitsToRecPoints(TTree *treeD,TTree *treeR,Int_t lastentry,Option_t *opt, Int_t optCluFind){
833 // cluster finding and reconstruction of space points
834 // the condition below will disappear when the geom class will be
835 // initialized for all versions - for the moment it is only for v5 !
836 // 7 is the SDD beam test version
838 // TTree *treeD Digits tree
839 // TTree *treeR Clusters tree
840 // Int_t lastentry Offset for module when not all of the modules
842 // Option_t *opt String indicating which ITS sub-detectors should
843 // be processed. If ="All" then all of the ITS
844 // sub detectors are processed.
846 const char *all = strstr(opt,"All");
847 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
850 SetDefaultClusterFindersV2();
851 AliInfo("V2 cluster finder has been selected \n");
853 SetDefaultClusterFindersV2();
854 AliInfo("Cluster Finder Option not implemented, V2 cluster finder will be used \n");
858 // Reset Fast-OR fired map
859 ResetFastOrFiredMap();
861 if (all || det[0]) { // SPD present
862 // Get the FO signals for this event
863 AliRunLoader* runLoader = AliRunLoader::Instance();
864 AliITSLoader* itsLoader = (AliITSLoader*) runLoader->GetLoader("ITSLoader");
866 AliError("ITS loader is NULL.");
869 fFOSignals = (AliITSFOSignalsSPD*)itsLoader->TreeD()->GetUserInfo()->FindObject("AliITSFOSignalsSPD");
870 if(!fFOSignals) AliError("FO signals not retrieved");
876 AliITSClusterFinder *rec = 0;
877 Int_t id,module,first=0;
878 for(module=0;module<GetITSgeom()->GetIndexMax();module++){
879 id = GetITSgeom()->GetModuleType(module);
880 if (!all && !det[id]) continue;
881 if(det[id]) first = GetITSgeom()->GetStartDet(id);
882 rec = (AliITSClusterFinder*)GetReconstructionModel(id);
883 TClonesArray *itsDigits = DigitsAddress(id);
885 AliFatal("The reconstruction class was not instanciated!");
886 ResetDigits(); // MvL: Not sure we neeed this when rereading anyways
888 treeD->GetEvent(lastentry+module);
891 treeD->GetEvent(lastentry+(module-first));
893 Int_t ndigits = itsDigits->GetEntriesFast();
894 if (ndigits>0 || id==0) { // for SPD we always want to call FindRawClusters (to process FO signals)
895 rec->SetDetTypeRec(this);
896 rec->SetDigits(DigitsAddress(id));
897 // rec->SetClusters(ClustersAddress(id));
898 rec->FindRawClusters(module);
904 // Remove PIT in-active chips from Fast-OR fired map
905 if (all || det[0]) { // SPD present
906 RemoveFastOrFiredInActive();
909 //______________________________________________________________________
910 void AliITSDetTypeRec::DigitsToRecPoints(AliRawReader* rawReader,TTree *treeR,Option_t *opt){
911 // cluster finding and reconstruction of space points
912 // the condition below will disappear when the geom class will be
913 // initialized for all versions - for the moment it is only for v5 !
914 // 7 is the SDD beam test version
916 // AliRawReader *rawReader Pointer to the raw-data reader
917 // TTree *treeR Clusters tree
922 const char *all = strstr(opt,"All");
923 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
926 // Reset Fast-OR fired map
927 ResetFastOrFiredMap();
929 AliITSClusterFinder *rec = 0;
932 TClonesArray *array=new TClonesArray("AliITSRecPoint",1000);
933 TBranch *branch = treeR->Branch("ITSRecPoints",&array);
936 TClonesArray** clusters = new TClonesArray*[GetITSgeom()->GetIndexMax()];
937 for (Int_t iModule = 0; iModule < GetITSgeom()->GetIndexMax(); iModule++) {
938 clusters[iModule] = NULL;
941 if (!all && !det[id]) continue;
942 rec = (AliITSClusterFinder*)GetReconstructionModel(id);
944 AliFatal("The reconstruction class was not instantiated");
945 rec->SetDetTypeRec(this);
946 rec->RawdataToClusters(rawReader,clusters);
949 TClonesArray *emptyArray=new TClonesArray("AliITSRecPoint");
950 for(Int_t iModule=0;iModule<GetITSgeom()->GetIndexMax();iModule++){
951 id = GetITSgeom()->GetModuleType(iModule);
952 if (!all && !det[id]) continue;
953 array = clusters[iModule];
955 AliDebug(1,Form("data for module %d missing!",iModule));
958 branch->SetAddress(&array);
960 nClusters+=array->GetEntriesFast();
962 if (array != emptyArray) {
970 Info("DigitsToRecPoints", "total number of found recpoints in ITS: %d\n",
973 // Remove PIT in-active chips from Fast-OR fired map
974 if (all || det[0]) { // SPD present
975 RemoveFastOrFiredInActive();
979 //______________________________________________________________________
980 void AliITSDetTypeRec::ReadOldSSDNoise(const TObjArray *array,
981 AliITSNoiseSSDv2 *noiseSSD) {
982 //Reads the old SSD calibration object and converts it to the new format
983 const Int_t fgkSSDSTRIPSPERMODULE = 1536;
984 const Int_t fgkSSDPSIDESTRIPSPERMODULE = 768;
986 Int_t gNMod = array->GetEntries();
987 cout<<"Converting old calibration object for noise..."<<endl;
990 Double_t noise = 0.0;
991 for (Int_t iModule = 0; iModule < gNMod; iModule++) {
992 AliITSNoiseSSD *noiseModule = (AliITSNoiseSSD*) (array->At(iModule));
993 for(Int_t iStrip = 0; iStrip < fgkSSDSTRIPSPERMODULE; iStrip++) {
994 noise = (iStrip < fgkSSDPSIDESTRIPSPERMODULE) ? noiseModule->GetNoiseP(iStrip) : noiseModule->GetNoiseN(1535 - iStrip);
995 if(iStrip < fgkSSDPSIDESTRIPSPERMODULE)
996 noiseSSD->AddNoiseP(iModule,iStrip,noise);
997 if(iStrip >= fgkSSDPSIDESTRIPSPERMODULE)
998 noiseSSD->AddNoiseN(iModule,1535 - iStrip,noise);
1000 }//loop over modules
1003 //______________________________________________________________________
1004 void AliITSDetTypeRec::ReadOldSSDBadChannels(const TObjArray *array,
1005 AliITSBadChannelsSSDv2 *badChannelsSSD) {
1006 //Reads the old SSD calibration object and converts it to the new format
1007 Int_t gNMod = array->GetEntries();
1008 cout<<"Converting old calibration object for bad channels..."<<endl;
1009 for (Int_t iModule = 0; iModule < gNMod; iModule++) {
1010 //for (Int_t iModule = 0; iModule < 1; iModule++) {
1011 AliITSBadChannelsSSD *bad = (AliITSBadChannelsSSD*) (array->At(iModule));
1012 TArrayI arrayPSide = bad->GetBadPChannelsList();
1013 for(Int_t iPCounter = 0; iPCounter < arrayPSide.GetSize(); iPCounter++)
1014 badChannelsSSD->AddBadChannelP(iModule,
1016 (Char_t)arrayPSide.At(iPCounter));
1018 TArrayI arrayNSide = bad->GetBadNChannelsList();
1019 for(Int_t iNCounter = 0; iNCounter < arrayNSide.GetSize(); iNCounter++)
1020 badChannelsSSD->AddBadChannelN(iModule,
1022 (Char_t)arrayNSide.At(iNCounter));
1024 }//loop over modules
1027 //______________________________________________________________________
1028 void AliITSDetTypeRec::ReadOldSSDGain(const TObjArray *array,
1029 AliITSGainSSDv2 *gainSSD) {
1030 //Reads the old SSD calibration object and converts it to the new format
1032 Int_t gNMod = array->GetEntries();
1033 cout<<"Converting old calibration object for gain..."<<endl;
1036 for (Int_t iModule = 0; iModule < gNMod; iModule++) {
1037 AliITSGainSSD *gainModule = (AliITSGainSSD*) (array->At(iModule));
1038 TArrayF arrayPSide = gainModule->GetGainP();
1039 for(Int_t iPCounter = 0; iPCounter < arrayPSide.GetSize(); iPCounter++)
1040 gainSSD->AddGainP(iModule,
1042 arrayPSide.At(iPCounter));
1043 TArrayF arrayNSide = gainModule->GetGainN();
1044 for(Int_t iNCounter = 0; iNCounter < arrayNSide.GetSize(); iNCounter++)
1045 gainSSD->AddGainN(iModule,
1047 arrayNSide.At(iNCounter));
1048 }//loop over modules
1050 //______________________________________________________________________
1051 void AliITSDetTypeRec::RemoveFastOrFiredInActive() {
1052 // Removes the chips that were in-active in the pixel trigger (from fast-or fired map)
1053 if (fTriggerConditions==NULL) {
1054 AliError("Pixel trigger conditions are missing.");
1060 while (fTriggerConditions->GetNextInActiveChip(eq,hs,chip)) {
1061 UInt_t chipKey = AliITSRawStreamSPD::GetOfflineChipKeyFromOnline(eq,hs,chip);
1062 fFastOrFiredMap.SetBitNumber(chipKey,kFALSE);
1065 //______________________________________________________________________
1066 void AliITSDetTypeRec::SetFastOrFiredMapOnline(UInt_t eq, UInt_t hs, UInt_t chip) {
1067 // Set fast-or fired map for this chip
1068 Int_t chipKey = AliITSRawStreamSPD::GetOfflineChipKeyFromOnline(eq,hs,chip);
1069 return SetFastOrFiredMap(chipKey);