FairRoot/PandaRoot
Public Member Functions | Private Member Functions | Private Attributes | List of all members
PndLmdTrackFinderTask Class Reference

#include <PndLmdTrackFinderTask.h>

Inheritance diagram for PndLmdTrackFinderTask:

Public Member Functions

 PndLmdTrackFinderTask (Int_t inFinderMode=0, TString hitBranch="LMDHitsStrip", Int_t innSensPP=8)
 
virtual ~PndLmdTrackFinderTask ()
 
virtual void SetParContainers ()
 
virtual InitStatus Init ()
 
virtual InitStatus ReInit ()
 
virtual void Exec (Option_t *opt)
 
void SetVerbose (Int_t verbose)
 
void SetInaccuracy (Double_t accu)
 
void SetSensStripFlag (bool fS)
 
void SetSensPixelFlag (bool fS)
 

Private Member Functions

Double_t GetTrackDip (PndMCTrack *myTrack)
 
Double_t GetTrackCurvature (PndMCTrack *myTrack)
 
void Register ()
 
void Reset ()
 
void ProduceHits ()
 
bool SortHitsByZ (std::vector< std::vector< std::pair< Int_t, bool > > > &hitsd, Int_t nStripHits)
 
bool SortHitsByDet (std::vector< std::vector< std::pair< Int_t, bool > > > &hitsd, Int_t nStripHits)
 
bool SortHitsByDet2 (std::vector< std::vector< std::pair< Int_t, bool > > > &hitsd, Int_t nStripHits)
 
void FindHitsI (std::vector< PndTrackCand > &tofill, std::vector< std::vector< std::pair< Int_t, bool > > > &hitsd, Int_t nStripHits)
 
void FindHitsII (std::vector< PndTrackCand > &tofill, std::vector< std::vector< std::pair< Int_t, bool > > > &hitsd, Int_t nStripHits)
 
void FindHitsIII (std::vector< PndTrackCand > &tofill, std::vector< std::vector< std::pair< Int_t, bool > > > &hitsd, Int_t nStripHits)
 
 ClassDef (PndLmdTrackFinderTask, 2)
 

Private Attributes

bool flagStipSens
 
bool flagPixelSens
 
Double_t dXY
 
Int_t fFinderMode
 
Int_t nSensPP
 
TString fHitBranchStrip
 
TClonesArray * fStripHitArray
 
TClonesArray * fTrackCandArray
 

Detailed Description

Definition at line 23 of file PndLmdTrackFinderTask.h.

Constructor & Destructor Documentation

PndLmdTrackFinderTask::PndLmdTrackFinderTask ( Int_t  inFinderMode = 0,
TString  hitBranch = "LMDHitsStrip",
Int_t  innSensPP = 8 
)

Default constructor

Definition at line 21 of file PndLmdTrackFinderTask.cxx.

References dXY, fFinderMode, fHitBranchStrip, flagPixelSens, and flagStipSens.

23  : FairTask("LMD Track Finding Task"),
24  nSensPP(innSensPP) {
25  fFinderMode = inFinderMode;
26  fHitBranchStrip = hitBranch;
27  dXY = 0.01;
28  flagStipSens = false;
29  flagPixelSens = false;
30  // dXY = 0.01;//TEST
31 }
PndLmdTrackFinderTask::~PndLmdTrackFinderTask ( )
virtual

Destructor

Definition at line 36 of file PndLmdTrackFinderTask.cxx.

36 {}

Member Function Documentation

PndLmdTrackFinderTask::ClassDef ( PndLmdTrackFinderTask  ,
 
)
private
void PndLmdTrackFinderTask::Exec ( Option_t *  opt)
virtual

Virtual method Exec

strip sensors

pixel sensors

Definition at line 680 of file PndLmdTrackFinderTask.cxx.

References Double_t, fFinderMode, FindHitsI(), FindHitsII(), FindHitsIII(), flagPixelSens, flagStipSens, fStripHitArray, fTrackCandArray, fVerbose, hit, i, SortHitsByDet(), SortHitsByDet2(), and t.

