]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliITSVertexerZ.cxx
warning on gcc 3.2.3 removed
[u/mrichter/AliRoot.git] / ITS / AliITSVertexerZ.cxx
CommitLineData
0c6af5c9 1/**************************************************************************
2 * Copyright(c) 1998-2003, 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#include <AliITSVertexerZ.h>
0c6af5c9 16#include <TString.h>
32e449be 17#include <AliITSLoader.h>
0c6af5c9 18#include<TBranch.h>
0c6af5c9 19#include<TH1.h>
0c6af5c9 20#include<TTree.h>
21#include <AliITSgeom.h>
7d62fb64 22#include "AliITSDetTypeRec.h"
0c6af5c9 23#include <AliITSRecPoint.h>
24
25/////////////////////////////////////////////////////////////////
26// this class implements a fast method to determine
27// the Z coordinate of the primary vertex
28// for p-p collisions it seems to give comparable or better results
29// with respect to what obtained with AliITSVertexerPPZ
30// It can be used successfully with Pb-Pb collisions
31////////////////////////////////////////////////////////////////
32
33ClassImp(AliITSVertexerZ)
34
35
36
37//______________________________________________________________________
38AliITSVertexerZ::AliITSVertexerZ():AliITSVertexer() {
39 // Default constructor
40 SetDiffPhiMax(0);
41 fX0 = 0.;
42 fY0 = 0.;
43 SetFirstLayerModules(0);
44 SetSecondLayerModules(0);
45 fZFound = 0;
46 fZsig = 0.;
0c6af5c9 47 fZCombc = 0;
48 fZCombf = 0;
ecc64c3f 49 fZCombv = 0;
0c6af5c9 50 SetLowLimit(0.);
51 SetHighLimit(0.);
52 SetBinWidthCoarse(0.);
53 SetBinWidthFine(0.);
54 SetTolerance(0.);
ecc64c3f 55 SetPPsetting(0.,0.);
6fd990e3 56 ConfigIterations();
0c6af5c9 57}
58
59//______________________________________________________________________
60AliITSVertexerZ::AliITSVertexerZ(TString fn, Float_t x0, Float_t y0):AliITSVertexer(fn) {
61 // Standard Constructor
62 SetDiffPhiMax();
63 fX0 = x0;
64 fY0 = y0;
65 SetFirstLayerModules();
66 SetSecondLayerModules();
67 fZFound = 0;
68 fZsig = 0.;
0c6af5c9 69 fZCombc = 0;
70 fZCombf = 0;
ecc64c3f 71 fZCombv = 0;
0c6af5c9 72 SetLowLimit();
73 SetHighLimit();
74 SetBinWidthCoarse();
75 SetBinWidthFine();
76 SetTolerance();
ecc64c3f 77 SetPPsetting();
6fd990e3 78 ConfigIterations();
0c6af5c9 79
80}
81
82//______________________________________________________________________
83AliITSVertexerZ::AliITSVertexerZ(const AliITSVertexerZ &vtxr) : AliITSVertexer(vtxr) {
84 // Copy constructor
85 // Copies are not allowed. The method is protected to avoid misuse.
86 Error("AliITSVertexerZ","Copy constructor not allowed\n");
87}
88
89//______________________________________________________________________
90AliITSVertexerZ& AliITSVertexerZ::operator=(const AliITSVertexerZ& /* vtxr */){
91 // Assignment operator
92 // Assignment is not allowed. The method is protected to avoid misuse.
93 Error("= operator","Assignment operator not allowed\n");
94 return *this;
95}
96
0c6af5c9 97//______________________________________________________________________
98AliITSVertexerZ::~AliITSVertexerZ() {
ecc64c3f 99 // Destructor
100 delete fZCombc;
101 delete fZCombf;
102 delete fZCombv;
0c6af5c9 103}
104
6fd990e3 105//______________________________________________________________________
106void AliITSVertexerZ::ConfigIterations(Int_t noiter,Float_t *ptr){
107 // configure the iterative procedure to gain efficiency for
108 // pp events with very low multiplicity
109 Float_t defaults[5]={0.05,0.1,0.2,0.3,0.5};
110 fMaxIter=noiter;
111 if(noiter>5){
112 Error("ConfigIterations","Maximum number of iterations is 5\n");
113 fMaxIter=5;
114 }
115 for(Int_t j=0;j<5;j++)fPhiDiffIter[j]=defaults[j];
116 if(ptr)for(Int_t j=0;j<fMaxIter;j++)fPhiDiffIter[j]=ptr[j];
117}
118
0c6af5c9 119//______________________________________________________________________
d681bb2d 120AliESDVertex* AliITSVertexerZ::FindVertexForCurrentEvent(Int_t evnumber){
121 // Defines the AliESDVertex for the current event
6fd990e3 122 VertexZFinder(evnumber);
123 Int_t ntrackl=0;
124 for(Int_t iteraz=0;iteraz<fMaxIter;iteraz++){
125 if(fCurrentVertex) ntrackl=fCurrentVertex->GetNContributors();
126 if(!fCurrentVertex || ntrackl==0 || ntrackl==-1){
127 Float_t diffPhiMaxOrig=fDiffPhiMax;
128 fDiffPhiMax=GetPhiMaxIter(iteraz);
129 VertexZFinder(evnumber);
130 fDiffPhiMax=diffPhiMaxOrig;
131 }
132 }
32e449be 133 FindMultiplicity(evnumber);
6fd990e3 134 return fCurrentVertex;
135}
136
32e449be 137
138
139
6fd990e3 140//______________________________________________________________________
141void AliITSVertexerZ::VertexZFinder(Int_t evnumber){
142 // Defines the AliESDVertex for the current event
0c6af5c9 143 fCurrentVertex = 0;
144 AliRunLoader *rl =AliRunLoader::GetRunLoader();
7d62fb64 145 AliITSLoader* itsLoader = (AliITSLoader*)rl->GetLoader("ITSLoader");
32e449be 146 AliITSgeom* geom = itsLoader->GetITSgeom();
2257f27e 147 itsLoader->LoadRecPoints();
0c6af5c9 148 rl->GetEvent(evnumber);
149
7d62fb64 150 AliITSDetTypeRec detTypeRec;
0c6af5c9 151
152 TTree *tR = itsLoader->TreeR();
7d62fb64 153 detTypeRec.SetTreeAddressR(tR);
0c6af5c9 154 TClonesArray *itsRec = 0;
ecc64c3f 155 // lc and gc are local and global coordinates for layer 1
0c6af5c9 156 Float_t lc[3]; for(Int_t ii=0; ii<3; ii++) lc[ii]=0.;
157 Float_t gc[3]; for(Int_t ii=0; ii<3; ii++) gc[ii]=0.;
ecc64c3f 158 // lc2 and gc2 are local and global coordinates for layer 2
0c6af5c9 159 Float_t lc2[3]; for(Int_t ii=0; ii<3; ii++) lc2[ii]=0.;
160 Float_t gc2[3]; for(Int_t ii=0; ii<3; ii++) gc2[ii]=0.;
161
7d62fb64 162 itsRec = detTypeRec.RecPoints();
2257f27e 163 TBranch *branch;
00a7cc50 164 branch = tR->GetBranch("ITSRecPoints");
0c6af5c9 165
0c6af5c9 166 Int_t nrpL1 = 0;
167 Int_t nrpL2 = 0;
ecc64c3f 168 // By default fFirstL1=0 and fLastL1=79
169 // This loop counts the number of recpoints on layer1 (central modules)
0c6af5c9 170 for(Int_t module= fFirstL1; module<=fLastL1;module++){
ecc64c3f 171 // Keep only central modules
0c6af5c9 172 if(module%4==0 || module%4==3)continue;
ecc64c3f 173 // cout<<"Procesing module "<<module<<" ";
0c6af5c9 174 branch->GetEvent(module);
ecc64c3f 175 // cout<<"Number of clusters "<<clusters->GetEntries()<<endl;
0c6af5c9 176 nrpL1+= itsRec->GetEntries();
7d62fb64 177 detTypeRec.ResetRecPoints();
0c6af5c9 178 }
ecc64c3f 179 //By default fFirstL2=80 and fLastL2=239
180 //This loop counts the number of RP on layer 2
0c6af5c9 181 for(Int_t module= fFirstL2; module<=fLastL2;module++){
182 branch->GetEvent(module);
183 nrpL2+= itsRec->GetEntries();
7d62fb64 184 detTypeRec.ResetRecPoints();
0c6af5c9 185 }
0c6af5c9 186 if(nrpL1 == 0 || nrpL2 == 0){
0c6af5c9 187 ResetHistograms();
6fd990e3 188 itsLoader->UnloadRecPoints();
189 fCurrentVertex = new AliESDVertex(0.,5.3,-2);
190 return;
0c6af5c9 191 }
ecc64c3f 192 // The vertex finding is attempted only if the number of RP is !=0 on
193 // both layers
194 Float_t *xc1 = new Float_t [nrpL1]; // coordinates of the L1 Recpoints
0c6af5c9 195 Float_t *yc1 = new Float_t [nrpL1];
196 Float_t *zc1 = new Float_t [nrpL1];
197 Float_t *phi1 = new Float_t [nrpL1];
ecc64c3f 198 Float_t *xc2 = new Float_t [nrpL2]; // coordinates of the L1 Recpoints
0c6af5c9 199 Float_t *yc2 = new Float_t [nrpL2];
200 Float_t *zc2 = new Float_t [nrpL2];
201 Float_t *phi2 = new Float_t [nrpL2];
ecc64c3f 202 Int_t ind = 0;// running index for RP
203 // Force a coarse bin size of 200 microns if the number of clusters on layer 2
204 // is low
205 if(nrpL2<fPPsetting[0])SetBinWidthCoarse(fPPsetting[1]);
206 // By default nbinfine = (10+10)/0.0005=40000
207 Int_t nbinfine = static_cast<Int_t>((fHighLim-fLowLim)/fStepFine);
208 // By default nbincoarse=(10+10)/0.01=2000
209 Int_t nbincoarse = static_cast<Int_t>((fHighLim-fLowLim)/fStepCoarse);
210 // Set stepverycoarse = 3*fStepCoarse
211 Int_t nbinvcoarse = static_cast<Int_t>((fHighLim-fLowLim)/(3.*fStepCoarse));
212 if(fZCombc)delete fZCombc;
213 fZCombc = new TH1F("fZCombc","Z",nbincoarse,fLowLim,fLowLim+nbincoarse*fStepCoarse);
214 if(fZCombv)delete fZCombv;
215 fZCombv = new TH1F("fZCombv","Z",nbinvcoarse,fLowLim,fLowLim+nbinvcoarse*3.*fStepCoarse);
216 if(fZCombf)delete fZCombf;
217 fZCombf = new TH1F("fZCombf","Z",nbinfine,fLowLim,fLowLim+nbinfine*fStepFine);
218 // Loop on modules of layer 1 (restricted to central modules)
0c6af5c9 219 for(Int_t module= fFirstL1; module<=fLastL1;module++){
220 if(module%4==0 || module%4==3)continue;
221 branch->GetEvent(module);
222 Int_t nrecp1 = itsRec->GetEntries();
223 for(Int_t j=0;j<nrecp1;j++){
224 AliITSRecPoint *recp = (AliITSRecPoint*)itsRec->At(j);
ecc64c3f 225 // Local coordinates of this recpoint
00a7cc50 226 lc[0]=recp->GetDetLocalX();
227 lc[2]=recp->GetDetLocalZ();
0c6af5c9 228 geom->LtoG(module,lc,gc);
ecc64c3f 229 // Global coordinates of this recpoints
230 gc[0]-=fX0; // Possible beam offset in the bending plane
231 gc[1]-=fY0; // " "
0c6af5c9 232 xc1[ind]=gc[0];
233 yc1[ind]=gc[1];
234 zc1[ind]=gc[2];
ecc64c3f 235 // azimuthal angle is computed in the interval 0 --> 2*pi
0c6af5c9 236 phi1[ind] = TMath::ATan2(gc[1],gc[0]);
237 if(phi1[ind]<0)phi1[ind]=2*TMath::Pi()+phi1[ind];
238 ind++;
239 }
7d62fb64 240 detTypeRec.ResetRecPoints();
0c6af5c9 241 }
ecc64c3f 242 ind = 0; // the running index is reset for Layer 2
0c6af5c9 243 for(Int_t module= fFirstL2; module<=fLastL2;module++){
244 branch->GetEvent(module);
245 Int_t nrecp2 = itsRec->GetEntries();
246 for(Int_t j=0;j<nrecp2;j++){
247 AliITSRecPoint *recp = (AliITSRecPoint*)itsRec->At(j);
00a7cc50 248 lc[0]=recp->GetDetLocalX();
249 lc[2]=recp->GetDetLocalZ();
0c6af5c9 250 geom->LtoG(module,lc,gc);
251 gc[0]-=fX0;
252 gc[1]-=fY0;
253 xc2[ind]=gc[0];
254 yc2[ind]=gc[1];
255 zc2[ind]=gc[2];
256 phi2[ind] = TMath::ATan2(gc[1],gc[0]);
257 if(phi2[ind]<0)phi2[ind]=2*TMath::Pi()+phi2[ind];
258 ind++;
259 }
7d62fb64 260 detTypeRec.ResetRecPoints();
0c6af5c9 261 }
32e449be 262
263 Int_t nolines=0;
ecc64c3f 264 for(Int_t i=0;i<nrpL1;i++){ // loop on L1 RP
265 Float_t r1=TMath::Sqrt(xc1[i]*xc1[i]+yc1[i]*yc1[i]); // radius L1 RP
266 for(Int_t j=0;j<nrpL2;j++){ // loop on L2 RP
267 Float_t diff = TMath::Abs(phi2[j]-phi1[i]); // diff in azimuth
268 if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff; //diff<pi
269 if(diff<fDiffPhiMax){ // cut on 10 milliradians by def.
270 Float_t r2=TMath::Sqrt(xc2[j]*xc2[j]+yc2[j]*yc2[j]); // radius L2 RP
271 Float_t zr0=(r2*zc1[i]-r1*zc2[j])/(r2-r1); //Z @ null radius
0c6af5c9 272 fZCombf->Fill(zr0);
273 fZCombc->Fill(zr0);
ecc64c3f 274 fZCombv->Fill(zr0);
32e449be 275 Double_t pA[3];
276 Double_t pB[3];
277 pA[0]=xc1[i];
278 pA[1]=yc1[i];
279 pA[2]=zc1[i];
280 pB[0]=xc2[j];
281 pB[1]=yc2[j];
282 pB[2]=zc2[j];
283 MakeTracklet(pA,pB,nolines);
0c6af5c9 284 }
285 }
286 }
287 delete [] xc1;
288 delete [] yc1;
289 delete [] zc1;
290 delete [] phi1;
291 delete [] xc2;
292 delete [] yc2;
293 delete [] zc2;
294 delete [] phi2;
6fd990e3 295 Double_t contents = fZCombc->GetEntries()- fZCombc->GetBinContent(0)-fZCombc->GetBinContent(nbincoarse+1);
296 if(contents<1.){
297 // Warning("FindVertexForCurrentEvent","Insufficient number of rec. points\n");
0c6af5c9 298 ResetHistograms();
6fd990e3 299 itsLoader->UnloadRecPoints();
300 fCurrentVertex = new AliESDVertex(0.,5.3,-1);
301 return;
0c6af5c9 302 }
ecc64c3f 303 /*
304 else {
305 if(fDebug>0)cout<<"Number of entries in hist. "<<fZCombc->GetEntries()<<endl;
306 }
307 */
308
309 TH1F *hc = fZCombc;
310 Bool_t goon = kFALSE;
311 Int_t cnt = 0;
312 Int_t bi;
313
314 do {
315 goon = kFALSE;
316 cnt++;
317 bi = hc->GetMaximumBin(); // bin with maximal content on coarse histogram
318 if(hc->GetBinContent(bi)<3){
319 if(cnt==1)goon = kTRUE;
320 hc = fZCombv;
321 }
322 } while(goon);
323
324
325 Float_t centre = hc->GetBinCenter(bi); // z value of the bin with maximal content
326
327 // evaluation of the centroid
328 Int_t ii1=TMath::Max(bi-3,1);
329 Int_t ii2=TMath::Min(bi+3,hc->GetNbinsX());
330 centre = 0.;
331 Int_t nn=0;
332 for(Int_t ii=ii1;ii<ii2;ii++){
333 centre+= hc->GetBinCenter(ii)*hc->GetBinContent(ii);
334 nn+=static_cast<Int_t>(hc->GetBinContent(ii));
335 }
3eb8b116 336 if (nn) centre/=nn;
ecc64c3f 337 /*
338 if(fDebug>0){
339 cout<<"Value of center "<<centre<<endl;
340 cout<<"Population of 3 central bins: "<<hc->GetBinContent(bi-1)<<", ";
341 cout<<hc->GetBinContent(bi)<<", ";
342 cout<<hc->GetBinContent(bi+1)<<endl;
343 }
344 */
345 // n1 is the bin number of fine histogram containing the point located 1 coarse bin less than "centre"
346 Int_t n1 = static_cast<Int_t>((centre-hc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fStepFine);
347 // n2 is the bin number of fine histogram containing the point located 1 coarse bin more than "centre"
348 Int_t n2 = static_cast<Int_t>((centre+hc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fStepFine);
349 if(n1<1)n1=1;
350 if(n2>nbinfine)n2=nbinfine;
0c6af5c9 351 Int_t niter = 0;
ecc64c3f 352 goon = kTRUE;
353 Int_t num=0;
354 Bool_t last = kFALSE;
355
356 while(goon || last){
0c6af5c9 357 fZFound = 0.;
358 fZsig = 0.;
359 num=0;
ecc64c3f 360 // at the end of the loop:
361 // fZFound = N*(Average Z) - where N is the number of tracklets
362 // num=N
363 // fZsig = N*Q - where Q is the average of Z*Z
0c6af5c9 364 for(Int_t n=n1;n<=n2;n++){
365 fZFound+=fZCombf->GetBinCenter(n)*fZCombf->GetBinContent(n);
366 num+=static_cast<Int_t>(fZCombf->GetBinContent(n));
367 fZsig+=fZCombf->GetBinCenter(n)*fZCombf->GetBinCenter(n)*fZCombf->GetBinContent(n);
368 }
369 if(num<2){
12e26c26 370 fZsig = 5.3; // Default error from the beam sigmoid
0c6af5c9 371 }
ecc64c3f 372 else{
979e3647 373 Float_t radi = fZsig/(num-1)-fZFound*fZFound/num/(num-1);
ecc64c3f 374 // radi = square root of sample variance of Z
979e3647 375 if(radi>0.)fZsig=TMath::Sqrt(radi);
12e26c26 376 else fZsig=5.3; // Default error from the beam sigmoid
ecc64c3f 377 // fZfound - Average Z
0c6af5c9 378 fZFound/=num;
379 }
ecc64c3f 380 if(!last){
381 // goon is true if the distance between the found Z and the lower bin differs from the distance between the found Z and
382 // the upper bin by more than tolerance (0.002)
383 goon = TMath::Abs(TMath::Abs(fZFound-fZCombf->GetBinCenter(n1))-TMath::Abs(fZFound-fZCombf->GetBinCenter(n2)))>fTolerance;
384 // a window in the fine grained histogram is centered aroung the found Z. The width is 2 bins of
385 // the coarse grained histogram
6fd990e3 386 if(num>0){
387 n1 = static_cast<Int_t>((fZFound-hc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fStepFine);
388 if(n1<1)n1=1;
389 n2 = static_cast<Int_t>((fZFound+hc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fStepFine);
390 if(n2>nbinfine)n2=nbinfine;
391 /*
392 if(fDebug>0){
393 cout<<"Search restricted to n1 = "<<n1<<", n2= "<<n2<<endl;
394 cout<<"z1= "<<fZCombf->GetBinCenter(n1)<<", n2 = "<<fZCombf->GetBinCenter(n2)<<endl;
395 }
396 */
397 }
398 else {
399 n1 = static_cast<Int_t>((centre-(niter+2)*hc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fStepFine);
400 n2 = static_cast<Int_t>((centre+(niter+2)*hc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fStepFine);
401 if(n1<1)n1=1;
402 if(n2>nbinfine)n2=nbinfine;
ecc64c3f 403 }
ecc64c3f 404 niter++;
405 // no more than 10 adjusting iterations
406 if(niter>=10){
407 goon = kFALSE;
6fd990e3 408 // Warning("FindVertexForCurrentEvent","The procedure does not converge\n");
ecc64c3f 409 }
410
411 if((fZsig> 0.0001) && !goon && num>5){
412 last = kTRUE;
413 n1 = static_cast<Int_t>((fZFound-fZsig-fZCombf->GetBinLowEdge(0))/fZCombf->GetBinWidth(0));
414 if(n1<1)n1=1;
415 n2 = static_cast<Int_t>((fZFound+fZsig-fZCombf->GetBinLowEdge(0))/fZCombf->GetBinWidth(0));
416 if(n2>nbinfine)n2=nbinfine;
417 /*
418 if(fDebug>0){
419 cout<<"FINAL Search restricted to n1 = "<<n1<<", n2= "<<n2<<endl;
420 cout<<"z1= "<<fZCombf->GetBinCenter(n1)<<", n2 = "<<fZCombf->GetBinCenter(n2)<<endl;
421 }
422 */
423 }
424 }
425 else {
426 last = kFALSE;
0c6af5c9 427 }
ecc64c3f 428
0c6af5c9 429 }
ecc64c3f 430 // if(fDebug>0)cout<<"Numer of Iterations "<<niter<<endl<<endl;
9679d5db 431 // if(num!=0)fZsig/=TMath::Sqrt(num);
12e26c26 432 if (fZsig<=0) fZsig=5.3; // Default error from the beam sigmoid
d681bb2d 433 fCurrentVertex = new AliESDVertex(fZFound,fZsig,num);
0c6af5c9 434 fCurrentVertex->SetTitle("vertexer: B");
435 ResetHistograms();
6fd990e3 436 itsLoader->UnloadRecPoints();
437 return;
0c6af5c9 438}
439
440//_____________________________________________________________________
441void AliITSVertexerZ::ResetHistograms(){
442 // delete TH1 data members
443 if(fZCombc)delete fZCombc;
444 if(fZCombf)delete fZCombf;
ecc64c3f 445 if(fZCombv)delete fZCombv;
0c6af5c9 446 fZCombc = 0;
447 fZCombf = 0;
ecc64c3f 448 fZCombv = 0;
0c6af5c9 449}
450
451//______________________________________________________________________
452void AliITSVertexerZ::FindVertices(){
453 // computes the vertices of the events in the range FirstEvent - LastEvent
454 AliRunLoader *rl = AliRunLoader::GetRunLoader();
455 AliITSLoader* itsLoader = (AliITSLoader*) rl->GetLoader("ITSLoader");
456 itsLoader->ReloadRecPoints();
457 for(Int_t i=fFirstEvent;i<=fLastEvent;i++){
6fd990e3 458 // cout<<"Processing event "<<i<<endl;
0c6af5c9 459 rl->GetEvent(i);
460 FindVertexForCurrentEvent(i);
461 if(fCurrentVertex){
462 WriteCurrentVertex();
463 }
ecc64c3f 464 /*
0c6af5c9 465 else {
466 if(fDebug>0){
467 cout<<"Vertex not found for event "<<i<<endl;
468 cout<<"fZFound = "<<fZFound<<", fZsig= "<<fZsig<<endl;
469 }
470 }
ecc64c3f 471 */
0c6af5c9 472 }
473}
474
475//________________________________________________________
476void AliITSVertexerZ::PrintStatus() const {
477 // Print current status
478 cout <<"=======================================================\n";
479 cout <<" First layer first and last modules: "<<fFirstL1<<", ";
480 cout <<fLastL1<<endl;
481 cout <<" Second layer first and last modules: "<<fFirstL2<<", ";
482 cout <<fLastL2<<endl;
483 cout <<" Max Phi difference: "<<fDiffPhiMax<<endl;
484 cout <<"Limits for Z histograms: "<<fLowLim<<"; "<<fHighLim<<endl;
485 cout <<"Bin sizes for coarse and fine z histos "<<fStepCoarse<<"; "<<fStepFine<<endl;
486 cout <<" Current Z "<<fZFound<<"; Z sig "<<fZsig<<endl;
487 cout <<" Debug flag: "<<fDebug<<endl;
488 cout <<"First event to be processed "<<fFirstEvent;
489 cout <<"\n Last event to be processed "<<fLastEvent<<endl;
ecc64c3f 490 if(fZCombc){
491 cout<<"fZCombc exists - entries="<<fZCombc->GetEntries()<<endl;
492 }
493 else{
494 cout<<"fZCombc does not exist\n";
495 }
496 if(fZCombv){
497 cout<<"fZCombv exists - entries="<<fZCombv->GetEntries()<<endl;
498 }
499 else{
500 cout<<"fZCombv does not exist\n";
501 }
502 if(fZCombf){
503 cout<<"fZCombf exists - entries="<<fZCombv->GetEntries()<<endl;
504 }
505 else{
506 cout<<"fZCombf does not exist\n";
507 }
0c6af5c9 508
509 cout <<"=======================================================\n";
510}
511