FairRoot/PandaRoot
Public Member Functions | List of all members
PndTrkCTFindTrackInXY2 Class Reference

#include <PndTrkCTFindTrackInXY2.h>

Inheritance diagram for PndTrkCTFindTrackInXY2:

Public Member Functions

 PndTrkCTFindTrackInXY2 ()
 
 ~PndTrkCTFindTrackInXY2 ()
 
void AddMvdHitsToSttTracks (Double_t delta, Double_t highqualitycut, Double_t FiRangeMvdLow, Double_t FiRangeMvdUp, Short_t maxmvdpixelhitsintrack, Short_t maxmvdstriphitsintrack, Short_t nMvdPixelHit, Short_t nMvdStripHit, Double_t Ox, Double_t Oy, Double_t R, Double_t *XMvdPixel, Double_t *XMvdStrip, Double_t *YMvdPixel, Double_t *YMvdStrip, Short_t &nPixelHitsinTrack, Short_t *ListPixelHitsinTrack, Short_t &nStripHitsinTrack, Short_t *ListStripHitsinTrack)
 
Short_t AssociateSciTilHit (Double_t dimensionscitil, Double_t *esse, bool *InclusionListSciTil, Short_t *List, Short_t maxscitilhitsintrack, Short_t nSciTilHits, Double_t Oxx, Double_t Oyy, Double_t posizSciTil[][3], Double_t Rr)
 
void DecideWhichAngularRangeAndCharge (Double_t fiCenter, Double_t fi_low_limit[2], Double_t fi_up_limit[2], Double_t(*info)[7], Short_t *ListHitsinTrack, Short_t nHitsinTrack, Double_t Oxx, Double_t Oyy, Short_t &charge, Double_t &FiRangeMvdLow, Double_t &FiRangeMvdUp, Double_t &Fi_low_limit, Double_t &Fi_up_limit)
 
bool FindTrackInXYProjection (struct FindTrackInXYProjection2_InputOutputData *InOut, int istampa, int IVOLTE)
 
void OrderingUsingConformal (Short_t Charge, Double_t info[][7], Int_t nHits, Double_t oX, Double_t oY, Short_t *ListHits)
 
void OrderingUsingFi (Short_t Charge, Double_t info[][7], Int_t nHits, Double_t oX, Double_t oY, Short_t *ListHits)
 
void OrderingUsingR (Double_t info[][7], Int_t nHits, Short_t *ListHits)
 
Short_t TrkAssociatedParallelHitsToHelix5 (Short_t *auxListHitsinTrack, bool *InclusionListStt, Double_t Fi_low, Double_t Fi_up, Double_t info[][7], Short_t *ListSttParHits, Int_t NhitsParallel, Double_t Oxx, Double_t Oyy, Double_t Rr, Double_t strawradius)
 
Short_t TrkAssociatedParallelHitsToHelix6 (Short_t *auxListHitsinTrack, bool *InclusionListStt, Double_t Fi_low, Double_t Fi_up, Double_t info[][7], Short_t *ListSttParHits, Int_t NhitsParallel, Double_t Oxx, Double_t Oyy, Double_t Rr, Double_t maximum_distance)
 
 ClassDef (PndTrkCTFindTrackInXY2, 1)
 

Detailed Description

Definition at line 108 of file PndTrkCTFindTrackInXY2.h.

Constructor & Destructor Documentation

PndTrkCTFindTrackInXY2::PndTrkCTFindTrackInXY2 ( )
inline

Default constructor

Definition at line 115 of file PndTrkCTFindTrackInXY2.h.

115 {};
PndTrkCTFindTrackInXY2::~PndTrkCTFindTrackInXY2 ( )
inline

Destructor

Definition at line 119 of file PndTrkCTFindTrackInXY2.h.

119 {};

Member Function Documentation

void PndTrkCTFindTrackInXY2::AddMvdHitsToSttTracks ( Double_t  delta,
Double_t  highqualitycut,
Double_t  FiRangeMvdLow,
Double_t  FiRangeMvdUp,
Short_t  maxmvdpixelhitsintrack,
Short_t  maxmvdstriphitsintrack,
Short_t  nMvdPixelHit,
Short_t  nMvdStripHit,
Double_t  Ox,
Double_t  Oy,
Double_t  R,
Double_t XMvdPixel,
Double_t XMvdStrip,
Double_t YMvdPixel,
Double_t YMvdStrip,
Short_t &  nPixelHitsinTrack,
Short_t *  ListPixelHitsinTrack,
Short_t &  nStripHitsinTrack,
Short_t *  ListStripHitsinTrack 
)

Definition at line 25 of file PndTrkCTFindTrackInXY2.cxx.

References angle, atan2(), Double_t, fabs(), PI, and sqrt().

47 {
48  //bool specialcase; //[R.K. 01/2017] unused variable?
49 
50  Short_t //i, //[R.K. 01/2017] unused variable?
51  jmvdhit;
52 
54  dist;
55 
56 // for(i=0; i<nSttTrackCand; i++){
57 
58 
59  // first try attach the Pixels;
60 
61  nMvdPixelHitsinTrack = 0;
62  for( jmvdhit=0; jmvdhit<nMvdPixelHit; jmvdhit++){
63  angle = atan2(YMvdPixel[jmvdhit]-Oy,XMvdPixel[jmvdhit]-Ox);
64  if(angle<0.) angle += 2.*PI;
65  if( angle>FiRangeMvdUp){
66  angle -= 2.*PI;
67  if( angle>FiRangeMvdUp) angle = FiRangeMvdUp;
68 
69  } else if (angle<FiRangeMvdLow){
70  angle += 2.*PI;
71  if (angle<FiRangeMvdLow) angle = FiRangeMvdLow;
72  }
73 
74  if(angle > FiRangeMvdLow && angle < FiRangeMvdUp)
75  {
76  dist=fabs( sqrt( (Ox-XMvdPixel[jmvdhit])* (Ox-XMvdPixel[jmvdhit])
77  +(Oy-YMvdPixel[jmvdhit])*(Oy-YMvdPixel[jmvdhit]))-R);
78  if(dist<delta)
79  {
80  ListMvdPixelHitsinTrack[nMvdPixelHitsinTrack]=jmvdhit;
81  nMvdPixelHitsinTrack++;
82  if( nMvdPixelHitsinTrack == maxmvdpixelhitsintrack ) break;
83  }
84  } // end of if(angle > FiRangeMvdLow)
85  } // end of for( jmvdhit=0; jmvdhit<nMvdPixelHits; jmvdhit++)
86 
87 
88  // then try attach the Strips;
89  nMvdStripHitsinTrack = 0;
90  for( jmvdhit=0; jmvdhit<nMvdStripHit; jmvdhit++){
91  angle = atan2(YMvdStrip[jmvdhit]-Oy,XMvdStrip[jmvdhit]-Ox);
92  if(angle<0.) angle += 2.*PI;
93  if( angle>FiRangeMvdUp){
94  angle -= 2.*PI;
95  if( angle>FiRangeMvdUp) angle = FiRangeMvdUp;
96  } else if (angle<FiRangeMvdLow){
97  angle += 2.*PI;
98  if (angle<FiRangeMvdLow) angle = FiRangeMvdLow;
99  }
100 
101  if(angle > FiRangeMvdLow && angle < FiRangeMvdUp)
102  {
103  dist=fabs( sqrt( (Ox-XMvdStrip[jmvdhit])* (Ox-XMvdStrip[jmvdhit])
104  +(Oy-YMvdStrip[jmvdhit])*(Oy-YMvdStrip[jmvdhit]))-R);
105  if(dist<delta)
106  {
107  ListMvdStripHitsinTrack[nMvdStripHitsinTrack]=jmvdhit;
108  nMvdStripHitsinTrack++;
109  if( nMvdStripHitsinTrack == maxmvdstriphitsintrack ) break;
110  }
111  } // end of if(angle > FiRangeMvdLow)
112  } // end of for( jmvdhit=0; jmvdhit<nMvdStripHits; jmvdhit++)
113 
114 // } // end of for(i=0; i<nSttTrackCand; i++)
115 
116  return;
117 }
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
Double_t
#define PI
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
friend F32vec4 atan2(const F32vec4 &y, const F32vec4 &x)
Definition: P4_F32vec4.h:117
Double_t angle
Double_t R
Definition: checkhelixhit.C:61
Short_t PndTrkCTFindTrackInXY2::AssociateSciTilHit ( Double_t  dimensionscitil,
Double_t esse,
bool *  InclusionListSciTil,
Short_t *  List,
Short_t  maxscitilhitsintrack,
Short_t  nSciTilHits,
Double_t  Oxx,
Double_t  Oyy,
Double_t  posizSciTil[][3],
Double_t  Rr 
)

Definition at line 127 of file PndTrkCTFindTrackInXY2.cxx.

References atan2(), Double_t, PndTrkCTGeometryCalculations::IntersectionSciTil_Circle(), CAMath::Nint(), and PI.