680  {
681  TStopwatch *timer_exec = new TStopwatch();
682  if (fVerbose > 2) timer_exec->Start();
683  if (fVerbose > 2) cout << "Evt started--------------" << endl << endl;
684 
685  // which combinations of planes to build pseudo-vectos
686  // Int_t planeFlag=0; //flag=... 0: planes 1-2, 1: 1-2 & 1-3, 2: 1-2 & 1-3
687  // & 2-3
688 
689  // Reset output array
690  if (!fTrackCandArray) Fatal("Exec", "No trackCandArray");
691  fTrackCandArray->Clear();
692 
693  Int_t nStripHits = fStripHitArray->GetEntriesFast();
694  if (fVerbose > 2) cout << "# Hits: \t" << nStripHits << endl << endl;
695  // bool usedFlag[nStripHits]; //for pseudo-vector building //[R.K.03/2017]
696  // unused variable
697  // for(Int_t seti=0; seti<nStripHits; seti++) //[R.K.03/2017] unused variable
698  // usedFlag[seti]=false; //[R.K.03/2017] unused variable
699 
700  if (nStripHits < 2) {
701  if (fVerbose > 2)
702  cout << "Evt finsihed: too less hits-----" << endl << endl;
703  return;
704  }
705 
706  //-----do some sorting first----------
707  std::vector<std::vector<std::pair<int, bool> > > hitsd(4);
708  bool resSortHits;
709  if (flagStipSens)
710  resSortHits = SortHitsByDet(hitsd, nStripHits);
711  else {
712  if (flagPixelSens)
713  resSortHits = SortHitsByDet2(hitsd, nStripHits);
714  else {
715  std::cout << "Algorithm is needed sensor type! Please, set it via "
716  "SetSensStripFlag(bool fS) or SetSensPixelFlag(bool fS)"
717  << std::endl;
718  return;
719  }
720  }
721  if (!resSortHits) {
722  if (fVerbose > 2)
723  cout << "Evt finsihed: too less planes-----" << endl << endl;
724  return;
725  }
726 
727  if (fVerbose > 2) {
728  cout << "HitMap size: " << hitsd.size() << endl;
729  if (hitsd.size() > 0) {
730  for (size_t i = 0; i < hitsd.at(0).size(); i++) {
731  PndSdsHit *hit =
732  (PndSdsHit *)fStripHitArray->At(hitsd.at(0).at(i).first);
733  cout << "Plane0 Hit=(" << hit->GetX() << ", " << hit->GetY() << ", "
734  << hit->GetZ() << ")" << endl;
735  }
736  }
737  if (hitsd.size() > 1) {
738  for (size_t i = 0; i < hitsd.at(1).size(); i++) {
739  PndSdsHit *hit =
740  (PndSdsHit *)fStripHitArray->At(hitsd.at(1).at(i).first);
741  cout << "Plane1 Hit=(" << hit->GetX() << ", " << hit->GetY() << ", "
742  << hit->GetZ() << ")" << endl;
743  }
744  }
745  if (hitsd.size() > 2) {
746  for (size_t i = 0; i < hitsd.at(2).size(); i++) {
747  PndSdsHit *hit =
748  (PndSdsHit *)fStripHitArray->At(hitsd.at(2).at(i).first);
749  cout << "Plane2 Hit=(" << hit->GetX() << ", " << hit->GetY() << ", "
750  << hit->GetZ() << ")" << endl;
751  }
752  }
753  if (hitsd.size() > 3) {
754  for (size_t i = 0; i < hitsd.at(3).size(); i++) {
755  PndSdsHit *hit =
756  (PndSdsHit *)fStripHitArray->At(hitsd.at(3).at(i).first);
757  cout << "Plane3 Hit=(" << hit->GetX() << ", " << hit->GetY() << ", "
758  << hit->GetZ() << ")" << endl;
759  }
760  }
761  }
762 
763  std::vector<PndTrackCand> theCands; // Track Candidates
764 
765  //--------find some tracks--------------
766  FindHitsI(theCands, hitsd, nStripHits);
767  if (fFinderMode) {
768  FindHitsII(theCands, hitsd, nStripHits);
769  FindHitsIII(theCands, hitsd, nStripHits);
770  }
771 
772  // fill tracklist für fitting
773  for (size_t t = 0; t < theCands.size(); t++)
774  new ((*fTrackCandArray)[t]) PndTrackCand(theCands.at(t));
775 
776  if (fVerbose > 2) cout << "Evt finsihed--------------" << endl << endl;
777  if (fVerbose > 2) {
778  timer_exec->Stop();
779  Double_t rtime_exec = timer_exec->RealTime();
780  Double_t ctime_exec = timer_exec->CpuTime();
781  cout << "Real time for Exec:" << rtime_exec << " s, CPU time " << ctime_exec
782  << " s" << endl;
783  cout << endl;
784  }
785 }
int fVerbose
Definition: poormantracks.C:24
Int_t i
Definition: run_full.C:25
Double_t
void FindHitsIII(std::vector< PndTrackCand > &tofill, std::vector< std::vector< std::pair< Int_t, bool > > > &hitsd, Int_t nStripHits)
bool SortHitsByDet2(std::vector< std::vector< std::pair< Int_t, bool > > > &hitsd, Int_t nStripHits)
void FindHitsI(std::vector< PndTrackCand > &tofill, std::vector< std::vector< std::pair< Int_t, bool > > > &hitsd, Int_t nStripHits)
TTree * t
Definition: bump_analys.C:13
void FindHitsII(std::vector< PndTrackCand > &tofill, std::vector< std::vector< std::pair< Int_t, bool > > > &hitsd, Int_t nStripHits)
PndSdsMCPoint * hit
Definition: anasim.C:70
bool SortHitsByDet(std::vector< std::vector< std::pair< Int_t, bool > > > &hitsd, Int_t nStripHits)
void PndLmdTrackFinderTask::FindHitsI ( std::vector< PndTrackCand > &  tofill,
std::vector< std::vector< std::pair< Int_t, bool > > > &  hitsd,
Int_t  nStripHits 
)
private

Add seed information to track---------—

Definition at line 509 of file PndLmdTrackFinderTask.cxx.

References PndTrackCand::AddHit(), Bool_t, Double_t, dXY, fStripHitArray, fVerbose, PndSdsHit::GetPosition(), hit, i, sqrt(), and vec.

Referenced by Exec().

