1 /**************************************************************************
2 * Copyright(c) 1998-2003, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ////////////////////////////////////////////////////
19 // Stand alone tracker class //
20 // Origin: Elisabetta Crescio //
21 // e-mail: crescio@to.infn.it //
22 // tracks are saved as AliITStrackV2 objects //
23 ////////////////////////////////////////////////////
29 #include <TObjArray.h>
33 #include "AliESDVertex.h"
34 #include "AliESDtrack.h"
35 #include "AliITSRiemannFit.h"
36 #include "AliITSVertexer.h"
37 #include "AliITSclusterTable.h"
38 #include "AliITSclusterV2.h"
39 #include "AliITSgeom.h"
40 #include "AliITStrackSA.h"
41 #include "AliITStrackerSA.h"
44 ClassImp(AliITStrackerSA)
46 //____________________________________________________________________________
47 AliITStrackerSA::AliITStrackerSA():AliITStrackerMI(){
48 // Default constructor
52 //____________________________________________________________________________
53 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom):AliITStrackerMI(geom)
55 // Standard constructor (Vertex is known and passed to this obj.)
62 //____________________________________________________________________________
63 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliESDVertex *vert):AliITStrackerMI(geom)
65 // Standard constructor (Vertex is known and passed to this obj.)
72 //____________________________________________________________________________
73 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliITSVertexer *vertexer):AliITStrackerMI(geom)
75 // Standard constructor (Vertex is unknown - vertexer is passed to this obj)
82 //____________________________________________________________________________
83 AliITStrackerSA::AliITStrackerSA(const AliITStrackerSA& tracker):AliITStrackerMI(){
85 fPhiEstimate = tracker.fPhiEstimate;
86 for(Int_t i=0;i<2;i++){
87 fPoint1[i]=tracker.fPoint1[i];
88 fPoint2[i]=tracker.fPoint2[i];
89 fPoint3[i]=tracker.fPoint3[i];
90 fPointc[i]=tracker.fPointc[i];
92 fLambdac = tracker.fLambdac;
93 fPhic = tracker.fPhic;
94 fCoef1 = tracker.fCoef1;
95 fCoef2 = tracker.fCoef2;
96 fCoef3 = tracker.fCoef3;
97 fNloop = tracker.fNloop;
98 fPhiWin = tracker.fPhiWin;
99 fLambdaWin = tracker.fLambdaWin;
100 if(tracker.fVertexer && tracker.fVert){
101 fVert = new AliESDVertex(*tracker.fVert);
104 fVert = tracker.fVert;
106 fVertexer = tracker.fVertexer;
107 fGeom = tracker.fGeom;
108 fTable = tracker.fTable;
109 fListOfTracks = tracker.fListOfTracks;
111 //______________________________________________________________________
112 AliITStrackerSA& AliITStrackerSA::operator=(const AliITStrackerSA& /*source*/){
113 // Assignment operator. This is a function which is not allowed to be
115 Error("operator=","Assignment operator not allowed\n");
119 //____________________________________________________________________________
120 AliITStrackerSA::~AliITStrackerSA(){
122 // if fVertexer is not null, the AliESDVertex obj. is owned by this class
123 // and is deleted here
125 if(fVert)delete fVert;
130 if(fPhiWin)delete []fPhiWin;
131 if(fLambdaWin)delete []fLambdaWin;
133 fListOfTracks->Delete();
136 //____________________________________________________________________________
137 void AliITStrackerSA::Init(){
138 // Reset all data members
140 for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
155 fListOfTracks=new TObjArray(0,0);
157 //_______________________________________________________________________
158 void AliITStrackerSA::ResetForFinding(){
159 // Reset data members used in all loops during track finding
161 for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
169 fListOfTracks->Delete();
171 //____________________________________________________________________________
172 void AliITStrackerSA::FindTracks(TTree *out,Int_t evnumber){
174 /**************************************************************************
175 * This function finds primary tracks.
178 * Example: to execute function with only the ITS (no combined tracking *
179 * with TPC+ITS) and requiring 5/6 points to define a good track *
180 * call SetSixPoinbts(kFALSE) in advance and then *
181 * use: FindTracks(treein,treeout,evnumber) *
182 * to execute combined tracking, before using FindTracks, use *
184 *************************************************************************/
187 Fatal("FindTracks","ITS cluster tree is not accessed - Abort!!!\n Please use method SetClusterTree to pass the pointer to the tree\n");
192 if(fVert)delete fVert;
193 fVert = fVertexer->FindVertexForCurrentEvent(evnumber);
196 gAlice->GetEvent(evnumber);
198 Fatal("FindTracks","Vertex is missing\n");
202 Double_t primaryVertex[3];
203 Double_t errorsprimvert[3];
204 fVert->GetXYZ(primaryVertex);
205 fVert->GetSigmaXYZ(errorsprimvert);
206 if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
207 Warning("FindTracks","Set errors on vertex positions x and y at 0.0001");
208 errorsprimvert[0]=0.0001;
209 errorsprimvert[1]=0.0001;
211 fVert->PrintStatus();
214 //Fill array with cluster indices for each module
216 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
217 fTable->FillArray(fITSclusters);
218 fTable->FillArrayCoorAngles();
222 //Fill tree for found tracks
223 AliITStrackV2* outrack=0;
224 TBranch* branch=out->Branch("tracks","AliITStrackV2",&outrack,32000,0);
225 if (!branch) out->Branch("tracks","AliITStrackV2",&outrack,32000,3);
226 else branch->SetAddress(&outrack);
229 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
230 for(Int_t i=0;i<fGeom->GetNlayers();i++){
231 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
233 // firstmod [i] number of the first module in the ITS layer i.
235 AliITSlayer &layer=fgLayers[0]; // first layer
237 Int_t dim=layer.GetNumberOfClusters();
238 //loop on the different windows
239 for(Int_t nloop=0;nloop<fNloop;nloop++){
240 for(Int_t ncl=0;ncl<dim;ncl++){
241 //loop starting from layer 0
244 AliITSclusterV2* cl = layer.GetCluster(ncl);
245 if(cl->IsUsed()==1) continue;
246 if(cl->TestBit(kSAflag)==kTRUE) continue;
248 fPhic = fTable->GetPhiCluster(0,ncl);
249 fLambdac = fTable->GetLambdaCluster(0,ncl);
250 fPhiEstimate = fPhic;
251 AliITStrackSA* trs = new AliITStrackSA();
252 fPoint1[0]=primaryVertex[0];
253 fPoint1[1]=primaryVertex[1];
254 fPoint2[0]=fTable->GetXCluster(0,ncl);
255 fPoint2[1]=fTable->GetYCluster(0,ncl);
257 Int_t * nn = new Int_t[fGeom->GetNlayers()];//counter for clusters on each layer
258 for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;}
259 nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
260 nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
264 fPoint3[0] = fPointc[0];
265 fPoint3[1] = fPointc[1];
267 nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
268 if(nn[1]==0 && nn[2]==0) pflag=0;
269 if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();}
270 if(nn[2]!=0 && nn[1]==0){
272 fPoint3[0]=fPointc[0];
273 fPoint3[1]=fPointc[1];
276 nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
278 if(nn[3]!=0) UpdatePoints();
279 nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
281 if(nn[4]!=0) UpdatePoints();
282 nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
286 Int_t numberofpoints;
287 if(fSixPoints) numberofpoints=6; //check of the candidate track
288 else numberofpoints=5; //if track is good (with the required number
289 for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){ //of points) it is written on file
290 if(nn[nnp]!=0) layOK+=1;
292 if(layOK>=numberofpoints){
293 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
295 Int_t nct = trs->GetNumberOfClustersSA();
297 Int_t index = trs->GetClusterIndexSA(nct);
298 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
299 if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
308 Int_t nct = tr2->GetNumberOfClusters();
311 Int_t index = tr2->GetClusterIndex(nct);
312 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
318 Int_t nct = trs->GetNumberOfClustersSA();
320 Int_t index = trs->GetClusterIndexSA(nct);
321 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
322 if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
328 }//end loop on clusters of layer1
332 //if 5/6 points are required, second loop starting
333 //from second layer, to find tracks with point of
337 // counter for clusters on each layer
338 Int_t * nn = new Int_t[fGeom->GetNlayers()-1];
339 for(Int_t nloop=0;nloop<fNloop;nloop++){
340 AliITSlayer &layer2=fgLayers[1]; //loop on layer 2
341 Int_t ncl2=layer2.GetNumberOfClusters();
342 while(ncl2--){ //loop starting from layer 2
345 AliITSclusterV2* cl = layer2.GetCluster(ncl2);
346 if(cl->IsUsed()==1) continue;
347 if(cl->TestBit(kSAflag)==kTRUE) continue;
348 fPhic = fTable->GetPhiCluster(1,ncl2);
349 fLambdac = fTable->GetLambdaCluster(1,ncl2);
350 fPhiEstimate = fPhic;
351 AliITStrackSA* trs = new AliITStrackSA();
353 fPoint1[0]=primaryVertex[0];
354 fPoint1[1]=primaryVertex[1];
355 fPoint2[0]=fTable->GetXCluster(1,ncl2);;
356 fPoint2[1]=fTable->GetYCluster(1,ncl2);;
358 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
359 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
360 nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
361 trs,primaryVertex[2],pflag);
369 for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
370 if(nn[nnp]!=0) fl+=1;
373 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
375 Int_t nct = trs->GetNumberOfClustersSA();
377 Int_t index = trs->GetClusterIndexSA(nct);
378 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
379 if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
386 Int_t nct = tr2->GetNumberOfClusters();
388 Int_t index = tr2->GetClusterIndex(nct);
389 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
395 Int_t nct = trs->GetNumberOfClustersSA();
397 Int_t index = trs->GetClusterIndexSA(nct);
398 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
400 if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
404 }//end loop on clusters of layer2
407 } // if(!fSixPoints....
410 delete fTable; fTable=0;
414 //______________________________________________________________________
415 Int_t AliITStrackerSA::FindTracks(AliESD* event){
417 // Track finder using the ESD object
421 Fatal("FindTracks","ITS cluster tree is not accessed - Abort!!!\n Please use method SetClusterTree to pass the pointer to the tree\n");
426 Double_t errorsprimvert[3];
427 Double_t primaryVertex[3];
428 event->GetVertex()->GetXYZ(primaryVertex);
429 event->GetVertex()->GetSigmaXYZ(errorsprimvert);
431 if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
432 // Warning("FindTracks","Set errors on vertex positions x and y at 0.005");
433 errorsprimvert[0]=0.005;
434 errorsprimvert[1]=0.005;
436 //Fill array with cluster indices for each module
438 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
439 fTable->FillArray(fITSclusters);
440 fTable->FillArrayCoorAngles();
442 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
443 for(Int_t i=0;i<fGeom->GetNlayers();i++){
444 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
446 // firstmod [i] number of the first module in the ITS layer i.
449 AliITSlayer &layer=fgLayers[0];
451 Int_t dim=layer.GetNumberOfClusters();
452 //loop on the different windows
453 for(Int_t nloop=0;nloop<fNloop;nloop++){
454 for(Int_t ncl=0;ncl<dim;ncl++){ //loop starting from layer 0
458 AliITSclusterV2* cl = (AliITSclusterV2*)layer.GetCluster(ncl);
460 if(cl->IsUsed()==1) continue;
461 if(cl->TestBit(kSAflag)==kTRUE) continue;
462 if (cl->GetQ()<=0) continue;
464 fPhic = fTable->GetPhiCluster(0,ncl);
465 fLambdac = fTable->GetLambdaCluster(0,ncl);
467 if (TMath::Abs(fLambdac)>0.26*TMath::Pi()) continue;
469 fPhiEstimate = fPhic;
470 AliITStrackSA* trs = new AliITStrackSA();
471 fPoint1[0]=primaryVertex[0];
472 fPoint1[1]=primaryVertex[1];
473 fPoint2[0]=fTable->GetXCluster(0,ncl);
474 fPoint2[1]=fTable->GetYCluster(0,ncl);
475 Int_t * nn = new Int_t[fGeom->GetNlayers()];//counter for clusters on each layer
476 for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;}
477 nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
479 nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
482 fPoint3[0] = fPointc[0];
483 fPoint3[1] = fPointc[1];
485 nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
486 if(nn[1]==0 && nn[2]==0) pflag=0;
487 if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();}
488 if(nn[2]!=0 && nn[1]==0){
490 fPoint3[0]=fPointc[0];
491 fPoint3[1]=fPointc[1];
494 nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
496 if(nn[3]!=0) UpdatePoints();
497 nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
499 if(nn[4]!=0) UpdatePoints();
500 nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
504 Int_t numberofpoints;
505 if(fSixPoints) numberofpoints=6; //check of the candidate track
506 else numberofpoints=5; //if track is good (with the required number
507 for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){ //of points) it is written on file
508 if(nn[nnp]!=0) layOK+=1;
510 if(layOK>=numberofpoints){
511 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
513 Int_t nct = trs->GetNumberOfClustersSA();
515 Int_t index = trs->GetClusterIndexSA(nct);
516 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
517 if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
523 AliESDtrack outtrack;
524 outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
525 event->AddTrack(&outtrack);
527 Int_t nct = tr2->GetNumberOfClusters();
529 Int_t index = tr2->GetClusterIndex(nct);
530 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
536 Int_t nct = trs->GetNumberOfClustersSA();
538 Int_t index = trs->GetClusterIndexSA(nct);
539 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
540 if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
547 }//end loop on clusters of layer1
553 //if 5/6 points are required, second loop starting
554 //from second layer, to find tracks with point of
558 // counter for clusters on each layer
559 Int_t * nn = new Int_t[fGeom->GetNlayers()-1];
560 for(Int_t nloop=0;nloop<fNloop;nloop++){
561 AliITSlayer &layer2=fgLayers[1];
562 Int_t ncl2=layer2.GetNumberOfClusters();
563 while(ncl2--){ //loop starting from layer 2
566 AliITSclusterV2* cl = layer2.GetCluster(ncl2);
567 if(cl->IsUsed()==1) continue;
568 if(cl->TestBit(kSAflag)==kTRUE) continue;
569 fPhic = fTable->GetPhiCluster(1,ncl2);
570 fLambdac = fTable->GetLambdaCluster(1,ncl2);
571 fPhiEstimate = fPhic;
572 AliITStrackSA* trs = new AliITStrackSA();
573 fPoint1[0]=primaryVertex[0];
574 fPoint1[1]=primaryVertex[1];
575 fPoint2[0]=fTable->GetXCluster(1,ncl2);
576 fPoint2[1]=fTable->GetYCluster(1,ncl2);
578 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
579 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
580 nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
581 trs,primaryVertex[2],pflag);
589 for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
590 if(nn[nnp]!=0) fl+=1;
593 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
595 Int_t nct = trs->GetNumberOfClustersSA();
597 Int_t index = trs->GetClusterIndexSA(nct);
598 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
599 if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
606 AliESDtrack outtrack;
607 outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
608 event->AddTrack(&outtrack);
610 Int_t nct = tr2->GetNumberOfClusters();
612 Int_t index = tr2->GetClusterIndex(nct);
613 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
619 Int_t nct = trs->GetNumberOfClustersSA();
621 Int_t index = trs->GetClusterIndexSA(nct);
622 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
624 if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
629 }//end loop on clusters of layer2
635 delete fTable;fTable=0;
636 Info("FindTracks","Number of found tracks: %d",event->GetNumberOfTracks());
642 //________________________________________________________________________
644 AliITStrackV2* AliITStrackerSA::FitTrack(AliITStrackSA* tr,Double_t *primaryVertex,Double_t *errorsprimvert){
645 //fit of the found track
648 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
649 for(Int_t i=0;i<fGeom->GetNlayers();i++){
650 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
653 Int_t nclusters = tr->GetNumberOfClustersSA();
654 TObjArray** listlayer = new TObjArray*[fGeom->GetNlayers()];
655 for(Int_t i=0;i<fGeom->GetNlayers();i++){
656 listlayer[i] = new TObjArray(0,0);
666 Int_t * nnn = new Int_t[fGeom->GetNlayers()];
667 for(Int_t i=0;i<fGeom->GetNlayers();i++)nnn[i]=0;
669 for(Int_t ncl=0;ncl<nclusters;ncl++){
670 Int_t index = tr->GetClusterIndexSA(ncl);
671 AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(index);
673 if(cl->TestBit(kSAflag)==kTRUE) cl->ResetBit(kSAflag);
674 Int_t lay = (index & 0xf0000000) >> 28;
675 if(lay==0) { listlayer[0]->AddLast(cl); clind0[nnn[0]]=index;nnn[0]++;}
676 if(lay==1) { listlayer[1]->AddLast(cl); clind1[nnn[1]]=index;nnn[1]++;}
677 if(lay==2) { listlayer[2]->AddLast(cl); clind2[nnn[2]]=index;nnn[2]++;}
678 if(lay==3) { listlayer[3]->AddLast(cl); clind3[nnn[3]]=index;nnn[3]++;}
679 if(lay==4) { listlayer[4]->AddLast(cl); clind4[nnn[4]]=index;nnn[4]++;}
680 if(lay==5) { listlayer[5]->AddLast(cl); clind5[nnn[5]]=index;nnn[5]++;}
685 Int_t * end = new Int_t[fGeom->GetNlayers()];
686 for(Int_t i=0;i<fGeom->GetNlayers();i++){
687 if(listlayer[i]->GetEntries()==0) end[i]=1;
688 else end[i]=listlayer[i]->GetEntries();
691 for(Int_t l1=0;l1<end[0];l1++){//loop on layer 1
692 AliITSclusterV2* cl0 = (AliITSclusterV2*)listlayer[0]->At(l1);
693 TVector3** recp = new TVector3*[3];
694 TVector3** errs = new TVector3*[3];
695 recp[0] = new TVector3(primaryVertex[0],primaryVertex[1],primaryVertex[2]);
696 errs[0] = new TVector3(errorsprimvert[0],errorsprimvert[1],errorsprimvert[2]);
697 Double_t x1,y1,z1,sx1,sy1,sz1;
698 Double_t x2,y2,z2,sx2,sy2,sz2;
699 AliITSclusterV2* p1=0;
700 AliITSclusterV2* p2=0;
701 Int_t index1=clind0[l1];
703 for(Int_t l2=0;l2<end[1];l2++){//loop on layer 2
704 AliITSclusterV2* cl1 = (AliITSclusterV2*)listlayer[1]->At(l2);
706 for(Int_t l3=0;l3<end[2];l3++){ //loop on layer 3
707 AliITSclusterV2* cl2 = (AliITSclusterV2*)listlayer[2]->At(l3);
709 if(cl0==0 && cl1!=0) {
710 p2 = cl2;index1=clind2[l3];
714 if(cl0!=0 && cl1==0){
716 p2=cl2;index2=clind2[l3];
718 if(cl0!=0 && cl1!=0){
722 Int_t lay1=(index1 & 0xf0000000) >> 28;
723 Int_t cln1=(index1 & 0x0fffffff) >> 00;
724 Int_t lay2=(index2 & 0xf0000000) >> 28;
725 Int_t cln2=(index2 & 0x0fffffff) >> 00;
726 x1 = fTable->GetXCluster(lay1,cln1);
727 x2 = fTable->GetXCluster(lay2,cln2);
728 y1 = fTable->GetYCluster(lay1,cln1);
729 y2 = fTable->GetYCluster(lay2,cln2);
730 z1 = fTable->GetZCluster(lay1,cln1);
731 z2 = fTable->GetZCluster(lay2,cln2);
732 sx1 = fTable->GetXClusterError(lay1,cln1);
733 sx2 = fTable->GetXClusterError(lay2,cln2);
734 sy1 = fTable->GetYClusterError(lay1,cln1);
735 sy2 = fTable->GetYClusterError(lay2,cln2);
736 sz1 = fTable->GetZClusterError(lay1,cln1);
737 sz2 = fTable->GetZClusterError(lay2,cln2);
738 Double_t phi1 = fTable->GetPhiCluster(lay1,cln1);
739 Int_t module1 = p1->GetDetectorIndex()+firstmod[0];
740 recp[1] = new TVector3(x1,y1,z1);
741 errs[1] = new TVector3(sx1,sy1,sz1);
742 recp[2] = new TVector3(x2,y2,z2);
743 errs[2] = new TVector3(sx2,sy2,sz2);
745 //fit on the Riemann sphere
746 Float_t seed1,seed2,seed3;
747 AliITSRiemannFit fit;
748 Int_t rf = fit.FitHelix(3,recp,errs,seed1,seed2,seed3); //this gives phi,tgl,curvature to start Kalman Filter
750 for(Int_t i=1;i<3;i++){
760 if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
761 phi=seed1+1.5*TMath::Pi();
764 if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
765 phi=seed1+0.5*TMath::Pi();
768 if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
769 phi=seed1-0.5*TMath::Pi();
774 if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
775 phi=seed1+0.5*TMath::Pi();
778 if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
779 phi=seed1-0.5*TMath::Pi();
782 if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
783 phi=seed1-1.5*TMath::Pi();
788 Int_t layer,ladder,detector;
789 fGeom->GetModuleId(module1,layer,ladder,detector);
790 Float_t yclu1 = p1->GetY();
791 Float_t zclu1 = p1->GetZ();
792 Double_t cv=Curvature(primaryVertex[0],primaryVertex[1],x1,y1,x2,y2);
794 for(Int_t l4=0;l4<end[3];l4++){ //loop on layer 4
795 AliITSclusterV2* cl3 = (AliITSclusterV2*)listlayer[3]->At(l4);
796 for(Int_t l5=0;l5<end[4];l5++){ //loop on layer 5
797 AliITSclusterV2* cl4 = (AliITSclusterV2*)listlayer[4]->At(l5);
798 for(Int_t l6=0;l6<end[5];l6++){ //loop on layer 6
799 AliITSclusterV2* cl5 = (AliITSclusterV2*)listlayer[5]->At(l6);
800 AliITStrackSA* trac = new AliITStrackSA(fGeom,layer,ladder,detector,yclu1,zclu1,phi,tgl,cv,1);
802 if(cl5!=0) trac->AddClusterV2(5,(clind5[l6] & 0x0fffffff)>>0);
803 if(cl4!=0) trac->AddClusterV2(4,(clind4[l5] & 0x0fffffff)>>0);
804 if(cl3!=0) trac->AddClusterV2(3,(clind3[l4] & 0x0fffffff)>>0);
805 if(cl2!=0) trac->AddClusterV2(2,(clind2[l3] & 0x0fffffff)>>0);
806 if(cl1!=0) trac->AddClusterV2(1,(clind1[l2] & 0x0fffffff)>>0);
807 if(cl0!=0) trac->AddClusterV2(0,(clind0[l1] & 0x0fffffff)>>0);
809 //fit with Kalman filter using AliITStrackerMI::RefitAt()
811 AliITStrackMI* ot = new AliITStrackSA(*trac);
813 ot->ResetCovariance();
816 if(RefitAt(49.,ot,trac)){ //fit from layer 1 to layer 6
818 AliITStrackMI *otrack2 = new AliITStrackMI(*ot);
819 otrack2->ResetCovariance();
820 otrack2->ResetClusters();
821 //fit from layer 6 to layer 1
822 if(RefitAt(3.7,otrack2,ot)) {
823 fListOfTracks->AddLast(otrack2);
836 for(Int_t i=1;i<3;i++){
850 Int_t dim=fListOfTracks->GetEntries();
852 for(Int_t i=0;i<fGeom->GetNlayers();i++){
860 AliITStrackV2* otrack =(AliITStrackV2*)FindTrackLowChiSquare(fListOfTracks,dim);
863 for(Int_t i=0;i<fGeom->GetNlayers();i++){
870 Int_t * indexc = new Int_t[fGeom->GetNlayers()];
871 for(Int_t i=0;i<fGeom->GetNlayers();i++) indexc[i]=0;
872 for(Int_t nind=0;nind<otrack->GetNumberOfClusters();nind++){
873 indexc[nind] = otrack->GetClusterIndex(nind);
875 AliITSclusterV2* cl0 = (AliITSclusterV2*)GetCluster(indexc[0]);
876 AliITSclusterV2* cl1 = (AliITSclusterV2*)GetCluster(indexc[1]);
877 AliITSclusterV2* cl2 = (AliITSclusterV2*)GetCluster(indexc[2]);
878 AliITSclusterV2* cl3 = (AliITSclusterV2*)GetCluster(indexc[3]);
879 AliITSclusterV2* cl4 = (AliITSclusterV2*)GetCluster(indexc[4]);
880 Int_t labl[3]={-1,-1,-1};
881 if(otrack->GetNumberOfClusters()==fGeom->GetNlayers()){
882 AliITSclusterV2* cl5 = (AliITSclusterV2*)GetCluster(indexc[5]);
883 labl[0]=cl5->GetLabel(0);
884 labl[1]=cl5->GetLabel(1);
885 labl[2]=cl5->GetLabel(2);
888 if(otrack->GetNumberOfClusters()==(fGeom->GetNlayers()-1)){
893 Int_t numberofpoints;
894 if(fSixPoints) numberofpoints=6;
895 else numberofpoints=5;
896 CookLabel(otrack,0.); //MI change - to see fake ratio
897 Int_t label = Label(cl0->GetLabel(0),cl1->GetLabel(0),
898 cl2->GetLabel(0),cl3->GetLabel(0),
899 cl4->GetLabel(0),labl[0],
900 cl0->GetLabel(1),cl1->GetLabel(1),
901 cl2->GetLabel(1),cl3->GetLabel(1),
902 cl4->GetLabel(1),labl[1],
903 cl0->GetLabel(2),cl1->GetLabel(2),
904 cl2->GetLabel(2),cl3->GetLabel(2),
905 cl4->GetLabel(2),labl[2],numberofpoints);
907 otrack->SetLabel(label);
908 for(Int_t i=0;i<fGeom->GetNlayers();i++){
917 //_______________________________________________________________________
918 void AliITStrackerSA::UseFoundTracksV2(Int_t evnum,TTree* treev2){
919 // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
920 //(or AliITStrackV2 tracks found with function FindTracks of this class)
925 if(fVert)delete fVert;
926 fVert = fVertexer->FindVertexForCurrentEvent(evnum);
929 gAlice->GetEvent(evnum);
931 Fatal("FindTracks","Vertex is missing\n");
935 Double_t primaryVertex[3];
936 fVert->GetXYZ(primaryVertex);
939 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
940 fTable->FillArray(fITSclusters);
941 fTable->FillArrayCoorAngles();
944 TBranch* bra = (TBranch*)treev2->GetBranch("tracks");
945 if(!bra) Warning("UseFoundTracksV2","No branch for track tree");
946 AliITStrackV2* ttrrt = new AliITStrackV2;
947 bra->SetAddress(&ttrrt);
949 for(Int_t nj=0;nj<treev2->GetEntries();nj++){
950 treev2->GetEvent(nj);
951 Int_t ncl = ttrrt->GetNumberOfClusters();
952 for(Int_t k=0;k<ncl;k++){
953 Int_t index = ttrrt->GetClusterIndex(k);
954 AliITSclusterV2* clui = (AliITSclusterV2*)GetCluster(index);
955 if(clui->IsUsed()==0) clui->Use();
963 //_______________________________________________________________________
964 void AliITStrackerSA::UseFoundTracksV2(AliESD *event){
965 // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
969 Double_t primaryVertex[3];
970 event->GetVertex()->GetXYZ(primaryVertex);
973 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
974 fTable->FillArray(fITSclusters);
975 fTable->FillArrayCoorAngles();
978 Int_t ntracks = event->GetNumberOfTracks();
980 AliESDtrack *esd=event->GetTrack(ntracks);
981 if ((esd->GetStatus()&
982 AliESDtrack::kITSin|AliESDtrack::kTPCin)==0) continue;
984 Int_t ncl = esd->GetITSclusters(idx);
985 for(Int_t clu=0; clu<ncl; clu++){
986 AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(idx[clu]);
987 if(cl->IsUsed()==0) cl->Use();
991 Info("UseFoundTracksV2","Clusters of tracks prolonged from TPC deleted");
997 //_______________________________________________________
998 Int_t AliITStrackerSA::SearchClusters(Int_t layer,Double_t phiwindow,Double_t lambdawindow, AliITStrackSA* trs,Double_t zvertex,Int_t pflag){
999 //function used to to find the clusters associated to the track
1001 AliITSlayer &lay = fgLayers[layer];
1002 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
1003 for(Int_t i=0;i<fGeom->GetNlayers();i++){
1004 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
1008 Float_t cx1,cx2,cy1,cy2;
1009 FindEquation(fPoint1[0],fPoint1[1],fPoint2[0],fPoint2[1],fPoint3[0],fPoint3[1],fCoef1,fCoef2,fCoef3);
1010 Int_t fun = FindIntersection(fCoef1,fCoef2,fCoef3,-(lay.GetR()*lay.GetR()),cx1,cy1,cx2,cy2);
1016 Double_t fi1 =TMath::ATan2(cy1,cx1);
1017 Double_t fi2 =TMath::ATan2(cy2,cx2);
1018 fPhiEstimate = ChoosePoint(fi1,fi2,fPhic);
1021 Double_t zed = TMath::Tan(fLambdac)*lay.GetR()+zvertex;
1022 Double_t zed1 = TMath::Tan(fLambdac+lambdawindow)*lay.GetR()+zvertex;
1023 Double_t zed2 = TMath::Tan(fLambdac-lambdawindow)*lay.GetR()+zvertex;
1025 Double_t fi = fPhiEstimate;
1026 Int_t nmod = lay.FindDetectorIndex(fi,zed);
1031 nmod += firstmod[layer];
1033 Int_t nm[8]={0,0,0,0,0,0,0,0};
1034 nm[0] = lay.FindDetectorIndex(fi+phiwindow,zed);
1035 nm[1] = lay.FindDetectorIndex(fi-phiwindow,zed);
1036 nm[2] = lay.FindDetectorIndex(fi,zed1);
1037 nm[3] = lay.FindDetectorIndex(fi,zed2);
1038 nm[4] = lay.FindDetectorIndex(fi+phiwindow,zed1);
1039 nm[5] = lay.FindDetectorIndex(fi-phiwindow,zed1);
1040 nm[6] = lay.FindDetectorIndex(fi+phiwindow,zed2);
1041 nm[7] = lay.FindDetectorIndex(fi-phiwindow,zed2);
1045 TArrayI* array =(TArrayI*)fTable->GetListOfClusters(nmod);
1046 TArrayI* listc = new TArrayI(array->GetSize());
1047 for(Int_t i=0;i<array->GetSize();i++){
1048 Int_t in=(Int_t)array->At(i);
1049 listc->AddAt(in,nn);
1056 for(Int_t h=k+1;h<8;h++){
1069 for(Int_t ii=0;ii<8;ii++){
1070 if(nm[ii]!=value && nm[ii]!=nmod && nm[ii]>=0){
1071 TArrayI* ar =(TArrayI*)fTable->GetListOfClusters(nm[ii]+firstmod[layer]);
1072 listc->Set(listc->GetSize()+ar->GetSize());
1073 for(Int_t j=0;j<ar->GetSize();j++){
1074 Int_t in=(Int_t)ar->At(j);
1075 listc->AddAt(in,nn);
1083 for(Int_t i=0;i<listc->GetSize();i++){
1084 Int_t index = (Int_t)listc->At(i);
1085 AliITSclusterV2* cllay = lay.GetCluster(index);
1086 if(cllay==0) continue;
1087 if(cllay->IsUsed()==1) continue;
1088 if(cllay->TestBit(kSAflag)==kTRUE) continue;
1089 Double_t phi = fTable->GetPhiCluster(layer,index);
1090 Double_t lambda= fTable->GetLambdaCluster(layer,index);
1092 if(TMath::Abs(fLambdac-lambda)<lambdawindow &&
1093 TMath::Abs(fPhiEstimate-phi)<phiwindow){
1096 if(trs->GetNumberOfClustersSA()==15){
1101 trs->AddClusterSA(layer,index);
1102 cllay->SetBit(kSAflag);
1104 fPointc[0]=fTable->GetXCluster(layer,index);
1105 fPointc[1]=fTable->GetYCluster(layer,index);
1117 //_______________________________________________________
1118 Int_t AliITStrackerSA::SearchClusters(Int_t layer,Double_t phiwindow,Double_t lambdawindow, AliITStrackSA* trs,Double_t zvertex,Int_t pflag){
1119 //function used to to find the clusters associated to the track
1121 AliITSlayer &lay = fgLayers[layer];
1122 Double_t r=lay.GetR(),tgl=TMath::Tan(fLambdac);
1125 Float_t cx1,cx2,cy1,cy2;
1126 FindEquation(fPoint1[0],fPoint1[1],fPoint2[0],fPoint2[1],fPoint3[0],fPoint3[1],fCoef1,fCoef2,fCoef3);
1127 if (FindIntersection(fCoef1,fCoef2,fCoef3,-r*r,cx1,cy1,cx2,cy2)==0)
1129 Double_t fi1=TMath::ATan2(cy1,cx1);
1130 Double_t fi2=TMath::ATan2(cy2,cx2);
1131 fPhiEstimate=ChoosePoint(fi1,fi2,fPhic);
1134 Double_t dz=r*lambdawindow*TMath::Sqrt(1+tgl*tgl) + 0.3*TMath::Abs(tgl);
1135 Double_t zmax=r*tgl + zvertex + dz;
1136 Double_t zmin=r*tgl + zvertex - dz;
1138 Int_t indmin = lay.FindClusterIndex(zmin);
1139 Int_t indmax = lay.FindClusterIndex(zmax);
1140 for (Int_t index=indmin; index<indmax; index++) {
1141 AliITSclusterV2 *c=lay.GetCluster(index);
1143 if (c->IsUsed()) continue;
1144 if (c->GetQ()<=0) continue;
1145 if (c->TestBit(kSAflag)==kTRUE) continue;
1147 Double_t phi =fTable->GetPhiCluster(layer,index);
1148 if (TMath::Abs(phi-fPhiEstimate)>phiwindow) continue;
1150 Double_t lambda=fTable->GetLambdaCluster(layer,index);
1151 if (TMath::Abs(lambda-fLambdac)>lambdawindow) continue;
1153 if(trs->GetNumberOfClustersSA()==15) return 0;
1155 trs->AddClusterSA(layer,index);
1159 fPointc[0]=fTable->GetXCluster(layer,index);
1160 fPointc[1]=fTable->GetYCluster(layer,index);
1168 //________________________________________________________________
1169 void AliITStrackerSA::UpdatePoints(){
1170 //update of points for the estimation of the curvature
1172 //fPoint1[0]=fPoint2[0];
1173 //fPoint1[1]=fPoint2[1];
1174 fPoint2[0]=fPoint3[0];
1175 fPoint2[1]=fPoint3[1];
1176 fPoint3[0]=fPointc[0];
1177 fPoint3[1]=fPointc[1];
1182 //___________________________________________________________________
1183 Int_t AliITStrackerSA::FindEquation(Float_t x1, Float_t y1, Float_t x2, Float_t y2, Float_t x3, Float_t y3,Float_t& a, Float_t& b, Float_t& c){
1185 //given (x,y) of three recpoints (in global coordinates)
1186 //returns the parameters a,b,c of circonference x*x + y*y +a*x + b*y +c
1188 Float_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1189 if(den==0) return 0;
1190 a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1191 b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1192 c = -x1*x1-y1*y1-a*x1-b*y1;
1195 //__________________________________________________________________________
1196 Int_t AliITStrackerSA::FindIntersection(Float_t a1, Float_t b1, Float_t c1, Float_t c2,Float_t& x1,Float_t& y1, Float_t& x2, Float_t& y2){
1198 //Finds the intersection between the circonference of the track and the circonference centered in (0,0) represented by one layer
1199 //c2 is -rlayer*rlayer
1203 Double_t aA = (b1*b1)/(a1*a1)+1;
1204 Double_t bB = (-2*m*b1/(a1*a1));
1205 Double_t cC = c2+(m*m)/(a1*a1);
1206 Double_t dD = bB*bB-4*aA*cC;
1209 y1 = (-bB+TMath::Sqrt(dD))/(2*aA);
1210 y2 = (-bB-TMath::Sqrt(dD))/(2*aA);
1211 x1 = (c2-c1-b1*y1)/a1;
1212 x2 = (c2-c1-b1*y2)/a1;
1216 //____________________________________________________________________
1217 Double_t AliITStrackerSA::Curvature(Double_t x1,Double_t y1,Double_t
1218 x2,Double_t y2,Double_t x3,Double_t y3){
1220 //calculates the curvature of track
1221 Double_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1222 if(den==0) return 0;
1223 Double_t a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1224 Double_t b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1225 Double_t c = -x1*x1-y1*y1-a*x1-b*y1;
1228 if((a*a+b*b-4*c)<0) return 0;
1229 Double_t rad = TMath::Sqrt(a*a+b*b-4*c)/2.;
1230 if(rad==0) return 0;
1232 if((x1>0 && y1>0 && x1<xc)) rad*=-1;
1233 if((x1<0 && y1>0 && x1<xc)) rad*=-1;
1234 // if((x1<0 && y1<0 && x1<xc)) rad*=-1;
1235 // if((x1>0 && y1<0 && x1<xc)) rad*=-1;
1240 //____________________________________________________________________
1241 Double_t AliITStrackerSA::ChoosePoint(Double_t p1, Double_t p2, Double_t pp){
1243 //Returns the point closest to pp
1245 Double_t diff1 = p1-pp;
1246 Double_t diff2 = p2-pp;
1248 if(TMath::Abs(diff1)<TMath::Abs(diff2)) fPhiEstimate=p1;
1249 else fPhiEstimate=p2;
1250 return fPhiEstimate;
1255 //_________________________________________________________________
1256 AliITStrackV2* AliITStrackerSA::FindTrackLowChiSquare(TObjArray* tracklist, Int_t dim) const {
1257 // returns track with lowes chi square
1259 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(0);
1262 if(dim==0) return 0;
1263 Double_t * chi2 = new Double_t[dim];
1264 Int_t * index = new Int_t[dim];
1265 for(Int_t i=0;i<dim;i++){
1266 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(i);
1267 chi2[i]=trk->GetChi2();
1271 Int_t w=0;Double_t value;
1274 for(Int_t j=w+1;j<dim;j++){
1275 if(chi2[w]<chi2[j]){
1287 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(index[dim-1]);
1294 //__________________________________________________________
1295 Int_t AliITStrackerSA::FindLabel(Int_t l1, Int_t l2, Int_t l3, Int_t l4, Int_t l5, Int_t l6){
1297 //function used to determine the track label
1299 Int_t lb[6] = {l1,l2,l3,l4,l5,l6};
1300 Int_t aa[6]={1,1,1,1,1,1};
1303 Int_t k=0;Int_t w=0;Int_t num=6;
1304 if(lb[5]==-1) num=5;
1308 for(Int_t i=k+1;i<num;i++){
1310 if(lb[k]==lb[i] && aa[k]!=0){
1321 for(Int_t j=0;j<6;j++){
1333 if(num==6) return lb[5];
1337 //_____________________________________________________________________________
1338 Int_t AliITStrackerSA::Label(Int_t gl1, Int_t gl2, Int_t gl3, Int_t gl4, Int_t gl5, Int_t gl6,Int_t gl7, Int_t gl8, Int_t gl9, Int_t gl10,Int_t gl11,
1339 Int_t gl12, Int_t gl13, Int_t gl14,Int_t gl15, Int_t gl16, Int_t gl17, Int_t gl18, Int_t numberofpoints){
1342 //function used to assign label to the found track. If track is fake, the label is negative
1344 Int_t lb0[6] = {gl1,gl2,gl3,gl4,gl5,gl6};
1345 Int_t lb1[6] = {gl7,gl8,gl9,gl10,gl11,gl12};
1346 Int_t lb2[6] = {gl13,gl14,gl15,gl16,gl17,gl18};
1347 Int_t ll=FindLabel(lb0[0],lb0[1],lb0[2],lb0[3],lb0[4],lb0[5]);
1348 Int_t lflag=0;Int_t num=6;
1349 if(lb0[5]==-1 && lb1[5]==-1 && lb2[5]==-1) num=5;
1351 for(Int_t i=0;i<num;i++){
1352 if(lb0[i]==ll || lb1[i]==ll || lb2[i]==ll) lflag+=1;
1355 if(lflag>=numberofpoints) return ll;
1361 //_____________________________________________________________________________
1362 void AliITStrackerSA::SetWindowSizes(Int_t n, Double_t *phi, Double_t *lam){
1363 // Set sizes of the phi and lambda windows used for track finding
1365 if(phi){ // user defined values
1366 fPhiWin = new Double_t[fNloop];
1367 fLambdaWin = new Double_t[fNloop];
1368 for(Int_t k=0;k<fNloop;k++){
1370 fLambdaWin[k]=lam[k];
1373 else { // default values
1375 Double_t phid[33] = {0.002,0.003,0.004,0.0045,0.0047,
1376 0.005,0.0053,0.0055,
1377 0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077,
1378 0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097,
1379 0.01,0.0105,0.011,0.0115,0.012,0.0125,0.013,0.0135,0.0140,0.0145};
1380 Double_t lambdad[33] = {0.003,0.004,0.005,0.005,0.005,
1382 0.006,0.006,0.006,0.007,0.007,0.007,0.007,
1383 0.007,0.007,0.007,0.007,0.007,0.007,0.007,
1384 0.008,0.008,0.008,0.008,0.008,0.008,0.008,0.008,0.008,0.008};
1391 fPhiWin = new Double_t[fNloop];
1392 fLambdaWin = new Double_t[fNloop];
1394 for(Int_t k=0;k<fNloop;k++){
1396 fLambdaWin[k]=lambdad[k];