139 {
140 
141  bool intersect;
142 
143  Short_t
144  igoodScit,
145  iScitHit,
146  Nint;
147 
148  Double_t
149  XintersectionList[2],
150  YintersectionList[2];
151 
152 
153  PndTrkCTGeometryCalculations GeomCalculator;
154 
155 
156  igoodScit=0;
157 
158  for(iScitHit=0; iScitHit<nSciTilHits; iScitHit++){
159  if(!InclusionListSciTil[iScitHit]) continue;
160 
161  intersect=GeomCalculator.IntersectionSciTil_Circle(
162 // 2.*dimensionscitil, // dimensionscitil = true width of SciTil tiles;
163  posizSciTil[iScitHit][0],
164  posizSciTil[iScitHit][1],
165  Oxx,
166  Oyy,
167  Rr,
168  &Nint, // output
169  XintersectionList, // output
170  YintersectionList // output
171  );
172 
173  if(intersect){
174  if( igoodScit == maxscitilhitsintrack) break;
175  List[igoodScit] = iScitHit;
176 
177  // calculate S on the lateral face of the Helix.
178  if ( Nint==1){ // the majority of the cases
179  esse[igoodScit] = atan2(YintersectionList[0]-Oyy,
180  XintersectionList[0]-Oxx);
181  } else { // in this case Nint=2 (it should be a very rare case).
182  // do an average of the two positions.
183  esse[igoodScit] = atan2( 0.5*(YintersectionList[0]+
184  YintersectionList[1])-Oyy,
185  0.5*(XintersectionList[0]+
186  XintersectionList[1])-Oxx);
187  } // end of if ( Nint==1)
188  if ( esse[igoodScit]<0.) esse[igoodScit] += 2.*PI;
189  igoodScit++;
190 
191  } // end of if(intersect)
192  } // end of for(iScitHit=0; iScitHit<nScitHits; iScitHit++)
193 
194  return igoodScit;
195 
196 }
Double_t
#define PI
friend F32vec4 atan2(const F32vec4 &y, const F32vec4 &x)
Definition: P4_F32vec4.h:117
bool IntersectionSciTil_Circle(Double_t posizSciTilx, Double_t posizSciTily, Double_t Oxx, Double_t Oyy, Double_t Rr, Short_t *Nintersections, Double_t XintersectionList[2], Double_t YintersectionList[2])
int Nint(float x)
Definition: PndCAMath.h:117
PndTrkCTFindTrackInXY2::ClassDef ( PndTrkCTFindTrackInXY2  ,
 
)
void PndTrkCTFindTrackInXY2::DecideWhichAngularRangeAndCharge ( Double_t  fiCenter,
Double_t  fi_low_limit[2],
Double_t  fi_up_limit[2],
Double_t(*)  info[7],
Short_t *  ListHitsinTrack,
Short_t  nHitsinTrack,
Double_t  Oxx,
Double_t  Oyy,
Short_t &  charge,
Double_t FiRangeMvdLow,
Double_t FiRangeMvdUp,
Double_t Fi_low_limit,
Double_t Fi_up_limit 
)

Definition at line 204 of file PndTrkCTFindTrackInXY2.cxx.

References atan2(), Double_t, fi, sin(), and two_pi.

221  {
222 
223  Short_t
224  ihit,
225  Nleft,
226  Nright
227  ;
228 
229  Double_t
230  fi,
231  sinus;
232 
233 //-------------------
234 
235  Nleft=0;
236  Nright=0;
237  for(ihit=0;ihit<nSttParHitsinTrack;ihit++){
238  // using the SIGN of the cross product [in the Helix reference frame
239  // in the XY projection] (Hit_Straw_Center) CROSS (center of the Laboratory Frame)
240  // to determine if
241  // a hit is on the left side [looking into the beam] or in the right
242  // side of the trajectory in the XY projection; if the sign is positive
243  // then it is on the right side with respect to the segment joining
244  // the (0,0) to the center of the Helix;
245  // it turns out the the sign of such cross product is simply the
246  // sign of sin( fi - fiCenter);
247 
248  // fi of the Straw hit :
249  fi = atan2( info[ ListSttParHitsinTrack[ihit] ][1]- Oyy,
250  info[ ListSttParHitsinTrack[ihit] ][0]- Oxx);
251 
252  // sign of the cross product :
253  sinus = sin(fi - fiCenter);
254  if (sinus > 0. ) Nright ++ ; else Nleft ++;
255 
256  } // end of for(ihit=0;ihit<nSttParHitsinTrack;ihit++)
257 
258  if( Nright > Nleft) {
259  // in this case the particle is travelling counterclockwise --> it is a
260  // negative particle; the Stt detector range is between fi_low_limit[0]
261  // and fi_up_limit[0] in ANY case (even when there is ONLY ONE entrance and exit);
262  // the Mvd range is between (0,0) and Fi_low_limit;
263  Fi_low_limit = fi_low_limit[0];
264  Fi_up_limit = fi_up_limit[0];
265  charge =-1;
266  FiRangeMvdLow = fiCenter; // fiCenter is between 0. and 2PI;
267  FiRangeMvdUp = Fi_low_limit;
268 
269  // put FiRangeMvdUp between 0 and 2PI;
270  FiRangeMvdUp = fmod(FiRangeMvdUp,two_pi);
271  if( FiRangeMvdUp < FiRangeMvdLow) {
272  FiRangeMvdUp += two_pi;
273  if( FiRangeMvdUp < FiRangeMvdLow) FiRangeMvdUp = FiRangeMvdLow;
274  } // end if( FiRangeMvdUp < FiRangeMvdLow)
275 
276  } else { // continuation of if( Nright > Nleft)
277  // in this case the particle is travelling clockwise --> it is a
278  // positive particle; the Mvd range is between and Fi_up_limit and (0,0) ;
279 
280  if( fi_low_limit[1] < -99.) {
281  // there are only 2 intersections;
282  Fi_low_limit = fi_low_limit[0];
283  Fi_up_limit = fi_up_limit[0];
284  } else {
285  // there are 4 intersections;
286  Fi_low_limit = fi_low_limit[1];
287  Fi_up_limit = fi_up_limit[1];
288  }
289  charge =1;
290  FiRangeMvdUp = fiCenter; // fiCenter is between 0. and 2PI;
291  FiRangeMvdLow = Fi_up_limit;
292  // put FiRangeMvdLow between 0 and 2PI;
293  FiRangeMvdLow = fmod(FiRangeMvdLow,two_pi);
294  if( FiRangeMvdUp < FiRangeMvdLow) {
295  FiRangeMvdUp += two_pi;
296  if( FiRangeMvdUp < FiRangeMvdLow) FiRangeMvdUp = FiRangeMvdLow;
297  }
298  } // end if( Nright > Nleft)
299 
300  // calculate the range for possible Mvd hits attached to this track;
301 
302 
303 
304  return;
305  }
friend F32vec4 sin(const F32vec4 &a)
Definition: P4_F32vec4.h:111
#define two_pi
TFile * fi
Double_t
friend F32vec4 atan2(const F32vec4 &y, const F32vec4 &x)
Definition: P4_F32vec4.h:117
bool PndTrkCTFindTrackInXY2::FindTrackInXYProjection ( struct FindTrackInXYProjection2_InputOutputData InOut,
int  istampa,
int  IVOLTE 
)

Definition at line 312 of file PndTrkCTFindTrackInXY2.cxx.

References FindTrackInXYProjection2_InputOutputData::ALFA, FindTrackInXYProjection2_InputOutputData::apotemastrawdetectormin, atan2(), FindTrackInXYProjection2_InputOutputData::BETA, FindTrackInXYProjection2_InputOutputData::Charge, Cleanup(), FindTrackInXYProjection2_InputOutputData::Cosine, FindTrackInXYProjection2_InputOutputData::dimensionscitil, Double_t, fabs(), FindTrackInXYProjection2_InputOutputData::Fi_initial_helix_referenceframe, FindTrackInXYProjection2_InputOutputData::Fi_low_limit, FindTrackInXYProjection2_InputOutputData::Fi_up_limit, PndTrkChi2Fits::FitHelixCylinder(), PndTrkLegendreFits::FitHelixCylinder2(), FindTrackInXYProjection2_InputOutputData::GAMMA, i, FindTrackInXYProjection2_InputOutputData::icounter, FindTrackInXYProjection2_InputOutputData::InclusionListSciTil, FindTrackInXYProjection2_InputOutputData::InclusionListStt, FindTrackInXYProjection2_InputOutputData::info, FindTrackInXYProjection2_InputOutputData::infoparalConformal, FindTrackInXYProjection2_InputOutputData::legiandre_nradiusdiv, FindTrackInXYProjection2_InputOutputData::legiandre_nthetadiv, FindTrackInXYProjection2_InputOutputData::ListHitsinTrack, FindTrackInXYProjection2_InputOutputData::ListMvdPixelHitsinTrack, FindTrackInXYProjection2_InputOutputData::ListMvdStripHitsinTrack, FindTrackInXYProjection2_InputOutputData::ListSciTilHitsinTrack, FindTrackInXYProjection2_InputOutputData::ListSttParHits, m, FindTrackInXYProjection2_InputOutputData::maxhitsinfit, FindTrackInXYProjection2_InputOutputData::maxmvdpixelhitsintrack, FindTrackInXYProjection2_InputOutputData::maxmvdstriphitsintrack, FindTrackInXYProjection2_InputOutputData::maxscitilhitsintrack, FindTrackInXYProjection2_InputOutputData::maxstthits, FindTrackInXYProjection2_InputOutputData::maxstthitsintrack, FindTrackInXYProjection2_InputOutputData::Mvdhits, FindTrackInXYProjection2_InputOutputData::nHitsinTrack, FindTrackInXYProjection2_InputOutputData::nMvdPixelHit, FindTrackInXYProjection2_InputOutputData::nMvdPixelHitsinTrack, FindTrackInXYProjection2_InputOutputData::nMvdStripHit, FindTrackInXYProjection2_InputOutputData::nMvdStripHitsinTrack, NN, FindTrackInXYProjection2_InputOutputData::nSciTilHits, FindTrackInXYProjection2_InputOutputData::nSciTilHitsinTrack, FindTrackInXYProjection2_InputOutputData::nsttparhit, FindTrackInXYProjection2_InputOutputData::Oxx, FindTrackInXYProjection2_InputOutputData::Oyy, FindTrackInXYProjection2_InputOutputData::posizSciT, R, r2, FindTrackInXYProjection2_InputOutputData::Rr, FindTrackInXYProjection2_InputOutputData::rstrawdetectormax, FindTrackInXYProjection2_InputOutputData::S_SciTilHitsinTrack, FindTrackInXYProjection2_InputOutputData::Sinus, sqrt(), status, FindTrackInXYProjection2_InputOutputData::strawradius, FindTrackInXYProjection2_InputOutputData::thetamax, FindTrackInXYProjection2_InputOutputData::thetamin, FindTrackInXYProjection2_InputOutputData::trajectory_vertex, two_pi, FindTrackInXYProjection2_InputOutputData::TypeConf, FindTrackInXYProjection2_InputOutputData::XMvdPixel, FindTrackInXYProjection2_InputOutputData::XMvdStrip, FindTrackInXYProjection2_InputOutputData::YMvdPixel, and FindTrackInXYProjection2_InputOutputData::YMvdStrip.