513 {
514  std::vector<TVector3> trackStart, trackVec; // pseudo tracks
515  std::vector<TVector3> trackStartd, trackVecd; // save errors
516  std::vector<Int_t> trackID1, trackID2; // for TrackCand
517 
518  // iterate first discs-hits with all seconds, save pseudo-tracks
519  TVector3 start, tmp, vec, dstart, dvec; // temp-vars
520 
521  if (hitsd.size() < 2) return;
522  for (size_t i = 0; i < hitsd.at(0).size(); i++) {
523  PndSdsHit *hit1 = (PndSdsHit *)fStripHitArray->At(hitsd.at(0).at(i).first);
524  start.SetXYZ(hit1->GetX(), hit1->GetY(), hit1->GetZ());
525  dstart.SetXYZ(hit1->GetDx(), hit1->GetDy(), hit1->GetDz());
526  for (size_t k = 0; k < hitsd.at(1).size(); k++) {
527  PndSdsHit *hit2 =
528  (PndSdsHit *)fStripHitArray->At(hitsd.at(1).at(k).first);
529  tmp.SetXYZ(hit2->GetX(), hit2->GetY(), hit2->GetZ());
530  vec = tmp - start; // calc direction vector for FINDING
531  dvec.SetXYZ(hit2->GetDx(), hit2->GetDy(), hit2->GetDz());
532  if (vec.Theta() > 0.03 && vec.Theta() < 0.05 && vec.Phi() > -0.25 &&
533  vec.Phi() < 0.25) { // ignore vectors with theta outside 2-9 mrad
534  // if(vec.Theta()<0.01){ //ignore vectors with theta outside 10
535  // mrad
536  trackStart.push_back(start);
537  trackStartd.push_back(dstart);
538  trackVec.push_back(vec); // save vector from start to second
539  trackVecd.push_back(dvec); // save error of second point for FIT, NOT
540  // error of direction vector
541  // trackID1.push_back(hitsd.at(0).at(i).first); //save Hit-Id's for
542  // TrackCand
543  // trackID2.push_back(hitsd.at(1).at(k).first);
544  trackID1.push_back(i); // save Hit-Id's for TrackCand
545  trackID2.push_back(k);
546  }
547  } // end of disc 2 (or 3 if none in disc 2) hits
548  } // end of disc 1 hits
549 
550  if (fVerbose > 1) cout << "Pseudos L1L2: " << trackStart.size() << endl;
551 
552  std::vector<Int_t> ids;
553 
554  // check if other discs have hits in track, add points for fitting
555  // ---------------
556  for (size_t i = 0; i < trackStart.size(); i++) // pseudo-loop
557  {
558  ids.clear();
559  Int_t pntcnt = 2;
560 
561  start = trackStart.at(i);
562  dstart = trackStartd.at(i);
563  vec = trackVec.at(i);
564  dvec = trackVecd.at(i);
565  ids.push_back(hitsd.at(0).at(trackID1.at(i)).first);
566  ids.push_back(hitsd.at(1).at(trackID2.at(i)).first);
567  std::vector<std::pair<Int_t, Int_t> > otherIDs;
568 
569  for (Int_t idet = 2; idet < 4; idet++) {
570  Double_t distClosest = 2 * idet * dXY; // just bigger as possible
571  Bool_t firstp = true;
572  for (size_t ihit = 0; ihit < hitsd.at(idet).size(); ihit++) {
573  PndSdsHit *hit =
574  (PndSdsHit *)fStripHitArray->At(hitsd.at(idet).at(ihit).first);
575  Double_t scale = (hit->GetZ() - start.z()) / vec.z();
576  tmp = start + scale * vec; // extend search-vector to hit-plane
577  Double_t distTmp =
578  sqrt((tmp.x() - hit->GetX()) * (tmp.x() - hit->GetX()) +
579  (tmp.y() - hit->GetY()) * (tmp.y() - hit->GetY()));
580  if (distTmp < (idet * dXY)) { // if in diameter
581  if (firstp) {
582  pntcnt++;
583  ids.push_back(hitsd.at(idet).at(ihit).first);
584  otherIDs.push_back(make_pair(idet, ihit));
585  distClosest = distTmp;
586  firstp = false;
587  } else {
588  if (distTmp < distClosest) { // change if nearer
589  ids.pop_back();
590  otherIDs.pop_back();
591  ids.push_back(hitsd.at(idet).at(ihit).first);
592  otherIDs.push_back(make_pair(idet, ihit));
593  }
594  }
595  } // endif in diameter
596  } // end of disc n hits
597  } // end of discs
598 
599  if (fVerbose > 2)
600  cout << " Track L1L2: " << i << "#Planes: " << ids.size() << endl;
601 
602  // create track für fitting
603  if (ids.size() > 2) { // third hit found <=> !track found!
604  // if(ids.size()>3){ //4 hits are needed because 6 parameters are used
605  // in trk-fit!
606  // // //third hit found <=> !track found!
607  // // //check direction for reduction of ghost tracks
608  // if(ids.size()>2 &&
609  // fabs(vec.Phi())<0.26 && vec.Theta()>3e-2 && vec.Theta()<5e-2){
610  if (fVerbose > 2)
611  cout << " Track: " << i << "#Planes: " << ids.size() << endl;
612  PndTrackCand *myTCand = new PndTrackCand();
613 
614  for (size_t id = 0; id < ids.size(); id++) {
615  PndSdsHit *myHit = (PndSdsHit *)(fStripHitArray->At(ids.at(id)));
616  myTCand->AddHit(myHit->GetDetectorID(), ids.at(id),
617  myHit->GetPosition().Mag());
618  }
619 
620  // mark used hits---------
621  hitsd.at(0).at(trackID1.at(i)).second = true;
622  hitsd.at(1).at(trackID2.at(i)).second = true;
623  for (size_t id = 0; id < otherIDs.size(); id++) {
624  hitsd.at(otherIDs.at(id).first).at(otherIDs.at(id).second).second =
625  true;
626  }
627 
629  // Double_t Z0 = 1099.;
630  // // Double_t Z0 = 0;
631  // PndSdsHit* myHit0 = (PndSdsHit*)(fStripHitArray->At(ids.at(0)));
632  // PndSdsHit* myHit1 = (PndSdsHit*)(fStripHitArray->At(ids.at(1)));
633  // TVector3 hit0 = myHit0->GetPosition(); TVector3 hit1 =
634  // myHit1->GetPosition();
635  // double p1seed = (hit0.X()-hit1.X())/(hit0.Z()-hit1.Z());
636  // double p0seed = hit0.X() - p1seed*(hit0.Z()-Z0);
637  // // double p0seed =
638  // 0.5*(hit0.X()+hit1.X()-p1seed*(hit0.Z()+hit1.Z()-2*1099.)); //TO DO:
639  // don't use const
640  // double p3seed = (hit0.Y()-hit1.Y())/(hit0.Z()-hit1.Z());
641  // double p2seed = hit0.Y() - p3seed*(hit0.Z()-Z0);
642  // //double p2seed =
643  // 0.5*(hit0.Y()+hit1.Y()-p1seed*(hit0.Z()+hit1.Z()-2*1099.)); //TO DO:
644  // don't use const
645  // TVector3 posSeed(p0seed,p2seed,Z0);
646  // TVector3 dirSeed(p1seed,p3seed,1.);
647  // cout<<"posSeed:"<<endl;
648  // posSeed.Print();
649  // cout<<"dirSeed:"<<endl;
650  // dirSeed.Print();
651 
652  // // Double_t Z0 = 1099.;
653  // // TVector3 posSeed(start.X(),start.Y(),Z0);
654  // TVector3 posSeed(start.X(),start.Y(),start.Z());
655  // // cout<<"Trk cand: pos"<<endl;
656  // // posSeed.Print();
657  // vec*=1./vec.Mag();
658  // //shift trk out of plane [needed for correct treatment in Kalman
659  // Fillter and GEANE]
660  // double sh_z = -0.035; //350 mkm
661  // double sh_x = vec.X()*sh_z;
662  // double sh_y = vec.Y()*sh_y;
663  // TVector3 sh_point(sh_x,sh_y,sh_z);
664  // posSeed +=sh_point;
665 
666  // // cout<<"Trk cand: vec"<<endl;
667  // // vec.Print();
668  // myTCand->setTrackSeed(posSeed,vec,-1);
669  // ///-------------------------------------
670 
671  tofill.push_back(*(myTCand)); // save Track Candidate
672  // new((*fTrackCandArray)[trackCnt]) PndTrackCand(*(myTCand));
673  delete myTCand;
674  } // Track Cand build
675  } // end of pseudo-tracks
676 }
int fVerbose
Definition: poormantracks.C:24
Int_t i
Definition: run_full.C:25
TVector3 GetPosition() const
Definition: PndSdsHit.h:93
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
void AddHit(UInt_t detId, UInt_t hitId, Double_t rho)
Double_t
PndSdsMCPoint * hit
Definition: anasim.C:70
dble_vec_t vec[12]
Definition: ranlxd.cxx:380
void PndLmdTrackFinderTask::FindHitsII ( std::vector< PndTrackCand > &  tofill,
std::vector< std::vector< std::pair< Int_t, bool > > > &  hitsd,
Int_t  nStripHits 
)
private

