]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HMPID/AliHMPIDCluster.cxx
Corrected index (icc)
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDCluster.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// **************************************************************************
15
16#include "AliHMPIDCluster.h" //class header
17#include <TMinuit.h> //Solve()
18#include <TClonesArray.h> //Solve()
d1bf51e1 19#include <TMarker.h> //Draw()
d3da6dc4 20ClassImp(AliHMPIDCluster)
21//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
22void AliHMPIDCluster::CoG()
23{
24// Calculates naive cluster position as a center of gravity of its digits.
25// Arguments: none
d1bf51e1 26// Returns: none
27
28// if(fDigs==0) return; //no digits in this cluster
29 fX=fY=fQ=0; //set cluster position to (0,0) to start to collect contributions
30 Int_t maxQpad=-1,maxQ=-1; //to calculate the pad with the highest charge
31 AliHMPIDDigit *pDig;
d3da6dc4 32 for(Int_t iDig=0;iDig<fDigs->GetEntriesFast();iDig++){//digits loop
d1bf51e1 33 pDig=(AliHMPIDDigit*)fDigs->At(iDig); //get pointer to next digit
d3da6dc4 34 Float_t q=pDig->Q(); //get QDC
35 fX += pDig->LorsX()*q;fY +=pDig->LorsY()*q; //add digit center weighted by QDC
d1bf51e1 36 fQ+=q; //increment total charge
37 if(q>maxQ) {maxQpad = pDig->Pad();maxQ=(Int_t)q;} // to find pad with highest charge
d3da6dc4 38 }//digits loop
d1bf51e1 39 if ( fQ != 0 ) fX/=fQ;fY/=fQ; //final center of gravity
d3da6dc4 40
d1bf51e1 41
42 CorrSin(); //correct it by sinoid
43 fCh=pDig->Ch(); //initialize chamber number
44 fMaxQpad = maxQpad; fMaxQ=maxQ; //store max charge pad to the field
45 fXi=fX+99; fYi=fY+99; fQi=fQ+99; //initial local max position is to be shifted artificially
d3da6dc4 46 fSt=kCoG;
47}//CoG()
48//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
49void AliHMPIDCluster::CorrSin()
50{
51// Correction of cluster x position due to sinoid, see HMPID TDR page 30
52// Arguments: none
53// Returns: none
3c6274c1 54 AliHMPIDDigit dig;dig.Manual1(Ch(),fX,fY); //tmp digit to get it center
d3da6dc4 55 Float_t x=fX-dig.LorsX();
56 fX+=3.31267e-2*TMath::Sin(2*TMath::Pi()/0.8*x)-2.66575e-3*TMath::Sin(4*TMath::Pi()/0.8*x)+2.80553e-3*TMath::Sin(6*TMath::Pi()/0.8*x)+0.0070;
57}
58//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
d1bf51e1 59void AliHMPIDCluster::Draw(Option_t*)
60{
61 TMarker *pMark=new TMarker(X(),Y(),5); pMark->SetMarkerColor(kBlue); pMark->Draw();
62}
63//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
d3da6dc4 64void AliHMPIDCluster::FitFunc(Int_t &iNpars, Double_t *, Double_t &chi2, Double_t *par, Int_t )
65{
66// Cluster fit function
67// par[0]=x par[1]=y par[2]=q for the first Mathieson shape
68// par[3]=x par[4]=y par[5]=q for the second Mathieson shape and so on up to iNpars/3 Mathieson shapes
69// For each pad of the cluster calculates the difference between actual pad charge and the charge induced to this pad by all Mathieson distributions
70// Then the chi2 is calculated as the sum of this value squared for all pad in the cluster.
71// Arguments: iNpars - number of parameters which is number of local maxima of cluster * 3
72// chi2 - function result to be minimised
73// par - parameters array of size iNpars
74// Returns: none
75 AliHMPIDCluster *pClu=(AliHMPIDCluster*)gMinuit->GetObjectFit();
76 Int_t iNshape = iNpars/3;
77
78 chi2 = 0;
79 for(Int_t i=0;i<pClu->Size();i++){ //loop on all pads of the cluster
80 Double_t dQpadMath = 0; //pad charge collector
81 for(Int_t j=0;j<iNshape;j++){ //Mathiesons loop as all of them may contribute to this pad
82 dQpadMath+=par[3*j+2]*pClu->Dig(i)->Mathieson(par[3*j],par[3*j+1]); // par[3*j+2] is charge par[3*j] is x par[3*j+1] is y of current Mathieson
83 }
84 chi2 +=TMath::Power((pClu->Dig(i)->Q()-dQpadMath),2); //
85 } //loop on all pads of the cluster
86}//FitFunction()
87//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
88void AliHMPIDCluster::Print(Option_t* opt)const
89{
90//Print current cluster
91 const char *status=0;
92 switch(fSt){
d1bf51e1 93 case kFrm : status="formed " ;break;
94 case kUnf : status="unfolded (fit)" ;break;
95 case kCoG : status="coged " ;break;
96 case kLo1 : status="locmax 1 (fit)" ;break;
97 case kAbn : status="abnorm (fit)" ;break;
98 case kMax : status="exceeded (cog)" ;break;
99 case kNot : status="not done (cog)" ;break;
100 case kEmp : status="empty " ;break;
101 case kEdg : status="edge (fit)" ;break;
102 case kSi1 : status="size 1 (cog)" ;break;
103 case kNoLoc: status="no LocMax(fit)" ;break;
104
105 default: status="??????" ;break;
d3da6dc4 106 }
d1bf51e1 107 Printf("%sCLU:(%7.3f,%7.3f) Q=%8.3f ch=%i, FormedSize=%2i N loc. max. %i Box %i Chi2 %7.3f %s",
108 opt, X(), Y(), Q(), Ch(), Size(), fNlocMax, fBox, fChi2, status);
109 if(fDigs) fDigs->Print();
d3da6dc4 110}//Print()
111//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
112Int_t AliHMPIDCluster::Solve(TClonesArray *pCluLst,Bool_t isTryUnfold)
113{
114//This methode is invoked when the cluster is formed to solve it. Solve the cluster means to try to unfold the cluster
115//into the local maxima number of clusters. This methode is invoked by AliHMPIDRconstructor::Dig2Clu() on cluster by cluster basis.
116//At this point, cluster contains a list of digits, cluster charge and size is precalculated in AddDigit(), position is preset to (-1,-1) in ctor,
117//status is preset to kFormed in AddDigit(), chamber-sector info is preseted to actual values in AddDigit()
118//Method first finds number of local maxima and if it's more then one tries to unfold this cluster into local maxima number of clusters
119//Arguments: pCluLst - cluster list pointer where to add new cluster(s)
120// isTryUnfold - flag to switch on/off unfolding
121// Returns: number of local maxima of original cluster
d1bf51e1 122 CoG();
123 // Printf("1 - fStatus: %d",fSt);
124 Int_t iCluCnt=pCluLst->GetEntriesFast(); //get current number of clusters already stored in the list by previous operations
125 if(isTryUnfold==kFALSE || Size()==1) { //if cluster contains single pad there is no way to improve the knowledge
126 (isTryUnfold)?fSt=kSi1:fSt=kNot;
127 new ((*pCluLst)[iCluCnt++]) AliHMPIDCluster(*this); //add this raw cluster
128 return 1;
129 }
130 // Printf("2 - fStatus: %d",fSt);
d3da6dc4 131//Phase 0. Initialise TMinuit
132 const Int_t kMaxLocMax=6; //max allowed number of loc max for fitting
133 TMinuit *pMinuit = new TMinuit(3*kMaxLocMax); //init MINUIT with this number of parameters (3 params per mathieson)
d1bf51e1 134 pMinuit->SetObjectFit((TObject*)this); pMinuit->SetFCN(AliHMPIDCluster::FitFunc); //set fit function
135 Double_t aArg=-1; Int_t iErrFlg; //tmp vars for TMinuit
d3da6dc4 136 pMinuit->mnexcm("SET PRI",&aArg,1,iErrFlg); //suspend all printout from TMinuit
137 pMinuit->mnexcm("SET NOW",&aArg,0,iErrFlg); //suspend all warning printout from TMinuit
d1bf51e1 138//Phase 1. Find number of local maxima. Strategy is to check if the current pad has QDC more then all neigbours. Also find the box contaning the cluster
139 fNlocMax=0;
140 Int_t minPadX=999,minPadY=999,maxPadX=-1,maxPadY=-1,pc=-1; //for box finding
141 //Double_t lowX,highX,lowY,highY;
142
143 // Printf("3 - fStatus: %d",fSt);
144 for(Int_t iDig1=0;iDig1<Size();iDig1++) { //first digits loop
d3da6dc4 145 AliHMPIDDigit *pDig1 = Dig(iDig1); //take next digit
d1bf51e1 146 pc=pDig1->Pc(); //finding the box
147
148 if(pDig1->PadPcX() > maxPadX) maxPadX = pDig1->PadPcX();
149 if(pDig1->PadPcY() > maxPadY) maxPadY = pDig1->PadPcY();
150 if(pDig1->PadPcX() < minPadX) minPadX = pDig1->PadPcX();
151 if(pDig1->PadPcY() < minPadY) minPadY = pDig1->PadPcY();
152
153 fBox=(maxPadX-minPadX+1)*100+maxPadY-minPadY+1;
154
155 Int_t iHowManyMoreCnt = 0; //counts how many neighbouring pads has QDC more then current one
156 for(Int_t iDig2=0;iDig2<Size();iDig2++) { //loop on all digits again
157 if(iDig1==iDig2) continue; //the same digit, no need to compare
d3da6dc4 158 AliHMPIDDigit *pDig2 = Dig(iDig2); //take second digit to compare with the first one
da08475b 159 Int_t dist = TMath::Sign(Int_t(pDig1->PadChX()-pDig2->PadChX()),1)+TMath::Sign(Int_t(pDig1->PadChY()-pDig2->PadChY()),1);//distance between pads
d1bf51e1 160 if(dist==1) //means dig2 is a neighbour of dig1
161 if(pDig2->Q()>=pDig1->Q()) iHowManyMoreCnt++; //count number of pads with Q more then Q of current pad
d3da6dc4 162 }//second digits loop
d1bf51e1 163 if(iHowManyMoreCnt==0&&fNlocMax<kMaxLocMax){ //this pad has Q more then any neighbour so it's local maximum
164
165 /*
166 lowX = AliHMPIDDigit::LorsX(pc,minPadX) - 0.5 *AliHMPIDDigit::SizePadX();
167 highX = AliHMPIDDigit::LorsX(pc,maxPadX) + 0.5 *AliHMPIDDigit::SizePadX();
168 lowY = AliHMPIDDigit::LorsY(pc,minPadY) - 0.5 *AliHMPIDDigit::SizePadY();
169 highY = AliHMPIDDigit::LorsY(pc,maxPadY) + 0.5 *AliHMPIDDigit::SizePadY();
170 */
171 //Double_t lowQ=0,highQ=30000;
172
173 fQi=pDig1->Q(); fXi=pDig1->LorsX(); fYi=pDig1->LorsY(); //initial position of this Mathieson is to be in the center of loc max pad
174 /*
175 pMinuit->mnparm(3*fNlocMax ,Form("x%i",fNlocMax),fXi,0.01,lowX,highX,iErrFlg);
176 pMinuit->mnparm(3*fNlocMax+1,Form("y%i",fNlocMax),fYi,0.01,lowY,highY,iErrFlg);
177 pMinuit->mnparm(3*fNlocMax+2,Form("q%i",fNlocMax),fQi,0.01,lowQ,highQ,iErrFlg);
178 */
179 pMinuit->mnparm(3*fNlocMax ,Form("x%i",fNlocMax),fXi,0.01,0,0,iErrFlg);
180 pMinuit->mnparm(3*fNlocMax+1,Form("y%i",fNlocMax),fYi,0.01,0,0,iErrFlg);
181 pMinuit->mnparm(3*fNlocMax+2,Form("q%i",fNlocMax),fQi,0.01,0,100000,iErrFlg);
182
183 fNlocMax++;
d3da6dc4 184 }//if this pad is local maximum
185 }//first digits loop
d1bf51e1 186
187 //Int_t fitChk=0;
188
d3da6dc4 189//Phase 2. Fit loc max number of Mathiesons or add this current cluster to the list
d1bf51e1 190// Printf("4 - fStatus: %d",fSt);
191 if ( fNlocMax == 0) { // case of no local maxima found: pads with same charge...
192 pMinuit->mnparm(3*fNlocMax ,Form("x%i",fNlocMax),fX,0.01,0,0,iErrFlg);
193 pMinuit->mnparm(3*fNlocMax+1,Form("y%i",fNlocMax),fY,0.01,0,0,iErrFlg);
194 pMinuit->mnparm(3*fNlocMax+2,Form("q%i",fNlocMax),fQ,0.01,0,100000,iErrFlg);
195 fNlocMax = 1;
196 fSt=kNoLoc;
197 }
198
199 if ( fNlocMax >= kMaxLocMax)
200 {
201 fSt = kMax; new ((*pCluLst)[iCluCnt++]) AliHMPIDCluster(*this); //add this raw cluster
202 }
203 else{ //resonable number of local maxima to fit and user requested it
204 Double_t arglist[10]; arglist[0] = 10000; arglist[1] = 1.; //number of steps and sigma on pads charges
205 pMinuit->mnexcm("MIGRAD" ,arglist,0,iErrFlg); //start fitting
206
207 if (iErrFlg)
208 {
209 fSt = kAbn; //fit fails, MINUIT returns error flag
210 new ((*pCluLst)[iCluCnt++]) AliHMPIDCluster(*this); //add this raw cluster
211 }
212 else
213 { //Only if MIGRAD converged normally
214 Double_t d2,d3; TString sName; //vars to get results from TMinuit
215 for(Int_t i=0;i<fNlocMax;i++){ //local maxima loop
216 pMinuit->mnpout(3*i ,sName, fX, fXe , d2, d3, iErrFlg);
217 pMinuit->mnpout(3*i+1 ,sName, fY, fYe , d2, d3, iErrFlg);
218 pMinuit->mnpout(3*i+2 ,sName, fQ, fQe , d2, d3, iErrFlg);
219 pMinuit->mnstat(fChi2,d2,d2,iErrFlg,iErrFlg,iErrFlg);
220
221 if(fNlocMax!=1)fSt=kUnf;
222 if(fNlocMax==1&&fSt!=kNoLoc) fSt=kLo1;
223 if ( !IsInPc()) fSt = kEdg;
224 if(fSt==kNoLoc) fNlocMax=0;
225 new ((*pCluLst)[iCluCnt++]) AliHMPIDCluster(*this); //add new unfolded cluster
226 }
227 }
228 }
229
230
231
232
233
d3da6dc4 234 delete pMinuit;
d1bf51e1 235 return fNlocMax;
d3da6dc4 236}//Solve()
237//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++