Referenced by PndTrkTracking2::Exec().

317 {
318 
319  // load the struct members in local variables;
320 
321  Double_t (*info)[7]=InOut->info;
322  Double_t (*infoparalConformal)[5] = InOut->infoparalConformal ;
323  Short_t &nSttParHitsinTrack = (*(InOut->nHitsinTrack)),
324  nMvdPixelHit = InOut->nMvdPixelHit,
325  &nMvdPixelHitsinTrack = *(InOut->nMvdPixelHitsinTrack),
326  nMvdStripHit = InOut->nMvdStripHit,
327  &nMvdStripHitsinTrack = *(InOut->nMvdStripHitsinTrack),
328  *ListSttParHitsinTrack = InOut->ListHitsinTrack ,
329  *ListMvdPixelHitsinTrack = InOut->ListMvdPixelHitsinTrack ,
330  *ListMvdStripHitsinTrack = InOut->ListMvdStripHitsinTrack ;
331 
332  Double_t &Ox = *(InOut->Oxx),
333  &Oy =*(InOut->Oyy),
334  &R = *(InOut->Rr),
335  *XMvdPixel = InOut->XMvdPixel,
336  *XMvdStrip = InOut->XMvdStrip,
337  *YMvdPixel = InOut->YMvdPixel,
338  *YMvdStrip = InOut->YMvdStrip;
339 
340 
341  // make a POINTER to an ARRAY[][3] of
342  // Double_t and assign value present in the calling sequence of this method;
343 
344  Double_t (*posizSciTil)[3]=(Double_t (*)[3])InOut->posizSciT;
345 
346 //---------------
347 
348  //bool
349  //Type; //[R.K. 01/2017] unused variable?
350 //int istampa = 2;
351 
352  Short_t
353  auxListHitsinTrack[InOut->maxstthits],
354  flagStt,
355  i,
356  //iexcl, //[R.K. 01/2017] unused variable?
357  j,
358  //ListHitsinTrackinWhichToSearch[InOut->maxstthits], //[R.K. 01/2017] unused variable?
359  //Naux, //[R.K. 01/2017] unused variable?
360  //Nbaux, //[R.K. 01/2017] unused variable?
361  nFitPoints,
362  //Nint, //[R.K. 01/2017] unused variable?
363  NN,
364  //Nouter, //[R.K. 01/2017] unused variable?
365  //OutputListHitsinTrack[InOut->maxstthits], //[R.K. 01/2017] unused variable?
366  //OutputList2HitsinTrack[InOut->maxstthits], //[R.K. 01/2017] unused variable?
367  status;
368 
369 
370  Double_t
371  aaa,
372  //d, //[R.K. 01/2017] unused variable?
373  //diff, //[R.K. 01/2017] unused variable?
374  fiCenter,
375  fi_low_limit[2],
376  fi_up_limit[2],
377  FiRangeMvdLow,
378  FiRangeMvdUp,
379  gamma,
380  m,
381  q,
382  r2;
383  //rotationangle, //[R.K. 01/2017] unused variable?
384  //rotationcos, //[R.K. 01/2017] unused variable?
385  //rotationsin; //[R.K. 01/2017] unused variable?
386 
387 
388 //---------------------
389 /*
390  Double_t
391  tDriftRadiusconformal[InOut->maxhitsinfit],
392  tErrorDriftRadiusconformal[InOut->maxhitsinfit],
393  tXconformal[InOut->maxhitsinfit],
394  tYconformal[InOut->maxhitsinfit];
395  Vec <Double_t>
396  DriftRadiusconformal ( tDriftRadiusconformal, InOut->maxhitsinfit, "DriftRadiusconformal"),
397  ErrorDriftRadiusconformal ( tErrorDriftRadiusconformal, InOut->maxhitsinfit, "ErrorDriftRadiusconformal"),
398  Xconformal(tXconformal, InOut->maxhitsinfit, "Xconformal"),
399  Yconformal(tYconformal, InOut->maxhitsinfit, "Yconformal");
400 */
401  Double_t
402  DriftRadiusconformal[InOut->maxhitsinfit],
403  ErrorDriftRadiusconformal[InOut->maxhitsinfit],
404  Xconformal[InOut->maxhitsinfit],
405  Yconformal[InOut->maxhitsinfit];
406 
407 //--------------------------
408 
409 
410 
412  PndTrkCTGeometryCalculations GeomCalculator;
413  PndTrkPrintouts Print2;
414 
415 
416  nFitPoints = nSttParHitsinTrack;
417  if( nFitPoints > InOut->maxhitsinfit ) nFitPoints = InOut->maxhitsinfit;
418 
419  for(j=0; j< nFitPoints ; j++){
420  Xconformal[j] =infoparalConformal[ListSttParHitsinTrack[j]][0];
421  Yconformal[j] =infoparalConformal[ListSttParHitsinTrack[j]][1];
422  ErrorDriftRadiusconformal[j]=
423  infoparalConformal[ListSttParHitsinTrack[j]][2];
424  DriftRadiusconformal[j]=
425  infoparalConformal[ListSttParHitsinTrack[j]][2];
426  }
427 
428 
429 // PndTrkGlpkFits fit;
430 
431  PndTrkLegendreFits fitLegiandre;
432  PndTrkChi2Fits fitChi2;
433 
434 
435 
436  status = fitLegiandre.FitHelixCylinder2(
437  InOut->Cosine,
438  InOut->legiandre_nthetadiv,
439  InOut->legiandre_nradiusdiv,
440  nFitPoints,
441  Xconformal,
442  Yconformal,
443  DriftRadiusconformal,
444  ErrorDriftRadiusconformal,
445 
446 
447  0., // rotationangle, input;
448  InOut->Sinus,
449  InOut->thetamax, // input;
450  InOut->thetamin, // input;
451  InOut->trajectory_vertex, // vertex in (X,Y) of this trajectory
452  InOut->maxhitsinfit, // maximum n. of hits allowed in fast fit
453  &m,
454  &q,
455  InOut->ALFA,
456  InOut->BETA,
457  InOut->GAMMA,
458  InOut->TypeConf,
459  2, // istampa
460  InOut->icounter // IVOLTE
461  );
462 
463  if(status < 0 ) return false;
464 
465 
466 // this trasformation is valid even if the equation is a straight line from the fit
467  Ox= -0.5*(*(InOut->ALFA));
468  Oy= -0.5*(*(InOut->BETA));
469  R= Ox * Ox + Oy * Oy - (*(InOut->GAMMA));
470 
471 
472 // some obvious preliminary cuts
473  if( R < 0. ) return false;
474  R= sqrt( R );
475  aaa = sqrt( Ox * Ox + Oy * Oy );
476 
477 
478  // the following is because the circumference is supposed to come from (0,0);
479  // here the factor 0.9 is used in order to be conservative.
480  if(aaa< 0.9*InOut->apotemastrawdetectormin/2.) return false;
481 
482 // here the factor 0.9 is used in order to be conservative.
483  if ( R + aaa < InOut->apotemastrawdetectormin *0.9 ) return false;
484 
485 //-------------- debug printout
486 if (istampa>=2){
487 bool tkeepit[10];
488 tkeepit[0] = true;
489 Short_t nSttSkewHitsinTrack[10],
490  nSciTilHitsinTrack[10];
491 Short_t ListSttSkewHitsinTrack[100];
492 Double_t KAPPA[10] ;
493 KAPPA[0] = 0.;
494 nSttSkewHitsinTrack[0] = 0;
495 nSciTilHitsinTrack[0] = 0;
496 nMvdPixelHitsinTrack = 0;
497 nMvdStripHitsinTrack = 0;
498 
499 cout<<"\tstampa in CTFind..2.cxx, dopo FitHelixCylinder :" <<endl;
500 Print2.stampetta2(tkeepit,ListMvdPixelHitsinTrack,
501 ListMvdStripHitsinTrack,ListSttParHitsinTrack,
502 ListSttSkewHitsinTrack,InOut->ListSciTilHitsinTrack,
503 &nMvdPixelHitsinTrack,&nMvdStripHitsinTrack,&nSttParHitsinTrack,
504 
505 nSttSkewHitsinTrack,nSciTilHitsinTrack,1, // questo e' nTotCand, cioe' 1
506 -1, // print all candidates;
508 InOut->Fi_initial_helix_referenceframe,KAPPA);
509 }
510 //-------------- end of debug printout
511 //--------------------- find the angular range for candidate track in the Stt detector (find
512 // the two possibilities in general);
513 
514 // fi_low_limit[0], fi_up_limit[0] is the solution (radians) corresponding to the intersection
515 // on the RIGHT side, in the XY projection of the trajectory [looking into the beam] with respect
516 // to the segment joining the Center of the Helix with the origin (0,0);
517 // fi_low_limit[1], fi_up_limit[1] is the solution corresponding to the LEFT intersection; in case
518 // there is no second solution it is set at -100.;
519  GeomCalculator.FindingParallelTrackAngularRange2(
520  Ox,
521  Oy,
522  InOut->rstrawdetectormax,
524  R,
525  fi_low_limit,
526  fi_up_limit,
527  &flagStt
528  );
529  if(flagStt < 0) return false;
530 
531 
532 // decide which fi_up_limit, fi_low_limit are good by looking at where it is one hit of the cluster;
533 
534 // the fi of the origin [= (0,0,0) ] of the trajectory;
535  fiCenter = atan2(-Oy,-Ox);
536  if( fiCenter < 0. ) fiCenter += two_pi;
537  *(InOut->Fi_initial_helix_referenceframe) = fiCenter;
538 
539 
541  fiCenter, // input, fi of the (0,0) in the Helix reference frame;
542  fi_low_limit, // input; fi low limit of the Stt detector in the Helix frame;
543  fi_up_limit, // input; fi up limit of the Stt detector in the Helix frame;
544  info, // input
545  InOut->ListHitsinTrack, // input
546  nSttParHitsinTrack, // input
547  Ox, // input
548  Oy, // input
549 
550  *(InOut->Charge), // output; charge of the particle calculated with the
551  // trajectory present parameters;
552  FiRangeMvdLow, // output; Fi range of possible Mvd hits; FiRangeMvdLow always
553  // between 0 and 2PI; and always FiRangeMvdLow<FiRangeMvdUp;
554  FiRangeMvdUp, // output; Fi range of possible Mvd hits;
555  *InOut->Fi_low_limit, // output; inside DecideWhichAngulaRange this is an alias;
556  *InOut->Fi_up_limit // output; inside DecideWhichAngulaRange this is an alias;
557  );
558 
559 
560 // finds the intersection of the trajectory with the outer and inner STT detector radius
561 // and calculates the corresponding fi_up_limit and fi_low_limit in the Helix reference
562 // frame both for a positive and for a negative particle;
563 // fi_up_limit[0], fi_low_limit[0] --> for a positive track (it moves clockwise when looking
564 // at the beam);
565 // fi_up_limit[1], fi_low_limit[1] --> for a negative track (it moves anticlockwise when looking
566 // at the beam); or they are -100. when there is only one entrance and exit for the trajectory;
567 
568 //----
569 
570 
572  auxListHitsinTrack, // this is the output
573  InOut->InclusionListStt,
574  *(InOut->Fi_low_limit),
575  *(InOut->Fi_up_limit),
576  info,
577  InOut->ListSttParHits,
578  InOut->nsttparhit,
579  Ox,
580  Oy,
581  R,
582  InOut->strawradius
583  );
584 
585  // check that there are at least 2 axial STT hits (minimumhitspertrack is usually 2);
586  if( NN < InOut->minimumhitspertrack ) return false;
587 
588 
589  if( NN>InOut->maxstthitsintrack) {
590  nSttParHitsinTrack = InOut->maxstthitsintrack;
591  } else {
592  nSttParHitsinTrack=NN;
593  }
594 
595 
596  for(i=0; i< nSttParHitsinTrack;i++){
597  ListSttParHitsinTrack[i]=auxListHitsinTrack[i];
598  }
599 
600 
601 //-------------- debug printout
602 
603 if (istampa>=2){
604 bool tkeepit[10];
605 tkeepit[0] = true;
606 Short_t nSttSkewHitsinTrack[10],
607  nSciTilHitsinTrack[10];
608 Short_t ListSttSkewHitsinTrack[100];
609 Double_t KAPPA[10] ;
610 KAPPA[0] = 0.;
611 nSttSkewHitsinTrack[0] = 0;
612 nSciTilHitsinTrack[0] = 0;
613 
614 cout<<"\tstampa in CTFind..2.cxx, dopo TrkAssociatedParallelHitsToHelix5 :"<<endl;
615 Print2.stampetta2(tkeepit,ListMvdPixelHitsinTrack,
616 ListMvdStripHitsinTrack,ListSttParHitsinTrack,
617 ListSttSkewHitsinTrack,InOut->ListSciTilHitsinTrack,
618 &nMvdPixelHitsinTrack,&nMvdStripHitsinTrack,&nSttParHitsinTrack,
619 
620 nSttSkewHitsinTrack,nSciTilHitsinTrack,1, // questo e' nTotCand, cioe' 1
621 -1, // print all candidates;
623 InOut->Fi_initial_helix_referenceframe,KAPPA);
624 }
625 //-------------- end of debug printout
626 
627 
628 // adding the Mvd hits;
629 
630  Double_t
631  delta = 0.5, // parameter of proximity for associating Mvd hits to Stt tracks;
632  highqualitycut = 0.2; // parameter of proximity for associating Mvd hits to Stt tracks; at the moment it is not used
633  // actually ;
634 
636  delta, // input;
637  highqualitycut, // input;
638  FiRangeMvdLow, // input;
639  FiRangeMvdUp, // input;
640  InOut->maxmvdpixelhitsintrack, // input;
641  InOut->maxmvdstriphitsintrack, // input;
642  nMvdPixelHit, // input;
643  nMvdStripHit, // input;
644  Ox, // input;
645  Oy, // input;
646  R, // input;
647  XMvdPixel, // input;
648  XMvdStrip, // input;
649  YMvdPixel, // input;
650  YMvdStrip, // input;
651 
652  nMvdPixelHitsinTrack, // output
653  ListMvdPixelHitsinTrack, // output; dimensionality : [MAXMVDPIXELHITSINTRACK];
654  nMvdStripHitsinTrack, // output
655  ListMvdStripHitsinTrack // output; dimensionality : [MAXMVDSTRIPHITSINTRACK];
656  );
657 
658 
659 //-------------- debug printout
660 
661 if (istampa>=2){
662 bool tkeepit[10];
663 tkeepit[0] = true;
664 Short_t nSttSkewHitsinTrack[10],
665  nSciTilHitsinTrack[10];
666 Short_t ListSttSkewHitsinTrack[100];
667 Double_t KAPPA[10] ;
668 KAPPA[0] = 0.;
669 nSttSkewHitsinTrack[0] = 0;
670 nSciTilHitsinTrack[0] = 0;
671 
672 cout<<"\tstampa in CTFind..2.cxx, dopo AddMvdHitsToSttTracks (m = "<<m<<", q = "<<q<<", atan(m) = "<<atan(m)<<"):"<<endl;
673 Print2.stampetta2(tkeepit,ListMvdPixelHitsinTrack,
674 ListMvdStripHitsinTrack,ListSttParHitsinTrack,
675 ListSttSkewHitsinTrack,InOut->ListSciTilHitsinTrack,
676 &nMvdPixelHitsinTrack,&nMvdStripHitsinTrack,&nSttParHitsinTrack,
677 
678 nSttSkewHitsinTrack,nSciTilHitsinTrack,1, // questo e' nTotCand, cioe' 1
679 -1, // print all candidates;
681 InOut->Fi_initial_helix_referenceframe,KAPPA);
682 }
683 //-------------- end of debug printout
684 
685 
686 
687 
688 //-----------------------------------------------------------------------
689 
690  // if there are new Mvd hits attached, redo the fit;
691 
692  if( nMvdPixelHitsinTrack + nMvdStripHitsinTrack == 0 ){
693  // in this case no Mvd hits were added, skip and go to the SciTil section;
694  *(InOut->Mvdhits) = false;
695  *(InOut->ALFA) = -2.* Ox ;
696  *(InOut->BETA) = -2.* Oy ;
697  *(InOut->GAMMA) = Ox * Ox + Oy * Oy -R*R;
698  }
699 
700  // in this case new Mvd hits were added;
701  else {
702 
703  *(InOut->Mvdhits) = true;
704 
705  // fit again in XY projection, this time with the Chi2 fit;
706  // redo the Xconformal and Yconformal arrays since the
707  // hit composition of this track candidate may have changed; use first of all
708  // the Mvd hits and also all the Stt hits (provided the total is <= InOut->maxstthitsintrack);
709  // in the latter case the Mvd points are used first;
710 
711 
712  nFitPoints = 0;
713  // Mvd Pixels;
714  for(j=0; j<nMvdPixelHitsinTrack; j++){
715  if(nFitPoints==InOut->maxhitsinfit) break;
716  nFitPoints++;
717  r2 = XMvdPixel[ListMvdPixelHitsinTrack[j]]*XMvdPixel[ListMvdPixelHitsinTrack[j]]+
718  YMvdPixel[ListMvdPixelHitsinTrack[j]]*YMvdPixel[ListMvdPixelHitsinTrack[j]];
719  Xconformal[j] = XMvdPixel[ListMvdPixelHitsinTrack[j]]/r2 ;
720  Yconformal[j] = YMvdPixel[ListMvdPixelHitsinTrack[j]]/r2 ;
721 
722  // I assume the 'drift radius' to be the max dimension of the Pixel;
723  gamma = r2 - 0.01 * 0.01 ;
724  // 0.01 is the dimension of the pixel;
725  ErrorDriftRadiusconformal[j]= 3.*delta /fabs(gamma); // 3 is an arbitrary loose safety factor;
726 
727 
728  DriftRadiusconformal[j]= -1; // // only to signal later this is a Mvd hit;
729 
730 
731 
732  } // end of for(j=0; j<nMvdPixelHitsinTrack; j++)
733 
734  // Mvd Strips;
735  for(j=0; j<nMvdStripHitsinTrack; j++){
736  if(nFitPoints==InOut->maxhitsinfit) break;
737  r2 = XMvdStrip[ListMvdStripHitsinTrack[j]]*XMvdStrip[ListMvdStripHitsinTrack[j]]+
738  YMvdStrip[ListMvdStripHitsinTrack[j]]*YMvdStrip[ListMvdStripHitsinTrack[j]];
739  Xconformal[nFitPoints] = XMvdStrip[ListMvdStripHitsinTrack[j]]/r2 ;
740  Yconformal[nFitPoints] = YMvdStrip[ListMvdStripHitsinTrack[j]]/r2 ;
741  // I assume the 'drift radius' to be the max dimension of the Strip;
742  gamma = r2 - 0.01 * 0.01 ;
743  ErrorDriftRadiusconformal[nFitPoints]= 3. * 0.01 /fabs(gamma); // 3 is an arbitrary loose safety factor;
744  // 0.01 is the dimension of the strip;
745  DriftRadiusconformal[nFitPoints]= -1; // // only to signal later this is a Mvd hit;
746 
747  nFitPoints++;
748 
749 
750  }
751 
752  // axial Stt;
753  for(j=0; j<nSttParHitsinTrack; j++){
754  if(nFitPoints==InOut->maxhitsinfit) break;
755  Xconformal[nFitPoints] =infoparalConformal[(InOut->ListHitsinTrack)[j]][0];
756  Yconformal[nFitPoints] =infoparalConformal[(InOut->ListHitsinTrack)[j]][1];
757 
758  // errors on a Stt hit are assumed to be 0.5 cm in cartesian XY variables;
759  r2 = info [(InOut->ListHitsinTrack)[j] ][0]*info [(InOut->ListHitsinTrack)[j] ][0]+
760  info [(InOut->ListHitsinTrack)[j] ][1]*info [(InOut->ListHitsinTrack)[j] ][1];
761  gamma = r2 - 0.5 * 0.5 ;
762  ErrorDriftRadiusconformal[nFitPoints]=3.* 0.01/fabs(gamma); // 3 is an arbitrary loose safety factor;
763 
764  DriftRadiusconformal[nFitPoints]=
765  infoparalConformal[(InOut->ListHitsinTrack)[j]][2];
766 
767 
768 
769  nFitPoints++;
770 
771 
772  }
773 
774  // the following fit in XY with the Chi2 methods requires already an existing fit (in order
775  // to decide a priori on the left-right ambiguity of the Stt axial straws);
776  // The calculation of a good rotationangle is very, very useful to the success of the fitting
777  // in those cases when the trajectory is nearly perpendicular in the conformal space;
778 
779  status = fitChi2.FitHelixCylinder(
780  nFitPoints,
781  Xconformal,
782  Yconformal,
783  DriftRadiusconformal,
784  ErrorDriftRadiusconformal,
785 
786 
787  atan(m) , // rotationangle, in radians; m was found by the Hough transform fit;
788 // 0. , // rotationangle, in radians; m was found by the Hough transform fit;
789  InOut->trajectory_vertex, // vertex in (X,Y) of this trajectory
790  InOut->maxhitsinfit, // maximum n. of hits allowed in fast fit
791  &m,
792  &q,
793  InOut->ALFA, // input and output;
794  InOut->BETA, // input and output;
795  InOut->GAMMA, // input and output;
796  InOut->TypeConf,
797  2, // istampa
798  InOut->icounter // IVOLTE
799  );
800 
801 
802  if(status> 0 ) { // in this case the previous fit was successful, so proceed with further (better)
803  // association of Mvd hits and STT hits;
804  // otherwise go directly to the association of the SciTil hits;
805  Ox= -0.5*(*(InOut->ALFA));
806  Oy= -0.5*(*(InOut->BETA));
807  R= sqrt( Ox * Ox + Oy * Oy - (*(InOut->GAMMA)) );
808 
809 
810 
811 //-------------- debug printout
812 
813 if (istampa>=2){
814 bool tkeepit[10];
815 tkeepit[0] = true;
816 Short_t nSttSkewHitsinTrack[10],
817  nSciTilHitsinTrack[10];
818 Short_t ListSttSkewHitsinTrack[100];
819 Double_t KAPPA[10] ;
820 KAPPA[0] = 0.;
821 nSttSkewHitsinTrack[0] = 0;
822 nSciTilHitsinTrack[0] = 0;
823 
824 cout<<"\tstampa in CTFind..2.cxx, dopo fitChi2.FitHelixCylinder(m = "<<m<<", q = "<<q<<", atan(m) = "<<atan(m)<<"):"<<endl;
825 Print2.stampetta2(tkeepit,ListMvdPixelHitsinTrack,
826 ListMvdStripHitsinTrack,ListSttParHitsinTrack,
827 ListSttSkewHitsinTrack,InOut->ListSciTilHitsinTrack,
828 &nMvdPixelHitsinTrack,&nMvdStripHitsinTrack,&nSttParHitsinTrack,
829 
830 nSttSkewHitsinTrack,nSciTilHitsinTrack,1, // questo e' nTotCand, cioe' 1
831 -1, // print all candidates;
833 InOut->Fi_initial_helix_referenceframe,KAPPA);
834 }
835 //-------------- end of the debug printout
836 
837 
838  // find again the angular range for candidate track in the Stt detector (find
839  // the two possibilities in general);
840 
841  // fi_low_limit[0], fi_up_limit[0] is the solution (radians) corresponding to the intersection
842  // on the RIGHT side, in the XY projection of the trajectory [looking into the beam] with respect
843  // to the segment joining the Center of the Helix with the origin (0,0);
844  // fi_low_limit[1], fi_up_limit[1] is the solution corresponding to the LEFT intersection; in case
845  // there is no second solution it is set at -100.;
846  GeomCalculator.FindingParallelTrackAngularRange2(
847  Ox,
848  Oy,
849  InOut->rstrawdetectormax,
851  R,
852  fi_low_limit, // output;
853  fi_up_limit, // output;
854  &flagStt // output;
855  );
856 
857 
858  if(flagStt < 0) return false;
859 
860  // decide which fi_up_limit, fi_low_limit are good by looking at where it is one hit of the cluster;
861 
862  // the fi of the origin [= (0,0,0) ] of the trajectory;
863  fiCenter = atan2(-Oy,-Ox);
864  if( fiCenter < 0. ) fiCenter += two_pi;
865  *(InOut->Fi_initial_helix_referenceframe) = fiCenter;
866 
868  fiCenter, // input, fi of the (0,0) in the Helix reference frame;
869  fi_low_limit, // input; fi low limit of the Stt detector in the Helix frame;
870  fi_up_limit, // input; fi up limit of the Stt detector in the Helix frame;
871  info, // input
872  InOut->ListHitsinTrack, // input
873  nSttParHitsinTrack, // input
874  Ox, // input
875  Oy, // input
876 
877  *(InOut->Charge), // output; charge of the particle calculated with the
878  // trajectory present parameters;
879  FiRangeMvdLow, // output; Fi range of possible Mvd hits; FiRangeMvdLow always
880  // between 0 and 2PI; and always FiRangeMvdLow<FiRangeMvdUp;
881  // inside DecideWhichAngulaRange this is an alias;
882  FiRangeMvdUp, // output; Fi range of possible Mvd hits;
883  // inside DecideWhichAngulaRange this is an alias;
884  *InOut->Fi_low_limit, // output; inside DecideWhichAngulaRange this is an alias;
885  *InOut->Fi_up_limit // output; inside DecideWhichAngulaRange this is an alias;
886  );
887 
888 
889 
890 
891  // ------------------------ now redo the selection of the Mvd hits belonging to this track;
892 
893  delta=0.5; // parameter of proximity for associating Mvd hits to Stt tracks
894  // highqualitycut=0.3; // parameter of proximity for associating Mvd hits to Stt tracks
895  highqualitycut=0.5; // parameter of proximity for associating Mvd hits to Stt tracks
896 
898  delta, // input;
899  highqualitycut, // input;
900  FiRangeMvdLow, // input;
901  FiRangeMvdUp, // input;
902  InOut->maxmvdpixelhitsintrack, // input;
903  InOut->maxmvdstriphitsintrack, // input;
904  nMvdPixelHit, // input;
905  nMvdStripHit, // input;
906  Ox, // input;
907  Oy, // input;
908  R, // input;
909  XMvdPixel, // input;
910  XMvdStrip, // input;
911  YMvdPixel, // input;
912  YMvdStrip, // input;
913 
914  nMvdPixelHitsinTrack, // output
915  ListMvdPixelHitsinTrack, // output; dimensionality : [MAXMVDPIXELHITSINTRACK];
916  nMvdStripHitsinTrack, // output
917  ListMvdStripHitsinTrack // output; dimensionality : [MAXMVDSTRIPHITSINTRACK];
918  );
919 
920  // ------------------------ now redo the selection of the STT Axial hits belonging to this track;
921 
922 
923  // try TrkAssociatedParallelHitsToHelix6 that is the same as TrkAssociatedParallelHitsToHelix5 except that
924  // it takes as input the maximum distance allowed for an associated hit;
926  auxListHitsinTrack, // this is the output
927  InOut->InclusionListStt,
928  *(InOut->Fi_low_limit),
929  *(InOut->Fi_up_limit),
930  info,
931  InOut->ListSttParHits,
932  InOut->nsttparhit,
933  Ox,
934  Oy,
935  R,
936  2.*InOut->strawradius // this is the maximum allowed distance;
937  );
938 
939  if( NN < InOut->minimumhitspertrack ) return false;
940  if( NN>InOut->maxstthitsintrack) {
941  nSttParHitsinTrack = InOut->maxstthitsintrack;
942  } else {
943  nSttParHitsinTrack=NN;
944  }
945 
946 
947  for(i=0; i< nSttParHitsinTrack;i++){
948  ListSttParHitsinTrack[i]=auxListHitsinTrack[i];
949  }
950 
951  } // end of if(status> 0 )
952 
953 
954 
955  } // end of if( nMvdPixelHitsinTrack + nMvdStripHitsinTrack == 0 )
956 
957 //---------------------------
958 
959  // now the section that associates the SciTil hits to this track candidate;
960 
961 // equation of the SciTil segment : y0 * y + x0 * x - x0**2 - y0**2 = 0
962 // where (x0,y0) = position of center of the SciTil.
963 
964 // delimiting points of the SciTil segment : define L = length of the SciTil,
965 // and RR = sqrt(x0**2+y0**2), SIGN = the sign of (-x0*y0) or SIGN=1 when y0=0,
966 // SIGN=irrelevant when x0=0; then :
967 // P1 = [ x0- abs{(L/2)*y0/RR}; y0-SIGN*abs{(L/2)*x0/RR} ],
968 // P2 = [ x0+abs{(L/2)*y0/RR}; y0+SIGN*abs{(L/2)*x0/RR} ].
969 
970 
971 
973  InOut->dimensionscitil,
974  InOut->S_SciTilHitsinTrack,// output; S on the lateral face of the Helix
975  // of the SciTil hit (if present).
976  InOut->InclusionListSciTil,
977  InOut->ListSciTilHitsinTrack,
978  InOut->maxscitilhitsintrack,
979  InOut->nSciTilHits,
980  *(InOut->Oxx),
981  *(InOut->Oyy),
982  posizSciTil,
983  *(InOut->Rr)
984  );
985 
986  // even though it should be impossible in principle, EXCLUDE the possibility of having more
987  // than TWO SciTil hits belonging to a track;
988 
989  if( *(InOut->nSciTilHitsinTrack) >0 ){
990  // even though it should be impossible in principle, EXCLUDE
991  // the possibility of having more
992  // than TWO SciTil hits belonging to a track;
993  if( *(InOut->nSciTilHitsinTrack) > 2 ) *(InOut->nSciTilHitsinTrack)=2;
994 
995 // for(j=0;j<*(InOut->nSciTilHitsinTrack);j++){
996 // (InOut->InclusionListSciTil)[(InOut->ListSciTilHitsinTrack)[j]]
997 // =false;
998 // }
999 
1000 
1001 
1002  }
1003 
1004 
1005  // orderig the hits in this track cand; for a trajectory Radius not too large
1006  // better the ordering with conformal.
1007 
1008 
1009  if( R< InOut->rstrawdetectormax){
1010 
1011  // the origin of the tack is assumed to be (0,0);
1012 
1013  // ordering the Stt axial using the conformal coordinates;
1015  *InOut->Charge, // input;
1016  info, // input;
1017  nSttParHitsinTrack, // input;
1018  Ox, // input;
1019  Oy, // input;
1020 
1021  ListSttParHitsinTrack // input and output;
1022  );
1023 
1024 
1025  } else { // otherwise it is better distance from (0,0) method.
1027  info, // input;
1028  nSttParHitsinTrack, // input;
1029  ListSttParHitsinTrack // input and output;
1030  );
1031  }
1032 
1033 
1034 //-------------- debug printout
1035 if (istampa>=2){
1036 bool tkeepit[10];
1037 tkeepit[0] = true;
1038 Short_t nSttSkewHitsinTrack[10],
1039  nSciTilHitsinTrack[10];
1040 Short_t ListSttSkewHitsinTrack[100];
1041 Double_t KAPPA[10] ;
1042 KAPPA[0] = 0.;
1043 nSttSkewHitsinTrack[0] = 0;
1044 nSciTilHitsinTrack[0] = *(InOut->nSciTilHitsinTrack);
1045 
1046 cout<<"\tstampa in CTFind..2.cxx, dopo ordering :"<<endl;
1047 Print2.stampetta2(tkeepit,ListMvdPixelHitsinTrack,
1048 ListMvdStripHitsinTrack,ListSttParHitsinTrack,
1049 ListSttSkewHitsinTrack,InOut->ListSciTilHitsinTrack,
1050 &nMvdPixelHitsinTrack,&nMvdStripHitsinTrack,&nSttParHitsinTrack,
1051 
1052 nSttSkewHitsinTrack,nSciTilHitsinTrack,1, // questo e' nTotCand, cioe' 1
1053 -1, // print all candidates;
1055 InOut->Fi_initial_helix_referenceframe,KAPPA);
1056 
1057 
1058 }
1059 //-------------- end debug printout
1060 
1061 
1062  return true;
1063 
1064 
1065 };
Int_t i
Definition: run_full.C:25
__m128 m
Definition: P4_F32vec4.h:28
Short_t FitHelixCylinder2(Double_t *Cosine, Short_t LEGIANDRE_NTHETADIV, Short_t LEGIANDRE_NRADIUSDIV, Short_t nHitsinTrack, Double_t *Xconformal, Double_t *Yconformal, Double_t *DriftRadiusconformal, Double_t *ErrorDriftRadiusconformal, Double_t rotationangle, Double_t *Sinus, Double_t THETAMAX, Double_t THETAMIN, Double_t trajectory_vertex[2], Short_t NMAX, Double_t *m, Double_t *q, Double_t *pAlfa, Double_t *pBeta, Double_t *pGamma, bool *Type, int istampa, int IVOLTE)
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
void AddMvdHitsToSttTracks(Double_t delta, Double_t highqualitycut, Double_t FiRangeMvdLow, Double_t FiRangeMvdUp, Short_t maxmvdpixelhitsintrack, Short_t maxmvdstriphitsintrack, Short_t nMvdPixelHit, Short_t nMvdStripHit, Double_t Ox, Double_t Oy, Double_t R, Double_t *XMvdPixel, Double_t *XMvdStrip, Double_t *YMvdPixel, Double_t *YMvdStrip, Short_t &nPixelHitsinTrack, Short_t *ListPixelHitsinTrack, Short_t &nStripHitsinTrack, Short_t *ListStripHitsinTrack)
Short_t FitHelixCylinder(Short_t nHitsinTrack, Double_t *Xconformal, Double_t *Yconformal, Double_t *DriftRadiusconformal, Double_t *ErrorDriftRadiusconformal, Double_t rotationangle, Double_t trajectory_vertex[2], Short_t NMAX, Double_t *m, Double_t *q, Double_t *pAlfa, Double_t *pBeta, Double_t *pGamma, bool *Type, int istampa, int IVOLTE)
Short_t AssociateSciTilHit(Double_t dimensionscitil, Double_t *esse, bool *InclusionListSciTil, Short_t *List, Short_t maxscitilhitsintrack, Short_t nSciTilHits, Double_t Oxx, Double_t Oyy, Double_t posizSciTil[][3], Double_t Rr)
Short_t TrkAssociatedParallelHitsToHelix5(Short_t *auxListHitsinTrack, bool *InclusionListStt, Double_t Fi_low, Double_t Fi_up, Double_t info[][7], Short_t *ListSttParHits, Int_t NhitsParallel, Double_t Oxx, Double_t Oyy, Double_t Rr, Double_t strawradius)
Short_t TrkAssociatedParallelHitsToHelix6(Short_t *auxListHitsinTrack, bool *InclusionListStt, Double_t Fi_low, Double_t Fi_up, Double_t info[][7], Short_t *ListSttParHits, Int_t NhitsParallel, Double_t Oxx, Double_t Oyy, Double_t Rr, Double_t maximum_distance)
#define two_pi
Double_t
void OrderingUsingR(Double_t info[][7], Int_t nHits, Short_t *ListHits)
void DecideWhichAngularRangeAndCharge(Double_t fiCenter, Double_t fi_low_limit[2], Double_t fi_up_limit[2], Double_t(*info)[7], Short_t *ListHitsinTrack, Short_t nHitsinTrack, Double_t Oxx, Double_t Oyy, Short_t &charge, Double_t &FiRangeMvdLow, Double_t &FiRangeMvdUp, Double_t &Fi_low_limit, Double_t &Fi_up_limit)
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
friend F32vec4 atan2(const F32vec4 &y, const F32vec4 &x)
Definition: P4_F32vec4.h:117
SttMvdTracking Cleanup()
void OrderingUsingFi(Short_t Charge, Double_t info[][7], Int_t nHits, Double_t oX, Double_t oY, Short_t *ListHits)
Double_t R
Definition: checkhelixhit.C:61
double r2
int status[10]
Definition: f_Init.h:28
void PndTrkCTFindTrackInXY2::OrderingUsingConformal ( Short_t  Charge,
Double_t  info[][7],
Int_t  nHits,
Double_t  oX,
Double_t  oY,
Short_t *  ListHits 
)

