]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HMPID/AliHMPIDCluster.cxx
Typo corrected.
[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()
27311693 20
21Bool_t AliHMPIDCluster::fDoCorrSin=kTRUE;
22
d3da6dc4 23ClassImp(AliHMPIDCluster)
24//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
25void AliHMPIDCluster::CoG()
26{
27// Calculates naive cluster position as a center of gravity of its digits.
28// Arguments: none
d1bf51e1 29// Returns: none
e4a3eae8 30 Int_t minPadX=999,minPadY=999,maxPadX=-1,maxPadY=-1; //for box finding
31 if(fDigs==0) return; //no digits in this cluster
32 fX=fY=fQRaw=0; //init summable parameters
33 Int_t maxQpad=-1,maxQ=-1; //to calculate the pad with the highest charge
d1bf51e1 34 AliHMPIDDigit *pDig;
e4a3eae8 35 for(Int_t iDig=0;iDig<fDigs->GetEntriesFast();iDig++){ //digits loop
36 pDig=(AliHMPIDDigit*)fDigs->At(iDig); //get pointer to next digit
37
38 if(pDig->PadPcX() > maxPadX) maxPadX = pDig->PadPcX(); // find the minimum box that contain the cluster MaxX
39 if(pDig->PadPcY() > maxPadY) maxPadY = pDig->PadPcY(); // MaxY
40 if(pDig->PadPcX() < minPadX) minPadX = pDig->PadPcX(); // MinX
41 if(pDig->PadPcY() < minPadY) minPadY = pDig->PadPcY(); // MinY
42
43 Float_t q=pDig->Q(); //get QDC
44 fX += pDig->LorsX()*q;fY +=pDig->LorsY()*q; //add digit center weighted by QDC
45 fQRaw+=q; //increment total charge
46 if(q>maxQ) {maxQpad = pDig->Pad();maxQ=(Int_t)q;} // to find pad with highest charge
d3da6dc4 47 }//digits loop
e4a3eae8 48
49 fBox=(maxPadX-minPadX+1)*100+maxPadY-minPadY+1; // dimension of the box: format Xdim*100+Ydim
50
51 if ( fQRaw != 0 ) fX/=fQRaw;fY/=fQRaw; //final center of gravity
27311693 52
53 if(fDigs->GetEntriesFast()>1&&fDoCorrSin)CorrSin(); //correct it by sinoid
e4a3eae8 54
c5c19d6a 55 fQ = fQRaw; // Before starting fit procedure, Q and QRaw must be equal
56 fCh=pDig->Ch(); //initialize chamber number
57 fMaxQpad = maxQpad; fMaxQ=maxQ; //store max charge pad to the field
58 fChi2=0; // no Chi2 to find
59 fNlocMax=0; // proper status from this method
d3da6dc4 60 fSt=kCoG;
61}//CoG()
62//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
63void AliHMPIDCluster::CorrSin()
64{
65// Correction of cluster x position due to sinoid, see HMPID TDR page 30
66// Arguments: none
67// Returns: none
1d4857c5 68 Int_t pc,px,py;
69 AliHMPIDDigit::Lors2Pad(fX,fY,pc,px,py); //tmp digit to get it center
70 Float_t x=fX-AliHMPIDDigit::LorsX(pc,px); //diff between cluster x and center of the pad contaning this cluster
d3da6dc4 71 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;
72}
73//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
d1bf51e1 74void AliHMPIDCluster::Draw(Option_t*)
75{
a1d55ff3 76 TMarker *pMark=new TMarker(X(),Y(),5); pMark->SetUniqueID(fSt);pMark->SetMarkerColor(kBlue); pMark->Draw();
d1bf51e1 77}
78//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
d3da6dc4 79void AliHMPIDCluster::FitFunc(Int_t &iNpars, Double_t *, Double_t &chi2, Double_t *par, Int_t )
80{
81// Cluster fit function
82// par[0]=x par[1]=y par[2]=q for the first Mathieson shape
83// par[3]=x par[4]=y par[5]=q for the second Mathieson shape and so on up to iNpars/3 Mathieson shapes
84// For each pad of the cluster calculates the difference between actual pad charge and the charge induced to this pad by all Mathieson distributions
85// Then the chi2 is calculated as the sum of this value squared for all pad in the cluster.
86// Arguments: iNpars - number of parameters which is number of local maxima of cluster * 3
87// chi2 - function result to be minimised
88// par - parameters array of size iNpars
89// Returns: none
90 AliHMPIDCluster *pClu=(AliHMPIDCluster*)gMinuit->GetObjectFit();
91 Int_t iNshape = iNpars/3;
92
93 chi2 = 0;
e4a3eae8 94 for(Int_t i=0;i<pClu->Size();i++){ //loop on all pads of the cluster
95 Double_t dQpadMath = 0; //pad charge collector
96 for(Int_t j=0;j<iNshape;j++){ //Mathiesons loop as all of them may contribute to this pad
c5c19d6a 97 dQpadMath+=par[3*j+2]*pClu->Dig(i)->IntMathieson(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
d3da6dc4 98 }
c5c19d6a 99// if(dQpadMath>0)chi2 +=TMath::Power((pClu->Dig(i)->Q()-dQpadMath),2)/dQpadMath; //
7c147ea1 100 if(dQpadMath>0)chi2 +=TMath::Power((pClu->Dig(i)->Q()-dQpadMath),2)/pClu->Dig(i)->Q(); //
e4a3eae8 101 } //loop on all pads of the cluster
d3da6dc4 102}//FitFunction()
103//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
104void AliHMPIDCluster::Print(Option_t* opt)const
105{
106//Print current cluster
107 const char *status=0;
108 switch(fSt){
d1bf51e1 109 case kFrm : status="formed " ;break;
110 case kUnf : status="unfolded (fit)" ;break;
111 case kCoG : status="coged " ;break;
112 case kLo1 : status="locmax 1 (fit)" ;break;
d1bf51e1 113 case kMax : status="exceeded (cog)" ;break;
114 case kNot : status="not done (cog)" ;break;
115 case kEmp : status="empty " ;break;
116 case kEdg : status="edge (fit)" ;break;
117 case kSi1 : status="size 1 (cog)" ;break;
118 case kNoLoc: status="no LocMax(fit)" ;break;
c5c19d6a 119 case kAbn : status="Abnormal fit " ;break;
d1bf51e1 120
121 default: status="??????" ;break;
d3da6dc4 122 }
e4a3eae8 123 Double_t ratio=0;
124 if(Q()>0&&QRaw()>0) ratio = Q()/QRaw()*100;
1d4857c5 125 Printf("%sCLU: ch=%i (%7.3f,%7.3f) Q=%8.3f Qraw=%8.3f(%3.0f%%) Size=%2i DimBox=%i LocMax=%i Chi2=%7.3f %s",
1006986d 126 opt,Ch(),X(),Y(),Q(),QRaw(),ratio,Size(),fBox,fNlocMax,fChi2,status);
d1bf51e1 127 if(fDigs) fDigs->Print();
d3da6dc4 128}//Print()
129//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
130Int_t AliHMPIDCluster::Solve(TClonesArray *pCluLst,Bool_t isTryUnfold)
131{
132//This methode is invoked when the cluster is formed to solve it. Solve the cluster means to try to unfold the cluster
133//into the local maxima number of clusters. This methode is invoked by AliHMPIDRconstructor::Dig2Clu() on cluster by cluster basis.
134//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,
135//status is preset to kFormed in AddDigit(), chamber-sector info is preseted to actual values in AddDigit()
136//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
137//Arguments: pCluLst - cluster list pointer where to add new cluster(s)
138// isTryUnfold - flag to switch on/off unfolding
139// Returns: number of local maxima of original cluster
d1bf51e1 140 CoG();
e4a3eae8 141 Int_t iCluCnt=pCluLst->GetEntriesFast(); //get current number of clusters already stored in the list by previous operations
142 if(isTryUnfold==kFALSE || Size()==1) { //if cluster contains single pad there is no way to improve the knowledge
d1bf51e1 143 (isTryUnfold)?fSt=kSi1:fSt=kNot;
144 new ((*pCluLst)[iCluCnt++]) AliHMPIDCluster(*this); //add this raw cluster
145 return 1;
146 }
d3da6dc4 147//Phase 0. Initialise TMinuit
e4a3eae8 148 const Int_t kMaxLocMax=6; //max allowed number of loc max for fitting
7c147ea1 149 if(!gMinuit) gMinuit = new TMinuit(100); //init MINUIT with this number of parameters (3 params per mathieson)
150 gMinuit->mncler(); // reset Minuit list of paramters
151 gMinuit->SetObjectFit((TObject*)this); gMinuit->SetFCN(AliHMPIDCluster::FitFunc); //set fit function
152 Double_t aArg=-1; //tmp vars for TMinuit
153 Int_t iErrFlg;
154 gMinuit->mnexcm("SET PRI",&aArg,1,iErrFlg); //suspend all printout from TMinuit
155 gMinuit->mnexcm("SET NOW",&aArg,0,iErrFlg); //suspend all warning printout from TMinuit
d1bf51e1 156//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
157 fNlocMax=0;
e4a3eae8 158
159 for(Int_t iDig1=0;iDig1<Size();iDig1++) { //first digits loop
160 AliHMPIDDigit *pDig1 = Dig(iDig1); //take next digit
c5c19d6a 161 Int_t iCnt = 0; //counts how many neighbouring pads has QDC more then current one
d1bf51e1 162 for(Int_t iDig2=0;iDig2<Size();iDig2++) { //loop on all digits again
163 if(iDig1==iDig2) continue; //the same digit, no need to compare
d3da6dc4 164 AliHMPIDDigit *pDig2 = Dig(iDig2); //take second digit to compare with the first one
da08475b 165 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 166 if(dist==1) //means dig2 is a neighbour of dig1
c5c19d6a 167 if(pDig2->Q()>=pDig1->Q()) iCnt++; //count number of pads with Q more then Q of current pad
d3da6dc4 168 }//second digits loop
c5c19d6a 169 if(iCnt==0&&fNlocMax<kMaxLocMax){ //this pad has Q more then any neighbour so it's local maximum
170 Double_t xStart=pDig1->LorsX();Double_t yStart=pDig1->LorsY();
171 Double_t xMin=xStart-AliHMPIDDigit::SizePadX();
172 Double_t xMax=xStart+AliHMPIDDigit::SizePadX();
173 Double_t yMin=yStart-AliHMPIDDigit::SizePadY();
174 Double_t yMax=yStart+AliHMPIDDigit::SizePadY();
7c147ea1 175 gMinuit->mnparm(3*fNlocMax ,Form("x%i",fNlocMax),xStart,0.1,xMin,xMax,iErrFlg); // X,Y,Q initial values of the loc max pad
176 gMinuit->mnparm(3*fNlocMax+1,Form("y%i",fNlocMax),yStart,0.1,yMin,yMax,iErrFlg); // X, Y constrained to be near the loc max
177 gMinuit->mnparm(3*fNlocMax+2,Form("q%i",fNlocMax),pDig1->Q(),0.1,0,100000,iErrFlg);// Q constrained to be positive
d1bf51e1 178 fNlocMax++;
d3da6dc4 179 }//if this pad is local maximum
180 }//first digits loop
d1bf51e1 181
d3da6dc4 182//Phase 2. Fit loc max number of Mathiesons or add this current cluster to the list
e4a3eae8 183// case 1 -> no loc max found
184 if ( fNlocMax == 0) { // case of no local maxima found: pads with same charge...
7c147ea1 185 gMinuit->mnparm(3*fNlocMax ,Form("x%i",fNlocMax),fX,0.1,0,0,iErrFlg); // Init values taken from CoG() -> fX,fY,fQRaw
186 gMinuit->mnparm(3*fNlocMax+1,Form("y%i",fNlocMax),fY,0.1,0,0,iErrFlg); //
187 gMinuit->mnparm(3*fNlocMax+2,Form("q%i",fNlocMax),fQRaw,0.1,0,100000,iErrFlg); //
d1bf51e1 188 fNlocMax = 1;
189 fSt=kNoLoc;
190 }
e4a3eae8 191
192// case 2 -> loc max found. Check # of loc maxima
193 if ( fNlocMax >= kMaxLocMax) // if # of local maxima exceeds kMaxLocMax...
d1bf51e1 194 {
e4a3eae8 195 fSt = kMax; new ((*pCluLst)[iCluCnt++]) AliHMPIDCluster(*this); //...add this raw cluster
196 } //or...
197 else{ //...resonable number of local maxima to fit and user requested it
198 Double_t arglist[10];arglist[0] = 10000;arglist[1] = 1.; //number of steps and sigma on pads charges
7c147ea1 199 gMinuit->mnexcm("SIMPLEX" ,arglist,2,iErrFlg); //start fitting with Simplex
200 gMinuit->mnexcm("MIGRAD" ,arglist,2,iErrFlg); //fitting improved by Migrad
c5c19d6a 201 if(iErrFlg) {
202 Double_t strategy=2;
7c147ea1 203 gMinuit->mnexcm("SET STR",&strategy,1,iErrFlg); //change level of strategy
c5c19d6a 204 if(!iErrFlg) {
7c147ea1 205 gMinuit->mnexcm("SIMPLEX" ,arglist,2,iErrFlg);
206 gMinuit->mnexcm("MIGRAD" ,arglist,2,iErrFlg); //fitting improved by Migrad
07dd8038 207// Printf("Try to improve fit --> err %d",iErrFlg);
c5c19d6a 208 }
209 }
210 if(iErrFlg) fSt=kAbn; //no convergence of the fit...
e4a3eae8 211 Double_t dummy; TString sName; //vars to get results from Minuit
212 for(Int_t i=0;i<fNlocMax;i++){ //store the local maxima parameters
7c147ea1 213 gMinuit->mnpout(3*i ,sName, fX, fErrX , dummy, dummy, iErrFlg); // X
214 gMinuit->mnpout(3*i+1 ,sName, fY, fErrY , dummy, dummy, iErrFlg); // Y
215 gMinuit->mnpout(3*i+2 ,sName, fQ, fErrQ , dummy, dummy, iErrFlg); // Q
216 gMinuit->mnstat(fChi2,dummy,dummy,iErrFlg,iErrFlg,iErrFlg); // Chi2 of the fit
c5c19d6a 217 if(fSt!=kAbn) {
218 if(fNlocMax!=1)fSt=kUnf; // if unfolded
219 if(fNlocMax==1&&fSt!=kNoLoc) fSt=kLo1; // if only 1 loc max
220 if ( !IsInPc()) fSt = kEdg; // if Out of Pc
221 if(fSt==kNoLoc) fNlocMax=0; // if with no loc max (pads with same charge..)
222 }
e4a3eae8 223 new ((*pCluLst)[iCluCnt++]) AliHMPIDCluster(*this); //add new unfolded cluster
224 }
225 }
d1bf51e1 226
e4a3eae8 227 return fNlocMax;
228
d3da6dc4 229}//Solve()
230//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++