Definition at line 367 of file PndLmdTrackFinderTask.cxx.

References PndTrackCand::AddHit(), Bool_t, Double_t, dXY, fStripHitArray, fVerbose, PndSdsHit::GetPosition(), hit, i, sqrt(), and vec.

Referenced by Exec().

371 {
372  std::vector<TVector3> trackStart, trackVec; // pseudo tracks
373  std::vector<TVector3> trackStartd, trackVecd; // save errors
374  std::vector<Int_t> trackID1, trackID3; // for TrackCand
375 
376  // iterate first discs-hits with all seconds, save pseudo-tracks
377  TVector3 start, tmp, vec, dstart, dvec; // temp-vars
378 
379  if (hitsd.size() < 3) return;
380  for (size_t i = 0; i < hitsd.at(0).size(); i++) {
381  if (hitsd.at(0).at(i).second) continue; // if already used
382  PndSdsHit *hit1 = (PndSdsHit *)fStripHitArray->At(hitsd.at(0).at(i).first);
383  start.SetXYZ(hit1->GetX(), hit1->GetY(), hit1->GetZ());
384  dstart.SetXYZ(hit1->GetDx(), hit1->GetDy(), hit1->GetDz());
385  for (size_t k = 0; k < hitsd.at(2).size(); k++) {
386  if (hitsd.at(2).at(k).second) continue; // if already used
387  PndSdsHit *hit2 =
388  (PndSdsHit *)fStripHitArray->At(hitsd.at(2).at(k).first);
389  tmp.SetXYZ(hit2->GetX(), hit2->GetY(), hit2->GetZ());
390  vec = tmp - start; // calc direction vector for FINDING
391  dvec.SetXYZ(hit2->GetDx(), hit2->GetDy(), hit2->GetDz());
392  if (vec.Theta() > 0.03 && vec.Theta() < 0.05 && vec.Phi() > -0.25 &&
393  vec.Phi() < 0.25) { // ignore vectors with theta outside 2-9 mrad
394  // if(0<1){
395  // if(vec.Theta()<0.01){ //ignore vectors with theta outside 10
396  //mrad
397  trackStart.push_back(start);
398  trackStartd.push_back(dstart);
399  trackVec.push_back(vec); // save vector from start to second
400  trackVecd.push_back(dvec); // save error of second point for FIT, NOT
401  // error of direction vector
402  // trackID1.push_back(hitsd.at(0).at(i).first); //save Hit-Id's for
403  // TrackCand
404  // trackID2.push_back(hitsd.at(1).at(k).first);
405  trackID1.push_back(i); // save Hit-Id's for TrackCand
406  trackID3.push_back(k);
407  }
408  } // end of disc 2 (or 3 if none in disc 2) hits
409  } // end of disc 1 hits
410 
411  if (fVerbose > 1) cout << "Pseudos L1L3: " << trackStart.size() << endl;
412 
413  std::vector<Int_t> ids;
414 
415  // check if other discs have hits in track, add points for fitting
416  // ---------------
417  for (size_t i = 0; i < trackStart.size(); i++) // pseudo-loop
418  {
419  ids.clear();
420  Int_t pntcnt = 2;
421 
422  start = trackStart.at(i);
423  dstart = trackStartd.at(i);
424  vec = trackVec.at(i);
425  dvec = trackVecd.at(i);
426  ids.push_back(hitsd.at(0).at(trackID1.at(i)).first);
427  ids.push_back(hitsd.at(2).at(trackID3.at(i)).first);
428  std::vector<std::pair<Int_t, Int_t> > otherIDs;
429 
430  for (Int_t idet = 3; idet < 4;
431  idet++) { // i know this is not a loop, but in case we someday will
432  // have more planes
433  Double_t distClosest = 2 * idet * dXY; // just bigger as possible
434  Bool_t firstp = true;
435  for (size_t ihit = 0; ihit < hitsd.at(idet).size(); ihit++) {
436  PndSdsHit *hit =
437  (PndSdsHit *)fStripHitArray->At(hitsd.at(idet).at(ihit).first);
438  Double_t scale = (hit->GetZ() - start.z()) / vec.z();
439  tmp = start + scale * vec; // extend search-vector to hit-plane
440  Double_t distTmp =
441  sqrt((tmp.x() - hit->GetX()) * (tmp.x() - hit->GetX()) +
442  (tmp.y() - hit->GetY()) * (tmp.y() - hit->GetY()));
443  if (distTmp < (idet * dXY)) { // if in diameter
444  if (firstp) {
445  pntcnt++;
446  ids.push_back(hitsd.at(idet).at(ihit).first);
447  otherIDs.push_back(make_pair(idet, ihit));
448  distClosest = distTmp;
449  firstp = false;
450  } else {
451  if (distTmp < distClosest) { // change if nearer
452  ids.pop_back();
453  otherIDs.pop_back();
454  ids.push_back(hitsd.at(idet).at(ihit).first);
455  otherIDs.push_back(make_pair(idet, ihit));
456  }
457  }
458  } // endif in diameter
459  } // end of disc n hits
460  } // end of discs
461 
462  if (fVerbose > 2)
463  cout << " Track L1L3: " << i << "#Planes: " << ids.size() << endl;
464 
465  // create track für fitting
466  if (ids.size() > 2) { // third hit found <=> !track found!
467  PndTrackCand *myTCand = new PndTrackCand();
468 
469  for (size_t id = 0; id < ids.size(); id++) {
470  PndSdsHit *myHit = (PndSdsHit *)(fStripHitArray->At(ids.at(id)));
471  myTCand->AddHit(myHit->GetDetectorID(), ids.at(id),
472  myHit->GetPosition().Mag());
473  }
474 
475  // mark used hits---------
476  hitsd.at(0).at(trackID1.at(i)).second = true;
477  hitsd.at(2).at(trackID3.at(i)).second = true;
478  for (size_t id = 0; id < otherIDs.size(); id++) {
479  hitsd.at(otherIDs.at(id).first).at(otherIDs.at(id).second).second =
480  true;
481  }
482 
483  // ///Add seed information to track------------
484  // PndSdsHit* myHit0 = (PndSdsHit*)(fStripHitArray->At(ids.at(0)));
485  // PndSdsHit* myHit1 = (PndSdsHit*)(fStripHitArray->At(ids.at(1)));
486  // TVector3 hit0 = myHit0->GetPosition(); TVector3 hit1 =
487  // myHit1->GetPosition();
488  // TVector3 posSeed(hit0.X(),hit0.Y(),hit0.Z());
489  // vec*=1./vec.Mag();
490  // //shift trk out of plane [needed for correct treatment in Kalman
491  // Fillter and GEANE]
492  // double sh_z = -0.035; //350 mkm
493  // double sh_x = vec.X()*sh_z;
494  // double sh_y = vec.Y()*sh_y;
495  // TVector3 sh_point(sh_x,sh_y,sh_z);
496  // posSeed +=sh_point;
497  // myTCand->setTrackSeed(posSeed,vec,-1);
498  // ///-------------------------------------
499 
500  tofill.push_back(*(myTCand)); // save Track Candidate
501  // new((*fTrackCandArray)[trackCnt]) PndTrackCand(*(myTCand));
502  delete myTCand;
503  } // Track Cand build
504  } // end of pseudo-tracks
505 }
int fVerbose
Definition: poormantracks.C:24
Int_t i
Definition: run_full.C:25
TVector3 GetPosition() const
Definition: PndSdsHit.h:93
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
void AddHit(UInt_t detId, UInt_t hitId, Double_t rho)
Double_t
PndSdsMCPoint * hit
Definition: anasim.C:70
dble_vec_t vec[12]
Definition: ranlxd.cxx:380
void PndLmdTrackFinderTask::FindHitsIII ( std::vector< PndTrackCand > &  tofill,
std::vector< std::vector< std::pair< Int_t, bool > > > &  hitsd,
Int_t  nStripHits 
)
private