Definition at line 1074 of file PndTrkCTFindTrackInXY2.cxx.

References atan2(), b1, Double_t, i, PndTrkMergeSort::Merge_Sort2(), nHits, and PI.

1082 {
1083 
1084  Short_t i,j,
1085  tmp[nHits];
1086  Double_t aaa,
1087  b1,
1088  U[nHits],
1089  V[nHits];
1090 
1091 
1092  PndTrkMergeSort MergeSort;
1093 
1094 
1095 
1096 // here there is the ordering of the hits, NOT under the assumption that the circumference
1097 // in XY goes through Trajectory_Start.
1098 // Moreover, the code before is supposed to have selected trajectories in XY with (fOx,fOy)
1099 // farther from (0,0) by > 0.9 * RminStrawDetector/2 and consequently fOx and fOy are not both 0.
1100 // The scheme for the ordering of the hit is as follows :
1101 // 1) order hits by increasing U or V of the conformal mapping; see Gianluigi's Logbook page 283;
1102 // 2) find the charge of the track by checking if it is closest to the center in XY
1103 // the first or the last of the ordered hits.
1104 // 3) in case, invert the ordering of U, V and ListHits such that the first hits in the
1105 // list are always those closer to the Trajectory_Start.
1106 
1107 
1108 // ordering of the hits
1109 
1110  aaa = atan2( oY, oX); // atan2 defined between -PI and PI.
1111 
1112  // the following statement is necessary since for unknown reason the root interpreter
1113  // gives a weird error when using PI directly in the if statement below!!!!!!! I lost
1114  // 2 hours trying to figure this out!
1115  b1 = PI/4.;
1116 
1117  if((aaa>b1&&aaa<3.*b1) || (aaa>-3.*b1&&aaa<-b1)){//use U as ordering variable;
1118  //[case 1 or 3 Gianluigi's Logbook page 285].
1119  for (j = 0; j< nHits; j++){
1120  U[j]= info[ListHits[j]][0]/(info[ListHits[j]][0]*info[ListHits[j]][0]+
1121  info[ListHits[j]][1]*info[ListHits[j]][1]);
1122  }
1123  MergeSort.Merge_Sort2( nHits, U, ListHits);
1124 
1125  if((aaa>b1&&aaa<3.*b1)){ // case #1;
1126  if( Charge == -1){
1127  // inverting the order of the hits.
1128  for(i=0;i<nHits;i++){
1129  tmp[i]=ListHits[nHits-1-i];
1130  }
1131  for(i=0;i<nHits;i++){
1132  ListHits[i]=tmp[i];
1133  }
1134  }
1135  } else{ // case # 3.
1136  if(Charge == 1){
1137  // inverting the order of the hits.
1138  for(i=0;i<nHits;i++){
1139  tmp[i]=ListHits[nHits-1-i];
1140  }
1141  for(i=0;i<nHits;i++){
1142  ListHits[i]=tmp[i];
1143  }
1144  }// end of if( Charge ==1)
1145  }// end of if((aaa>b1&&aaa<3.*b1))
1146 
1147  } else { // use V as ordering variable [case 2 or 4 Gianluigi's Logbook page 285].
1148  for (j = 0; j< nHits; j++){
1149  V[j]= info[ListHits[j]][1]/(info[ListHits[j]][0]*info[ListHits[j]][0]+
1150  info[ListHits[j]][1]*info[ListHits[j]][1]);
1151  }
1152  MergeSort.Merge_Sort2( nHits, V, ListHits);
1153 
1154  if((aaa<=-3.*b1 || aaa>=3.*b1)){ // case #2;
1155  if( Charge == -1){
1156  // inverting the order of the hits.
1157  for(i=0;i<nHits;i++){
1158  tmp[i]=ListHits[nHits-1-i];
1159  }
1160  for(i=0;i<nHits;i++){
1161  ListHits[i]=tmp[i];
1162  }
1163  }
1164  } else{ // case # 4.
1165  if( Charge == 1){
1166  // inverting the order of the hits.
1167  for(i=0;i<nHits;i++){
1168  tmp[i]=ListHits[nHits-1-i];
1169  }
1170  for(i=0;i<nHits;i++){
1171  ListHits[i]=tmp[i];
1172  }
1173  }
1174  }
1175 
1176  } // end of if((aaa>b1&& ....
1177 
1178 
1179 
1180  return;
1181 
1182 
1183 }
Int_t i
Definition: run_full.C:25
int nHits
Definition: RiemannTest.C:16
Double_t
void Merge_Sort2(Short_t n_ele, Double_t *array, Short_t *ind)
#define PI
friend F32vec4 atan2(const F32vec4 &y, const F32vec4 &x)
Definition: P4_F32vec4.h:117
void PndTrkCTFindTrackInXY2::OrderingUsingFi ( Short_t  Charge,
Double_t  info[][7],
Int_t  nHits,
Double_t  oX,
Double_t  oY,
Short_t *  ListHits 
)

