StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFtpcConfMapper.cc
1 // $Id: StFtpcConfMapper.cc,v 1.38 2012/11/07 23:30:18 fisyak Exp $
2 // $Log: StFtpcConfMapper.cc,v $
3 // Revision 1.38 2012/11/07 23:30:18 fisyak
4 // Supress warnings
5 //
6 // Revision 1.37 2010/06/03 11:09:39 jcs
7 // correction for Bug #1939 - TObjArray track was still used after it was removed
8 // removed inactive code
9 //
10 // Revision 1.36 2010/06/02 12:12:01 jcs
11 // undo correction for Bug#1939 since it was incorrect
12 //
13 // Revision 1.35 2010/05/19 14:55:29 jcs
14 // Correction for Bug #1939 - variables with the same name were defined twice;
15 // once in StFtpcTrackingParams.hh and again in StFtpcConfMapper.hh
16 // The names of the StFtpcConfMapper variables have been changed to make them unique
17 //
18 // Revision 1.34 2007/02/06 11:42:14 jcs
19 // move unessential output messages from INFO to DEBUG
20 //
21 // Revision 1.33 2007/01/15 08:23:01 jcs
22 // replace printf, cout and gMesMgr with Logger commands
23 //
24 // Revision 1.32 2004/09/03 20:36:22 perev
25 // Big LeakOff + mem optimisation
26 //
27 // Revision 1.31 2004/02/12 19:37:09 oldi
28 // *** empty log message ***
29 //
30 // Revision 1.30 2004/01/28 01:41:32 jeromel
31 // *** empty log message ***
32 //
33 // Revision 1.29 2003/11/20 03:14:37 perev
34 // LeakOff
35 //
36 // Revision 1.28 2003/10/12 04:15:20 perev
37 // few bugs fixed. division by zero disappeared
38 //
39 // Revision 1.27 2003/10/08 12:48:25 jcs
40 // initialize mNumMainVertexTracks which I overlooked the first time
41 //
42 // Revision 1.26 2003/10/02 15:07:34 jcs
43 // initialize mNumMainVertexTracks
44 //
45 // Revision 1.25 2003/09/16 15:27:01 jcs
46 // removed inline as it would leave a few undefined reference
47 //
48 // Revision 1.24 2003/01/21 10:04:13 jcs
49 // initialize variables to eliminate compiler warnings for NODEBUG=yes
50 //
51 // Revision 1.23 2003/01/20 13:16:23 oldi
52 // Additional volume segment added as garbage container. Hits which give a
53 // segment index which is out of range (esp. those ones sitting exactly on the
54 // beam line) are put in here.
55 // Handling of function GetSegm() simplified.
56 //
57 // Revision 1.22 2002/11/06 14:05:35 oldi
58 // Minor (minimal) clean up.
59 //
60 // Revision 1.21 2002/11/06 13:44:47 oldi
61 // Vertex handling simplifed.
62 // Flag for clusters not to be used for tracking introduced.
63 // Code clean ups.
64 //
65 // Revision 1.20 2002/10/11 15:44:58 oldi
66 // Get FTPC geometry and dimensions from database.
67 // No field fit activated: Returns momentum = 0 but fits a helix.
68 // Bug in TrackMaker fixed (events with z_vertex > outer_ftpc_radius were cut).
69 // QA histograms corrected (0 was supressed).
70 // Code cleanup (several lines of code changed due to *params -> Instance()).
71 // cout -> gMessMgr.
72 //
73 // Revision 1.19 2002/04/29 15:49:33 oldi
74 // All tracking parameters moved to StFtpcTrackingParameters.cc/hh.
75 // In a future version the actual values should be moved to an .idl file (the
76 // interface should be kept as is in StFtpcTrackingParameters.cc/hh).
77 //
78 // Revision 1.18 2002/04/08 15:37:43 oldi
79 // Switch for magnetic field factor installed.
80 // Minor corrections/improvements.
81 //
82 // Revision 1.17 2002/04/05 16:50:04 oldi
83 // Cleanup of MomentumFit (StFtpcMomentumFit is now part of StFtpcTrack).
84 // Each Track inherits from StHelix, now.
85 // Therefore it is possible to calculate, now:
86 // - residuals
87 // - vertex estimations obtained by back extrapolations of FTPC tracks
88 // Chi2 was fixed.
89 // Many additional minor (and major) changes.
90 //
91 // Revision 1.16 2002/02/21 22:57:56 oldi
92 // Fixes to avoid warnings during optimized compilation.
93 //
94 // Revision 1.15 2001/09/19 21:01:28 jcs
95 // change track finding and fitting settings
96 //
97 // Revision 1.14 2001/07/12 08:29:02 oldi
98 // New constructor introduced to be able to use only cluster which were
99 // found to be good in a previous run.
100 // New function NoFieldTracking() introduced.
101 // Tracking parameters for LaserTracking() changed.
102 //
103 // Revision 1.13 2001/04/25 17:53:06 perev
104 // HPcorrs
105 //
106 // Revision 1.12 2001/03/30 21:30:27 jeromel
107 // Constructor #2 argument 1 had an implicit/hidden copy of object.
108 // Removed const in declaration.
109 //
110 // Revision 1.11 2001/01/30 13:31:27 oldi
111 // New variable mTime introduced to count total time consumption.
112 //
113 // Revision 1.10 2001/01/25 15:21:33 oldi
114 // Review of the complete code.
115 // Fix of several bugs which caused memory leaks:
116 // - Tracks were not allocated properly.
117 // - Tracks (especially split tracks) were not deleted properly.
118 // - TClonesArray seems to have a problem (it could be that I used it in a
119 // wrong way). I changed all occurences to TObjArray which makes the
120 // program slightly slower but much more save (in terms of memory usage).
121 // Speed up of HandleSplitTracks() which is now 12.5 times faster than before.
122 // Cleanup.
123 //
124 // Revision 1.9 2000/11/10 18:37:01 oldi
125 // New functions introduced to be able to extend tracks after main vertex tracking.
126 // This implied many changes in other parts of the class (loop over clusters can run backwards now).
127 // TBenchmark object ('mBench') moved to StFtpcTracker.
128 // Cleanup.
129 //
130 // Revision 1.8 2000/07/24 02:42:47 oldi
131 // Problem of memory exhaustion solved. This was introduced during the last changes.
132 //
133 // Revision 1.7 2000/07/18 21:22:15 oldi
134 // Changes due to be able to find laser tracks.
135 // Cleanup: - new functions in StFtpcConfMapper, StFtpcTrack, and StFtpcPoint
136 // to bundle often called functions
137 // - short functions inlined
138 // - formulas of StFormulary made static
139 // - avoid streaming of objects of unknown size
140 // (removes the bunch of CINT warnings during compile time)
141 // - two or three minor bugs cured
142 //
143 // Revision 1.6 2000/07/03 12:39:21 jcs
144 // correct comment
145 //
146 // Revision 1.5 2000/06/13 14:34:25 oldi
147 // Excluded function call to ExtendTracks().
148 // Changed cout to gMessMgr->Message().
149 // Limits for pseudorapidity eta changed. They are now calculated from the
150 // z-position of the main vertex directly. This solves bug #574.
151 //
152 // Revision 1.4 2000/06/07 10:04:18 oldi
153 // Changed 0 pointers to NULL pointers.
154 // Constructor changed. Setting of mNumRowSegment is not possible anymore to
155 // avoid that it is set to a value != 20. This was necessary to avoid the
156 // exit(-20) function call.
157 // New funtion CalcChiSquared(StFtpcTrack *track, StFtpcConfMapPoint *point,
158 // Double_t *chi2) introduced.
159 // StraightLineFit(StFtpcTrack *track, Double_t *a, Int_t n) now calculates the
160 // chi squared. Therefore the pointer to 'a' has to be a pointer to an 6-dim array
161 // (a[6]). Also the code of this function was modified to take care of the first
162 // point which is not allowed to be used in the circle fit (to avoid infinities).
163 // Introduction of new function CreateTrack() to cleanup ClusterLoop().
164 // Introduction of new function ExtendTracks() and TrackExtension() to get cleaner
165 // tracks.
166 // Cleanup.
167 //
168 // Revision 1.3 2000/05/12 12:59:12 oldi
169 // removed delete operator for mSegment in StFtpcConfMapper (mSegment was deleted twice),
170 // add two new constructors for StFtpcTracker to be able to refit already existing tracks,
171 // minor cosmetics
172 //
173 // Revision 1.2 2000/05/11 15:14:41 oldi
174 // Changed class names *Hit.* due to already existing class StFtpcHit.cxx in StEvent
175 //
176 // Revision 1.1 2000/05/10 13:39:09 oldi
177 // Initial version of StFtpcTrackMaker
178 //
179 
180 //----------Author: Markus D. Oldenburg
181 //----------Last Modified: 10.11.2000
182 //----------Copyright: &copy MDO Production 1999
183 
184 #include "StFtpcConfMapper.hh"
185 #include "StFtpcTrackingParams.hh"
186 #include "StFtpcConfMapPoint.hh"
187 #include "StFtpcTrack.hh"
188 #include "StFormulary.hh"
189 #include "MIntArray.h"
190 
191 #include "TMath.h"
192 #include "TBenchmark.h"
193 #include "TCanvas.h"
194 #include "TH1.h"
195 #include "TH2.h"
196 #include "TH3.h"
197 #include "TPolyMarker.h"
198 #include "TFile.h"
199 #include "TLine.h"
200 #include "TGraph.h"
201 #include "TMarker.h"
202 #include "TPolyLine3D.h"
203 #include "TPolyMarker3D.h"
204 #include "TTUBE.h"
205 #include "TBRIK.h"
206 
207 #include "StMessMgr.h"
208 
210 // //
211 // StFtpcConfMapper class - tracking class to do tracking with conformal mapping. //
212 // //
214 
215 
216 ClassImp(StFtpcConfMapper)
217 
218 
220 {
221  // Default constructor.
222 
223  mLaser = (Bool_t)kFALSE;
224 
225  mVolume = 0;
226 
227  mVertexConstraint = (Bool_t)kTRUE;
228 }
229 
230 
231 StFtpcConfMapper::StFtpcConfMapper(TObjArray *inputHits, MIntArray *good_hits, StFtpcVertex *vertex, Bool_t bench,
232  Int_t phi_segments, Int_t eta_segments)
233  : StFtpcTracker(inputHits, vertex, bench)
234 {
235  // Constructor to fill in evaluated hits.
236 
237  if (mBench) {
238  mBench->Start("init");
239  }
240 
241  mLaser = (Bool_t)kFALSE;
242 
243  mNumRowSegment = StFtpcTrackingParams::Instance()->RowSegments();
244  mNumPhiSegment = StFtpcTrackingParams::Instance()->PhiSegments();
245  mNumEtaSegment = StFtpcTrackingParams::Instance()->EtaSegments();
246  mBounds = mNumRowSegment * mNumPhiSegment * mNumEtaSegment;
247  mMaxFtpcRow = StFtpcTrackingParams::Instance()->NumberOfPadRowsPerSide();
248 
249  CalcEtaMinMax();
250 
251  mMainVertexTracks = 0;
252  mMergedTracks = 0;
253  mMergedTracklets = 0;
254  mMergedSplits = 0;
255  mExtendedTracks = 0;
256  mDiffHits = 0;
257  mDiffHitsStill = 0;
258  mLengthFitNaN = 0;
259 
260  Int_t n_clusters = inputHits->GetEntriesFast(); // number of clusters
261  Int_t n_goodclusters = good_hits->CountAppearance(1); // number of good clusters
262  mClustersUnused = n_goodclusters;
263  mBadClusters = 0;
264 
265 //VP Dangerous play. mHit already set in base class and here reseting again
266  mHit = new TObjArray(n_goodclusters); // create TObjArray
267  mHitsCreated = (Bool_t)kTRUE;
268  Int_t good_hits_counted = 0;
269 
270  {for (Int_t i = 0; i < n_clusters; i++) {
271  StFtpcConfMapPoint *point = (StFtpcConfMapPoint*)inputHits->At(i);
272 
273  if (good_hits->At(i) == 1) {
274  StFtpcConfMapPoint *hit = new StFtpcConfMapPoint(point, mVertex);
275  mHit->AddAt(hit, good_hits_counted);
276  hit->SetHitNumber(i);
277  if (hit->GetUnusableForTrackingFlag()) mBadClusters++;
278  good_hits_counted++;
279  }
280  }}
281 
282  mVolume = new TObjArray(mBounds+1); // create ObjArray for volume cells (of size mBounds + one cell for garbage)
283 
284  {for (Int_t i = 0; i <= mBounds; i++) {
285  mVolume->AddAt(new TObjArray(0), i); // Fill ObjArray with empty ObjArrays
286  }}
287 
288  StFtpcConfMapPoint *hit;
289 
290  {for (Int_t i = 0; i < mHit->GetEntriesFast(); i++) {
291  hit = (StFtpcConfMapPoint *)mHit->At(i);
292  hit->Setup(mVertex);
293  ((TObjArray *)mVolume->At(GetSegm(hit)))->AddLast(hit);
294  }}
295 
296  if (mBench) {
297  mBench->Stop("init");
298  LOG_DEBUG << "Setup finished (" << mBench->GetCpuTime("init") << " s)." << endm;
299  mTime += mBench->GetCpuTime("init");
300  }
301 }
302 
303 
304 StFtpcConfMapper::StFtpcConfMapper(TObjArray *hits, StFtpcVertex *vertex, Bool_t bench,
305  Int_t phi_segments, Int_t eta_segments)
306  : StFtpcTracker(hits, vertex, bench)
307 {
308  // Constructor which needs a ClonesArray of hits as an input.
309 
310  if (mBench) {
311  mBench->Start("init");
312  }
313 
314  mLaser = (Bool_t)kFALSE;
315 
316  mNumRowSegment = StFtpcTrackingParams::Instance()->RowSegments();
317  mNumPhiSegment = StFtpcTrackingParams::Instance()->PhiSegments();
318  mNumEtaSegment = StFtpcTrackingParams::Instance()->EtaSegments();
319  mBounds = mNumRowSegment * mNumPhiSegment * mNumEtaSegment;
320  mMaxFtpcRow = StFtpcTrackingParams::Instance()->NumberOfPadRowsPerSide();
321 
322  CalcEtaMinMax();
323 
324  mMainVertexTracks = 0;
325  mMergedTracks = 0;
326  mMergedTracklets = 0;
327  mMergedSplits = 0;
328  mExtendedTracks = 0;
329  mDiffHits = 0;
330  mDiffHitsStill = 0;
331  mLengthFitNaN = 0;
332 
333  mClustersUnused = mHit->GetEntriesFast();
334  mBadClusters = 0;
335 
336  mVolume = new TObjArray(mBounds+1); // create ObjArray for volume cells (of size mBounds + one cell for garbage)
337 
338  {for (Int_t i = 0; i <= mBounds; i++) {
339  mVolume->AddAt(new TObjArray(0), i); // Fill ObjArray with empty ObjArrays
340  }}
341 
342  StFtpcConfMapPoint *hit;
343 
344  {for (Int_t i = 0; i < mHit->GetEntriesFast(); i++) {
345  hit = (StFtpcConfMapPoint *)mHit->At(i);
346  hit->SetHitNumber(i);
347  if (hit->GetUnusableForTrackingFlag()) mBadClusters++;
348  hit->Setup(mVertex);
349  ((TObjArray *)mVolume->At(GetSegm(hit)))->AddLast(hit);
350  }}
351 
352  if (mBench) {
353  mBench->Stop("init");
354  LOG_DEBUG << "Setup finished (" << mBench->GetCpuTime("init") << " s)." << endm;
355  mTime += mBench->GetCpuTime("init");
356  }
357 }
358 
359 
360 StFtpcConfMapper::~StFtpcConfMapper()
361 {
362  // Destructor.
363 
364  if (mVolume) {
365  mVolume->Delete();
366  delete mVolume; mVolume=0;
367  }
368 }
369 
370 
371 void StFtpcConfMapper::CompleteTrack(StFtpcTrack *track)
372 {
373  // Do everything which is still left to be done for the given track.
374  // Updates the statistics.
375 
376  track->SetPointDependencies();
377  track->ComesFromMainVertex(mVertexConstraint);
378  track->CalculateNMax();
379 
380  mClustersUnused -= track->GetNumberOfPoints();
381  if (mVertexConstraint) mMainVertexTracks++;
382 
383  return;
384 }
385 
386 
387 void StFtpcConfMapper::MainVertexTracking()
388 {
389  // Tracking with vertex constraint.
390 
391  if (mBench) {
392  mBench->Start("main_vertex");
393  }
394 
395  Settings("main_vertex");
396  Cuts("main_vertex");
397 
398  ClusterLoop();
399 
400  if (mBench) {
401  mBench->Stop("main_vertex");
402  LOG_DEBUG << "Main vertex tracking finished (" << mBench->GetCpuTime("main_vertex") << " s)." << endm;
403  mTime += mBench->GetCpuTime("main_vertex");
404 
405  mBench->Start("extend");
406  }
407 
408  ExtendTracks();
409 
410  if (mBench) {
411  mBench->Stop("extend");
412  LOG_DEBUG << "Track extension finished (" << mBench->GetCpuTime("extend") << " s)." << endm;
413  mTime += mBench->GetCpuTime("extend");
414 
415  mBench->Start("splits");
416  }
417 
418  HandleSplitTracks();
419 
420  if (mBench) {
421  mBench->Stop("splits");
422  LOG_DEBUG << "Split track merging finished (" << mBench->GetCpuTime("splits") << " s)." << endm;
423  mTime += mBench->GetCpuTime("splits");
424  }
425 
426  return;
427 }
428 
429 
430 void StFtpcConfMapper::FreeTracking()
431 {
432  // Tracking without vertex constraint.
433 
434  if (mBench) {
435  mBench->Start("non_vertex");
436  }
437 
438  Settings("non_vertex");
439  Cuts("non_vertex");
440 
441  ClusterLoop();
442 
443  if (mBench) {
444  mBench->Stop("non_vertex");
445  LOG_DEBUG << "Non vertex tracking finished (" << mBench->GetCpuTime("non_vertex") << " s)." << endm;
446  mTime += mBench->GetCpuTime("non_vertex");
447  }
448 
449  return;
450 }
451 
452 
453 void StFtpcConfMapper::LaserTracking()
454 {
455  // Tracking of straight tracks in arbitrary directions.
456  // Cuts and settings are optimized to find straight tracks.
457 
458  if (mBench) {
459  mBench->Start("laser");
460  }
461 
462  Settings("laser");
463  Cuts("laser");
464 
465  ClusterLoop();
466 
467  if (mBench) {
468  mBench->Stop("laser");
469  LOG_DEBUG << "Laser tracking finished (" << mBench->GetCpuTime("laser") << " s)." << endm;
470  mTime += mBench->GetCpuTime("laser");
471  mBench->GetCpuTime("nofield");
472 
473  mBench->Start("extend");
474  }
475 
476  ExtendTracks();
477 
478  if (mBench) {
479  mBench->Stop("extend");
480  LOG_DEBUG << "Track extension finished (" << mBench->GetCpuTime("extend") << " s)." << endm;
481  mTime += mBench->GetCpuTime("extend");
482 
483  mBench->Start("splits");
484  }
485 
486  HandleSplitTracks();
487 
488  if (mBench) {
489  mBench->Stop("splits");
490  LOG_DEBUG << "Split track merging finished (" << mBench->GetCpuTime("splits") << " s)." << endm;
491  mTime += mBench->GetCpuTime("splits");
492  }
493 
494  return;
495 }
496 
497 
498 void StFtpcConfMapper::NoFieldTracking()
499 {
500  // Tracking of straight tracks originating from main vertex.
501  // Cuts and settings are optimized to find straight tracks.
502 
503  if (mBench) {
504  mBench->Start("no_field");
505  }
506 
507  Settings("no_field");
508  Cuts("no_field");
509 
510  ClusterLoop();
511 
512  if (mBench) {
513  mBench->Stop("no_field");
514  LOG_DEBUG << "No field tracking finished (" << mBench->GetCpuTime("no_field") << " s)." << endm;
515  mBench->GetCpuTime("no_field");
516 
517  mBench->Start("extend");
518  }
519 
520  ExtendTracks();
521 
522  if (mBench) {
523  mBench->Stop("extend");
524  LOG_DEBUG << "Track extension finished (" << mBench->GetCpuTime("extend") << " s)." << endm;
525  mTime += mBench->GetCpuTime("extend");
526 
527  mBench->Start("splits");
528  }
529 
530  HandleSplitTracks();
531 
532  if (mBench) {
533  mBench->Stop("splits");
534  LOG_DEBUG << "Split track merging finished (" << mBench->GetCpuTime("splits") << " s)." << endm;
535  mTime += mBench->GetCpuTime("splits");
536  }
537 
538  return;
539 }
540 
541 void StFtpcConfMapper::Settings(TString method) {
542  // Sets settings with default values by given tracking method.
543 
544  Int_t method_id = 0;
545 
546  if (method == "main_vertex") {
547  method_id = 0;
548  }
549 
550  else if (method == "non_vertex") {
551  method_id = 1;
552  }
553 
554  else if (method == "no_field") {
555  method_id = 2;
556  }
557 
558  else if (method == "laser") {
559  method_id = 3;
560  }
561 
562  Settings(method_id);
563 
564  return;
565 }
566 
567 
568 void StFtpcConfMapper::Settings(Int_t method_id) {
569  // Sets settings with default values by given tracking method id.
570 
571  SetLaser(StFtpcTrackingParams::Instance()->Laser(method_id));
572  SetMaxDca(StFtpcTrackingParams::Instance()->MaxDca(method_id));
573  SetVertexConstraint(StFtpcTrackingParams::Instance()->VertexConstraint(method_id));
574  SetTrackletLength(StFtpcTrackingParams::Instance()->MaxTrackletLength(method_id));
575  SetRowScopeTracklet(StFtpcTrackingParams::Instance()->RowScopeTracklet(method_id));
576  SetMinPoints(StFtpcTrackingParams::Instance()->MinTrackLength(method_id));
577  SetRowScopeTrack(StFtpcTrackingParams::Instance()->RowScopeTrack(method_id));
578  SetPhiScope(StFtpcTrackingParams::Instance()->PhiScope(method_id));
579  SetEtaScope(StFtpcTrackingParams::Instance()->EtaScope(method_id));
580 
581  return;
582 }
583 
584 
585 void StFtpcConfMapper::Cuts(TString method) {
586  // Sets cuts with default values by given tracking method.
587 
588  Int_t method_id = 0;
589 
590  if (method == "main_vertex") {
591  method_id = 0;
592  }
593 
594  else if (method == "non_vertex") {
595  method_id = 1;
596  }
597 
598  else if (method == "no_field") {
599  method_id = 2;
600  }
601 
602  else if (method == "laser") {
603  method_id = 3;
604  }
605 
606  Cuts(method_id);
607 
608  return;
609 }
610 
611 
612 void StFtpcConfMapper::Cuts(Int_t method_id) {
613  // Sets cuts with default values by given tracking method id.
614 
615  SetMaxAngleTracklet(StFtpcTrackingParams::Instance()->MaxAngleTracklet(method_id));
616  SetMaxAngleTrack(StFtpcTrackingParams::Instance()->MaxAngleTrack(method_id));
617  SetMaxCircleDistTrack(StFtpcTrackingParams::Instance()->MaxCircleDist(method_id));
618  SetMaxLengthDistTrack(StFtpcTrackingParams::Instance()->MaxLengthDist(method_id));
619 
620  return;
621 }
622 
623 
624 void StFtpcConfMapper::Settings(Int_t trackletlength1, Int_t trackletlength2, Int_t tracklength1, Int_t tracklength2, Int_t rowscopetracklet1, Int_t rowscopetracklet2, Int_t rowscopetrack1, Int_t rowscopetrack2, Int_t phiscope1, Int_t phiscope2, Int_t etascope1, Int_t etascope2)
625 {
626  // Sets all settings of the tracker.
627  //
628  // This is the order of settings given to this function:
629  //
630  // - number of points to perform 'tracklet' search for main vertex tracks
631  // - number of points to perform 'tracklet' search for non vertex tracks
632  // These two mean no fitting but just looking for the nearest point
633  // in the direction to the main vertex (also for non vertex tracks!).
634  //
635  // - minimum number of points on a main vertex track
636  // - minimum number of points on a non vertex track
637  // These remove tracks and release points again already during the tracking
638  // if the tracker hasn't found enough hits on the track.
639  //
640  // - number of row segments to look in inward direction for main vertex tracklets
641  // - number of row segments to look in inward direction for non vertex tracklets
642  // - number of row segments to look in inward direction for main vertex tracks
643  // - number of row segments to look in inward direction for non vertex tracks
644  // These should be set to 1 for tracklets and to a value not less than 1 for track.
645  // Otherwise (if the value is set to 0) the tracker looks for the next point
646  // only in the same row. It should be set to a value higher than 1 if you
647  // want to be able to miss a cluster (because it is not there) but still
648  // to extend the track. This is (may be) not a good idea for tracklets.
649  // The tracker will extend the search in the next to the following padrow only
650  // if it has not found a cluster in the following row.
651  //
652  // - number of phi segments to look on both sides for main vertex tracks
653  // - number of phi segments to look on both sides for non vertex tracks
654  // These values have the same meaning as the values above for the row segemnts
655  // but now for the azimuth angle phi.
656  // The diffences are that the search is performed in any case over all ofthese
657  // segemnts (not only if in the nearest segments nothing was found) and that the
658  // reason to be able to set these values differential for main and non vertex tracks
659  // is the fact that non vertex tracks may be bent more than the high momentum
660  // main vertex tracks.
661  //
662  // - number of eta segments to look on both sides for main vertex tracks
663  // - number of eta segments to look on both sides for non vertex tracks
664  // Same as for the phi segments (above but now in eta (pseudorapidity) space).
665 
666  SetTrackletLength(trackletlength1, trackletlength2);
667  SetRowScopeTracklet(rowscopetracklet1, rowscopetracklet2);
668  SetRowScopeTrack(rowscopetrack1, rowscopetrack2);
669  SetPhiScope(phiscope1, phiscope2);
670  SetEtaScope(etascope1, etascope2);
671  SetMinPoints(tracklength1, tracklength2);
672 }
673 
674 
675 void StFtpcConfMapper::Cuts(Double_t maxangletracklet1, Double_t maxangletracklet2, Double_t maxangletrack1, Double_t maxangletrack2, Double_t maxcircletrack1, Double_t maxcircletrack2, Double_t maxlengthtrack1, Double_t maxlengthtrack2)
676 {
677  // Sets all cuts for the tracking.
678  //
679  // This is the order of settings given to this function:
680  //
681  // - maximum angle of main vertex tracklets
682  // - maximum angle of non vertex tracklets
683  // - maximum angle of main vertex tracks
684  // - maximum angle of non vertex tracks
685  // - maximal distance from circle fit for main vertex tracks
686  // - maximal distance from circle fit for non vertex tracks
687  // - maximal distance from length fit for main vertex tracks
688  // - maximal distance from length fit for non vertex tracks
689 
690  SetMaxAngleTracklet(maxangletracklet1, maxangletracklet2);
691  SetMaxAngleTrack(maxangletrack1, maxangletrack2);
692  SetMaxCircleDistTrack(maxcircletrack1, maxcircletrack2);
693  SetMaxLengthDistTrack(maxlengthtrack1, maxlengthtrack2);
694 }
695 
696 
697 void StFtpcConfMapper::LoopUpdate(Int_t *sub_row_segm, Bool_t backward)
698 {
699  // Increments or decrements *sub_row_segm, depending on backward.
700 
701  if (backward) {
702  (*sub_row_segm)--;
703  }
704 
705  else { // forward
706  (*sub_row_segm)++;
707  }
708 }
709 
710 
711 Bool_t const StFtpcConfMapper::TestExpression(Int_t sub_row_segm, Int_t end_row, Bool_t backward)
712 {
713  // Tests if loop should be continued or not.
714 
715  if (backward) {
716 
717  if (sub_row_segm >= end_row) {
718  return (Bool_t)kTRUE;
719  }
720 
721  else {
722  return (Bool_t) kFALSE;
723  }
724  }
725 
726  else { // forward
727 
728  if (sub_row_segm <= end_row) {
729  return (Bool_t)kTRUE;
730  }
731 
732  else {
733  return (Bool_t) kFALSE;
734  }
735  }
736 
737  return (Bool_t) kFALSE;
738 }
739 
740 
741 Double_t const StFtpcConfMapper::CalcDistance(const StFtpcConfMapPoint *hit1, const StFtpcConfMapPoint *hit2, Bool_t laser)
742 {
743  // Returns the distance of two given clusters. The distance in this respect (conformal mapping)
744  // is defined in the paper "A Fast track pattern recognition" by Pablo Yepes, NIM A 380 (1996) 585-585.
745  // If laser == true I expect straight tracks and the returned value will be the ordinary 3d distance.
746 
747  if (!laser) {
748 
749  Double_t phi_diff = TMath::Abs( hit1->GetPhi() - hit2->GetPhi() );
750  if (phi_diff > TMath::Pi()) phi_diff = 2*TMath::Pi() - phi_diff;
751 
752  return TMath::Abs(hit1->GetPadRow() - hit2->GetPadRow()) * (phi_diff + TMath::Abs( hit1->GetEta() - hit2->GetEta() ));
753  }
754 
755  else {
756  return TMath::Sqrt(TMath::Power(hit1->GetX() - hit2->GetX(), 2.) + TMath::Power(hit1->GetY() - hit2->GetY(), 2.) + TMath::Power(hit1->GetZ() - hit2->GetZ(), 2.));
757  }
758 }
759 
760 
761 Double_t const StFtpcConfMapper::CalcDistance(const StFtpcConfMapPoint *hit, Double_t *coeff)
762 {
763  // Returns the distance of a point to a straight line.
764  // The point is given by the to conformal coordinates of a cluster and the
765  // straight line is given by its to coefficients: y = coeff[0]*x + coeff[1].
766 
767  Double_t x = (coeff[0] / (1 + coeff[0]*coeff[0])) * (1/coeff[0] * hit->GetXprime() + hit->GetYprime() - coeff[1]);
768 
769  return TMath::Sqrt(TMath::Power(x - hit->GetXprime(), 2) + TMath::Power(coeff[0]*x + coeff[1] - hit->GetYprime(), 2));
770 }
771 
772 
773 Bool_t const StFtpcConfMapper::VerifyCuts(const StFtpcConfMapPoint *lasttrackhit, const StFtpcConfMapPoint *newhit, Bool_t backward)
774 {
775  // Returns true if circle, length, and angle cut holds.
776 
777  if (newhit->GetCircleDist() < mMaxCircleDist[mVertexConstraint] &&
778  newhit->GetLengthDist() < mMaxLengthDist[mVertexConstraint] &&
779  TrackAngle(lasttrackhit, newhit, backward) < mMaxAngleTrack[mVertexConstraint]) {
780  return kTRUE;
781  }
782 
783  else {
784  return kFALSE;
785  }
786 }
787 
788 
789 Double_t const StFtpcConfMapper::TrackAngle(const StFtpcPoint *lasthitoftrack, const StFtpcPoint *hit, Bool_t backward = (Bool_t)kTRUE)
790 {
791  // Returns the 'angle' between the last two points on the track (of which the last point is
792  // given as input) and the second given point.
793 
794  Double_t x1[3];
795  Double_t x2[3];
796 
797  StFtpcTrack *track = lasthitoftrack->GetTrack(mTrack);
798  TObjArray *hits = track->GetHits();
799  Int_t n = track->GetNumberOfPoints();
800 
801  if (n<2) {
802  LOG_ERROR << "StFtpcConfMapper::TrackAngle(StFtpcPoint *lasthitoftrack, StFtpcPoint *hit)" << endm;
803  LOG_ERROR << " - Call this function only if you are sure to have at least two points on the track already!" << endm;
804 
805  return 0.;
806  }
807 
808  if (backward) {
809  x1[0] = ((StFtpcPoint *)hits->At(n-1))->GetX() - ((StFtpcPoint *)hits->At(n-2))->GetX();
810  x1[1] = ((StFtpcPoint *)hits->At(n-1))->GetY() - ((StFtpcPoint *)hits->At(n-2))->GetY();
811  x1[2] = ((StFtpcPoint *)hits->At(n-1))->GetZ() - ((StFtpcPoint *)hits->At(n-2))->GetZ();
812 
813  x2[0] = hit->GetX() - ((StFtpcPoint *)hits->At(n-1))->GetX();
814  x2[1] = hit->GetY() - ((StFtpcPoint *)hits->At(n-1))->GetY();
815  x2[2] = hit->GetZ() - ((StFtpcPoint *)hits->At(n-1))->GetZ();
816  }
817 
818  else { // forward
819  x1[0] = ((StFtpcPoint *)hits->At(0))->GetX() - ((StFtpcPoint *)hits->At(1))->GetX();
820  x1[1] = ((StFtpcPoint *)hits->At(0))->GetY() - ((StFtpcPoint *)hits->At(1))->GetY();
821  x1[2] = ((StFtpcPoint *)hits->At(0))->GetZ() - ((StFtpcPoint *)hits->At(1))->GetZ();
822 
823  x2[0] = hit->GetX() - ((StFtpcPoint *)hits->At(0))->GetX();
824  x2[1] = hit->GetY() - ((StFtpcPoint *)hits->At(0))->GetY();
825  x2[2] = hit->GetZ() - ((StFtpcPoint *)hits->At(0))->GetZ();
826  }
827 
828  return StFormulary::Angle(x1, x2, 3);
829 }
830 
831 
832 Double_t const StFtpcConfMapper::TrackletAngle(StFtpcTrack *track, Int_t n)
833 {
834  // Returns the angle 'between' the last three points (started at point number n) on this track.
835 
836  Double_t x1[3];
837  Double_t x2[3];
838 
839  TObjArray *hits = track->GetHits();
840  if (n > track->GetNumberOfPoints() || n == 0) {
841  n = track->GetNumberOfPoints();
842  }
843 
844  if (n<3) {
845  LOG_ERROR << "StFtpcConfMapper::TrackletAngle(StFtpcTrack *track, Int_t n)" << endm;
846  LOG_ERROR << " - Call this function only if you are sure to have at least three points on this track already!" << endm;
847 
848  return 0.;
849  }
850 
851  x1[0] = ((StFtpcPoint *)hits->At(n-2))->GetX() - ((StFtpcPoint *)hits->At(n-3))->GetX();
852  x1[1] = ((StFtpcPoint *)hits->At(n-2))->GetY() - ((StFtpcPoint *)hits->At(n-3))->GetY();
853  x1[2] = ((StFtpcPoint *)hits->At(n-2))->GetZ() - ((StFtpcPoint *)hits->At(n-3))->GetZ();
854 
855  x2[0] = ((StFtpcPoint *)hits->At(n-1))->GetX() - ((StFtpcPoint *)hits->At(n-2))->GetX();
856  x2[1] = ((StFtpcPoint *)hits->At(n-1))->GetY() - ((StFtpcPoint *)hits->At(n-2))->GetY();
857  x2[2] = ((StFtpcPoint *)hits->At(n-1))->GetZ() - ((StFtpcPoint *)hits->At(n-2))->GetZ();
858 
859  return StFormulary::Angle(x1, x2, 3);
860 }
861 
862 
863 void StFtpcConfMapper::CalcChiSquared(StFtpcTrack *track, StFtpcConfMapPoint *point, Double_t *chi2)
864 {
865  // Calculates the chi squared for the circle and the length fit for a given point
866  // as an extension for the also given track.
867 
868  TObjArray *trackpoint = track->GetHits();
869  StFtpcConfMapPoint *last_point = (StFtpcConfMapPoint *)trackpoint->Last();
870 
871  if (!mVertexConstraint) {
872  point->SetAllCoord(last_point);
873  }
874 
875  trackpoint->AddLast(point);
876 
877  Double_t a[6];
878  StraightLineFit(track, a);
879  chi2[0] = a[4];
880  chi2[1] = a[5];
881 
882  trackpoint->Remove(point);
883 
884  return;
885 }
886 
887 /*
888 Double_t StFtpcConfMapper::CalcLength(StFtpcTrack *track, StFtpcTrackPoint* point)
889 {
890 
891  Double_t s;
892 
893  asin_arg = (point->GetY() - track->GetCenterY()) / track->GetRadius();
894 
895  // The following lines were inserted because ~1% of all tracks produce arguments of arcsin
896  // which are above the |1| limit. But they were differing only in the 5th digit after the point.
897  if (TMath::Abs(asin_arg) > 1.) {
898  asin_arg = (asin_arg >= 0) ? +1. : -1.;
899  //mLengthFitNaN++;
900  }
901 
902  if (point->GetX() >= track->GetCenterX() && point->GetY() > track->GetCenterY()) {
903  angle = TMath::ASin(asin_arg);
904  }
905 
906  else if (point->GetX() < track->GetCenterX() && point->GetY() >= track->GetCenterY()) {
907  angle = TMath::ASin(-asin_arg) + TMath::Pi();
908  }
909 
910  else if (point->GetX() <= track->GetCenterX() && point->GetY() < track->GetCenterY()) {
911  angle = TMath::ASin(-asin_arg) + TMath::Pi();
912  }
913 
914  else if (point->GetX() > track->GetCenterX() && point->GetY() <= track->GetCenterY()) {
915  angle = TMath::ASin(asin_arg) + 2 * TMath::Pi();
916  }
917 
918  s = (angle - track->GetAlpha0()) *
919  TMath::Sqrt(TMath::Power(track->GetRadius(), 2.) + TMath::Power(point->GetZv() , 2.));
920 
921  return s;
922 }
923 
924 
925 void StFtcpConfMapper::CalcAngle(StFtpcTrack *track, StFtpcTrackPoint* point)
926 {
927  // Calculates the the angle of the given point of a track with respect to the vertex of this point.
928 
929  Double_t asin_arg;
930 
931  asin_arg = (point->GetY() - track->GetCenterY()) / track->GetRadius();
932 
933  // The following lines were inserted because ~1% of all tracks produce arguments of arcsin
934  // which are above the |1| limit. But they were differing only in the 5th digit after the point.
935  if (TMath::Abs(asin_arg) > 1.) {
936  asin_arg = (asin_arg >= 0) ? +1. : -1.;
937  //mLengthFitNaN++;
938  }
939 
940  if (point->GetX() >= track->GetCenterX() && point->GetY() > track->GetCenterY()) {
941  point->SetPhiv(TMath::ASin(asin_arg));
942  }
943 
944  else if (point->GetX() < track->GetCenterX() && point->GetY() >= track->GetCenterY()) {
945  point->SetPhiv(-TMath::ASin(asin_arg) + TMath::Pi());
946  }
947 
948  else if (point->GetX() <= track->GetCenterX() && point->GetY() < track->GetCenterY()) {
949  point->SetPhiv(TMath::ASin(-asin_arg) + TMath::Pi());
950  }
951 
952  else if (point->GetX() > track->GetCenterX() && point->GetY() <= track->GetCenterY()) {
953  point->SetPhiv(-TMath::ASin(-asin_arg) + 2 * TMath::Pi());
954  }
955 
956  return;
957 }
958 */
959 
960 void StFtpcConfMapper::StraightLineFit(StFtpcTrack *track, Double_t *a, Int_t n)
961 {
962  // Calculates two straight line fits with the given clusters
963  //
964  // The first calculation is performed in the conformal mapping space (Xprime, Yprime):
965  // Yprime(Xprime) = a[0]*Xprime + a[1].
966  // The second calculates the fit for length s vs. z:
967  // s(z) = a[2]*z +a[3].
968  // a[4] is the chi squared per degrees of freedom for the circle fit.
969  // a[5] is the chi squared per degrees of freedom for the length fit.
970 
971  TObjArray *trackpoints = track->GetHits();
972 
973  if (n>0) {
974 
975  if (n > trackpoints->GetEntriesFast()) {
976  n = trackpoints->GetEntriesFast();
977  }
978  }
979 
980  else {
981  n = trackpoints->GetEntriesFast();
982  }
983 
984  Double_t L11 = 0.;
985  Double_t L12 = 0.;
986  Double_t L22 = 0.;
987  Double_t g1 = 0.;
988  Double_t g2 = 0.;
989 
990  Double_t asin_arg, asin_arg2;
991  Int_t start_counter = 0;
992 
993  if (!mVertexConstraint) {
994  start_counter = 1;
995  }
996 
997  // Circle Fit
998  StFtpcConfMapPoint *trackpoint = 0;
999 
1000  {for (Int_t i = start_counter; i < n; i++) {
1001  trackpoint = (StFtpcConfMapPoint *)trackpoints->At(i);
1002 
1003  L11 += 1./* / (trackpoint->GetYprimeerr() * trackpoint->GetYprimeerr())*/;
1004  L12 += trackpoint->GetXprime()/* / (trackpoint->GetYprimeerr() * trackpoint->GetYprimeerr())*/;
1005  L22 += (trackpoint->GetXprime() * trackpoint->GetXprime())/* / (trackpoint->GetYprimeerr() * trackpoint->GetYprimeerr())*/;
1006  g1 += trackpoint->GetYprime()/* / (trackpoint->GetYprimeerr() * trackpoint->GetYprimeerr())*/;
1007  g2 += (trackpoint->GetXprime() * trackpoint->GetYprime())/* / (trackpoint->GetYprimeerr() * trackpoint->GetYprimeerr())*/;
1008  }}
1009 
1010  Double_t D = L11*L22 - L12*L12;
1011 
1012  a[0] = (g2*L11 - g1*L12)/D;
1013  a[1] = (g1*L22 - g2*L12)/D;
1014 
1015  // Set circle parameters
1016  track->SetCenterX(- a[0] / (2. * a[1]) + trackpoint->GetXt());
1017  track->SetCenterY(- 1. / (2. * a[1]) + trackpoint->GetYt());
1018  track->SetRadius(TMath::Sqrt(a[0]*a[0] + 1.) / (2. * TMath::Abs(a[1])));
1019  track->CalcAndSetAlpha0();
1020 
1021  //LOG_DEBUG << track->GetRadius() << " " << track->GetAlpha0() << endm;
1022 
1023  // Tracklength Fit
1024  Double_t s = 0.;
1025  Double_t angle = 0.;
1026  Double_t angle_diff = 0.;
1027 
1028  // Set variables again
1029  // first track point is main vertex or first track point (both at (0, 0, 0) [shifted coordinates])
1030  L11 = 1.;
1031  L12 = L22 = g1 = g2 = 0.;
1032 
1033  {for (Int_t i = start_counter; i < n; i++ ) {
1034 
1035  trackpoint= (StFtpcConfMapPoint *)trackpoints->At(i);
1036 
1037  // calculate angle
1038 
1039  // The following lines were inserted because ~1% of all tracks produce arguments of arcsin
1040  // which are above the |1| limit. But they were differing only in the 5th digit after the point.
1041  asin_arg = StFormulary::CheckASinArg(asin_arg2 = (trackpoint->GetY() - track->GetCenterY()) / track->GetRadius());
1042  if (asin_arg != asin_arg2) mLengthFitNaN++;
1043 
1044  if (trackpoint->GetX() >= track->GetCenterX() && trackpoint->GetY() > track->GetCenterY()) {
1045  angle = TMath::ASin(asin_arg);
1046  }
1047 
1048  else if (trackpoint->GetX() < track->GetCenterX() && trackpoint->GetY() >= track->GetCenterY()) {
1049  angle = TMath::ASin(-asin_arg) + TMath::Pi();
1050  }
1051 
1052  else if (trackpoint->GetX() <= track->GetCenterX() && trackpoint->GetY() < track->GetCenterY()) {
1053  angle = TMath::ASin(-asin_arg) + TMath::Pi();
1054  }
1055 
1056  else if (trackpoint->GetX() > track->GetCenterX() && trackpoint->GetY() <= track->GetCenterY()) {
1057  angle = TMath::ASin(asin_arg) + 2 * TMath::Pi();
1058  }
1059 
1060  angle_diff = angle - track->GetAlpha0();
1061 
1062  if (TMath::Abs(angle_diff) > TMath::Pi()) {
1063  if (angle_diff > 0.) {
1064  angle_diff -= 2*TMath::Pi();
1065  }
1066 
1067  else {
1068  angle_diff += 2*TMath::Pi();
1069  }
1070  }
1071 
1072  s = TMath::Sqrt(TMath::Power(track->GetRadius() * angle_diff, 2.)
1073  + TMath::Power(trackpoint->GetZv() , 2.));
1074 
1075  //LOG_DEBUG << angle << " " << angle - track->GetAlpha0() << " " << s << endm;
1076 
1077  L11 += 1;
1078  L12 += trackpoint->GetZv();
1079  L22 += (trackpoint->GetZv() * trackpoint->GetZv());
1080  g1 += s;
1081  g2 += (s * trackpoint->GetZv());
1082  }}
1083 
1084  D = L11*L22 - L12*L12;
1085 
1086  //LOG_DEBUG << endm;
1087 
1088  a[2] = (g2*L11 - g1*L12)/D;
1089  a[3] = (g1*L22 - g2*L12)/D;
1090 
1091  // Set chi squared
1092  Double_t chi2circle = 0.;
1093  Double_t chi2length = 0.;
1094 
1095  {for (Int_t i = start_counter; i < n; i++) {
1096  trackpoint = (StFtpcConfMapPoint *)trackpoints->At(i);
1097  asin_arg = StFormulary::CheckASinArg((trackpoint->GetYv() - track->GetCenterY()) / track->GetRadius());
1098 
1099 //VP s = TMath::Sqrt(TMath::Power(track->GetRadius() * angle_diff, 2.)
1100 //VP + TMath::Power(trackpoint->GetZv() , 2.));
1101 //VP
1102 //VP chi2circle += TMath::Power(trackpoint->GetYprime() - a[0]*trackpoint->GetXprime()+a[1], 2.) / (a[0]*trackpoint->GetXprime()+a[1]);
1103 //VP chi2length += TMath::Power(s - a[2]*trackpoint->GetZv()+a[3], 2.) / (a[2]*trackpoint->GetZv()+a[3]);
1104  double x = trackpoint->GetXprime();
1105  double y = trackpoint->GetYprime();
1106  double z = trackpoint->GetZv();
1107  double yf = a[0]*x + a[1];
1108  double sf = a[2]*z + a[3];
1109  chi2circle += TMath::Power((y - yf)/yf,2);
1110  chi2length += TMath::Power((s - sf)/sf,2);
1111  }}
1112 
1113 
1114  if (mVertexConstraint) {
1115  n++;
1116  }
1117 
1118 //VP a[4] = chi2circle/(n-3.);//-3 is wrong. only 2 values were fitted a[0] & a[1]
1119  a[4] = chi2circle/(n-2.);
1120  a[5] = chi2length/(n-2.);
1121 
1122  //track->SetChi2Circle(chi2circle/(n-2.));
1123  //track->SetChi2Length(chi2length/(n-2.));
1124 
1125  return;
1126 }
1127 
1128 
1129 void StFtpcConfMapper::HandleSplitTracks() {
1130  // Handle split tracks with default values.
1131 
1132  HandleSplitTracks(StFtpcTrackingParams::Instance()->MaxDist(),
1133  StFtpcTrackingParams::Instance()->MinPointRatio(),
1134  StFtpcTrackingParams::Instance()->MaxPointRatio());
1135 
1136  return;
1137 }
1138 
1139 
1140 
1141 void StFtpcConfMapper::HandleSplitTracks(Double_t max_dist, Double_t ratio_min, Double_t ratio_max)
1142 {
1143  // Looks for split tracks and passes them to the track merger.
1144  Double_t ratio;
1145  Double_t dist;
1146 
1147  Double_t r1;
1148  Double_t r2;
1149  Double_t R1;
1150  Double_t R2;
1151  Double_t phi1;
1152  Double_t phi2;
1153  Double_t Phi1;
1154  Double_t Phi2;
1155 
1156  Int_t first_split = -1;
1157 
1158  StFtpcTrack *t1;
1159  StFtpcTrack *t2;
1160 
1161  Int_t entries = mTrack->GetEntriesFast();
1162 
1163  for (Int_t i = 0; i < entries; i++) {
1164  t1 = (StFtpcTrack *)mTrack->At(i);
1165 
1166  if (!t1) {
1167  // track was removed before (already merged)
1168  continue;
1169  }
1170 
1171  for (Int_t j = i+1; j < entries; j++) {
1172  t2 = (StFtpcTrack *)mTrack->At(j);
1173 
1174  if (!t2) {
1175  // track was removed before (already merged)
1176  continue;
1177  }
1178 
1179  if (!t1) {
1180  // track t1 was removed before (already merged) - leave inner loop
1181  break;
1182  }
1183 
1184  if (!(t1->GetRowsWithPoints() & t2->GetRowsWithPoints()) &&
1185  (t1->GetHemisphere() == t2->GetHemisphere())) {
1186 
1187  r1 = t1->GetRLast();
1188  r2 = t2->GetRLast();
1189  phi1 = t1->GetAlphaLast();
1190  phi2 = t2->GetAlphaLast();
1191  R1 = t1->GetRFirst();
1192  R2 = t2->GetRFirst();
1193  Phi1 = t1->GetAlphaFirst();
1194  Phi2 = t2->GetAlphaFirst();
1195 
1196  dist = (TMath::Sqrt(r2*r2+r1*r1-2*r1*r2*(TMath::Cos(phi1)*TMath::Cos(phi2)+TMath::Sin(phi1)*TMath::Sin(phi2))) +
1197  TMath::Sqrt(R2*R2+R1*R1-2*R1*R2*(TMath::Cos(Phi1)*TMath::Cos(Phi2)+TMath::Sin(Phi1)*TMath::Sin(Phi2)))) / 2.;
1198  ratio = (Double_t)(t1->GetNumberOfPoints() + t2->GetNumberOfPoints()) / (Double_t)(t1->GetNMax() + t2->GetNMax());
1199 
1200  if (dist <= max_dist && ratio <= ratio_max && ratio >= ratio_min) {
1201 
1202  if (first_split == -1) {
1203  first_split = i;
1204  }
1205 
1206  MergeSplitTracks(t1, t2);
1207  t1 = t2 = 0;
1208  }
1209  }
1210  }
1211  }
1212 
1213  if (first_split != -1) {
1214  // adjust track numbers only if something has changed
1215  AdjustTrackNumbers(first_split);
1216  }
1217 
1218  return;
1219 }
1220 
1221 
1222 void StFtpcConfMapper::MergeSplitTracks(StFtpcTrack *t1, StFtpcTrack *t2)
1223 {
1224  // Merges two tracks which are split.
1225 
1226  Int_t new_track_number = mTrack->GetEntriesFast();
1227  Int_t num_t1_points = t1->GetNumberOfPoints();
1228  Int_t num_t2_points = t2->GetNumberOfPoints();
1229  TObjArray *t1_points = t1->GetHits();
1230  TObjArray *t2_points = t2->GetHits();
1231 
1232  StFtpcTrack *track = new StFtpcTrack(new_track_number);
1233  mTrack->AddAt(track, new_track_number);
1234 
1235  TObjArray *trackpoint = track->GetHits();
1236  trackpoint->Expand(mMaxFtpcRow);
1237 
1238  MIntArray *trackhitnumber = track->GetHitNumbers();
1239 
1240  {for (Int_t i = 0; i < mMaxFtpcRow; i++) {
1241 
1242  if (i < num_t1_points) {
1243  trackpoint->AddAt(t1_points->At(i), mMaxFtpcRow - 1 -(((StFtpcPoint *)t1_points->At(i))->GetPadRow()-1)%mMaxFtpcRow);
1244  }
1245 
1246  if (i < num_t2_points) {
1247  trackpoint->AddAt(t2_points->At(i), mMaxFtpcRow - 1 -(((StFtpcPoint *)t2_points->At(i))->GetPadRow()-1)%mMaxFtpcRow);
1248  }
1249  }}
1250 
1251  trackpoint->Compress();
1252  trackpoint->Expand(trackpoint->GetLast()+1);
1253 
1254  {for (Int_t i = 0; i < trackpoint->GetEntriesFast(); i++) {
1255  trackhitnumber->AddLast(((StFtpcPoint *)trackpoint->At(i))->GetHitNumber());
1256  }}
1257 
1258  track->SetProperties(kTRUE, new_track_number);
1259  track->ComesFromMainVertex(mVertexConstraint);
1260  track->CalculateNMax();
1261 
1262  if (mVertexConstraint) mMainVertexTracks++;
1263  if (t1->ComesFromMainVertex()) mMainVertexTracks--;
1264  if (t2->ComesFromMainVertex()) mMainVertexTracks--;
1265 
1266 
1267  mTrack->Remove(t1);
1268  mTrack->Remove(t2);
1269  delete t1;
1270  delete t2;
1271 
1272  mMergedSplits++;
1273 
1274  return;
1275 }
1276 
1277 
1278 void StFtpcConfMapper::ClusterLoop()
1279 {
1280  // This function loops over all clusters to do the tracking.
1281  // These clusters act as starting points for new tracks.
1282 
1283  Int_t row_segm;
1284  Int_t phi_segm;
1285  Int_t eta_segm;
1286  Int_t entries;
1287  Int_t hit_num;
1288 
1289  TObjArray *segment;
1290  StFtpcConfMapPoint *hit;
1291 
1292  // loop over two Ftpcs
1293  for (mFtpc = 1; mFtpc <= 2; mFtpc++) {
1294 
1295  // loop over the respective 10 RowSegments ("layers") per Ftpc
1296  // loop only so far to where you can still put a track in the remaining padrows (due to length)
1297  for (row_segm = mFtpc * mMaxFtpcRow - 1; row_segm >= (mFtpc-1) * mMaxFtpcRow + mMinPoints[mVertexConstraint] - 1; row_segm--) {
1298 
1299  // loop over phi segments
1300  for (Int_t phi_segm_counter = 0; phi_segm_counter < mNumPhiSegment; phi_segm_counter++) {
1301 
1302  // go over phi in two directions, one segment in each direction alternately
1303  if(phi_segm_counter%2) {
1304  phi_segm = mNumPhiSegment - phi_segm_counter - mNumPhiSegment%2;
1305  }
1306 
1307  else {
1308  phi_segm = phi_segm_counter;
1309  }
1310 
1311  // loop over eta segments
1312  for(eta_segm = 0; eta_segm < mNumEtaSegment; eta_segm++) {
1313 
1314  // loop over entries in one segment
1315  if ((entries = (segment = (TObjArray *)mVolume->At(GetSegm(row_segm, phi_segm, eta_segm)))->GetEntriesFast())) {
1316 
1317  for (hit_num = 0; hit_num < entries; hit_num++) {
1318  hit = (StFtpcConfMapPoint *)segment->At(hit_num);
1319 
1320  if (hit->IsUnusable()) { // start hit was used before or flagged as not to be used for tracking
1321  continue;
1322  }
1323 
1324  else { // start hit was not used before
1325  CreateTrack(hit);
1326  }
1327  }
1328  }
1329 
1330  else continue; // no entries in this segment
1331  }
1332  }
1333  }
1334  }
1335 
1336  return; // end of track finding
1337 }
1338 
1339 
1340 void StFtpcConfMapper::CreateTrack(StFtpcConfMapPoint *hit)
1341 {
1342  // This function takes as input a point in the Ftpc which acts as starting point for a new track.
1343  // It forms tracklets then extends them to tracks.
1344 
1345  Double_t *coeff = 0,coeffT[6];
1346  //Double_t chi2[2];
1347 
1348  StFtpcConfMapPoint *closest_hit = 0;
1349  StFtpcTrack *track = 0;
1350 
1351  Int_t point;
1352  Int_t tracks = GetNumberOfTracks();
1353  if (tracks >= mTrack->GetSize()) mTrack->Expand(mTrack->GetSize()+1000);
1354  track = new StFtpcTrack(tracks);
1355  mTrack->AddAt(track, tracks);
1356  TObjArray *trackpoint = track->GetHits();
1357 
1358  track->AddPoint(hit);
1359 
1360  // set conformal mapping coordinates if looking for non vertex tracks
1361  if (!mVertexConstraint) {
1362  hit->SetAllCoord(hit);
1363  }
1364 
1365  // create tracklets
1366  for (point = 1; point < mTrackletLength[mVertexConstraint]; point++) {
1367 
1368  if ((closest_hit = GetNextNeighbor(hit, coeff, (Bool_t)kTRUE))) {
1369  // closest_hit for hit exists
1370  track->AddPoint(closest_hit);
1371  hit = closest_hit;
1372 
1373 
1374  if (GetLaser() && track->GetNumberOfPoints() == mTrackletLength[mVertexConstraint]) {
1375  tracks++;
1376  CompleteTrack(track);
1377 
1378  return; // no additional track search
1379  }
1380  }
1381 
1382  else { // closest hit does not exist
1383 
1384  if (!GetLaser()) { // this is not a laser event
1385  RemoveTrack(track);
1386  return; // continue with next hit in segment
1387  }
1388 
1389  else { // this is a laser event
1390 
1391  if (track->GetNumberOfPoints() < mMinPoints[mVertexConstraint]) { // not enough points
1392  RemoveTrack(track);
1393  return; // continue with next hit in segment
1394  }
1395 
1396  else { // enough points
1397  tracks++;
1398  CompleteTrack(track);
1399 
1400  return; // no additional track search
1401  }
1402  }
1403  }
1404  }
1405 
1406  // tracklet is long enough to be extended to a track
1407  if (trackpoint->GetEntriesFast() == mTrackletLength[mVertexConstraint]) {
1408 
1409  if (TrackletAngle(track) > mMaxAngleTracklet[mVertexConstraint]) { // proof if the first points seem to be a beginning of a track
1410  RemoveTrack(track);
1411  }
1412 
1413  else { // good tracklet -> proceed
1414  tracks++;
1415 
1416  // create tracks
1417  for (point = mTrackletLength[mVertexConstraint]; point < mMaxFtpcRow; point++) {
1418 
1419  if (!coeff) coeff = coeffT;
1420  //Double_t chi_circ = track->GetChi2Circle();
1421  //Double_t chi_len = track->GetChi2Length();
1422 
1423  StraightLineFit(track, coeff);
1424  closest_hit = GetNextNeighbor((StFtpcConfMapPoint *)trackpoint->Last(), coeff, (Bool_t)kTRUE);
1425 
1426  if (closest_hit) {
1427  // add closest hit to track
1428  track->AddPoint(closest_hit);
1429  }
1430 
1431  else {
1432  // closest hit does not exist
1433  point = mMaxFtpcRow; // continue with next hit in segment
1434  }
1435  }
1436 
1437  // remove tracks with not enough points already now
1438  if (track->GetNumberOfPoints() < mMinPoints[mVertexConstraint]) {
1439  RemoveTrack(track);
1440  tracks--;
1441  }
1442 
1443  else {
1444  //LOG_DEBUG << coeff[2] << endm;
1445  CompleteTrack(track);
1446  }
1447 
1448  // cleanup
1449  coeff = 0;
1450  }
1451  }
1452 
1453  return;
1454 }
1455 
1456 
1457 StFtpcConfMapPoint *StFtpcConfMapper::GetNextNeighbor(StFtpcConfMapPoint *start_hit, Double_t *coeff = 0, Bool_t backward = (Bool_t)kTRUE)
1458 {
1459  // Returns the nearest cluster to a given start_hit.
1460 
1461  Double_t dist, closest_dist = 1.e7;
1462  Double_t closest_circle_dist = 1.e7;
1463  Double_t closest_length_dist = 1.e7;
1464 
1465  StFtpcConfMapPoint *hit = 0;
1466  StFtpcConfMapPoint *closest_hit = 0;
1467  StFtpcConfMapPoint *closest_circle_hit = 0;
1468  StFtpcConfMapPoint *closest_length_hit = 0;
1469 
1470  TObjArray *sub_segment;
1471  Int_t sub_entries;
1472 
1473  Int_t sub_row_segm;
1474  Int_t sub_phi_segm;
1475  Int_t sub_eta_segm;
1476  Int_t sub_hit_num;
1477 
1478  Int_t start_row;
1479  Int_t end_row;
1480 
1481  if (backward) {
1482  start_row = GetRowSegm(start_hit) - 1;
1483 
1484  if (coeff) {
1485  end_row = GetRowSegm(start_hit) - mRowScopeTrack[mVertexConstraint];
1486  }
1487 
1488  else {
1489  end_row = GetRowSegm(start_hit) - mRowScopeTracklet[mVertexConstraint];
1490  }
1491 
1492  while (end_row < (mFtpc-1) * mMaxFtpcRow) {
1493  end_row++;
1494  }
1495 
1496  if (start_row < end_row) return 0;
1497  }
1498 
1499  else { // forward
1500  start_row = GetRowSegm(start_hit) + 1;
1501 
1502  if (coeff) {
1503  end_row = GetRowSegm(start_hit) + mRowScopeTrack[mVertexConstraint];
1504  }
1505 
1506  else {
1507  end_row = GetRowSegm(start_hit) + mRowScopeTracklet[mVertexConstraint];
1508  }
1509 
1510  while (end_row > mFtpc * mMaxFtpcRow - 1) {
1511  end_row--;
1512  }
1513 
1514  if (start_row > end_row) return 0;
1515  }
1516 
1517  // loop over sub rows
1518  for (sub_row_segm = start_row; TestExpression(sub_row_segm, end_row, backward); LoopUpdate(&sub_row_segm, backward)) {
1519 
1520  // loop over sub phi segments
1521  for (Int_t i = -(mPhiScope[mVertexConstraint]); i <= mPhiScope[mVertexConstraint]; i++) {
1522  sub_phi_segm = GetPhiSegm(start_hit) + i; // neighboring phi segment
1523 
1524  if (sub_phi_segm < 0) { // find neighboring segment if #segment < 0
1525  sub_phi_segm += mNumPhiSegment;
1526  }
1527 
1528  else if (sub_phi_segm >= mNumPhiSegment) { // find neighboring segment if #segment > fNum_phi_segm
1529  sub_phi_segm -= mNumPhiSegment;
1530  }
1531 
1532  // loop over sub eta segments
1533  for (Int_t j = -(mEtaScope[mVertexConstraint]); j <= mEtaScope[mVertexConstraint]; j++) {
1534  sub_eta_segm = GetEtaSegm(start_hit) + j; // neighboring eta segment
1535 
1536  if (sub_eta_segm < 0 || sub_eta_segm >= mNumEtaSegment) {
1537  continue; // #segment exceeds bounds -> skip
1538  }
1539 
1540  // loop over entries in one sub segment
1541  if ((sub_entries = ((sub_segment = (TObjArray *)mVolume->At(GetSegm(sub_row_segm, sub_phi_segm, sub_eta_segm)))->GetEntriesFast()))) {
1542 
1543  for (sub_hit_num = 0; sub_hit_num < sub_entries; sub_hit_num++) {
1544  hit = (StFtpcConfMapPoint *)sub_segment->At(sub_hit_num);
1545 
1546  if ((hit = (StFtpcConfMapPoint *)sub_segment->At(sub_hit_num))->IsUsable()) {
1547  // hit was not used before and is flagged ok for tracking
1548 
1549  // set conformal mapping coordinates if looking for non vertex tracks
1550  if (!mVertexConstraint) {
1551  hit->SetAllCoord(start_hit);
1552  }
1553 
1554  if (coeff) { // track search - look for nearest neighbor to extrapolated track
1555 
1556  // test distance
1557  hit->SetDist(CalcDistance(hit, coeff+0), CalcDistance(hit, coeff+2));
1558 
1559  if (hit->GetCircleDist() < closest_circle_dist) {
1560  closest_circle_dist = hit->GetCircleDist();
1561  closest_circle_hit = hit;
1562  }
1563 
1564  if (hit->GetLengthDist() < closest_length_dist) {
1565  closest_length_dist = hit->GetLengthDist();
1566  closest_length_hit = hit;
1567  }
1568  }
1569 
1570  else {
1571  // tracklet search - just look for the nearest neighbor (distance defined by Pablo Jepes)
1572 
1573  // test distance
1574  if ((dist = CalcDistance(start_hit, hit, GetLaser())) < closest_dist) { // hit found that is closer than the hits before
1575 
1576  if (GetLaser() && ((StFtpcTrack*)start_hit->GetTrack(mTrack))->GetNumberOfPoints() > 1) { // laser tracking mode and already 2 hits on this track
1577 
1578  if (TrackAngle(start_hit, hit, backward) <= mMaxAngleTracklet[mVertexConstraint]) { // angle between last two hits and new within limits
1579 
1580  closest_dist = dist;
1581  closest_hit = hit;
1582  }
1583 
1584  else { // angle not in limits
1585  continue;
1586  }
1587  }
1588 
1589  else { // no laser mode or laser mode but less the two points on the track up to now
1590  closest_dist = dist;
1591  closest_hit = hit;
1592  }
1593  }
1594 
1595  else { // sub hit was farther away than a hit before
1596  continue;
1597  }
1598  }
1599  }
1600 
1601  else continue; // sub hit was used before
1602  }
1603  }
1604 
1605  else continue; // no sub hits
1606  }
1607  }
1608 
1609 
1610  if ((coeff && (closest_circle_hit || closest_length_hit)) || ((!coeff) && closest_hit)) {
1611 
1612  if ((start_row - sub_row_segm) >= 1) {
1613 
1614  if (coeff) {
1615  mMergedTracks++;
1616  }
1617 
1618  else {
1619  mMergedTracklets++;
1620  }
1621  }
1622 
1623  // found a hit in a sub layer - don't look in other sub layers
1624  break;
1625  }
1626 
1627  else {
1628  // didn't find a hit in this sub layer - try next sub layer
1629  continue;
1630  }
1631  }
1632 
1633  if (coeff) {
1634 
1635  if (closest_circle_hit && closest_length_hit) { // hits are not zero
1636 
1637  if (closest_circle_hit == closest_length_hit) { // both found hits are identical
1638 
1639  if (VerifyCuts(start_hit, closest_circle_hit, backward)) {
1640 
1641  // closest hit within limits found
1642  return closest_circle_hit;
1643  }
1644 
1645  else { // limits exceeded
1646  return 0;
1647  }
1648  }
1649 
1650  else { // found hits are different
1651  mDiffHits++;
1652 
1653  Bool_t cut_circle = VerifyCuts(start_hit, closest_circle_hit, backward);
1654  Bool_t cut_length = VerifyCuts(start_hit, closest_length_hit, backward);
1655 
1656  if (cut_circle && cut_length) { // both hits are within the limit
1657 
1658  if (GetDistanceFromFit(closest_circle_hit) < GetDistanceFromFit(closest_length_hit)) { // circle_hit is closer
1659  return closest_circle_hit;
1660  }
1661 
1662  else { // length_hit is closer
1663  return closest_length_hit;
1664  }
1665  }
1666 
1667  else if (!(cut_circle || cut_length)) { // both hits exceed limits
1668  return 0;
1669  }
1670 
1671  else if (cut_circle) { // closest_circle_hit is the only one within limits
1672  return closest_circle_hit;
1673  }
1674 
1675  else { // closest_length_hit is the only one within limits
1676  return closest_length_hit;
1677  }
1678  }
1679  }
1680 
1681  else { // no hits found
1682  return 0;
1683  }
1684  }
1685 
1686 
1687  else { // closest hit for tracklet found
1688 
1689  if (closest_hit) { // closest hit exists
1690  return closest_hit;
1691  }
1692 
1693  else { // hit does not exist
1694  return 0;
1695  }
1696  }
1697 }
1698 
1699 
1700 void StFtpcConfMapper::ExtendTracks()
1701 {
1702  // Loops over all found tracks and passes them to the part of the program where each track is tried to be extended.
1703 
1704  for (Int_t t = 0; t < mTrack->GetEntriesFast(); t++) {
1705 
1706  StFtpcTrack *track = (StFtpcTrack *)mTrack->At(t);
1707 
1708  if (track->GetHemisphere() == 1) {
1709  mFtpc = (Char_t)1;
1710  }
1711 
1712  else {
1713  mFtpc = (Char_t)2;
1714  }
1715 
1716  if (TrackExtension(track)) {
1717  mExtendedTracks++;
1718  }
1719  }
1720 
1721  return;
1722 }
1723 
1724 
1725 Bool_t StFtpcConfMapper::TrackExtension(StFtpcTrack *track)
1726 {
1727  // Trys to extend a given track.
1728 
1729  Int_t point;
1730  Int_t number_of_points = track->GetNumberOfPoints();
1731  Double_t *coeff = 0,coeffT[6];
1732 
1733  StFtpcConfMapPoint *closest_hit;
1734  StFtpcConfMapPoint *hit;
1735  TObjArray *trackpoint = track->GetHits();
1736 
1737  for (Int_t direction = 0; direction < 2; direction ++) {
1738 
1739  for (point = number_of_points; point < mMaxFtpcRow; point++) {
1740 
1741  if (direction == 1) { // backward
1742  hit = (StFtpcConfMapPoint *)trackpoint->Last();
1743  }
1744 
1745  else { // forward
1746  hit = (StFtpcConfMapPoint *)trackpoint->First();
1747  }
1748 
1749  Int_t padrow = hit->GetPadRow()%StFtpcTrackingParams::Instance()->NumberOfPadRowsPerSide();
1750 
1751  if ((padrow != 1 && direction == 1) ||
1752  (padrow != 0 && direction == 0)) {
1753 
1754  if (GetLaser()) { // Laser event
1755 
1756  if ((closest_hit = GetNextNeighbor(hit, coeff, (Bool_t)direction))) {
1757 
1758  // add closest hit to track
1759  if(direction == 1) { // backward
1760  track->AddPoint(closest_hit);
1761  }
1762 
1763  else { // forward
1764  track->AddForwardPoint(closest_hit);
1765  }
1766  }
1767  }
1768 
1769  else { // usual event
1770 
1771  if (!coeff) coeff = coeffT;
1772  StraightLineFit(track, coeff);
1773 
1774  if ((closest_hit = GetNextNeighbor(hit, coeff, (Bool_t)direction))) {
1775 
1776  // add closest hit to track
1777  if(direction == 1) { // backward
1778  track->AddPoint(closest_hit);
1779  }
1780 
1781  else { // forward
1782  track->AddForwardPoint(closest_hit);
1783  }
1784  }
1785 
1786  else {
1787  point = mMaxFtpcRow; // continue with next hit in segment
1788  }
1789  }
1790  }
1791  }
1792  }
1793 
1794  // cleanup
1795  coeff = 0;
1796 
1797  if (track->GetNumberOfPoints() - number_of_points) {
1798 
1799  track->SetPointDependencies();
1800  track->ComesFromMainVertex(mVertexConstraint);
1801  track->CalculateNMax();
1802 
1803  mClustersUnused -= (track->GetNumberOfPoints() - number_of_points);
1804 
1805  return (Bool_t)kTRUE;
1806  }
1807 
1808  else {
1809  return (Bool_t)kFALSE;
1810  }
1811 }
1812 
1813 
1814 void StFtpcConfMapper::CalcEtaMinMax()
1815 {
1816  // Calculates the min. and max. value of eta (pseudorapidity) with the given main vertex.
1817  // The FTPCs are placed in a distance to the point of origin between 162.75 and 256.45 cm.
1818  // Their inner radius is 8, the outer one 30 cm. This means they are seen from (0, 0, 0)
1819  // under an angle between 1.79 and 10.62 degrees which translates directly into max. eta and min. eta.
1820  // Due to the fact that the main vertex is shifted, the values of min./max. eta is calculated for each
1821  // event. To be save, 0.01 is substracted/added.
1822 
1823  Double_t vertex_rad_squared = TMath::Power(mVertex->GetRadius2(), 2);
1824  Double_t ftpc_inner_rad_squared = TMath::Power(StFtpcTrackingParams::Instance()->InnerRadius(), 2);
1825  Double_t ftpc_outer_rad_squared = TMath::Power(StFtpcTrackingParams::Instance()->OuterRadius(), 2);
1826 
1827  mEtaMin = -TMath::Log(TMath::Tan( TMath::ASin(TMath::Sqrt(ftpc_outer_rad_squared + vertex_rad_squared)/ (StFtpcTrackingParams::Instance()->PadRowPosZ(0) - TMath::Abs(mVertex->GetZ()) ) ) /2.) ) - 0.01;
1828  mEtaMax = -TMath::Log(TMath::Tan( TMath::ASin(TMath::Sqrt(ftpc_inner_rad_squared - vertex_rad_squared)/ (StFtpcTrackingParams::Instance()->PadRowPosZ(9) + TMath::Abs(mVertex->GetZ()) ) ) /2.) ) + 0.01;
1829 
1830  return;
1831 }
1832 
1833 
1834 void StFtpcConfMapper::TrackingInfo()
1835 {
1836  // Information about the tracking process.
1837 
1838  LOG_INFO << endm;
1839  LOG_INFO << "Tracking information" << endm;
1840  LOG_INFO << "--------------------" << endm;
1841 
1842  LOG_INFO << Form("%5d (%5d/%5d) tracks (main vertex/non vertex) found.",
1843  GetNumberOfTracks(), GetNumMainVertexTracks(), GetNumberOfTracks() - GetNumMainVertexTracks()) << endm;
1844 
1845  LOG_INFO << Form("%5d (%5d/%5d) clusters (used/unused, of which %d were flagged 'bad').",
1846  GetNumberOfClusters(), GetNumberOfClusters() - GetNumClustersUnused(),
1847  GetNumClustersUnused(), GetNumBadClusters()) << endm;
1848 
1849  LOG_INFO << Form(" %5d/%5d tracks/tracklets merged.",
1850  GetNumMergedTracks(), GetNumMergedTracklets()) << endm;
1851 
1852  LOG_INFO << Form("%18d times different hits for circle and length fit found.",
1853  GetNumDiffHits()) << endm;
1854 
1855  LOG_INFO << Form("%18d times argument of arcsin set to +/-1.",
1856  GetNumLengthFitNaN()) << endm;
1857 
1858  LOG_INFO << Form("%18d tracks extended.", GetNumExtendedTracks()) << endm;
1859 
1860  LOG_INFO << Form("%18d split tracks merged.", GetNumMergedSplits()) << endm;
1861 
1862  return;
1863 }
1864 
1865 
1866 void StFtpcConfMapper::CutInfo()
1867 {
1868  // Information about cuts.
1869 
1870  LOG_INFO << endm;
1871  LOG_INFO << "Cuts for main vertex constraint on / off" << endm;
1872  LOG_INFO << "----------------------------------------" << endm;
1873  LOG_INFO << Form("Max. angle between last three points of tracklets: %6f / %6f",mMaxAngleTracklet[1],mMaxAngleTracklet[0]) << endm;
1874  LOG_INFO << Form("Max. angle between last three points of tracks: %6f / %6f",mMaxAngleTrack[1],mMaxAngleTrack[0]) << endm;
1875  LOG_INFO << Form("Max. distance between circle fit and trackpoint: %6f / %6f",mMaxCircleDist[1],mMaxCircleDist[0]) << endm;
1876  LOG_INFO << Form("Max. distance between length fit and trackpoint: %6f / %6f",mMaxLengthDist[1],mMaxLengthDist[0])<< endm;
1877 
1878  return;
1879 }
1880 
1881 
1882 void StFtpcConfMapper::SettingInfo()
1883 {
1884  // Information about settings.
1885 
1886  LOG_INFO << endm;
1887  LOG_INFO << "Settings for main vertex constraint on / off" << endm;
1888  LOG_INFO << "--------------------------------------------" << endm;
1889  LOG_INFO << "Points required to create a tracklet: " << mTrackletLength[1] << " / " << mTrackletLength[0] << endm;
1890  LOG_INFO << "Points required for a track: " << mMinPoints[1] << " / " << mMinPoints[0] << endm;
1891  LOG_INFO << "Subsequent padrows to look for next tracklet point: " << mRowScopeTracklet[1] << " / " << mRowScopeTracklet[0] << endm;
1892  LOG_INFO << "Subsequent padrows to look for next track point: " << mRowScopeTrack[1] << " / " << mRowScopeTrack[0] << endm;
1893  LOG_INFO << "Adjacent phi segments to look for next point: " << mPhiScope[1] << " / " << mPhiScope[0] << endm;
1894  LOG_INFO << "Adjacent eta segments to look for next point: " << mEtaScope[1] << " / " << mEtaScope[0] << endm;
1895 
1896  return;
1897 }
1898 
1899 
1900 void StFtpcConfMapper::TwoCycleTracking()
1901 {
1902  // Tracking in 2 cycles:
1903  // 1st cycle: tracking with vertex constraint
1904  // 2nd cycle: without vertex constraint (of remaining clusters)Begin_Html<a name="settings"></a>End_Html
1905 
1906  MainVertexTracking();
1907  FreeTracking();
1908  return;
1909 }
1910 
1911 
1912 void StFtpcConfMapper::NonVertexTracking()
1913 {
1914  // "Shortcut" to FreeTracking().
1915 
1916  FreeTracking();
1917 }
1918 
1919 
1920 void StFtpcConfMapper::RemoveTrack(StFtpcTrack *track)
1921 {
1922  // Removes track from ObjArry and takes care that the points are released again.
1923 
1924  track->SetProperties(kFALSE, -1); // release points
1925  mTrack->Remove(track); // actual removement of TObjArray
1926  delete track; // delete track
1927 
1928  return;
1929 }
1930 
1931 
1932 Int_t StFtpcConfMapper::GetRowSegm(StFtpcConfMapPoint *hit)
1933 {
1934  // Returns number of pad segment of a specific hit.
1935 
1936  return hit->GetPadRow() - 1; // fPadRow (1-20) already segmented, only offset substraction
1937 }
1938 
1939 
1940 Int_t StFtpcConfMapper::GetRow(Int_t segm)
1941 {
1942  // Returns the row number of a given row segment.
1943 
1944  return segm + 1;
1945 }
1946 
1947 
1948 Int_t StFtpcConfMapper::GetPhiSegm(StFtpcConfMapPoint *hit)
1949 {
1950  // Returns number of phi segment of a specific hit.
1951 
1952  return (Int_t)(hit->GetPhi() * mNumPhiSegment / (2.*TMath::Pi())); // fPhi has no offset but needs to be segmented (this is done by type conversion to Int_t)
1953 }
1954 
1955 
1956 Double_t StFtpcConfMapper::GetPhi(Int_t segm)
1957 {
1958  // Returns the angle phi of a given segment.
1959 
1960  return 2 * TMath::Pi() * segm / (Double_t)mNumPhiSegment;
1961 }
1962 
1963 
1964 Int_t StFtpcConfMapper::GetEtaSegm(StFtpcConfMapPoint *hit)
1965 {
1966  // Returns number of eta segment of a specific hit.
1967 
1968 
1969  Double_t eta;
1970  Int_t eta_segm;
1971 
1972  if ((eta = hit->GetEta()) > 0.) { // positive values
1973  eta_segm = (Int_t)((eta - mEtaMin) * mNumEtaSegment/(mEtaMax - mEtaMin) /2.); // Only use n_eta_segm/2. bins because of negative eta values.
1974  }
1975 
1976  else { // negative eta values
1977  eta_segm = (Int_t)((-eta - mEtaMin) * mNumEtaSegment/(mEtaMax - mEtaMin) /2. + mNumEtaSegment/2.);
1978  }
1979 
1980  return eta_segm;
1981 }
1982 
1983 
1984 Double_t StFtpcConfMapper::GetEta(Int_t segm)
1985 {
1986  // Returns the pseudorapidity eta of the given segment.
1987 
1988  Bool_t minus_sign = (Bool_t)kFALSE;
1989 
1990  if (segm >= mNumEtaSegment/2.) {
1991  minus_sign = (Bool_t)kTRUE;
1992  segm -= mNumEtaSegment/2;
1993  }
1994 
1995  if (minus_sign) {
1996  return -segm * (mEtaMax - mEtaMin)/ (2. * mNumEtaSegment) + mEtaMin;
1997  }
1998 
1999  else {
2000  return +segm * (mEtaMax - mEtaMin)/ (2. * mNumEtaSegment) + mEtaMin;
2001  }
2002 }
2003 
2004 
2005 Int_t StFtpcConfMapper::GetSegm(StFtpcConfMapPoint *hit)
2006 {
2007  // Calculates segment of a given hit.
2008 
2009  return GetSegm(GetRowSegm(hit), GetPhiSegm(hit), GetEtaSegm(hit));
2010 }
2011 
2012 
2013 Int_t StFtpcConfMapper::GetSegm(Int_t row_segm, Int_t phi_segm, Int_t eta_segm)
2014 {
2015  // Calculates the volume segment number from the segmented volumes (segm = segm(pad,phi,eta)).
2016 
2017  Int_t segm = row_segm * (mNumPhiSegment * mNumEtaSegment) + phi_segm * (mNumEtaSegment) + eta_segm;
2018 
2019  if (segm < 0 || segm >= mBounds) {
2020  LOG_WARN << "Segment calculation out of bounds (row = " << GetRow(row_segm) << ", phi = " << GetPhi(phi_segm) << ", eta = " << GetEta(eta_segm) << ")! Garbage segment returned." << endm;
2021  return mBounds;
2022  }
2023 
2024  else return segm;
2025 }
2026 
2027 
2028 Int_t StFtpcConfMapper::GetRowSegm(Int_t segm)
2029 {
2030  // Returns number of pad segment of a specifiv segment.
2031  return (segm - GetEtaSegm(segm) - GetPhiSegm(segm)) / (mNumPhiSegment * mNumEtaSegment);
2032 }
2033 
2034 
2035 Int_t StFtpcConfMapper::GetPhiSegm(Int_t segm)
2036 {
2037  // Returns number of phi segment of a specifiv segment.
2038  return (segm - GetEtaSegm(segm)) % (mNumPhiSegment * mNumEtaSegment) / (mNumEtaSegment);
2039 }
2040 
2041 
2042 Int_t StFtpcConfMapper::GetEtaSegm(Int_t segm)
2043 {
2044  // Returns number of eta segment of a specifiv segment.
2045 
2046  return (segm % (mNumPhiSegment * mNumEtaSegment)) % (mNumEtaSegment);
2047 }
2048 
2049 
2050 Double_t const StFtpcConfMapper::GetPhiDiff(const StFtpcConfMapPoint *hit1, const StFtpcConfMapPoint *hit2)
2051 {
2052  // Returns the difference in angle phi of the two given clusters.
2053  // Normalizes the result to the arbitrary angle between two subsequent padrows.
2054 
2055  Double_t angle = TMath::Abs(hit1->GetPhi() - hit2->GetPhi());
2056 
2057  if (angle > TMath::Pi()) {
2058  angle = 2*TMath::Pi() - angle;
2059  }
2060 
2061  return angle / (TMath::Abs(hit1->GetPadRow() - hit2->GetPadRow()));
2062 }
2063 
2064 
2065 Double_t const StFtpcConfMapper::GetEtaDiff(const StFtpcConfMapPoint *hit1, const StFtpcConfMapPoint *hit2)
2066 {
2067  // Returns the difference in pseudrapidity eta of the two given clusters.
2068  // Normalizes the result to the arbitrary pseudorapidity between two subsequent padrows.
2069 
2070  return (TMath::Abs(hit1->GetEta() - hit2->GetEta())) / (TMath::Abs(hit1->GetPadRow() - hit2->GetPadRow()));
2071 }
2072 
2073 
2074 Double_t const StFtpcConfMapper::GetClusterDistance(const StFtpcConfMapPoint *hit1, const StFtpcConfMapPoint *hit2)
2075 {
2076  // Returns the distance of two clusters measured in terms of angle phi and pseudorapidity eta weighted by the
2077  // maximal allowed values for phi and eta.
2078 
2079  return TMath::Sqrt(TMath::Power(GetPhiDiff(hit1, hit2)/mMaxCircleDist[mVertexConstraint], 2) + TMath::Power(GetEtaDiff(hit1, hit2)/mMaxLengthDist[mVertexConstraint], 2));
2080 }
2081 
2082 
2083 Double_t const StFtpcConfMapper::GetDistanceFromFit(const StFtpcConfMapPoint *hit)
2084 {
2085  // Returns the distance of the given cluster to the track to which it probably belongs.
2086  // The distances to the circle and length fit are weighted by the cuts on these values.
2087  // Make sure that the variables mCircleDist and mLengthDist for the hit are set already.
2088 
2089  return TMath::Sqrt(TMath::Power((hit->GetCircleDist() / mMaxCircleDist[mVertexConstraint]), 2) + TMath::Power((hit->GetLengthDist() / mMaxLengthDist[mVertexConstraint]), 2));
2090 }
2091 
2092 
2093 void StFtpcConfMapper::AdjustTrackNumbers(Int_t first_split)
2094 {
2095  // After split tracks are merged the size of the ObjArry holding the tracks needs to be changed.
2096  // Afterwards the tracknumbers need to be changed to be the number of the slot of the array again.
2097 
2098  mTrack->Compress();
2099  mTrack->Expand(mTrack->GetLast()+1);
2100 
2101  for (Int_t i = first_split; i < mTrack->GetEntriesFast(); ((StFtpcTrack*)mTrack->At(i))->SetTrackNumber(i), i++);
2102 
2103  return;
2104 }
2105