Definition at line 224 of file PndLmdTrackFinderTask.cxx.

References PndTrackCand::AddHit(), Bool_t, Double_t, dXY, fStripHitArray, fVerbose, PndSdsHit::GetPosition(), hit, i, sqrt(), and vec.

Referenced by Exec().

228 {
229  std::vector<TVector3> trackStart, trackVec; // pseudo tracks
230  std::vector<TVector3> trackStartd, trackVecd; // save errors
231  std::vector<Int_t> trackID2, trackID3; // for TrackCand
232 
233  // iterate first discs-hits with all seconds, save pseudo-tracks
234  TVector3 start, tmp, vec, dstart, dvec; // temp-vars
235 
236  if (hitsd.size() < 3) return;
237  for (size_t i = 0; i < hitsd.at(1).size(); i++) {
238  if (hitsd.at(1).at(i).second) continue; // if already used
239  PndSdsHit *hit1 = (PndSdsHit *)fStripHitArray->At(hitsd.at(1).at(i).first);
240  start.SetXYZ(hit1->GetX(), hit1->GetY(), hit1->GetZ());
241  dstart.SetXYZ(hit1->GetDx(), hit1->GetDy(), hit1->GetDz());
242  for (size_t k = 0; k < hitsd.at(2).size(); k++) {
243  if (hitsd.at(2).at(k).second) continue; // if already used
244  PndSdsHit *hit2 =
245  (PndSdsHit *)fStripHitArray->At(hitsd.at(2).at(k).first);
246  tmp.SetXYZ(hit2->GetX(), hit2->GetY(), hit2->GetZ());
247  vec = tmp - start; // calc direction vector for FINDING
248  dvec.SetXYZ(hit2->GetDx(), hit2->GetDy(), hit2->GetDz());
249  if (vec.Theta() > 0.03 && vec.Theta() < 0.05 && vec.Phi() > -0.25 &&
250  vec.Phi() < 0.25) { // ignore vectors with theta outside 2-9 mrad
251  // if(0<1){
252  // if(vec.Theta()<0.01){ //ignore vectors with theta outside 10
253  // mrad
254  trackStart.push_back(start);
255  trackStartd.push_back(dstart);
256  trackVec.push_back(vec); // save vector from start to second
257  trackVecd.push_back(dvec); // save error of second point for FIT, NOT
258  // error of direction vector
259  // trackID1.push_back(hitsd.at(0).at(i).first); //save Hit-Id's for
260  // TrackCand
261  // trackID2.push_back(hitsd.at(1).at(k).first);
262  trackID2.push_back(i); // save Hit-Id's for TrackCand
263  trackID3.push_back(k);
264  }
265  } // end of disc 2 (or 3 if none in disc 2) hits
266  } // end of disc 1 hits
267 
268  if (fVerbose > 1) cout << "Pseudos L2L3: " << trackStart.size() << endl;
269 
270  std::vector<Int_t> ids;
271 
272  // check if other discs have hits in track, add points for fitting
273  // ---------------
274  for (size_t i = 0; i < trackStart.size(); i++) // pseudo-loop
275  {
276  ids.clear();
277  Int_t pntcnt = 2;
278 
279  start = trackStart.at(i);
280  dstart = trackStartd.at(i);
281  vec = trackVec.at(i);
282  dvec = trackVecd.at(i);
283  ids.push_back(hitsd.at(1).at(trackID2.at(i)).first);
284  ids.push_back(hitsd.at(2).at(trackID3.at(i)).first);
285  std::vector<std::pair<Int_t, Int_t> > otherIDs;
286 
287  for (Int_t idet = 3; idet < 4;
288  idet++) { // i know this is not a loop, but in case we someday will
289  // have more planes
290  Double_t distClosest = 2 * idet * dXY; // just bigger as possible
291  Bool_t firstp = true;
292  for (size_t ihit = 0; ihit < hitsd.at(idet).size(); ihit++) {
293  PndSdsHit *hit =
294  (PndSdsHit *)fStripHitArray->At(hitsd.at(idet).at(ihit).first);
295  Double_t scale = (hit->GetZ() - start.z()) / vec.z();
296  tmp = start + scale * vec; // extend search-vector to hit-plane
297  Double_t distTmp =
298  sqrt((tmp.x() - hit->GetX()) * (tmp.x() - hit->GetX()) +
299  (tmp.y() - hit->GetY()) * (tmp.y() - hit->GetY()));
300  if (distTmp < (idet * dXY)) { // if in diameter
301  if (firstp) {
302  pntcnt++;
303  ids.push_back(hitsd.at(idet).at(ihit).first);
304  otherIDs.push_back(make_pair(idet, ihit));
305  distClosest = distTmp;
306  firstp = false;
307  } else {
308  if (distTmp < distClosest) { // change if nearer
309  ids.pop_back();
310  otherIDs.pop_back();
311  ids.push_back(hitsd.at(idet).at(ihit).first);
312  otherIDs.push_back(make_pair(idet, ihit));
313  }
314  }
315  } // endif in diameter
316  } // end of disc n hits
317  } // end of discs
318 
319  if (fVerbose > 2)
320  cout << " Track L2L3: " << i << "#Planes: " << ids.size() << endl;
321 
322  // create track for fitting
323  if (ids.size() > 2) { // third hit found <=> !track found!
324  PndTrackCand *myTCand = new PndTrackCand();
325 
326  for (size_t id = 0; id < ids.size(); id++) {
327  PndSdsHit *myHit = (PndSdsHit *)(fStripHitArray->At(ids.at(id)));
328  myTCand->AddHit(myHit->GetDetectorID(), ids.at(id),
329  myHit->GetPosition().Mag());
330  }
331 
332  // mark used hits---------
333  hitsd.at(1).at(trackID2.at(i)).second = true;
334  hitsd.at(2).at(trackID3.at(i)).second = true;
335  for (size_t id = 0; id < otherIDs.size(); id++) {
336  hitsd.at(otherIDs.at(id).first).at(otherIDs.at(id).second).second =
337  true;
338  }
339 
340  // ///Add seed information to track------------
341  // PndSdsHit* myHit0 = (PndSdsHit*)(fStripHitArray->At(ids.at(0)));
342  // PndSdsHit* myHit1 = (PndSdsHit*)(fStripHitArray->At(ids.at(1)));
343  // TVector3 hit0 = myHit0->GetPosition(); TVector3 hit1 =
344  // myHit1->GetPosition();
345  // TVector3 posSeed(hit0.X(),hit0.Y(),hit0.Z());
346  // vec*=1./vec.Mag();
347  // //shift trk out of plane [needed for correct treatment in Kalman
348  // Fillter and GEANE]
349  // double sh_z = -0.035; //350 mkm
350  // double sh_x = vec.X()*sh_z;
351  // double sh_y = vec.Y()*sh_y;
352  // TVector3 sh_point(sh_x,sh_y,sh_z);
353  // posSeed +=sh_point;
354  // myTCand->setTrackSeed(posSeed,vec,-1);
355  // ///-------------------------------------
356 
357  tofill.push_back(*(myTCand)); // save Track Candidate
358  // new((*fTrackCandArray)[trackCnt]) PndTrackCand(*(myTCand));
359  delete myTCand;
360  } // Track Cand build
361  } // end of pseudo-tracks
362 }
int fVerbose
Definition: poormantracks.C:24
Int_t i
Definition: run_full.C:25
TVector3 GetPosition() const
Definition: PndSdsHit.h:93
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
void AddHit(UInt_t detId, UInt_t hitId, Double_t rho)
Double_t
PndSdsMCPoint * hit
Definition: anasim.C:70
dble_vec_t vec[12]
Definition: ranlxd.cxx:380
Double_t PndLmdTrackFinderTask::GetTrackCurvature ( PndMCTrack myTrack)
private