Definition at line 1192 of file PndTrkCTFindTrackInXY2.cxx.

References atan2(), Double_t, fi, PndTrkMergeSort::Merge_Sort2(), nHits, and PI.

1200 {
1201 
1202  Short_t iaux[nHits],
1203  j;
1204  Double_t fi[nHits],
1205  Fi0;
1206 
1207  PndTrkMergeSort MergeSort;
1208 
1209 
1210 
1211 // here there is the ordering of the hits under the assumption that the circumference
1212 // in XY goes through (0,0).
1213 
1214 // ordering of the hits
1215 
1216  Fi0 = atan2( -oY, -oX); // atan2 defined between -PI and PI.
1217 
1218  // the following statement is necessary since for unknown reason the root interpreter
1219  if( Charge <0 ){ // particle rotates counterclockwise (beam direction along Z);
1220  for (j = 0; j< nHits; j++){
1221  fi[j]= atan2( info[ListHits[j]][1]-oY, info[ListHits[j]][0]-oX);
1222  if( fi[j] < Fi0 ) fi[j] += 2.* PI;
1223  if( fi[j] < Fi0 ) fi[j] = Fi0;
1224  }
1225  MergeSort.Merge_Sort2( nHits, fi, ListHits);
1226 
1227 
1228  } else { // particle rotates clockwise;
1229  for (j = 0; j< nHits; j++){
1230  fi[j]= atan2( info[ListHits[j]][1]-oY, info[ListHits[j]][0]-oX);
1231  if( fi[j] > Fi0 ) fi[j] -= 2.* PI;
1232  if( fi[j] > Fi0 ) fi[j] = Fi0;
1233  }
1234  MergeSort.Merge_Sort2( nHits, fi, ListHits);
1235  // it is necessary now to invert the order;
1236  for (j = 0; j< nHits; j++){
1237  iaux[j] = ListHits[nHits-j-1];
1238  }
1239 
1240  for (j = 0; j< nHits; j++){
1241  ListHits[j]=iaux[j];
1242  }
1243  } // end of if( Charge <0 )
1244 
1245 
1246 
1247  return;
1248 
1249 
1250 }
int nHits
Definition: RiemannTest.C:16
TFile * fi
Double_t
void Merge_Sort2(Short_t n_ele, Double_t *array, Short_t *ind)
#define PI
friend F32vec4 atan2(const F32vec4 &y, const F32vec4 &x)
Definition: P4_F32vec4.h:117
void PndTrkCTFindTrackInXY2::OrderingUsingR ( Double_t  info[][7],
Int_t  nHits,
Short_t *  ListHits 
)

