3 /**************************************************************************
4 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
6 * Author: The ALICE Off-line Project. *
7 * Contributors are mentioned in the code where appropriate. *
9 * Permission to use, copy, modify and distribute this software and its *
10 * documentation strictly for non-commercial purposes is hereby granted *
11 * without fee, provided that the above copyright notice appears in all *
12 * copies and that both the copyright notice and this permission notice *
13 * appear in the supporting documentation. The authors make no claims *
14 * about the suitability of this software for any purpose. It is *
15 * provided "as is" without express or implied warranty. *
16 **************************************************************************/
17 #include "AliHLTITSVertexerZ.h"
23 #include "AliITSLoader.h"
24 #include <AliITSgeom.h>
25 #include <AliITSRecPoint.h>
26 #include <AliITSclusterV2.h>
28 //-------------------------------------------------------------------------
29 // Implementation of the HLT ITS vertexer class
31 // Origin: Cvetan Cheshkov, CERN, Cvetan.Cheshkov@cern.ch
32 //-------------------------------------------------------------------------
34 ClassImp(AliHLTITSVertexerZ)
36 AliHLTITSVertexerZ::AliHLTITSVertexerZ():
41 // Constructor in case that there is no runloader
45 AliHLTITSVertexerZ::AliHLTITSVertexerZ(Float_t x0, Float_t y0):
46 AliITSVertexerZ(x0,y0),
50 // Standard Constructor
54 AliHLTITSVertexerZ::~AliHLTITSVertexerZ()
57 if (fZCombf) delete fZCombf;
60 //______________________________________________________________________
61 AliESDVertex* AliHLTITSVertexerZ::FindVertexForCurrentEvent(AliITSgeom *geom,TTree *tR){
62 // Defines the AliESDVertex for the current event
66 Float_t lc[3]; for(Int_t ii=0; ii<3; ii++) lc[ii]=0.;
67 Float_t gc[3]; for(Int_t ii=0; ii<3; ii++) gc[ii]=0.;
68 Float_t lc2[3]; for(Int_t ii=0; ii<3; ii++) lc2[ii]=0.;
69 Float_t gc2[3]; for(Int_t ii=0; ii<3; ii++) gc2[ii]=0.;
71 TClonesArray dummy("AliITSclusterV2",10000), *clusters=&dummy;
73 branch = tR->GetBranch("Clusters");
74 branch->SetAddress(&clusters);
78 for(Int_t module= fFirstL1; module<=fLastL1;module++){
79 if(module%4==0 || module%4==3)continue;
80 // cout<<"Procesing module "<<module<<" ";
82 // cout<<"Number of clusters "<<clusters->GetEntries()<<endl;
83 nrpL1+= clusters->GetEntriesFast();
85 for(Int_t module= fFirstL2; module<=fLastL2;module++){
87 nrpL2+= clusters->GetEntriesFast();
89 if(nrpL1 == 0 || nrpL2 == 0){
91 return fCurrentVertex;
94 Int_t nPhiBins = (Int_t)(TMath::Pi()/fDiffPhiMax);
95 Float_t phiBinSize = 2*TMath::Pi()/(Float_t)nPhiBins;
97 Int_t maxind1 = 2*nrpL1/nPhiBins;
98 Float_t **zc1 = new Float_t *[nPhiBins];
99 Float_t **phi1 = new Float_t *[nPhiBins];
100 Float_t **r1 = new Float_t *[nPhiBins];
101 Int_t *ind1 = new Int_t [nPhiBins];
102 Int_t maxind2 = 2*nrpL2/nPhiBins;
103 Float_t **zc2 = new Float_t *[nPhiBins];
104 Float_t **phi2 = new Float_t *[nPhiBins];
105 Float_t **r2 = new Float_t *[nPhiBins];
106 Int_t *ind2 = new Int_t [nPhiBins];
107 for(Int_t i=0;i<nPhiBins;i++) {
108 zc1[i] = new Float_t [maxind1];
109 phi1[i] = new Float_t [maxind1];
110 r1[i] = new Float_t [maxind1];
111 zc2[i] = new Float_t [maxind2];
112 phi2[i] = new Float_t [maxind2];
113 r2[i] = new Float_t [maxind2];
117 Float_t zshift[4] = {-10.708000, -3.536000, 3.536000, 10.708000};
120 memset(ind1,0,nPhiBins*sizeof(Int_t));
121 for(Int_t module= fFirstL1; module<=fLastL1;module++){
122 if(module%4==0 || module%4==3)continue;
123 tR->GetEvent(module);
124 Int_t nrecp1 = clusters->GetEntriesFast();
125 for(Int_t j=0;j<nrecp1;j++){
126 AliITSclusterV2 *recp = (AliITSclusterV2*)clusters->UncheckedAt(j);
127 lc[0]=-recp->GetY()+yshift;
128 lc[2]=-recp->GetZ()+zshift[module%4];
129 geom->LtoG(module,lc,gc);
130 gc[0]-=GetNominalPos()[0];
131 gc[1]-=GetNominalPos()[1];
135 Float_t phi = TMath::ATan2(gc[1],gc[0]);
136 if(phi<0)phi=2*TMath::Pi()+phi;
137 Int_t bin = (Int_t)(phi/phiBinSize);
138 if(bin>=nPhiBins || bin<0) bin = 0;
139 Int_t ind = ind1[bin];
141 phi1[bin][ind] = phi;
142 zc1[bin][ind]=gc[2]/fStepFine;
143 r1[bin][ind]=sqrt(xc1*xc1+yc1*yc1);
150 memset(ind2,0,nPhiBins*sizeof(Int_t));
151 for(Int_t module= fFirstL2; module<=fLastL2;module++){
152 tR->GetEvent(module);
153 Int_t nrecp2 = clusters->GetEntriesFast();
154 for(Int_t j=0;j<nrecp2;j++){
155 AliITSclusterV2 *recp = (AliITSclusterV2*)clusters->UncheckedAt(j);
156 lc[0]=recp->GetY()+yshift;
157 lc[2]=-recp->GetZ()+zshift[module%4];
158 geom->LtoG(module,lc,gc);
159 gc[0]-=GetNominalPos()[0];
160 gc[1]-=GetNominalPos()[1];
164 Float_t phi = TMath::ATan2(gc[1],gc[0]);
165 if(phi<0)phi=2*TMath::Pi()+phi;
166 Int_t bin = (Int_t)(phi/phiBinSize+0.5);
167 if(bin>=nPhiBins || bin<0) bin = 0;
168 Int_t ind = ind2[bin];
170 phi2[bin][ind] = phi;
171 zc2[bin][ind]=gc[2]/fStepFine;
172 r2[bin][ind]=sqrt(xc2*xc2+yc2*yc2);
178 Int_t nbinfine = static_cast<Int_t>((fHighLim-fLowLim)/fStepFine);
179 Float_t lowz = fLowLim/fStepFine;
180 Int_t *harray = new Int_t[nbinfine];
181 memset(harray,0,nbinfine*sizeof(Int_t));
182 for(Int_t ibin=0;ibin<nPhiBins;ibin++) {
183 Float_t *pphi1 = phi1[ibin];
184 Float_t *pr1 = r1[ibin];
185 Float_t *pzc1 = zc1[ibin];
186 for(Int_t i=0;i<ind1[ibin];i++){
187 for(Int_t jbin=ibin;jbin<=(ibin+1);jbin++) {
189 if(jbin==nPhiBins) jbin2 = 0;
190 Float_t *pphi2 = phi2[jbin2];
191 Float_t *pr2 = r2[jbin2];
192 Float_t *pzc2 = zc2[jbin2];
193 for(Int_t j=0;j<ind2[jbin2];j++){
194 Float_t diff = TMath::Abs(pphi2[j]-pphi1[i]);
195 if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff;
196 if(diff<fDiffPhiMax){
197 Float_t zr0=(pr2[j]*pzc1[i]-pr1[i]*pzc2[j])/(pr2[j]-pr1[i]);
198 Int_t bin = (Int_t)(zr0-lowz);
199 if(bin>=0 && bin<nbinfine){
207 for(Int_t i=0;i<nPhiBins;i++) {
224 Int_t nbinratio = (Int_t)(fStepCoarse/fStepFine+0.5);
225 Int_t nbincoarse = nbinfine/nbinratio;
227 if(fZCombc)delete fZCombc;
228 fZCombc = new TH1F("fZCombc","Z",nbincoarse,fLowLim,fLowLim+nbincoarse*fStepCoarse);
229 if(fZCombf)delete fZCombf;
230 fZCombf = new TH1F("fZCombf","Z",nbinfine,fLowLim,fLowLim+nbinfine*fStepFine);
233 for(Int_t bin=0;bin<nbinfine;bin++) {
234 fZCombf->SetBinContent(bin+1,(Stat_t)harray[bin]);
235 fZCombf->SetBinError(bin+1,TMath::Sqrt((Stat_t)harray[bin]));
236 contents+=(Stat_t)harray[bin];
238 if(counter==nbinratio) {
239 Int_t binc=bin/nbinratio;
240 fZCombc->SetBinContent(binc+1,contents);
241 fZCombc->SetBinError(binc+1,TMath::Sqrt(contents));
249 if(fZCombc->GetEntries()==0){
250 Warning("FindVertexForCurrentEvent","Insufficient number of rec. points\n");
252 return fCurrentVertex;
255 // cout<<"Number of entries in hist. "<<fZCombc->GetEntries()<<endl;
257 Int_t bi = fZCombc->GetMaximumBin();
258 Float_t centre = fZCombc->GetBinCenter(bi);
259 Int_t n1 = static_cast<Int_t>((centre-fZCombc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fZCombf->GetBinWidth(0));
260 Int_t n2 = static_cast<Int_t>((centre+fZCombc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fZCombf->GetBinWidth(0));
268 for(Int_t n=n1;n<=n2;n++){
269 fZFound+=fZCombf->GetBinCenter(n)*fZCombf->GetBinContent(n);
270 num+=static_cast<Int_t>(fZCombf->GetBinContent(n));
271 fZsig+=fZCombf->GetBinCenter(n)*fZCombf->GetBinCenter(n)*fZCombf->GetBinContent(n);
277 Float_t radi = fZsig/(num-1)-fZFound*fZFound/num/(num-1);
278 if(radi>0.)fZsig=TMath::Sqrt(radi);
282 goon = TMath::Abs(TMath::Abs(fZFound-fZCombf->GetBinCenter(n1))-TMath::Abs(fZFound-fZCombf->GetBinCenter(n2)))>fTolerance;
283 n1 = static_cast<Int_t>((fZFound-fZCombc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fZCombf->GetBinWidth(0));
284 n2 = static_cast<Int_t>((fZFound+fZCombc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fZCombf->GetBinWidth(0));
288 Warning("FindVertexForCurrentEvent","The procedure dows not converge\n");
291 // cout<<"Numer of Iterations "<<niter<<endl<<endl;
292 fCurrentVertex = new AliESDVertex(fZFound,fZsig,num);
293 fCurrentVertex->SetTitle("vertexer: HLT");
296 return fCurrentVertex;