Definition at line 788 of file PndLmdTrackFinderTask.cxx.

References PndMCTrack::GetMomentum(), p, and CAMath::Sqrt().

788  {
789  TVector3 p = myTrack->GetMomentum();
790  return (2 / TMath::Sqrt(p.Px() * p.Px() + p.Py() * p.Py()));
791 }
Double_t p
Definition: anasim.C:58
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
TVector3 GetMomentum() const
Definition: PndMCTrack.h:78
Double_t PndLmdTrackFinderTask::GetTrackDip ( PndMCTrack myTrack)
private

Definition at line 794 of file PndLmdTrackFinderTask.cxx.

References PndMCTrack::GetMomentum(), p, and CAMath::Sqrt().

794  {
795  TVector3 p = myTrack->GetMomentum();
796  return (p.Mag() / TMath::Sqrt(p.Px() * p.Px() + p.Py() * p.Py()));
797 }
Double_t p
Definition: anasim.C:58
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
TVector3 GetMomentum() const
Definition: PndMCTrack.h:78
InitStatus PndLmdTrackFinderTask::Init ( )
virtual

Definition at line 48 of file PndLmdTrackFinderTask.cxx.

References fHitBranchStrip, fStripHitArray, and fTrackCandArray.

48  {
49  FairRootManager *ioman = FairRootManager::Instance();
50 
51  if (!ioman) {
52  std::cout << "-E- PndLmdTrackFinderTask::Init: "
53  << "RootManager not instantiated!" << std::endl;
54  return kFATAL;
55  }
56 
57  // Get input array
58  fStripHitArray = (TClonesArray *)ioman->GetObject(fHitBranchStrip);
59  if (!fStripHitArray) {
60  std::cout << "-W- PndLmdTrackFinderTask::Init: "
61  << "No fStripHitArray!" << std::endl;
62  return kERROR;
63  }
64 
65  fTrackCandArray = new TClonesArray("PndTrackCand");
66  ioman->Register("LMDTrackCand", "PndLmd", fTrackCandArray, kTRUE);
67 
68  // fTruePointArray = new TClonesArray("PndSdsMCPoint");
69 
70  std::cout << "-I- PndLmdTrackFinderTask: Initialisation successfull"
71  << std::endl;
72  return kSUCCESS;
73 }
void PndLmdTrackFinderTask::ProduceHits ( )
private
void PndLmdTrackFinderTask::Register ( )
private
InitStatus PndLmdTrackFinderTask::ReInit ( )
virtual