Definition at line 1258 of file PndTrkCTFindTrackInXY2.cxx.

References Double_t, fabs(), i, PndTrkMergeSort::Merge_Sort2(), and nHits.

1263 {
1264  Short_t i,
1265  iaux;
1266  //j, //[R.K. 01/2017] unused variable?
1267  //ipar, //[R.K. 01/2017] unused variable?
1268  //iskew; //[R.K. 01/2017] unused variable?
1269  Double_t aux,
1270  auxR2[nHits],
1271  distq1,
1272  distq2;
1273 
1274  PndTrkMergeSort MergeSort;
1275 
1276 // ordering all the hits belonging to the candidate track, by increasing fR;
1277 // forming the new track with Mvd+Stt hits
1278 
1279 
1280  for(i=0; i<nHits; i++){
1281  auxR2[i] =
1282  info[ ListHits[i] ][0]*
1283  info[ ListHits[i] ][0]+
1284  info[ ListHits[i] ][1]*
1285  info[ ListHits[i] ][1];
1286  }
1287 
1288 
1289  // ordering the Stt Hits
1290  MergeSort.Merge_Sort2( nHits, auxR2, ListHits);
1291 
1292  // take care of the case when 2 hits ave the same distance from (0,0). this case can
1293  // happen (as a difference with the Ordering with the Conformal method) typically
1294  // when the last 2 hits lie at the boundary of the outer axial layer OR the first 2 hits lie at the boundary
1295  // of the inner axial layer ;
1296 
1297  // case of the two first hits at the same R : decide order based on proximity to the 3 hit;
1298  if( fabs( auxR2[0] - auxR2[1]) < 0.1) {
1299  distq1 = (info[ListHits[0]][0]-info[ListHits[2]][0])*(info[ListHits[0]][0]-info[ListHits[2]][0])+
1300  (info[ListHits[0]][1]-info[ListHits[2]][1])*(info[ListHits[0]][1]-info[ListHits[2]][1]);
1301  distq2 = (info[ListHits[1]][0]-info[ListHits[2]][0])*(info[ListHits[1]][0]-info[ListHits[2]][0])+
1302  (info[ListHits[1]][1]-info[ListHits[2]][1])*(info[ListHits[1]][1]-info[ListHits[2]][1]);
1303 
1304  if(distq1 < distq2 ){ // exchange place of hits;
1305  iaux = ListHits[0];
1306  ListHits[0] = ListHits[1];
1307  ListHits[1] = iaux;
1308  aux = auxR2[0];
1309  auxR2[0]=auxR2[1];
1310  auxR2[1]=aux;
1311  } // end of if(distq1 < distq2 )
1312  }
1313 
1314 
1315 
1316  for(i=2; i<nHits; i++){
1317  if( auxR2[i] == auxR2[i-1]){
1318  // ambiguity case; decide the order based on which of the two hits (hit # i and # i-1 )is closer the
1319  // hit # i-2;
1320  distq1 = (info[ListHits[i]][0]-info[ListHits[i-2]][0])*(info[ListHits[i]][0]-info[ListHits[i-2]][0])+
1321  (info[ListHits[i]][1]-info[ListHits[i-2]][1])*(info[ListHits[i]][1]-info[ListHits[i-2]][1]);
1322 
1323  distq2 = (info[ListHits[i-1]][0]-info[ListHits[i-2]][0])*(info[ListHits[i-1]][0]-info[ListHits[i-2]][0])+
1324  (info[ListHits[i-1]][1]-info[ListHits[i-2]][1])*(info[ListHits[i-1]][1]-info[ListHits[i-2]][1]);
1325 
1326  if(distq1<distq2){
1327  // exchange the order;
1328  iaux = ListHits[i];
1329  ListHits[i] = ListHits[i-1];
1330  ListHits[i-1] = iaux;
1331  aux = auxR2[i];
1332  auxR2[i]=auxR2[i-1];
1333  auxR2[i-1]=aux;
1334  }
1335 
1336  }
1337  } // end of for(i=0; i<nHits; i++)
1338 
1339  return;
1340 }
Int_t i
Definition: run_full.C:25
int nHits
Definition: RiemannTest.C:16
Double_t
void Merge_Sort2(Short_t n_ele, Double_t *array, Short_t *ind)
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
Short_t PndTrkCTFindTrackInXY2::TrkAssociatedParallelHitsToHelix5 ( Short_t *  auxListHitsinTrack,
bool *  InclusionListStt,
Double_t  Fi_low,
Double_t  Fi_up,
Double_t  info[][7],
Short_t *  ListSttParHits,
Int_t  NhitsParallel,
Double_t  Oxx,
Double_t  Oyy,
Double_t  Rr,
Double_t  strawradius 
)

