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 **************************************************************************/
16 ////////////////////////////////////////////////////
17 // Stand alone tracker class //
18 // Origin: Elisabetta Crescio //
19 // e-mail: crescio@to.infn.it //
20 // tracks are saved as AliITStrackV2 objects //
21 ////////////////////////////////////////////////////
27 #include <TObjArray.h>
30 #include "AliITSclusterTable.h"
31 #include "AliITSclusterV2.h"
32 #include "AliITSgeom.h"
33 #include "AliITSRiemannFit.h"
34 #include "AliITStrackerSA.h"
35 #include "AliITStrackSA.h"
36 #include "AliITSVertexer.h"
37 #include "AliESDVertex.h"
39 #include "AliESDtrack.h"
41 ClassImp(AliITStrackerSA)
43 //____________________________________________________________________________
44 AliITStrackerSA::AliITStrackerSA():AliITStrackerMI(){
45 // Default constructor
49 //____________________________________________________________________________
50 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom):AliITStrackerMI(geom)
52 // Standard constructor (Vertex is known and passed to this obj.)
59 //____________________________________________________________________________
60 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliESDVertex *vert):AliITStrackerMI(geom)
62 // Standard constructor (Vertex is known and passed to this obj.)
69 //______________________________________________________________________
70 AliITStrackerSA::AliITStrackerSA(const AliITStrackerSA &trkr) :
71 AliITStrackerMI(trkr) {
73 // Copies are not allowed. The method is protected to avoid misuse.
74 Error("AliITStrackerSA","Copy constructor not allowed\n");
77 //______________________________________________________________________
78 AliITStrackerSA& AliITStrackerSA::operator=(const
79 AliITStrackerSA& /* trkr */){
80 // Assignment operator
81 // Assignment is not allowed. The method is protected to avoid misuse.
82 Error("= operator","Assignment operator not allowed\n");
86 //____________________________________________________________________________
87 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliITSVertexer *vertexer):AliITStrackerMI(geom)
89 // Standard constructor (Vertex is unknown - vertexer is passed to this obj)
96 //____________________________________________________________________________
97 AliITStrackerSA::AliITStrackerSA(AliITStrackerSA& tracker):AliITStrackerMI(){
99 fPhiEstimate = tracker.fPhiEstimate;
100 for(Int_t i=0;i<2;i++){
101 fPoint1[i]=tracker.fPoint1[i];
102 fPoint2[i]=tracker.fPoint2[i];
103 fPoint3[i]=tracker.fPoint3[i];
104 fPointc[i]=tracker.fPointc[i];
106 fLambdac = tracker.fLambdac;
107 fPhic = tracker.fPhic;
108 fCoef1 = tracker.fCoef1;
109 fCoef2 = tracker.fCoef2;
110 fCoef3 = tracker.fCoef3;
111 fNloop = tracker.fNloop;
112 fPhiWin = tracker.fPhiWin;
113 fLambdaWin = tracker.fLambdaWin;
114 if(tracker.fVertexer && tracker.fVert){
115 fVert = new AliESDVertex(*tracker.fVert);
118 fVert = tracker.fVert;
120 fVertexer = tracker.fVertexer;
121 fGeom = tracker.fGeom;
122 fTable = tracker.fTable;
123 fListOfTracks = tracker.fListOfTracks;
126 //____________________________________________________________________________
127 AliITStrackerSA::~AliITStrackerSA(){
129 // if fVertexer is not null, the AliESDVertex obj. is owned by this class
130 // and is deleted here
132 if(fVert)delete fVert;
137 if(fPhiWin)delete []fPhiWin;
138 if(fLambdaWin)delete []fLambdaWin;
140 fListOfTracks->Delete();
143 //____________________________________________________________________________
144 void AliITStrackerSA::Init(){
145 // Reset all data members
147 for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
162 fListOfTracks=new TObjArray(0,0);
164 //_______________________________________________________________________
165 void AliITStrackerSA::ResetForFinding(){
166 // Reset data members used in all loops during track finding
168 for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
176 fListOfTracks->Delete();
178 //____________________________________________________________________________
179 void AliITStrackerSA::FindTracks(TTree *out,Int_t evnumber){
181 /**************************************************************************
182 * This function finds primary tracks.
185 * Example: to execute function with only the ITS (no combined tracking *
186 * with TPC+ITS) and requiring 5/6 points to define a good track *
187 * call SetSixPoinbts(kFALSE) in advance and then *
188 * use: FindTracks(treein,treeout,evnumber) *
189 * to execute combined tracking, before using FindTracks, use *
191 *************************************************************************/
194 Fatal("FindTracks","ITS cluster tree is not accessed - Abort!!!\n Please use method SetClusterTree to pass the pointer to the tree\n");
199 if(fVert)delete fVert;
200 fVert = fVertexer->FindVertexForCurrentEvent(evnumber);
203 gAlice->GetEvent(evnumber);
205 Fatal("FindTracks","Vertex is missing\n");
209 Double_t primaryVertex[3];
210 Double_t errorsprimvert[3];
211 fVert->GetXYZ(primaryVertex);
212 fVert->GetSigmaXYZ(errorsprimvert);
213 if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
214 Warning("FindTracks","Set errors on vertex positions x and y at 0.0001");
215 errorsprimvert[0]=0.0001;
216 errorsprimvert[1]=0.0001;
218 fVert->PrintStatus();
221 //Fill array with cluster indices for each module
223 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
224 fTable->FillArray(fITSclusters);
225 fTable->FillArrayCoorAngles();
229 //Fill tree for found tracks
230 AliITStrackV2* outrack=0;
231 TBranch* branch=out->Branch("tracks","AliITStrackV2",&outrack,32000,0);
232 if (!branch) out->Branch("tracks","AliITStrackV2",&outrack,32000,3);
233 else branch->SetAddress(&outrack);
236 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
237 for(Int_t i=0;i<fGeom->GetNlayers();i++){
238 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
240 // firstmod [i] number of the first module in the ITS layer i.
242 AliITSlayer &layer=fgLayers[0]; // first layer
244 Int_t dim=layer.GetNumberOfClusters();
245 //loop on the different windows
246 for(Int_t nloop=0;nloop<fNloop;nloop++){
247 for(Int_t ncl=0;ncl<dim;ncl++){
248 //loop starting from layer 0
251 AliITSclusterV2* cl = layer.GetCluster(ncl);
252 if(cl->IsUsed()==1) continue;
253 if(cl->TestBit(kSAflag)==kTRUE) continue;
255 fPhic = fTable->GetPhiCluster(0,ncl);
256 fLambdac = fTable->GetLambdaCluster(0,ncl);
257 fPhiEstimate = fPhic;
258 AliITStrackSA* trs = new AliITStrackSA();
259 fPoint1[0]=primaryVertex[0];
260 fPoint1[1]=primaryVertex[1];
261 fPoint2[0]=fTable->GetXCluster(0,ncl);
262 fPoint2[1]=fTable->GetYCluster(0,ncl);
264 Int_t * nn = new Int_t[fGeom->GetNlayers()];//counter for clusters on each layer
265 for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;}
266 nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
267 nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
271 fPoint3[0] = fPointc[0];
272 fPoint3[1] = fPointc[1];
274 nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
275 if(nn[1]==0 && nn[2]==0) pflag=0;
276 if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();}
277 if(nn[2]!=0 && nn[1]==0){
279 fPoint3[0]=fPointc[0];
280 fPoint3[1]=fPointc[1];
283 nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
285 if(nn[3]!=0) UpdatePoints();
286 nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
288 if(nn[4]!=0) UpdatePoints();
289 nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
293 Int_t numberofpoints;
294 if(fSixPoints) numberofpoints=6; //check of the candidate track
295 else numberofpoints=5; //if track is good (with the required number
296 for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){ //of points) it is written on file
297 if(nn[nnp]!=0) layOK+=1;
299 if(layOK>=numberofpoints){
300 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
302 Int_t nct = trs->GetNumberOfClustersSA();
304 Int_t index = trs->GetClusterIndexSA(nct);
305 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
306 if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
315 Int_t nct = tr2->GetNumberOfClusters();
318 Int_t index = tr2->GetClusterIndex(nct);
319 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
325 Int_t nct = trs->GetNumberOfClustersSA();
327 Int_t index = trs->GetClusterIndexSA(nct);
328 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
329 if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
335 }//end loop on clusters of layer1
339 //if 5/6 points are required, second loop starting
340 //from second layer, to find tracks with point of
344 // counter for clusters on each layer
345 Int_t * nn = new Int_t[fGeom->GetNlayers()-1];
346 for(Int_t nloop=0;nloop<fNloop;nloop++){
347 AliITSlayer &layer2=fgLayers[1]; //loop on layer 2
348 Int_t ncl2=layer2.GetNumberOfClusters();
349 while(ncl2--){ //loop starting from layer 2
352 AliITSclusterV2* cl = layer2.GetCluster(ncl2);
353 if(cl->IsUsed()==1) continue;
354 if(cl->TestBit(kSAflag)==kTRUE) continue;
355 fPhic = fTable->GetPhiCluster(1,ncl2);
356 fLambdac = fTable->GetLambdaCluster(1,ncl2);
357 fPhiEstimate = fPhic;
358 AliITStrackSA* trs = new AliITStrackSA();
360 fPoint1[0]=primaryVertex[0];
361 fPoint1[1]=primaryVertex[1];
362 fPoint2[0]=fTable->GetXCluster(1,ncl2);;
363 fPoint2[1]=fTable->GetYCluster(1,ncl2);;
365 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
366 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
367 nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
368 trs,primaryVertex[2],pflag);
376 for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
377 if(nn[nnp]!=0) fl+=1;
380 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
382 Int_t nct = trs->GetNumberOfClustersSA();
384 Int_t index = trs->GetClusterIndexSA(nct);
385 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
386 if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
393 Int_t nct = tr2->GetNumberOfClusters();
395 Int_t index = tr2->GetClusterIndex(nct);
396 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
402 Int_t nct = trs->GetNumberOfClustersSA();
404 Int_t index = trs->GetClusterIndexSA(nct);
405 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
407 if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
411 }//end loop on clusters of layer2
414 } // if(!fSixPoints....
417 delete fTable; fTable=0;
421 //______________________________________________________________________
422 Int_t AliITStrackerSA::FindTracks(AliESD* event){
424 // Track finder using the ESD object
428 Fatal("FindTracks","ITS cluster tree is not accessed - Abort!!!\n Please use method SetClusterTree to pass the pointer to the tree\n");
433 Double_t errorsprimvert[3];
434 Double_t primaryVertex[3];
435 event->GetVertex()->GetXYZ(primaryVertex);
436 event->GetVertex()->GetSigmaXYZ(errorsprimvert);
438 if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
439 // Warning("FindTracks","Set errors on vertex positions x and y at 0.005");
440 errorsprimvert[0]=0.005;
441 errorsprimvert[1]=0.005;
444 //Fill array with cluster indices for each module
446 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
447 fTable->FillArray(fITSclusters);
448 fTable->FillArrayCoorAngles();
451 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
452 for(Int_t i=0;i<fGeom->GetNlayers();i++){
453 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
455 // firstmod [i] number of the first module in the ITS layer i.
458 AliITSlayer &layer=fgLayers[0];
460 Int_t dim=layer.GetNumberOfClusters();
461 //loop on the different windows
462 for(Int_t nloop=0;nloop<fNloop;nloop++){
463 for(Int_t ncl=0;ncl<dim;ncl++){ //loop starting from layer 0
467 AliITSclusterV2* cl = (AliITSclusterV2*)layer.GetCluster(ncl);
469 if(cl->IsUsed()==1) continue;
470 if(cl->TestBit(kSAflag)==kTRUE) continue;
471 if (cl->GetQ()<=0) continue;
473 fPhic = fTable->GetPhiCluster(0,ncl);
474 fLambdac = fTable->GetLambdaCluster(0,ncl);
476 if (TMath::Abs(fLambdac)>0.26*TMath::Pi()) continue;
478 fPhiEstimate = fPhic;
479 AliITStrackSA* trs = new AliITStrackSA();
480 fPoint1[0]=primaryVertex[0];
481 fPoint1[1]=primaryVertex[1];
482 fPoint2[0]=fTable->GetXCluster(0,ncl);
483 fPoint2[1]=fTable->GetYCluster(0,ncl);
484 Int_t * nn = new Int_t[fGeom->GetNlayers()];//counter for clusters on each layer
485 for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;}
486 nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
488 nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
491 fPoint3[0] = fPointc[0];
492 fPoint3[1] = fPointc[1];
494 nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
495 if(nn[1]==0 && nn[2]==0) pflag=0;
496 if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();}
497 if(nn[2]!=0 && nn[1]==0){
499 fPoint3[0]=fPointc[0];
500 fPoint3[1]=fPointc[1];
503 nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
505 if(nn[3]!=0) UpdatePoints();
506 nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
508 if(nn[4]!=0) UpdatePoints();
509 nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
513 Int_t numberofpoints;
514 if(fSixPoints) numberofpoints=6; //check of the candidate track
515 else numberofpoints=5; //if track is good (with the required number
516 for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){ //of points) it is written on file
517 if(nn[nnp]!=0) layOK+=1;
519 if(layOK>=numberofpoints){
520 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
522 Int_t nct = trs->GetNumberOfClustersSA();
524 Int_t index = trs->GetClusterIndexSA(nct);
525 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
526 if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
532 AliESDtrack outtrack;
533 outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
534 event->AddTrack(&outtrack);
536 Int_t nct = tr2->GetNumberOfClusters();
538 Int_t index = tr2->GetClusterIndex(nct);
539 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
545 Int_t nct = trs->GetNumberOfClustersSA();
547 Int_t index = trs->GetClusterIndexSA(nct);
548 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
549 if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
556 }//end loop on clusters of layer1
562 //if 5/6 points are required, second loop starting
563 //from second layer, to find tracks with point of
567 // counter for clusters on each layer
568 Int_t * nn = new Int_t[fGeom->GetNlayers()-1];
569 for(Int_t nloop=0;nloop<fNloop;nloop++){
570 AliITSlayer &layer2=fgLayers[1];
571 Int_t ncl2=layer2.GetNumberOfClusters();
572 while(ncl2--){ //loop starting from layer 2
575 AliITSclusterV2* cl = layer2.GetCluster(ncl2);
576 if(cl->IsUsed()==1) continue;
577 if(cl->TestBit(kSAflag)==kTRUE) continue;
578 fPhic = fTable->GetPhiCluster(1,ncl2);
579 fLambdac = fTable->GetLambdaCluster(1,ncl2);
580 fPhiEstimate = fPhic;
581 AliITStrackSA* trs = new AliITStrackSA();
582 fPoint1[0]=primaryVertex[0];
583 fPoint1[1]=primaryVertex[1];
584 fPoint2[0]=fTable->GetXCluster(1,ncl2);
585 fPoint2[1]=fTable->GetYCluster(1,ncl2);
587 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
588 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
589 nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
590 trs,primaryVertex[2],pflag);
598 for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
599 if(nn[nnp]!=0) fl+=1;
602 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
604 Int_t nct = trs->GetNumberOfClustersSA();
606 Int_t index = trs->GetClusterIndexSA(nct);
607 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
608 if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
615 AliESDtrack outtrack;
616 outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
617 event->AddTrack(&outtrack);
619 Int_t nct = tr2->GetNumberOfClusters();
621 Int_t index = tr2->GetClusterIndex(nct);
622 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
628 Int_t nct = trs->GetNumberOfClustersSA();
630 Int_t index = trs->GetClusterIndexSA(nct);
631 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
633 if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
638 }//end loop on clusters of layer2
644 delete fTable;fTable=0;
645 Info("FindTracks","Number of found tracks: %d",event->GetNumberOfTracks());
651 //________________________________________________________________________
653 AliITStrackV2* AliITStrackerSA::FitTrack(AliITStrackSA* tr,Double_t *primaryVertex,Double_t *errorsprimvert){
654 //fit of the found track
657 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
658 for(Int_t i=0;i<fGeom->GetNlayers();i++){
659 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
662 Int_t nclusters = tr->GetNumberOfClustersSA();
663 TObjArray** listlayer = new TObjArray*[fGeom->GetNlayers()];
664 for(Int_t i=0;i<fGeom->GetNlayers();i++){
665 listlayer[i] = new TObjArray(0,0);
675 Int_t * nnn = new Int_t[fGeom->GetNlayers()];
676 for(Int_t i=0;i<fGeom->GetNlayers();i++)nnn[i]=0;
678 for(Int_t ncl=0;ncl<nclusters;ncl++){
679 Int_t index = tr->GetClusterIndexSA(ncl);
680 AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(index);
682 if(cl->TestBit(kSAflag)==kTRUE) cl->ResetBit(kSAflag);
683 Int_t lay = (index & 0xf0000000) >> 28;
684 if(lay==0) { listlayer[0]->AddLast(cl); clind0[nnn[0]]=index;nnn[0]++;}
685 if(lay==1) { listlayer[1]->AddLast(cl); clind1[nnn[1]]=index;nnn[1]++;}
686 if(lay==2) { listlayer[2]->AddLast(cl); clind2[nnn[2]]=index;nnn[2]++;}
687 if(lay==3) { listlayer[3]->AddLast(cl); clind3[nnn[3]]=index;nnn[3]++;}
688 if(lay==4) { listlayer[4]->AddLast(cl); clind4[nnn[4]]=index;nnn[4]++;}
689 if(lay==5) { listlayer[5]->AddLast(cl); clind5[nnn[5]]=index;nnn[5]++;}
694 Int_t * end = new Int_t[fGeom->GetNlayers()];
695 for(Int_t i=0;i<fGeom->GetNlayers();i++){
696 if(listlayer[i]->GetEntries()==0) end[i]=1;
697 else end[i]=listlayer[i]->GetEntries();
700 for(Int_t l1=0;l1<end[0];l1++){//loop on layer 1
701 AliITSclusterV2* cl0 = (AliITSclusterV2*)listlayer[0]->At(l1);
702 TVector3** recp = new TVector3*[3];
703 TVector3** errs = new TVector3*[3];
704 recp[0] = new TVector3(primaryVertex[0],primaryVertex[1],primaryVertex[2]);
705 errs[0] = new TVector3(errorsprimvert[0],errorsprimvert[1],errorsprimvert[2]);
706 Double_t x1,y1,z1,sx1,sy1,sz1;
707 Double_t x2,y2,z2,sx2,sy2,sz2;
708 AliITSclusterV2* p1=0;
709 AliITSclusterV2* p2=0;
710 Int_t index1=clind0[l1];
712 for(Int_t l2=0;l2<end[1];l2++){//loop on layer 2
713 AliITSclusterV2* cl1 = (AliITSclusterV2*)listlayer[1]->At(l2);
715 for(Int_t l3=0;l3<end[2];l3++){ //loop on layer 3
716 AliITSclusterV2* cl2 = (AliITSclusterV2*)listlayer[2]->At(l3);
718 if(cl0==0 && cl1!=0) {
719 p2 = cl2;index1=clind2[l3];
723 if(cl0!=0 && cl1==0){
725 p2=cl2;index2=clind2[l3];
727 if(cl0!=0 && cl1!=0){
731 Int_t lay1=(index1 & 0xf0000000) >> 28;
732 Int_t cln1=(index1 & 0x0fffffff) >> 00;
733 Int_t lay2=(index2 & 0xf0000000) >> 28;
734 Int_t cln2=(index2 & 0x0fffffff) >> 00;
735 x1 = fTable->GetXCluster(lay1,cln1);
736 x2 = fTable->GetXCluster(lay2,cln2);
737 y1 = fTable->GetYCluster(lay1,cln1);
738 y2 = fTable->GetYCluster(lay2,cln2);
739 z1 = fTable->GetZCluster(lay1,cln1);
740 z2 = fTable->GetZCluster(lay2,cln2);
741 sx1 = fTable->GetXClusterError(lay1,cln1);
742 sx2 = fTable->GetXClusterError(lay2,cln2);
743 sy1 = fTable->GetYClusterError(lay1,cln1);
744 sy2 = fTable->GetYClusterError(lay2,cln2);
745 sz1 = fTable->GetZClusterError(lay1,cln1);
746 sz2 = fTable->GetZClusterError(lay2,cln2);
747 Double_t phi1 = fTable->GetPhiCluster(lay1,cln1);
748 Int_t module1 = p1->GetDetectorIndex()+firstmod[0];
749 recp[1] = new TVector3(x1,y1,z1);
750 errs[1] = new TVector3(sx1,sy1,sz1);
751 recp[2] = new TVector3(x2,y2,z2);
752 errs[2] = new TVector3(sx2,sy2,sz2);
754 //fit on the Riemann sphere
755 Float_t seed1,seed2,seed3;
756 AliITSRiemannFit fit;
757 Int_t rf = fit.FitHelix(3,recp,errs,seed1,seed2,seed3); //this gives phi,tgl,curvature to start Kalman Filter
759 for(Int_t i=1;i<3;i++){
769 if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
770 phi=seed1+1.5*TMath::Pi();
773 if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
774 phi=seed1+0.5*TMath::Pi();
777 if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
778 phi=seed1-0.5*TMath::Pi();
783 if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
784 phi=seed1+0.5*TMath::Pi();
787 if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
788 phi=seed1-0.5*TMath::Pi();
791 if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
792 phi=seed1-1.5*TMath::Pi();
797 Int_t layer,ladder,detector;
798 fGeom->GetModuleId(module1,layer,ladder,detector);
799 Float_t yclu1 = p1->GetY();
800 Float_t zclu1 = p1->GetZ();
801 Double_t cv=Curvature(primaryVertex[0],primaryVertex[1],x1,y1,x2,y2);
803 for(Int_t l4=0;l4<end[3];l4++){ //loop on layer 4
804 AliITSclusterV2* cl3 = (AliITSclusterV2*)listlayer[3]->At(l4);
805 for(Int_t l5=0;l5<end[4];l5++){ //loop on layer 5
806 AliITSclusterV2* cl4 = (AliITSclusterV2*)listlayer[4]->At(l5);
807 for(Int_t l6=0;l6<end[5];l6++){ //loop on layer 6
808 AliITSclusterV2* cl5 = (AliITSclusterV2*)listlayer[5]->At(l6);
809 AliITStrackSA* trac = new AliITStrackSA(layer,ladder,detector,yclu1,zclu1,phi,tgl,cv,1);
811 if(cl5!=0) trac->AddClusterV2(5,(clind5[l6] & 0x0fffffff)>>0);
812 if(cl4!=0) trac->AddClusterV2(4,(clind4[l5] & 0x0fffffff)>>0);
813 if(cl3!=0) trac->AddClusterV2(3,(clind3[l4] & 0x0fffffff)>>0);
814 if(cl2!=0) trac->AddClusterV2(2,(clind2[l3] & 0x0fffffff)>>0);
815 if(cl1!=0) trac->AddClusterV2(1,(clind1[l2] & 0x0fffffff)>>0);
816 if(cl0!=0) trac->AddClusterV2(0,(clind0[l1] & 0x0fffffff)>>0);
818 //fit with Kalman filter using AliITStrackerMI::RefitAt()
820 AliITStrackMI* ot = new AliITStrackSA(*trac);
822 ot->ResetCovariance();
825 if(RefitAt(49.,ot,trac)){ //fit from layer 1 to layer 6
827 AliITStrackMI *otrack2 = new AliITStrackMI(*ot);
828 otrack2->ResetCovariance();
829 otrack2->ResetClusters();
830 //fit from layer 6 to layer 1
831 if(RefitAt(3.7,otrack2,ot)) {
832 fListOfTracks->AddLast(otrack2);
845 for(Int_t i=1;i<3;i++){
859 Int_t dim=fListOfTracks->GetEntries();
861 for(Int_t i=0;i<fGeom->GetNlayers();i++){
869 AliITStrackV2* otrack =(AliITStrackV2*)FindTrackLowChiSquare(fListOfTracks,dim);
872 for(Int_t i=0;i<fGeom->GetNlayers();i++){
879 Int_t * indexc = new Int_t[fGeom->GetNlayers()];
880 for(Int_t i=0;i<fGeom->GetNlayers();i++) indexc[i]=0;
881 for(Int_t nind=0;nind<otrack->GetNumberOfClusters();nind++){
882 indexc[nind] = otrack->GetClusterIndex(nind);
884 AliITSclusterV2* cl0 = (AliITSclusterV2*)GetCluster(indexc[0]);
885 AliITSclusterV2* cl1 = (AliITSclusterV2*)GetCluster(indexc[1]);
886 AliITSclusterV2* cl2 = (AliITSclusterV2*)GetCluster(indexc[2]);
887 AliITSclusterV2* cl3 = (AliITSclusterV2*)GetCluster(indexc[3]);
888 AliITSclusterV2* cl4 = (AliITSclusterV2*)GetCluster(indexc[4]);
889 Int_t labl[3]={-1,-1,-1};
890 if(otrack->GetNumberOfClusters()==fGeom->GetNlayers()){
891 AliITSclusterV2* cl5 = (AliITSclusterV2*)GetCluster(indexc[5]);
892 labl[0]=cl5->GetLabel(0);
893 labl[1]=cl5->GetLabel(1);
894 labl[2]=cl5->GetLabel(2);
897 if(otrack->GetNumberOfClusters()==(fGeom->GetNlayers()-1)){
902 Int_t numberofpoints;
903 if(fSixPoints) numberofpoints=6;
904 else numberofpoints=5;
905 CookLabel(otrack,0.); //MI change - to see fake ratio
906 Int_t label = Label(cl0->GetLabel(0),cl1->GetLabel(0),
907 cl2->GetLabel(0),cl3->GetLabel(0),
908 cl4->GetLabel(0),labl[0],
909 cl0->GetLabel(1),cl1->GetLabel(1),
910 cl2->GetLabel(1),cl3->GetLabel(1),
911 cl4->GetLabel(1),labl[1],
912 cl0->GetLabel(2),cl1->GetLabel(2),
913 cl2->GetLabel(2),cl3->GetLabel(2),
914 cl4->GetLabel(2),labl[2],numberofpoints);
916 otrack->SetLabel(label);
917 for(Int_t i=0;i<fGeom->GetNlayers();i++){
926 //_______________________________________________________________________
927 void AliITStrackerSA::UseFoundTracksV2(Int_t evnum,TTree* treev2){
928 // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
929 //(or AliITStrackV2 tracks found with function FindTracks of this class)
934 if(fVert)delete fVert;
935 fVert = fVertexer->FindVertexForCurrentEvent(evnum);
938 gAlice->GetEvent(evnum);
940 Fatal("FindTracks","Vertex is missing\n");
944 Double_t primaryVertex[3];
945 fVert->GetXYZ(primaryVertex);
948 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
949 fTable->FillArray(fITSclusters);
950 fTable->FillArrayCoorAngles();
953 TBranch* bra = (TBranch*)treev2->GetBranch("tracks");
954 if(!bra) Warning("UseFoundTracksV2","No branch for track tree");
955 AliITStrackV2* ttrrt = new AliITStrackV2;
956 bra->SetAddress(&ttrrt);
958 for(Int_t nj=0;nj<treev2->GetEntries();nj++){
959 treev2->GetEvent(nj);
960 Int_t ncl = ttrrt->GetNumberOfClusters();
961 for(Int_t k=0;k<ncl;k++){
962 Int_t index = ttrrt->GetClusterIndex(k);
963 AliITSclusterV2* clui = (AliITSclusterV2*)GetCluster(index);
964 if(clui->IsUsed()==0) clui->Use();
972 //_______________________________________________________________________
973 void AliITStrackerSA::UseFoundTracksV2(AliESD *event){
974 // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
978 Double_t primaryVertex[3];
979 event->GetVertex()->GetXYZ(primaryVertex);
982 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
983 fTable->FillArray(fITSclusters);
984 fTable->FillArrayCoorAngles();
987 Int_t ntracks = event->GetNumberOfTracks();
989 AliESDtrack *esd=event->GetTrack(ntracks);
990 if ((esd->GetStatus()&
991 AliESDtrack::kITSin|AliESDtrack::kTPCin)==0) continue;
993 Int_t ncl = esd->GetITSclusters(idx);
994 for(Int_t clu=0; clu<ncl; clu++){
995 AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(idx[clu]);
996 if(cl->IsUsed()==0) cl->Use();
1000 Info("UseFoundTracksV2","Clusters of tracks prolonged from TPC deleted");
1006 //_______________________________________________________
1007 Int_t AliITStrackerSA::SearchClusters(Int_t layer,Double_t phiwindow,Double_t lambdawindow, AliITStrackSA* trs,Double_t zvertex,Int_t pflag){
1008 //function used to to find the clusters associated to the track
1010 AliITSlayer &lay = fgLayers[layer];
1011 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
1012 for(Int_t i=0;i<fGeom->GetNlayers();i++){
1013 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
1017 Float_t cx1,cx2,cy1,cy2;
1018 FindEquation(fPoint1[0],fPoint1[1],fPoint2[0],fPoint2[1],fPoint3[0],fPoint3[1],fCoef1,fCoef2,fCoef3);
1019 Int_t fun = FindIntersection(fCoef1,fCoef2,fCoef3,-(lay.GetR()*lay.GetR()),cx1,cy1,cx2,cy2);
1025 Double_t fi1 =TMath::ATan2(cy1,cx1);
1026 Double_t fi2 =TMath::ATan2(cy2,cx2);
1027 fPhiEstimate = ChoosePoint(fi1,fi2,fPhic);
1030 Double_t zed = TMath::Tan(fLambdac)*lay.GetR()+zvertex;
1031 Double_t zed1 = TMath::Tan(fLambdac+lambdawindow)*lay.GetR()+zvertex;
1032 Double_t zed2 = TMath::Tan(fLambdac-lambdawindow)*lay.GetR()+zvertex;
1034 Double_t fi = fPhiEstimate;
1035 Int_t nmod = lay.FindDetectorIndex(fi,zed);
1040 nmod += firstmod[layer];
1042 Int_t nm[8]={0,0,0,0,0,0,0,0};
1043 nm[0] = lay.FindDetectorIndex(fi+phiwindow,zed);
1044 nm[1] = lay.FindDetectorIndex(fi-phiwindow,zed);
1045 nm[2] = lay.FindDetectorIndex(fi,zed1);
1046 nm[3] = lay.FindDetectorIndex(fi,zed2);
1047 nm[4] = lay.FindDetectorIndex(fi+phiwindow,zed1);
1048 nm[5] = lay.FindDetectorIndex(fi-phiwindow,zed1);
1049 nm[6] = lay.FindDetectorIndex(fi+phiwindow,zed2);
1050 nm[7] = lay.FindDetectorIndex(fi-phiwindow,zed2);
1054 TArrayI* array =(TArrayI*)fTable->GetListOfClusters(nmod);
1055 TArrayI* listc = new TArrayI(array->GetSize());
1056 for(Int_t i=0;i<array->GetSize();i++){
1057 Int_t in=(Int_t)array->At(i);
1058 listc->AddAt(in,nn);
1065 for(Int_t h=k+1;h<8;h++){
1078 for(Int_t ii=0;ii<8;ii++){
1079 if(nm[ii]!=value && nm[ii]!=nmod && nm[ii]>=0){
1080 TArrayI* ar =(TArrayI*)fTable->GetListOfClusters(nm[ii]+firstmod[layer]);
1081 listc->Set(listc->GetSize()+ar->GetSize());
1082 for(Int_t j=0;j<ar->GetSize();j++){
1083 Int_t in=(Int_t)ar->At(j);
1084 listc->AddAt(in,nn);
1092 for(Int_t i=0;i<listc->GetSize();i++){
1093 Int_t index = (Int_t)listc->At(i);
1094 AliITSclusterV2* cllay = lay.GetCluster(index);
1095 if(cllay==0) continue;
1096 if(cllay->IsUsed()==1) continue;
1097 if(cllay->TestBit(kSAflag)==kTRUE) continue;
1098 Double_t phi = fTable->GetPhiCluster(layer,index);
1099 Double_t lambda= fTable->GetLambdaCluster(layer,index);
1101 if(TMath::Abs(fLambdac-lambda)<lambdawindow &&
1102 TMath::Abs(fPhiEstimate-phi)<phiwindow){
1105 if(trs->GetNumberOfClustersSA()==15){
1110 trs->AddClusterSA(layer,index);
1111 cllay->SetBit(kSAflag);
1113 fPointc[0]=fTable->GetXCluster(layer,index);
1114 fPointc[1]=fTable->GetYCluster(layer,index);
1126 //_______________________________________________________
1127 Int_t AliITStrackerSA::SearchClusters(Int_t layer,Double_t phiwindow,Double_t lambdawindow, AliITStrackSA* trs,Double_t zvertex,Int_t pflag){
1128 //function used to to find the clusters associated to the track
1130 AliITSlayer &lay = fgLayers[layer];
1131 Double_t r=lay.GetR(),tgl=TMath::Tan(fLambdac);
1134 Float_t cx1,cx2,cy1,cy2;
1135 FindEquation(fPoint1[0],fPoint1[1],fPoint2[0],fPoint2[1],fPoint3[0],fPoint3[1],fCoef1,fCoef2,fCoef3);
1136 if (FindIntersection(fCoef1,fCoef2,fCoef3,-r*r,cx1,cy1,cx2,cy2)==0)
1138 Double_t fi1=TMath::ATan2(cy1,cx1);
1139 Double_t fi2=TMath::ATan2(cy2,cx2);
1140 fPhiEstimate=ChoosePoint(fi1,fi2,fPhic);
1143 Double_t dz=r*lambdawindow*TMath::Sqrt(1+tgl*tgl) + 0.3*TMath::Abs(tgl);
1144 Double_t zmax=r*tgl + zvertex + dz;
1145 Double_t zmin=r*tgl + zvertex - dz;
1147 Int_t indmin = lay.FindClusterIndex(zmin);
1148 Int_t indmax = lay.FindClusterIndex(zmax);
1149 for (Int_t index=indmin; index<indmax; index++) {
1150 AliITSclusterV2 *c=lay.GetCluster(index);
1151 if (c->IsUsed()) continue;
1152 if (c->GetQ()<=0) continue;
1153 if (c->TestBit(kSAflag)==kTRUE) continue;
1155 Double_t phi =fTable->GetPhiCluster(layer,index);
1156 if (TMath::Abs(phi-fPhiEstimate)>phiwindow) continue;
1158 Double_t lambda=fTable->GetLambdaCluster(layer,index);
1159 if (TMath::Abs(lambda-fLambdac)>lambdawindow) continue;
1161 if(trs->GetNumberOfClustersSA()==15) return 0;
1163 trs->AddClusterSA(layer,index);
1167 fPointc[0]=fTable->GetXCluster(layer,index);
1168 fPointc[1]=fTable->GetYCluster(layer,index);
1176 //________________________________________________________________
1177 void AliITStrackerSA::UpdatePoints(){
1178 //update of points for the estimation of the curvature
1180 //fPoint1[0]=fPoint2[0];
1181 //fPoint1[1]=fPoint2[1];
1182 fPoint2[0]=fPoint3[0];
1183 fPoint2[1]=fPoint3[1];
1184 fPoint3[0]=fPointc[0];
1185 fPoint3[1]=fPointc[1];
1190 //___________________________________________________________________
1191 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){
1193 //given (x,y) of three recpoints (in global coordinates)
1194 //returns the parameters a,b,c of circonference x*x + y*y +a*x + b*y +c
1196 Float_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1197 if(den==0) return 0;
1198 a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1199 b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1200 c = -x1*x1-y1*y1-a*x1-b*y1;
1203 //__________________________________________________________________________
1204 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){
1206 //Finds the intersection between the circonference of the track and the circonference centered in (0,0) represented by one layer
1207 //c2 is -rlayer*rlayer
1211 Double_t aA = (b1*b1)/(a1*a1)+1;
1212 Double_t bB = (-2*m*b1/(a1*a1));
1213 Double_t cC = c2+(m*m)/(a1*a1);
1214 Double_t dD = bB*bB-4*aA*cC;
1217 y1 = (-bB+TMath::Sqrt(dD))/(2*aA);
1218 y2 = (-bB-TMath::Sqrt(dD))/(2*aA);
1219 x1 = (c2-c1-b1*y1)/a1;
1220 x2 = (c2-c1-b1*y2)/a1;
1224 //____________________________________________________________________
1225 Double_t AliITStrackerSA::Curvature(Double_t x1,Double_t y1,Double_t
1226 x2,Double_t y2,Double_t x3,Double_t y3){
1228 //calculates the curvature of track
1229 Double_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1230 if(den==0) return 0;
1231 Double_t a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1232 Double_t b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1233 Double_t c = -x1*x1-y1*y1-a*x1-b*y1;
1236 if((a*a+b*b-4*c)<0) return 0;
1237 Double_t rad = TMath::Sqrt(a*a+b*b-4*c)/2.;
1238 if(rad==0) return 0;
1240 if((x1>0 && y1>0 && x1<xc)) rad*=-1;
1241 if((x1<0 && y1>0 && x1<xc)) rad*=-1;
1242 // if((x1<0 && y1<0 && x1<xc)) rad*=-1;
1243 // if((x1>0 && y1<0 && x1<xc)) rad*=-1;
1248 //____________________________________________________________________
1249 Double_t AliITStrackerSA::ChoosePoint(Double_t p1, Double_t p2, Double_t pp){
1251 //Returns the point closest to pp
1253 Double_t diff1 = p1-pp;
1254 Double_t diff2 = p2-pp;
1256 if(TMath::Abs(diff1)<TMath::Abs(diff2)) fPhiEstimate=p1;
1257 else fPhiEstimate=p2;
1258 return fPhiEstimate;
1263 //_________________________________________________________________
1264 AliITStrackV2* AliITStrackerSA::FindTrackLowChiSquare(TObjArray* tracklist, Int_t dim) const {
1265 // returns track with lowes chi square
1267 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(0);
1270 if(dim==0) return 0;
1271 Double_t * chi2 = new Double_t[dim];
1272 Int_t * index = new Int_t[dim];
1273 for(Int_t i=0;i<dim;i++){
1274 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(i);
1275 chi2[i]=trk->GetChi2();
1279 Int_t w=0;Double_t value;
1282 for(Int_t j=w+1;j<dim;j++){
1283 if(chi2[w]<chi2[j]){
1295 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(index[dim-1]);
1302 //__________________________________________________________
1303 Int_t AliITStrackerSA::FindLabel(Int_t l1, Int_t l2, Int_t l3, Int_t l4, Int_t l5, Int_t l6){
1305 //function used to determine the track label
1307 Int_t lb[6] = {l1,l2,l3,l4,l5,l6};
1308 Int_t aa[6]={1,1,1,1,1,1};
1311 Int_t k=0;Int_t w=0;Int_t num=6;
1312 if(lb[5]==-1) num=5;
1316 for(Int_t i=k+1;i<num;i++){
1318 if(lb[k]==lb[i] && aa[k]!=0){
1329 for(Int_t j=0;j<6;j++){
1341 if(num==6) return lb[5];
1345 //_____________________________________________________________________________
1346 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,
1347 Int_t gl12, Int_t gl13, Int_t gl14,Int_t gl15, Int_t gl16, Int_t gl17, Int_t gl18, Int_t numberofpoints){
1350 //function used to assign label to the found track. If track is fake, the label is negative
1352 Int_t lb0[6] = {gl1,gl2,gl3,gl4,gl5,gl6};
1353 Int_t lb1[6] = {gl7,gl8,gl9,gl10,gl11,gl12};
1354 Int_t lb2[6] = {gl13,gl14,gl15,gl16,gl17,gl18};
1355 Int_t ll=FindLabel(lb0[0],lb0[1],lb0[2],lb0[3],lb0[4],lb0[5]);
1356 Int_t lflag=0;Int_t num=6;
1357 if(lb0[5]==-1 && lb1[5]==-1 && lb2[5]==-1) num=5;
1359 for(Int_t i=0;i<num;i++){
1360 if(lb0[i]==ll || lb1[i]==ll || lb2[i]==ll) lflag+=1;
1363 if(lflag>=numberofpoints) return ll;
1369 //_____________________________________________________________________________
1370 void AliITStrackerSA::SetWindowSizes(Int_t n, Double_t *phi, Double_t *lam){
1371 // Set sizes of the phi and lambda windows used for track finding
1373 if(phi){ // user defined values
1374 fPhiWin = new Double_t[fNloop];
1375 fLambdaWin = new Double_t[fNloop];
1376 for(Int_t k=0;k<fNloop;k++){
1378 fLambdaWin[k]=lam[k];
1381 else { // default values
1383 Double_t phid[33] = {0.002,0.003,0.004,0.0045,0.0047,
1384 0.005,0.0053,0.0055,
1385 0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077,
1386 0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097,
1387 0.01,0.0105,0.011,0.0115,0.012,0.0125,0.013,0.0135,0.0140,0.0145};
1388 Double_t lambdad[33] = {0.003,0.004,0.005,0.005,0.005,
1390 0.006,0.006,0.006,0.007,0.007,0.007,0.007,
1391 0.007,0.007,0.007,0.007,0.007,0.007,0.007,
1392 0.008,0.008,0.008,0.008,0.008,0.008,0.008,0.008,0.008,0.008};
1399 fPhiWin = new Double_t[fNloop];
1400 fLambdaWin = new Double_t[fNloop];
1402 for(Int_t k=0;k<fNloop;k++){
1404 fLambdaWin[k]=lambdad[k];