Definition at line 42 of file PndLmdTrackFinderTask.cxx.

42  {
43  InitStatus stat = kERROR;
44  return stat;
45 }
void PndLmdTrackFinderTask::Reset ( )
private
void PndLmdTrackFinderTask::SetInaccuracy ( Double_t  accu)
inline

Definition at line 41 of file PndLmdTrackFinderTask.h.

References dXY.

Referenced by reco_LMD().

41 { dXY = accu; };
void PndLmdTrackFinderTask::SetParContainers ( )
virtual

Virtual method Init

Definition at line 40 of file PndLmdTrackFinderTask.cxx.

40 {}
void PndLmdTrackFinderTask::SetSensPixelFlag ( bool  fS)
inline

Definition at line 43 of file PndLmdTrackFinderTask.h.

References flagPixelSens.

Referenced by reco_LMD().

void PndLmdTrackFinderTask::SetSensStripFlag ( bool  fS)
inline

Definition at line 42 of file PndLmdTrackFinderTask.h.

References flagStipSens.

void PndLmdTrackFinderTask::SetVerbose ( Int_t  verbose)
inline

Definition at line 40 of file PndLmdTrackFinderTask.h.

References fVerbose, and verbose.

Referenced by reco_LMD().

40 { fVerbose = verbose; };
int fVerbose
Definition: poormantracks.C:24
#define verbose
bool PndLmdTrackFinderTask::SortHitsByDet ( std::vector< std::vector< std::pair< Int_t, bool > > > &  hitsd,
Int_t  nStripHits 
)
private

Definition at line 78 of file PndLmdTrackFinderTask.cxx.

References fStripHitArray, fVerbose, PndSdsHit::GetSensorID(), and nSensPP.

Referenced by Exec().

80  {
81  Int_t nPlanes = 0;
82 
83  // sort in plane's
84  for (Int_t iHit = 0; iHit < nStripHits; iHit++) {
85  PndSdsHit *myHit = (PndSdsHit *)(fStripHitArray->At(iHit));
86 
87  Int_t sensid = myHit->GetSensorID(); // Sensors: 1..32
88  // cout<<"sensid = "<<sensid<<endl;
89  Int_t planeid = floor(
90  (sensid) / (double)nSensPP); // nSensPP sensors/plane => Planes: 0..3
91  hitsd.at(planeid).push_back(make_pair(iHit, false));
92  }
93 
94  for (Int_t iPlane = 0; iPlane < 4; iPlane++) {
95  if (hitsd.at(iPlane).size() > 0) nPlanes++;
96  }
97 
98  // cout << "Hits: " << nStripHits << endl;
99  if (fVerbose > 2) {
100  cout << "Hits: " << nStripHits << " in " << nPlanes << " plane(s)." << endl;
101  for (Int_t idet = 0; idet < 4; idet++)
102  cout << "Plane: " << idet << " DiscHits: " << hitsd.at(idet).size()
103  << endl;
104  }
105 
106  if (nPlanes > 2) return true;
107  return false;
108 }
int fVerbose
Definition: poormantracks.C:24
Int_t GetSensorID() const
Definition: PndSdsHit.h:90
bool PndLmdTrackFinderTask::SortHitsByDet2 ( std::vector< std::vector< std::pair< Int_t, bool > > > &  hitsd,
Int_t  nStripHits 
)
private

Definition at line 113 of file PndLmdTrackFinderTask.cxx.

References fStripHitArray, fVerbose, PndLmdGeometryHelper::getInstance(), and PndSdsHit::GetSensorID().

Referenced by Exec().