Definition at line 1348 of file PndTrkCTFindTrackInXY2.cxx.

References angle, atan2(), Double_t, dx, dy, fabs(), i, PI, and sqrt().

1361 {
1362 
1363  Short_t i;
1364 
1365  Short_t nAssociatedHits;
1366 
1367  Double_t angle,
1368  dx,
1369  dy,
1370  distance,
1371 // NTIMES=1.5; // number of Straw radia allowed in association.
1372  NTIMES=2.; // number of Straw radia allowed in association.
1373 
1374 
1375  nAssociatedHits=0;
1376 // find the Hits belonging to this Track.
1377 
1378 
1379 
1380  for(i=0; i<NhitsParallel;i++){
1381  if( !InclusionListStt[ ListSttParHits[i] ] ) continue;
1382 // check if the hit position is near the circle of the Helix found by the fit
1383  dx = -Oxx+info[ListSttParHits[i]][0];
1384  dy = -Oyy+info[ListSttParHits[i]][1];
1385  angle=atan2(dy,dx);
1386  if(angle<0.) angle += 2.*PI;
1387  if(angle<0.) angle =0.;
1388  distance = sqrt(dx*dx+dy*dy);
1389 
1390  if ( fabs(Rr - distance ) > NTIMES*strawradius ) continue;
1391  if(angle<Fi_low) angle += 2.*PI;
1392  if(angle>Fi_up) continue;
1393  auxListHitsinTrack[nAssociatedHits]= ListSttParHits[i];
1394  nAssociatedHits++;
1395  } // end for(i=0; i<NhitsParallel;i++)
1396 
1397  return nAssociatedHits;
1398 
1399 };
double dy
Int_t i
Definition: run_full.C:25
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
Double_t
#define PI
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
double dx
friend F32vec4 atan2(const F32vec4 &y, const F32vec4 &x)
Definition: P4_F32vec4.h:117
Double_t angle
Short_t PndTrkCTFindTrackInXY2::TrkAssociatedParallelHitsToHelix6 ( Short_t *  auxListHitsinTrack,
bool *  InclusionListStt,
Double_t  Fi_low,
Double_t  Fi_up,
Double_t  info[][7],
Short_t *  ListSttParHits,
Int_t  NhitsParallel,
Double_t  Oxx,
Double_t  Oyy,
Double_t  Rr,
Double_t  maximum_distance 
)

