]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HMPID/AliHMPIDReconstructor.cxx
Script to create a random bad channel map.
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDReconstructor.cxx
CommitLineData
d3da6dc4 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 **************************************************************************/
423554a3 15//.
16// HMPID base class to reconstruct an event
17//.
18//.
19//.
d3da6dc4 20#include "AliHMPIDReconstructor.h" //class header
21#include "AliHMPID.h" //Reconstruct()
a1d55ff3 22#include "AliHMPIDCluster.h" //Dig2Clu()
3278403b 23#include "AliHMPIDPid.h" //FillEsd()
d3da6dc4 24#include "AliHMPIDParam.h" //FillEsd()
a0e5c1b9 25#include <AliCDBEntry.h> //ctor
26#include <AliCDBManager.h> //ctor
555a0dbf 27#include <AliESDEvent.h> //FillEsd()
a0e5c1b9 28#include <AliRawReader.h> //Reconstruct() for raw digits
65201c45 29#include <AliLog.h> //
d76c31f4 30#include "AliHMPIDRawStream.h" //ConvertDigits()
55a829a5 31#include "AliHMPIDRecoParam.h" //ctor
d3da6dc4 32ClassImp(AliHMPIDReconstructor)
33
55a829a5 34AliHMPIDRecoParam* AliHMPIDReconstructor::fgkRecoParam =0; //
94b1fbfa 35//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
a0e5c1b9 36AliHMPIDReconstructor::AliHMPIDReconstructor():AliReconstructor(),fUserCut(0),fDaqSig(0),fDig(0),fClu(0)
94b1fbfa 37{
38//
39//ctor
40//
909bfbdc 41 AliHMPIDParam::Instance(); //geometry loaded for reconstruction
14dc2e6c 42 fUserCut = new Int_t[7];
ae5a42aa 43 fClu=new TObjArray(AliHMPIDParam::kMaxCh+1); fClu->SetOwner(kTRUE);
44 fDig=new TObjArray(AliHMPIDParam::kMaxCh+1); fDig->SetOwner(kTRUE);
9cc1a411 45
ae5a42aa 46 for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++){
56c73976 47 fDig->AddAt(new TClonesArray("AliHMPIDDigit",24000),i);
48 TClonesArray *pClus = new TClonesArray("AliHMPIDCluster",24000);
611e810d 49 pClus->SetUniqueID(i);
50 fClu->AddAt(pClus,i);
94b1fbfa 51 }
55a829a5 52
65201c45 53 if(fgkRecoParam!=0x0 && fgkRecoParam->GetUserCutMode()==kTRUE) {
54 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) {
55a829a5 55 fUserCut[iCh] = fgkRecoParam->GetUserCut(iCh);
4f9a834e 56 AliDebug(1,Form("UserCut successfully loaded (from RecoParam) for chamber %i -> %i ",iCh,fUserCut[iCh]));
65201c45 57 }
55a829a5 58 }
59 else {
65201c45 60 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) {
61 fUserCut[iCh] = 3; // minimal requirement for sigma cut
62 AliDebug(1,Form("UserCut loaded from defaults for chamber %i -> %i ",iCh,fUserCut[iCh]));
85692374 63 }
65201c45 64 }
49881df7 65
66 AliCDBEntry *pQthreEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Qthre"); //contains TObjArray of 7 TF1
67 if(!pQthreEnt) AliFatal("No Qthre available");
68 TObjArray *pQthre = (TObjArray*)pQthreEnt->GetObject();
69 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) {
0d2b5802 70 if(pQthre->GetEntriesFast()==AliHMPIDParam::kMaxCh+1) { // for backward compatibility...
71 TF1 *pfQthre = (TF1*)pQthre->At(iCh);
72 Double_t tMin,tMax;
73 pfQthre->GetRange(tMin,tMax);
74 Double_t qthre=pfQthre->Eval(tMin);
75 AliDebug(1,Form("Qthre successfully loaded for chamber %i -> %f ",iCh,qthre));
76 } else {
77 for(Int_t isec=0;isec<=5;isec++) {
78 TF1 *pfQthre = (TF1*)pQthre->At(6*iCh+isec);
79 Double_t tMin,tMax;
80 pfQthre->GetRange(tMin,tMax);
81 Double_t qthre=pfQthre->Eval(tMin);
82 AliDebug(1,Form("Qthre successfully loaded for chamber %i sector %i -> %f ",iCh,isec,qthre));
83 }
49881df7 84 }
85 }
86
87
9cc1a411 88 AliCDBEntry *pDaqSigEnt =AliCDBManager::Instance()->Get("HMPID/Calib/DaqSig"); //contains TObjArray of TObjArray 14 TMatrixF sigmas values for pads
89 if(!pDaqSigEnt) AliFatal("No pedestals from DAQ!");
90 fDaqSig = (TObjArray*)pDaqSigEnt->GetObject();
91 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){ //chambers loop
4f9a834e 92 AliDebug(1,Form("DaqSigCut successfully loaded for chamber %i -> %i ",iCh,(Int_t)fDaqSig->At(iCh)->GetUniqueID()));
a0e5c1b9 93 }
94b1fbfa 94}//AliHMPIDReconstructor
d3da6dc4 95//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
54104a7c 96void AliHMPIDReconstructor::Dig2Clu(TObjArray *pDigAll,TObjArray *pCluAll,Int_t *pUserCut,Bool_t isTryUnfold)
d3da6dc4 97{
98// Finds all clusters for a given digits list provided not empty. Currently digits list is a list of all digits for a single chamber.
99// Puts all found clusters in separate lists, one per clusters.
3c6274c1 100// Arguments: pDigAll - list of digits for all chambers
101// pCluAll - list of clusters for all chambers
d3da6dc4 102// isTryUnfold - flag to choose between CoG and Mathieson fitting
103// Returns: none
ae5a42aa 104 TMatrixF padMap(AliHMPIDParam::kMinPx,AliHMPIDParam::kMaxPcx,AliHMPIDParam::kMinPy,AliHMPIDParam::kMaxPcy); //pads map for single chamber 0..159 x 0..143
3c6274c1 105
ae5a42aa 106 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){ //chambers loop
3c6274c1 107 TClonesArray *pDigCur=(TClonesArray*)pDigAll->At(iCh); //get list of digits for current chamber
108 if(pDigCur->GetEntriesFast()==0) continue; //no digits for current chamber
d3da6dc4 109
3c6274c1 110 padMap=(Float_t)-1; //reset map to -1 (means no digit for this pad)
111 TClonesArray *pCluCur=(TClonesArray*)pCluAll->At(iCh); //get list of clusters for current chamber
112
113 for(Int_t iDig=0;iDig<pDigCur->GetEntriesFast();iDig++){ //digits loop to fill pads map
114 AliHMPIDDigit *pDig= (AliHMPIDDigit*)pDigCur->At(iDig); //get current digit
115 padMap( pDig->PadChX(), pDig->PadChY() )=iDig; //fill the map, (padx,pady) cell takes digit index
116 }//digits loop to fill digits map
117
118 AliHMPIDCluster clu; //tmp cluster to be used as current
3c6274c1 119 for(Int_t iDig=0;iDig<pDigCur->GetEntriesFast();iDig++){ //digits loop for current chamber
120 AliHMPIDDigit *pDig=(AliHMPIDDigit*)pDigCur->At(iDig); //take current digit
121 if(!(pDig=UseDig(pDig->PadChX(),pDig->PadChY(),pDigCur,&padMap))) continue; //this digit is already taken in FormClu(), go after next digit
122 FormClu(&clu,pDig,pDigCur,&padMap); //form cluster starting from this digit by recursion
54104a7c 123 clu.Solve(pCluCur,pUserCut,isTryUnfold); //solve this cluster and add all unfolded clusters to provided list
3c6274c1 124 clu.Reset(); //empty current cluster
125 }//digits loop for current chamber
126 }//chambers loop
d3da6dc4 127}//Dig2Clu()
128//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
129void AliHMPIDReconstructor::FormClu(AliHMPIDCluster *pClu,AliHMPIDDigit *pDig,TClonesArray *pDigLst,TMatrixF *pDigMap)
130{
a1d55ff3 131// Forms the initial cluster as a combination of all adjascent digits. Starts from the given digit then calls itself recursevly for all neighbours.
d3da6dc4 132// Arguments: pClu - pointer to cluster being formed
a1d55ff3 133// Returns: none
d3da6dc4 134 pClu->DigAdd(pDig);//take this digit in cluster
135
3c6274c1 136 Int_t cnt=0,cx[4],cy[4];
d3da6dc4 137
ae5a42aa 138 if(pDig->PadPcX() != AliHMPIDParam::kMinPx){cx[cnt]=pDig->PadChX()-1; cy[cnt]=pDig->PadChY() ;cnt++;} //left
139 if(pDig->PadPcX() != AliHMPIDParam::kMaxPx){cx[cnt]=pDig->PadChX()+1; cy[cnt]=pDig->PadChY() ;cnt++;} //right
140 if(pDig->PadPcY() != AliHMPIDParam::kMinPy){cx[cnt]=pDig->PadChX() ; cy[cnt]=pDig->PadChY()-1;cnt++;} //down
141 if(pDig->PadPcY() != AliHMPIDParam::kMaxPy){cx[cnt]=pDig->PadChX() ; cy[cnt]=pDig->PadChY()+1;cnt++;} //up
d3da6dc4 142
a1d55ff3 143 for (Int_t i=0;i<cnt;i++){//neighbours loop
3c6274c1 144 if((pDig=UseDig(cx[i],cy[i],pDigLst,pDigMap))) FormClu(pClu,pDig,pDigLst,pDigMap); //check if this neighbour pad fired and mark it as taken
a1d55ff3 145 }//neighbours loop
d3da6dc4 146}//FormClu()
147//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
94b1fbfa 148void AliHMPIDReconstructor::Reconstruct(TTree *pDigTree,TTree *pCluTree)const
d3da6dc4 149{
150//Invoked by AliReconstruction to convert digits to clusters i.e. reconstruct simulated data
94b1fbfa 151//Arguments: pDigTree - pointer to Digit tree
152// pCluTree - poitner to Cluster tree
d3da6dc4 153// Returns: none
154 AliDebug(1,"Start.");
ae5a42aa 155 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) {
56c73976 156 pCluTree->Branch(Form("HMPID%d",iCh),&((*fClu)[iCh]),7);
94b1fbfa 157 pDigTree->SetBranchAddress(Form("HMPID%d",iCh),&((*fDig)[iCh]));
158 }
159 pDigTree->GetEntry(0);
54104a7c 160 Dig2Clu(fDig,fClu,fUserCut); //cluster finder
94b1fbfa 161 pCluTree->Fill(); //fill tree for current event
162
ae5a42aa 163 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
94b1fbfa 164 fDig->At(iCh)->Clear();
165 fClu->At(iCh)->Clear();
166 }
d3da6dc4 167
d3da6dc4 168 AliDebug(1,"Stop.");
169}//Reconstruct(for simulated digits)
a1d55ff3 170//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
94b1fbfa 171void AliHMPIDReconstructor::ConvertDigits(AliRawReader *pRR,TTree *pDigTree)const
172{
173//Invoked by AliReconstruction to convert raw digits from DDL files to digits
174//Arguments: pRR - ALICE raw reader pointer
175// pDigTree - pointer to Digit tree
176// Returns: none
177 AliDebug(1,"Start.");
21f61e25 178// Int_t digcnt=0;
179
56c73976 180 Int_t iDigCnt[7]={0,0,0,0,0,0,0};
21f61e25 181
56c73976 182 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
183 pDigTree->Branch(Form("HMPID%d",iCh),&((*fDig)[iCh]));
184 }
21f61e25 185
56c73976 186 AliHMPIDRawStream stream(pRR);
187
188 while(stream.Next())
189 {
190 Int_t ch = AliHMPIDParam::DDL2C(stream.GetDDLNumber());
191 for(Int_t iPad=0;iPad<stream.GetNPads();iPad++) {
192 AliHMPIDDigit dig(stream.GetPadArray()[iPad],stream.GetChargeArray()[iPad]);
193 if(!IsDigSurvive(&dig)) continue;
194 new((*((TClonesArray*)fDig->At(ch)))[iDigCnt[ch]++]) AliHMPIDDigit(dig); //add this digit to the tmp list
21f61e25 195 }
21f61e25 196 }
56c73976 197
94b1fbfa 198 pDigTree->Fill();
21f61e25 199
ae5a42aa 200 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++)fDig->At(iCh)->Clear();
21f61e25 201
94b1fbfa 202 AliDebug(1,"Stop.");
203}//Reconstruct digits from raw digits
204//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
d76c31f4 205void AliHMPIDReconstructor::FillESD(TTree */*digitsTree*/, TTree */*clustersTree*/, AliESDEvent *pESD) const
d3da6dc4 206{
3278403b 207// Fill ESD with all the infos from HMPID
208// Probability vector from AliHMPIDPid
209//...
9eb30f08 210 AliHMPIDPid pID;
dac53a45 211 Double_t prob[AliPID::kSPECIES];
9eb30f08 212
213 for(Int_t iTrk=0;iTrk<pESD->GetNumberOfTracks();iTrk++){//ESD tracks loop
214 AliESDtrack *pTrk = pESD->GetTrack(iTrk);// get next reconstructed track
dac53a45 215 pID.FindPid(pTrk,AliPID::kSPECIES,prob);
216 pTrk->SetHMPIDpid(prob);
9eb30f08 217 }//ESD tracks loop
3278403b 218
a1d55ff3 219}//FillESD()