115  {
116  Int_t nPlanes = 0;
117 
118  // sort in plane's
119  for (Int_t iHit = 0; iHit < nStripHits; iHit++) {
120  PndSdsHit *myHit = (PndSdsHit *)(fStripHitArray->At(iHit));
121 
122  Int_t sensid = myHit->GetSensorID();
123 
124  auto& lmd_geo_helper(PndLmdGeometryHelper::getInstance());
125  auto const &hit_info = lmd_geo_helper.getHitLocationInfo(sensid);
126  // hitsd.at(iplane).push_back(iHit);
127  hitsd.at(hit_info.plane).push_back(make_pair(iHit, false));
128  // cout<<"sensid = "<<sensid<<endl;
129  // Int_t planeid = floor((sensid)/(double)nSensPP); //nSensPP sensors/plane
130  // => Planes: 0..3
131  // hitsd.at(planeid).push_back( make_pair (iHit,false) );
132  }
133 
134  for (Int_t iPlane = 0; iPlane < 4; iPlane++) {
135  if (hitsd.at(iPlane).size() > 0) nPlanes++;
136  }
137 
138  // cout << "Hits: " << nStripHits << endl;
139  if (fVerbose > 2) {
140  cout << "Hits: " << nStripHits << " in " << nPlanes << " plane(s)." << endl;
141  for (Int_t idet = 0; idet < 4; idet++)
142  cout << "Plane: " << idet << " DiscHits: " << hitsd.at(idet).size()
143  << endl;
144  }
145 
146  if (nPlanes > 2) return true;
147  return false;
148 }
int fVerbose
Definition: poormantracks.C:24
static PndLmdGeometryHelper & getInstance()
Int_t GetSensorID() const
Definition: PndSdsHit.h:90
bool PndLmdTrackFinderTask::SortHitsByZ ( std::vector< std::vector< std::pair< Int_t, bool > > > &  hitsd,
Int_t  nStripHits 
)
private

Definition at line 153 of file PndLmdTrackFinderTask.cxx.

References Double_t, fabs(), fStripHitArray, fVerbose, i, pos, and z.

155  {
156  std::vector<Double_t> detZ;
157 
158  // find plane positions
159  for (Int_t iHit = 0; iHit < nStripHits; iHit++) {
160  Double_t tmp = ((PndSdsHit *)(fStripHitArray->At(iHit)))->GetZ();
161  bool newZ = true;
162  for (size_t idet = 0; idet < detZ.size(); idet++) {
163  // if(tmp == detZ.at(idet)){ //check if already found
164  // cout<<"tmp = "<<tmp<<" detZ.at(idet) = "<<detZ.at(idet)
165  // <<" fabs(tmp-detZ.at(idet))="<<fabs(tmp-detZ.at(idet))<<endl;
166  if (fabs(tmp - detZ.at(idet)) <
167  9.) { // check if already found [for using with Dipole]
168  newZ = false;
169  }
170  }
171  if (newZ) {
172  detZ.push_back(tmp);
173 
174  // sort positions
175  Int_t pos = -1;
176  for (Int_t idet = detZ.size() - 1; idet >= 0; idet--) {
177  if (tmp < detZ.at(idet)) pos = idet;
178  }
179  if (pos != -1) {
180  Double_t swap = detZ.at(pos);
181  detZ.at(pos) = tmp;
182  tmp = swap;
183  for (size_t i = pos + 1; i < detZ.size(); i++) {
184  swap = detZ.at(i);
185  detZ.at(i) = tmp;
186  tmp = swap;
187  }
188  }
189  }
190  }
191  // cout<<"Attention! detZ.size()="<<detZ.size()<<endl;
192 
193  // sort in plane's
194  for (Int_t iHit = 0; iHit < nStripHits; iHit++) {
195  PndSdsHit *myHit = (PndSdsHit *)(fStripHitArray->At(iHit));
196 
197  Double_t z = myHit->GetZ();
198  for (size_t idet = 0; idet < detZ.size(); idet++) { // planes
199 
200  // if( z == detZ.at(idet) ){
201  if (fabs(z - detZ.at(idet)) < 9.) { //[for using with Dipole]
202  hitsd.at(idet).push_back(make_pair(iHit, false));
203  // cout<<"detZ.at("<<idet<<")="<<detZ.at(idet)<<" z="<<z<<endl;
204  }
205  }
206  }
207 
208  // cout << "Hits: " << nStripHits << endl;
209  if (fVerbose > 2) {
210  cout << "Hits: " << nStripHits << " in " << detZ.size() << " plane(s)."
211  << endl;
212  for (size_t idet = 0; idet < detZ.size(); idet++)
213  cout << "Plane: " << idet << " DiscHits: " << hitsd.at(idet).size()
214  << endl;
215  }
216 
217  if (detZ.size() > 3) return true;
218  return false;
219 }
TVector3 pos
int fVerbose
Definition: poormantracks.C:24
Int_t i
Definition: run_full.C:25
Double_t
Double_t z
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47

Member Data Documentation

Double_t PndLmdTrackFinderTask::dXY
private
Int_t PndLmdTrackFinderTask::fFinderMode
private

Definition at line 49 of file PndLmdTrackFinderTask.h.

Referenced by Exec(), and PndLmdTrackFinderTask().

TString PndLmdTrackFinderTask::fHitBranchStrip
private

Definition at line 56 of file PndLmdTrackFinderTask.h.

Referenced by Init(), and PndLmdTrackFinderTask().

bool PndLmdTrackFinderTask::flagPixelSens
private

Definition at line 47 of file PndLmdTrackFinderTask.h.

Referenced by Exec(), PndLmdTrackFinderTask(), and SetSensPixelFlag().

bool PndLmdTrackFinderTask::flagStipSens
private

Definition at line 43 of file PndLmdTrackFinderTask.h.

Referenced by Exec(), PndLmdTrackFinderTask(), and SetSensStripFlag().

TClonesArray* PndLmdTrackFinderTask::fStripHitArray
private

Input array of PndSdsDigis

Definition at line 59 of file PndLmdTrackFinderTask.h.

Referenced by Exec(), FindHitsI(), FindHitsII(), FindHitsIII(), Init(), SortHitsByDet(), SortHitsByDet2(), and SortHitsByZ().

TClonesArray* PndLmdTrackFinderTask::fTrackCandArray
private

Output array of PndSdsHits

Definition at line 62 of file PndLmdTrackFinderTask.h.

Referenced by Exec(), and Init().

Int_t PndLmdTrackFinderTask::nSensPP
private

Definition at line 50 of file PndLmdTrackFinderTask.h.

Referenced by SortHitsByDet().


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