Definition at line 1408 of file PndTrkCTFindTrackInXY2.cxx.

References angle, atan2(), Double_t, dx, dy, fabs(), i, PI, and sqrt().

1421 {
1422 
1423  Short_t i;
1424 
1425  Short_t nAssociatedHits;
1426 
1427  Double_t angle,
1428  dx,
1429  dy,
1430  distance;
1431 
1432  nAssociatedHits=0;
1433 // find the Hits belonging to this Track.
1434 
1435 
1436 
1437  for(i=0; i<NhitsParallel;i++){
1438  if( !InclusionListStt[ ListSttParHits[i] ] ) continue;
1439 // check if the hit position is near the circle of the Helix found by the fit
1440  dx = -Oxx+info[ListSttParHits[i]][0];
1441  dy = -Oyy+info[ListSttParHits[i]][1];
1442  angle=atan2(dy,dx);
1443  if(angle<0.) angle += 2.*PI;
1444  if(angle<0.) angle =0.;
1445  distance = sqrt(dx*dx+dy*dy);
1446 
1447  if ( fabs(Rr - distance ) > maximum_distance ) continue;
1448  if(angle<Fi_low) angle += 2.*PI;
1449  if(angle>Fi_up) continue;
1450  auxListHitsinTrack[nAssociatedHits]= ListSttParHits[i];
1451  nAssociatedHits++;
1452  } // end for(i=0; i<NhitsParallel;i++)
1453 
1454  return nAssociatedHits;
1455 
1456 };
double dy
Int_t i
Definition: run_full.C:25
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
Double_t
#define PI
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
double dx
friend F32vec4 atan2(const F32vec4 &y, const F32vec4 &x)
Definition: P4_F32vec4.h:117
Double_t angle

The documentation for this class was generated from the following files: