CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Private Types | Private Member Functions | Private Attributes | Friends
CRackTrajectoryBuilder Class Reference

#include <CRackTrajectoryBuilder.h>

Classes

class  CompareDetByTraj
 

Public Member Functions

 CRackTrajectoryBuilder (const edm::ParameterSet &conf)
 
Trajectory createStartingTrajectory (const TrajectorySeed &seed) const
 
void init (const edm::EventSetup &es, bool)
 
void run (const TrajectorySeedCollection &collseed, const SiStripRecHit2DCollection &collstereo, const SiStripRecHit2DCollection &collrphi, const SiStripMatchedRecHit2DCollection &collmatched, const SiPixelRecHitCollection &collpixel, const edm::EventSetup &es, edm::Event &e, std::vector< Trajectory > &trajoutput)
 Runs the algorithm. More...
 
 ~CRackTrajectoryBuilder ()
 

Private Types

typedef std::pair
< TrackingRecHitRangeIterator,
TSOS
PairTrackingRecHitTsos
 
typedef TrajectoryMeasurement TM
 
typedef std::vector< const
TrackingRecHit * >::iterator 
TrackingRecHitIterator
 
typedef std::pair
< TrackingRecHitIterator,
TrackingRecHitIterator
TrackingRecHitRange
 
typedef std::vector
< TrackingRecHitRange >
::iterator 
TrackingRecHitRangeIterator
 
typedef TrajectoryStateOnSurface TSOS
 

Private Member Functions

void AddHit (Trajectory &traj, std::vector< const TrackingRecHit * >Hits, Propagator *currPropagator)
 
std::pair
< TrajectoryStateOnSurface,
const GeomDet * > 
innerState (const Trajectory &traj) const
 
bool isDifferentStripReHit2D (const SiStripRecHit2D &hitA, const SiStripRecHit2D &hitB)
 
bool qualityFilter (Trajectory traj)
 
std::vector
< TrajectoryMeasurement
seedMeasurements (const TrajectorySeed &seed) const
 
std::vector< const
TrackingRecHit * > 
SortHits (const SiStripRecHit2DCollection &collstereo, const SiStripRecHit2DCollection &collrphi, const SiStripMatchedRecHit2DCollection &collmatched, const SiPixelRecHitCollection &collpixel, const TrajectorySeed &seed, const bool bAddSeedHits)
 
TSOS startingTSOS (const TrajectorySeed &seed) const
 
void updateTrajectory (Trajectory &traj, const TM &tm, const TransientTrackingRecHit &hit) const
 

Private Attributes

double chi2cut
 
edm::ParameterSet conf_
 
bool debug_info
 
bool fastPropagation
 
std::string geometry
 
TransientTrackingRecHit::RecHitContainer hits
 
edm::ESHandle< MagneticFieldmagfield
 
const
TransientTrackingRecHitBuilder
RHBuilder
 
bool seed_plus
 
Chi2MeasurementEstimatortheEstimator
 
const KFTrajectoryFittertheFitter
 
int theMinHits
 
PropagatorWithMaterialthePropagator
 
PropagatorWithMaterialthePropagatorOp
 
const KFTrajectorySmoothertheSmoother
 
KFUpdatortheUpdator
 
edm::ESHandle< TrackerGeometrytracker
 
std::vector< TrajectorytrajFit
 
TrajectoryStateTransform tsTransform
 
bool useMatchedHits
 

Friends

class CompareDetByTraj
 

Detailed Description

Definition at line 128 of file CRackTrajectoryBuilder.h.

Member Typedef Documentation

Definition at line 141 of file CRackTrajectoryBuilder.h.

Definition at line 133 of file CRackTrajectoryBuilder.h.

typedef std::vector<const TrackingRecHit*>::iterator CRackTrajectoryBuilder::TrackingRecHitIterator
private

Definition at line 135 of file CRackTrajectoryBuilder.h.

Definition at line 137 of file CRackTrajectoryBuilder.h.

Definition at line 138 of file CRackTrajectoryBuilder.h.

Definition at line 132 of file CRackTrajectoryBuilder.h.

Constructor & Destructor Documentation

CRackTrajectoryBuilder::CRackTrajectoryBuilder ( const edm::ParameterSet conf)

Definition at line 32 of file CRackTrajectoryBuilder.cc.

References chi2cut, conf_, debug_info, fastPropagation, geometry, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), theMinHits, and useMatchedHits.

32  : conf_(conf) {
33  //minimum number of hits per tracks
34 
35  theMinHits=conf_.getParameter<int>("MinHits");
36  //cut on chi2
37  chi2cut=conf_.getParameter<double>("Chi2Cut");
38  edm::LogInfo("CosmicTrackFinder")<<"Minimum number of hits "<<theMinHits<<" Cut on Chi2= "<<chi2cut;
39 
40  debug_info=conf_.getUntrackedParameter<bool>("debug", false);
41  fastPropagation=conf_.getUntrackedParameter<bool>("fastPropagation", false);
42  useMatchedHits=conf_.getUntrackedParameter<bool>("useMatchedHits", true);
43 
44 
45  geometry=conf_.getUntrackedParameter<std::string>("GeometricStructure","STANDARD");
46 
47 
48 
49 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
CRackTrajectoryBuilder::~CRackTrajectoryBuilder ( )

Definition at line 52 of file CRackTrajectoryBuilder.cc.

52  {
53 // delete theInitialState;
54 }

Member Function Documentation

void CRackTrajectoryBuilder::AddHit ( Trajectory traj,
std::vector< const TrackingRecHit * >  Hits,
Propagator currPropagator 
)
private

do the old version ....

Definition at line 619 of file CRackTrajectoryBuilder.cc.

References TransientTrackingRecHitBuilder::build(), chi2cut, Trajectory::chiSquared(), CompareDetByTraj, gather_cfg::cout, debug_info, Chi2MeasurementEstimator::estimate(), fastPropagation, spr::find(), TrackingRecHit::geographicalId(), TrajectoryStateOnSurface::globalPosition(), hits, TrajectoryStateOnSurface::isValid(), Trajectory::lastMeasurement(), Propagator::propagate(), Trajectory::push(), DetId::rawId(), RHBuilder, theEstimator, theUpdator, tracker, KFUpdator::update(), and TrajectoryMeasurement::updatedState().

Referenced by run().

620  {
621 
622  if ( Hits.size() == 0 )
623  return;
624 
625  if (debug_info) cout << "CRackTrajectoryBuilder::AddHit" << endl;
626  if (debug_info) cout << "START " << traj.lastMeasurement().updatedState() << endl;
627 
628  vector <TrackingRecHitRange> hitRangeByDet;
629  TrackingRecHitIterator prevDet;
630 
631  prevDet = Hits.begin();
632  for( TrackingRecHitIterator iHit = Hits.begin(); iHit != Hits.end(); iHit++ )
633  {
634  if ( (*prevDet)->geographicalId() == (*iHit)->geographicalId() )
635  continue;
636 
637  hitRangeByDet.push_back( make_pair( prevDet, iHit ) );
638  prevDet = iHit;
639  }
640  hitRangeByDet.push_back( make_pair( prevDet, Hits.end() ) );
641 
643 
644  if (fastPropagation) {
645  for( TrackingRecHitRangeIterator iHitRange = hitRangeByDet.begin(); iHitRange != hitRangeByDet.end(); iHitRange++ )
646  {
647  const TrackingRecHit *currHit = *(iHitRange->first);
648  DetId currDet = currHit->geographicalId();
649 
650  TSOS prSt= currPropagator->propagate(traj.lastMeasurement().updatedState(),
651  tracker->idToDet(currDet)->surface());
652 
653  if ( !prSt.isValid())
654  {
655  if (debug_info) cout << "Not Valid: PRST" << prSt.globalPosition();
656  // if (debug_info) cout << "Not Valid: HIT" << *currHit;
657 
658 
659  continue;
660  }
661 
663  double chi2min = theEstimator->estimate( prSt, *bestHit).second;
664 
665  if (debug_info) cout << "Size " << iHitRange->first - (*iHitRange).second << endl;
666  for( TrackingRecHitIterator iHit = (*iHitRange).first+1; iHit != iHitRange->second; iHit++ )
667  {
668  if (debug_info) cout << "loop3 " <<" "<< Hits.end() - iHit << endl;
669 
671  double currChi2 = theEstimator->estimate(prSt, *tmpHit).second;
672  if ( currChi2 < chi2min )
673  {
674  currChi2 = chi2min;
675  bestHit = tmpHit;
676  }
677  }
678  //now we have check if the track can be added to the trajectory
679  if (debug_info) cout << chi2min << endl;
680  if (chi2min < chi2cut)
681  {
682  if (debug_info) cout << "chi2 fine : " << chi2min << endl;
683  TSOS UpdatedState= theUpdator->update( prSt, *bestHit );
684  if (UpdatedState.isValid()){
685  hits.push_back(&(*bestHit));
686  traj.push( TM(prSt,UpdatedState, bestHit, chi2min) );
687  if (debug_info) edm::LogInfo("CosmicTrackFinder") <<
688  "STATE UPDATED WITH HIT AT POSITION "
689 
690  << bestHit->globalPosition()
691  <<UpdatedState<<" "
692  <<traj.chiSquared();
693  if (debug_info) cout <<
694  "STATE UPDATED WITH HIT AT POSITION "
695 
696  << bestHit->globalPosition()
697  <<UpdatedState<<" "
698  <<traj.chiSquared();
699  if (debug_info) cout << "State is valid ..." << endl;
700  break; // now we need to
701  }
702  else
703  {
704  edm::LogWarning("CosmicTrackFinder")<<" State can not be updated with hit at position " << endl;
705  TSOS UpdatedState= theUpdator->update( prSt, *bestHit );
706  if (UpdatedState.isValid()){
707  cout <<
708  "NOT! UPDATED WITH HIT AT POSITION "
709 
710  << bestHit->globalPosition()
711  <<UpdatedState<<" "
712  <<traj.chiSquared();
713 
714  }
715  }
716  }
717  }
718  } //simple version end
719  else
720  {
721  //first sort the dets in the order they are traversed by the trajectory
722  // we need three loops:
723  // 1: loop as long as there can be an new hit added to the trajectory
724  // 2: loop over all dets that might be hit
725  // 3: loop over all hits on a certain det
726 
727 
728  std::vector < std::pair<TrackingRecHitRangeIterator, TSOS> > trackHitCandidates;
729  std::vector <std::pair<TrackingRecHitRangeIterator, TSOS> >::iterator iHitRange;
730  std::vector <uint32_t> processedDets;
731  do
732  {
733 
734  //create vector of possibly hit detectors...
735  trackHitCandidates.clear();
736  DetId currDet;
737  for( TrackingRecHitRangeIterator iHit = hitRangeByDet.begin(); iHit != hitRangeByDet.end(); iHit++ )
738  {
739  const TrackingRecHit *currHit = *(iHit->first);
740  currDet = currHit->geographicalId();
741 
742  if ( find(processedDets.begin(), processedDets.end(), currDet.rawId()) != processedDets.end() )
743  continue;
744 
745  TSOS prSt= currPropagator->propagate(traj.lastMeasurement().updatedState(),
746  tracker->idToDet(currDet)->surface());
747  if ( ( !prSt.isValid() ) || (theEstimator->Chi2MeasurementEstimatorBase::estimate(prSt,tracker->idToDet(currDet)->surface() ) == false) )
748  // if ( ( !prSt.isValid() ) || (theEstimator->estimate(prSt,tracker->idToDet(currDet)->surface() ) == false) )
749  continue;
750 
751  trackHitCandidates.push_back( make_pair(iHit, prSt) );
752  }
753 
754  if (!trackHitCandidates.size())
755  break;
756 
757  if (debug_info) cout << Hits.size() << " (int) trackHitCandidates.begin() " << trackHitCandidates.size() << endl;
758  if (debug_info) cout << "Before sorting ... " << endl;
759 
760  if (debug_info)
761  for( iHitRange = trackHitCandidates.begin(); iHitRange != trackHitCandidates.end(); iHitRange++ )
762  {
763  if (debug_info) cout << (tracker->idToDet((*(iHitRange->first->first))->geographicalId()))->position();
764  }
765  if (debug_info) cout << endl;
766 
767 
768  stable_sort( trackHitCandidates.begin(), trackHitCandidates.end(), CompareDetByTraj(traj.lastMeasurement().updatedState()) );
769 
770  if (debug_info) cout << "After sorting ... " << endl;
771  if (debug_info)
772  {
773  for( iHitRange = trackHitCandidates.begin(); iHitRange != trackHitCandidates.end(); iHitRange++ )
774  {
775  if (debug_info) cout << (tracker->idToDet((*(iHitRange->first->first))->geographicalId()))->position();
776  }
777  cout << endl;
778  }
779 
780  for( iHitRange = trackHitCandidates.begin(); iHitRange != trackHitCandidates.end(); iHitRange++ ) //loop over dets
781  {
782 
783  //now find the best hit of the detector
784  if (debug_info) cout << "loop2 " << trackHitCandidates.size() <<" " << trackHitCandidates.end() - iHitRange << endl;
785  const TrackingRecHit *currHit = *(iHitRange->first->first);
786 
788  TSOS currPrSt = (*iHitRange).second;
789 
790  if (debug_info) cout << "curr position" << bestHit->globalPosition();
791  for( TrackingRecHitIterator iHit = (*iHitRange).first->first+1; iHit != iHitRange->first->second; iHit++ )
792  {
794  if (debug_info) cout << "curr position" << tmpHit->globalPosition() ;
795 
796  }
797  }
798  if (debug_info) cout << "Cross check end ..." << endl;
799 
800 
801  //just a simple test if the same hit can be added twice ...
802  // for( iHitRange = trackHitCandidates.begin(); iHitRange != trackHitCandidates.end(); iHitRange++ ) //loop over all hits
803 
804  // break;
805 
806  for( iHitRange = trackHitCandidates.begin(); iHitRange != trackHitCandidates.end(); iHitRange++ ) //loop over detsall hits
807  {
808 
809  //now find the best hit of the detector
810  if (debug_info) cout << "loop2 " << trackHitCandidates.size() <<" " << trackHitCandidates.end() - iHitRange << endl;
811 
812  const TrackingRecHit *currHit = *(iHitRange->first->first);
813 
814  processedDets.push_back(currHit->geographicalId().rawId());
815 
816 
818 
819  if (debug_info) cout << "curr position A" << bestHit->globalPosition() << endl;
820  TSOS currPrSt = (*iHitRange).second;
821  double chi2min = theEstimator->estimate( currPrSt, *bestHit).second;
822 
823  if (debug_info) cout << "Size " << iHitRange->first->second - (*iHitRange).first->first << endl;
824  for( TrackingRecHitIterator iHit = (*iHitRange).first->first+1; iHit != iHitRange->first->second; iHit++ )
825  {
826  if (debug_info) cout << "loop3 " <<" "<< Hits.end() - iHit << endl;
827 
829  if (debug_info) cout << "curr position B" << tmpHit->globalPosition() << endl;
830  double currChi2 = theEstimator->estimate(currPrSt, *tmpHit).second;
831  if ( currChi2 < chi2min )
832  {
833  if (debug_info) cout << "Is best hit" << endl;
834  chi2min = currChi2;
835  bestHit = tmpHit;
836  }
837  }
838  //now we have checked the det and can remove the entry from the vector...
839 
840  //if (debug_info) cout << "before erase ..." << endl;
841  //this is to slow ...
842  // hitRangeByDet.erase( (*iHitRange).first,(*iHitRange).first+1 );
843  //if (debug_info) cout << "after erase ..." << endl;
844 
845  if (debug_info) cout << chi2min << endl;
846  //if the track can be added to the trajectory
847  if (chi2min < chi2cut)
848  {
849  if (debug_info) cout << "chi2 fine : " << chi2min << endl;
850 
851  // if (debug_info) cout << "previaous state " << traj.lastMeasurement().updatedState() <<endl;
852  TSOS UpdatedState= theUpdator->update( currPrSt, *bestHit );
853  if (UpdatedState.isValid()){
854 
855  hits.push_back(&(*bestHit));
856  traj.push( TM(currPrSt,UpdatedState, bestHit, chi2min) );
857  if (debug_info) edm::LogInfo("CosmicTrackFinder") <<
858  "STATE UPDATED WITH HIT AT POSITION "
859  // <<tmphitbestdet->globalPosition()
860  <<UpdatedState<<" "
861  <<traj.chiSquared();
862  if (debug_info) cout << "Added Hit" << bestHit->globalPosition() << endl;
863  if (debug_info) cout << "State is valid ..." << UpdatedState << endl;
864  //cout << "updated state " << traj.lastMeasurement().updatedState() <<endl;
865 
866  // return; //break;
867  //
868  // TSOS prSt= currPropagator->propagate(traj.lastMeasurement().updatedState(),
869  // tracker->idToDet( bestHit->geographicalId() )->surface());
870  //
871  // if ( prSt.isValid())
872  // cout << "the same hit can be added twice ..." << endl;
873  //
874 
875  break;
876  }
877  else
878  {
879  if (debug_info) edm::LogWarning("CosmicTrackFinder")<<" State can not be updated with hit at position "
880  << bestHit->globalPosition();
881  // cout << "State can not be updated with hit at " << bestHit->globalPosition() << endl;
882  }
883  continue;
884  }
885  else
886  {
887  // cout << "chi2 to big : " << chi2min << endl;
888  }
889  if (debug_info) cout << " continue 1 " << endl;
890  }
891  //now we remove all already processed dets from the list ...
892  // hitRangeByDet.erase( (*trackHitCandidates.begin()).first,(*iHitRange).first+1 );
893 
894  if (debug_info) cout << " continue 2 " << endl;
895  }
896  //if this was the last exit
897  while ( iHitRange != trackHitCandidates.end() );
898  }
899 
900 
901 }
std::vector< TrackingRecHitRange >::iterator TrackingRecHitRangeIterator
Chi2MeasurementEstimator * theEstimator
GlobalPoint globalPosition() const
const TransientTrackingRecHitBuilder * RHBuilder
TransientTrackingRecHit::RecHitContainer hits
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
std::vector< const TrackingRecHit * >::iterator TrackingRecHitIterator
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
edm::ESHandle< TrackerGeometry > tracker
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
TrajectoryMeasurement const & lastMeasurement() const
Definition: Trajectory.h:147
TrajectoryStateOnSurface updatedState() const
virtual std::pair< bool, double > estimate(const TrajectoryStateOnSurface &, const TransientTrackingRecHit &) const
virtual TrajectoryStateOnSurface propagate(const FreeTrajectoryState &, const Surface &) const
Definition: Propagator.cc:9
Definition: DetId.h:20
TrajectoryStateOnSurface update(const TrajectoryStateOnSurface &, const TransientTrackingRecHit &) const
Definition: KFUpdator.cc:10
tuple cout
Definition: gather_cfg.py:41
DetId geographicalId() const
void push(const TrajectoryMeasurement &tm)
Definition: Trajectory.cc:35
double chiSquared() const
Definition: Trajectory.h:208
TrajectoryMeasurement TM
Trajectory CRackTrajectoryBuilder::createStartingTrajectory ( const TrajectorySeed seed) const

Definition at line 255 of file CRackTrajectoryBuilder.cc.

References TrajectorySeed::direction(), i, query::result, and seedMeasurements().

Referenced by run(), and SortHits().

256 {
257  Trajectory result( seed, seed.direction());
258  std::vector<TM> seedMeas = seedMeasurements(seed);
259  if ( !seedMeas.empty()) {
260  for (std::vector<TM>::const_iterator i=seedMeas.begin(); i!=seedMeas.end(); i++){
261  result.push(*i);
262  }
263  }
264 
265  return result;
266 }
PropagationDirection direction() const
int i
Definition: DBlmapReader.cc:9
tuple result
Definition: query.py:137
std::vector< TrajectoryMeasurement > seedMeasurements(const TrajectorySeed &seed) const
void CRackTrajectoryBuilder::init ( const edm::EventSetup es,
bool  seedplus 
)

Definition at line 57 of file CRackTrajectoryBuilder.cc.

References alongMomentum, chi2cut, Chi2MeasurementEstimatorESProducer_cfi::Chi2MeasurementEstimator, conf_, edm::EventSetup::get(), edm::ParameterSet::getParameter(), KFTrajectoryFitterESProducer_cfi::KFTrajectoryFitter, KFTrajectorySmootherESProducer_cfi::KFTrajectorySmoother, magfield, oppositeToMomentum, edm::ESHandle< class >::product(), RHBuilder, seed_plus, theEstimator, theFitter, thePropagator, thePropagatorOp, theSmoother, theUpdator, and tracker.

Referenced by run().

57  {
58 
59 
60 // edm::ParameterSet tise_params = conf_.getParameter<edm::ParameterSet>("TransientInitialStateEstimatorParameters") ;
61 // theInitialState = new TransientInitialStateEstimator( es,tise_params);
62 
63  //services
66 
67 
68 
69  if (seedplus) {
70  seed_plus=true;
73  else {
74  seed_plus=false;
77  }
78 
79  theUpdator= new KFUpdator();
80 // theUpdator= new KFStripUpdator();
82 
83 
85  std::string builderName = conf_.getParameter<std::string>("TTRHBuilder");
86  es.get<TransientRecHitRecord>().get(builderName,theBuilder);
87 
88 
89  RHBuilder= theBuilder.product();
90 
91 
92 
93 
95  *theUpdator,
96  *theEstimator) ;
97 
98 
100  *theUpdator,
101  *theEstimator);
102 
103 }
T getParameter(std::string const &) const
edm::ESHandle< MagneticField > magfield
Chi2MeasurementEstimator * theEstimator
const TransientTrackingRecHitBuilder * RHBuilder
const KFTrajectorySmoother * theSmoother
edm::ESHandle< TrackerGeometry > tracker
PropagatorWithMaterial * thePropagatorOp
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
PropagatorWithMaterial * thePropagator
const KFTrajectoryFitter * theFitter
std::pair< TrajectoryStateOnSurface, const GeomDet * > CRackTrajectoryBuilder::innerState ( const Trajectory traj) const
private

Definition at line 950 of file CRackTrajectoryBuilder.cc.

References alongMomentum, funct::C, Chi2MeasurementEstimatorESProducer_cfi::Chi2MeasurementEstimator, Trajectory::direction(), Trajectory::foundHits(), i, TrajectoryStateOnSurface::localParameters(), PropagatorWithMaterial::magneticField(), Trajectory::measurements(), oppositeToMomentum, TrajectoryMeasurement::recHit(), TrajectoryStateOnSurface::surface(), thePropagator, and TrajectoryMeasurement::updatedState().

Referenced by run().

951 {
952  int lastFitted = 999;
953  int nhits = traj.foundHits();
954  if (nhits < lastFitted+1) lastFitted = nhits-1;
955 
956  std::vector<TrajectoryMeasurement> measvec = traj.measurements();
958 
959  bool foundLast = false;
960  int actualLast = -99;
961  for (int i=lastFitted; i >= 0; i--) {
962  if(measvec[i].recHit()->isValid()){
963  if(!foundLast){
964  actualLast = i;
965  foundLast = true;
966  }
967  firstHits.push_back( measvec[i].recHit());
968  }
969  }
970  TSOS unscaledState = measvec[actualLast].updatedState();
971  AlgebraicSymMatrix C(5,1);
972  // C *= 100.;
973 
974  TSOS startingState( unscaledState.localParameters(), LocalTrajectoryError(C),
975  unscaledState.surface(),
977 
978  // cout << endl << "FitTester starts with state " << startingState << endl;
979 
980  KFTrajectoryFitter backFitter( *thePropagator,
981  KFUpdator(),
982  Chi2MeasurementEstimator( 100., 3));
983 
985 
986  // only direction matters in this contest
989  backFitDirection);
990 
991  vector<Trajectory> fitres = backFitter.fit( fakeSeed, firstHits, startingState);
992 
993  if (fitres.size() != 1) {
994  // cout << "FitTester: first hits fit failed!" << endl;
995  return std::pair<TrajectoryStateOnSurface, const GeomDet*>();
996  }
997 
998  TrajectoryMeasurement firstMeas = fitres[0].lastMeasurement();
999  TSOS firstState = firstMeas.updatedState();
1000 
1001  // cout << "FitTester: Fitted first state " << firstState << endl;
1002  //cout << "FitTester: chi2 = " << fitres[0].chiSquared() << endl;
1003 
1004  TSOS initialState( firstState.localParameters(), LocalTrajectoryError(C),
1005  firstState.surface(),
1007 
1008  return std::pair<TrajectoryStateOnSurface, const GeomDet*>( initialState,
1009  firstMeas.recHit()->det());
1010 }
virtual const MagneticField * magneticField() const
int foundHits() const
Definition: Trajectory.h:190
int i
Definition: DBlmapReader.cc:9
const LocalTrajectoryParameters & localParameters() const
ConstRecHitPointer recHit() const
PropagationDirection
PropagationDirection const & direction() const
Definition: Trajectory.cc:195
DataContainer const & measurements() const
Definition: Trajectory.h:169
TrajectoryStateOnSurface updatedState() const
std::vector< ConstRecHitPointer > ConstRecHitContainer
PropagatorWithMaterial * thePropagator
const Surface & surface() const
CLHEP::HepSymMatrix AlgebraicSymMatrix
bool CRackTrajectoryBuilder::isDifferentStripReHit2D ( const SiStripRecHit2D hitA,
const SiStripRecHit2D hitB 
)
private

Definition at line 929 of file CRackTrajectoryBuilder.cc.

References TrackingRecHit::geographicalId(), BaseSiTrackerRecHit2DLocalPos::localPosition(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by SortHits().

930 {
931  if ( hitA.geographicalId() != hitB.geographicalId() )
932  return true;
933  if ( hitA.localPosition().x() != hitB.localPosition().x() )
934  return true;
935  if ( hitA.localPosition().y() != hitB.localPosition().y() )
936  return true;
937  if ( hitA.localPosition().z() != hitB.localPosition().z() )
938  return true;
939 
940  // if (debug_info) cout << hitA.localPosition() << endl;
941  // if (debug_info) cout << hitB << endl;
942 
943  return false;
944 }
T y() const
Definition: PV3DBase.h:57
virtual LocalPoint localPosition() const
T z() const
Definition: PV3DBase.h:58
DetId geographicalId() const
T x() const
Definition: PV3DBase.h:56
bool CRackTrajectoryBuilder::qualityFilter ( Trajectory  traj)
private

Definition at line 905 of file CRackTrajectoryBuilder.cc.

References Trajectory::foundHits(), geometry, hits, Trajectory::recHits(), and theMinHits.

Referenced by run().

905  {
906  int ngoodhits=0;
907  if(geometry=="MTCC"){
908  std::vector< ConstReferenceCountingPointer< TransientTrackingRecHit> > hits= traj.recHits();
909  std::vector< ConstReferenceCountingPointer< TransientTrackingRecHit> >::const_iterator hit;
910  for(hit=hits.begin();hit!=hits.end();hit++){
911  unsigned int iid=(*hit)->hit()->geographicalId().rawId();
912  //CHECK FOR 3 hits r-phi
913  if(((iid>>0)&0x3)!=1) ngoodhits++;
914  }
915  }
916  else ngoodhits=traj.foundHits();
917 
918  if ( ngoodhits >= theMinHits) {
919  return true;
920  }
921  else {
922  return false;
923  }
924 }
int foundHits() const
Definition: Trajectory.h:190
TransientTrackingRecHit::RecHitContainer hits
ConstRecHitContainer recHits(bool splitting=false) const
Definition: Trajectory.cc:67
void CRackTrajectoryBuilder::run ( const TrajectorySeedCollection collseed,
const SiStripRecHit2DCollection collstereo,
const SiStripRecHit2DCollection collrphi,
const SiStripMatchedRecHit2DCollection collmatched,
const SiPixelRecHitCollection collpixel,
const edm::EventSetup es,
edm::Event e,
std::vector< Trajectory > &  trajoutput 
)

Runs the algorithm.

Definition at line 106 of file CRackTrajectoryBuilder.cc.

References AddHit(), alongMomentum, TransientTrackingRecHitBuilder::build(), gather_cfg::cout, createStartingTrajectory(), debug_info, Trajectory::firstMeasurement(), KFTrajectoryFitter::fit(), Trajectory::foundHits(), hits, init(), innerState(), iseed, PropagatorWithMaterial::propagate(), qualityFilter(), Trajectory::recHits(), RHBuilder, Trajectory::seed(), seed_plus, SortHits(), theEstimator, theFitter, thePropagator, thePropagatorOp, theSmoother, theUpdator, tracker, KFTrajectorySmoother::trajectories(), trajFit, and TrajectoryMeasurement::updatedState().

Referenced by cms::CosmicTrackFinder::produce().

114 {
115 
116  std::vector<Trajectory> trajSmooth;
117  std::vector<Trajectory>::iterator trajIter;
118 
119  TrajectorySeedCollection::const_iterator iseed;
120  unsigned int IS=0;
121  for(iseed=collseed.begin();iseed!=collseed.end();iseed++){
122  bool seedplus=((*iseed).direction()==alongMomentum);
123  init(es,seedplus);
124  hits.clear();
125  trajFit.clear();
126  trajSmooth.clear();
127 
128  Trajectory startingTraj = createStartingTrajectory(*iseed);
129  Trajectory trajTmp; // used for the opposite direction
130  Trajectory traj = startingTraj;
131 
132 
133  //first propagate track in opposite direction ...
134  seed_plus = !seed_plus;
135  vector<const TrackingRecHit*> allHitsOppsite = SortHits(collstereo,collrphi,collmatched,collpixel, *iseed, true);
136  seed_plus = !seed_plus;
137  if (allHitsOppsite.size())
138  {
139  //there are hits which are above the seed,
140  //cout << "Number of hits higher than seed " <<allHitsOppsite.size() << endl;
141 
142  AddHit( traj, allHitsOppsite, thePropagatorOp);
143 
144 
145  if (debug_info)
146  {
147  cout << "Hits in opposite direction..." << endl;
148  TransientTrackingRecHit::RecHitContainer::const_iterator iHit;
149  for ( iHit = hits.begin(); iHit!=hits.end(); iHit++)
150  {
151  cout << (**iHit).globalPosition() << endl;
152  }
153  }
154 
155 
156  // now add hist opposite to seed
157  //now crate new starting trajectory
158  reverse(hits.begin(), hits.end());
160  tracker->idToDet((hits.front())->geographicalId())->surface()).isValid()){
161  TSOS startingStateNew = //TrajectoryStateWithArbitraryError()//
163  tracker->idToDet((hits.front())->geographicalId())->surface()));
164 
165  if (debug_info)
166  {
167  cout << "Hits in opposite direction reversed..." << endl;
168  TransientTrackingRecHit::RecHitContainer::const_iterator iHit;
169  for ( iHit = hits.begin(); iHit!=hits.end(); iHit++)
170  {
171  cout << (**iHit).globalPosition() << endl;
172  }
173  }
174 
175  const TrajectorySeed& tmpseed=traj.seed();
176  trajTmp = theFitter->fit(tmpseed,hits, startingStateNew ).front();
177 
178  if(debug_info){
179  cout << "Debugging show fitted hits" << endl;
180  std::vector< ConstReferenceCountingPointer< TransientTrackingRecHit> > hitsFit= trajTmp.recHits();
181  std::vector< ConstReferenceCountingPointer< TransientTrackingRecHit> >::const_iterator hit;
182  for(hit=hitsFit.begin();hit!=hitsFit.end();hit++){
183 
184  cout << RHBuilder->build( &(*(*hit)->hit()) )->globalPosition() << endl;
185  }
186  }
187 
188  }
189  }
190  else
191  {
192  if(debug_info) cout << "There are no hits in opposite direction ..." << endl;
193  }
194 
195  vector<const TrackingRecHit*> allHits;
196  if (trajTmp.foundHits())
197  {
198  traj = trajTmp;
199  allHits= SortHits(collstereo,collrphi,collmatched,collpixel,*iseed, false);
200  }
201  else
202  {
203  traj = startingTraj;
204  hits.clear();
205  allHits= SortHits(collstereo,collrphi,collmatched,collpixel,*iseed, true);
206  }
207 
208  AddHit( traj,allHits, thePropagator);
209 
210 
211  if(debug_info){
212  cout << "Debugging show All fitted hits" << endl;
213  std::vector< ConstReferenceCountingPointer< TransientTrackingRecHit> > hits= traj.recHits();
214  std::vector< ConstReferenceCountingPointer< TransientTrackingRecHit> >::const_iterator hit;
215  for(hit=hits.begin();hit!=hits.end();hit++){
216 
217  cout << RHBuilder->build( &(*(*hit)->hit()) )->globalPosition() << endl;
218  }
219 
220  cout << qualityFilter( traj) << " <- quality filter good?" << endl;
221  }
222 
223 
224  if (debug_info) cout << "now do quality checks" << endl;
225  if ( qualityFilter( traj) ) {
226  const TrajectorySeed& tmpseed=traj.seed();
227  std::pair<TrajectoryStateOnSurface, const GeomDet*> initState = innerState(traj); //theInitialState->innerState(traj);
228  if (initState.first.isValid())
229  trajFit = theFitter->fit(tmpseed,hits, initState.first);
230  }
231 
232  for (trajIter=trajFit.begin(); trajIter!=trajFit.end();trajIter++){
233  trajSmooth=theSmoother->trajectories((*trajIter));
234  }
235  for (trajIter= trajSmooth.begin(); trajIter!=trajSmooth.end();trajIter++){
236  if((*trajIter).isValid()){
237 
238  if (debug_info) cout << "adding track ... " << endl;
239  trajoutput.push_back((*trajIter));
240  }
241  }
242  delete theUpdator;
243  delete theEstimator;
244  delete thePropagator;
245  delete thePropagatorOp;
246  delete theFitter;
247  delete theSmoother;
248  //Only the first 30 seeds are considered
249  if (IS>30) return;
250  IS++;
251 
252  }
253 }
int foundHits() const
Definition: Trajectory.h:190
Trajectory createStartingTrajectory(const TrajectorySeed &seed) const
bool qualityFilter(Trajectory traj)
TrajectorySeed const & seed() const
Access to the seed used to reconstruct the Trajectory.
Definition: Trajectory.h:231
Chi2MeasurementEstimator * theEstimator
const TransientTrackingRecHitBuilder * RHBuilder
TransientTrackingRecHit::RecHitContainer hits
ConstRecHitContainer recHits(bool splitting=false) const
Definition: Trajectory.cc:67
const KFTrajectorySmoother * theSmoother
virtual std::vector< Trajectory > trajectories(const Trajectory &aTraj) const
std::pair< TrajectoryStateOnSurface, const GeomDet * > innerState(const Trajectory &traj) const
std::vector< Trajectory > trajFit
virtual std::vector< Trajectory > fit(const Trajectory &aTraj) const
edm::ESHandle< TrackerGeometry > tracker
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
PropagatorWithMaterial * thePropagatorOp
TrajectoryStateOnSurface updatedState() const
virtual TrajectoryStateOnSurface propagate(const TrajectoryStateOnSurface &tsos, const Plane &plane) const
void AddHit(Trajectory &traj, std::vector< const TrackingRecHit * >Hits, Propagator *currPropagator)
TrajectoryMeasurement const & firstMeasurement() const
Definition: Trajectory.h:160
int iseed
Definition: AMPTWrapper.h:124
void init(const edm::EventSetup &es, bool)
PropagatorWithMaterial * thePropagator
const KFTrajectoryFitter * theFitter
std::vector< const TrackingRecHit * > SortHits(const SiStripRecHit2DCollection &collstereo, const SiStripRecHit2DCollection &collrphi, const SiStripMatchedRecHit2DCollection &collmatched, const SiPixelRecHitCollection &collpixel, const TrajectorySeed &seed, const bool bAddSeedHits)
tuple cout
Definition: gather_cfg.py:41
std::vector< TrajectoryMeasurement > CRackTrajectoryBuilder::seedMeasurements ( const TrajectorySeed seed) const
private

Definition at line 270 of file CRackTrajectoryBuilder.cc.

References TransientTrackingRecHitBuilder::build(), TrajectorySeed::recHits(), query::result, RHBuilder, startingTSOS(), and GeomDet::surface().

Referenced by createStartingTrajectory().

271 {
272  std::vector<TrajectoryMeasurement> result;
273  TrajectorySeed::range hitRange = seed.recHits();
274  for (TrajectorySeed::const_iterator ihit = hitRange.first;
275  ihit != hitRange.second; ihit++) {
276  //RC TransientTrackingRecHit* recHit = RHBuilder->build(&(*ihit));
278  const GeomDet* hitGeomDet = (&(*tracker))->idToDet( ihit->geographicalId());
279  TSOS invalidState( new BasicSingleTrajectoryState( hitGeomDet->surface()));
280 
281  if (ihit == hitRange.second - 1) {
282  TSOS updatedState=startingTSOS(seed);
283  result.push_back(TM( invalidState, updatedState, recHit));
284 
285  }
286  else {
287  result.push_back(TM( invalidState, recHit));
288  }
289 
290  }
291 
292  return result;
293 }
const TransientTrackingRecHitBuilder * RHBuilder
recHitContainer::const_iterator const_iterator
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
tuple result
Definition: query.py:137
std::pair< const_iterator, const_iterator > range
TSOS startingTSOS(const TrajectorySeed &seed) const
range recHits() const
virtual const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
TrajectoryMeasurement TM
vector< const TrackingRecHit * > CRackTrajectoryBuilder::SortHits ( const SiStripRecHit2DCollection collstereo,
const SiStripRecHit2DCollection collrphi,
const SiStripMatchedRecHit2DCollection collmatched,
const SiPixelRecHitCollection collpixel,
const TrajectorySeed seed,
const bool  bAddSeedHits 
)
private

Definition at line 296 of file CRackTrajectoryBuilder.cc.

References TransientTrackingRecHitBuilder::build(), gather_cfg::cout, createStartingTrajectory(), edmNew::DetSetVector< T >::data(), debug_info, spr::find(), TrackingRecHit::geographicalId(), TrajectoryStateOnSurface::globalPosition(), hits, isDifferentStripReHit2D(), TrajectoryStateOnSurface::isValid(), Trajectory::lastMeasurement(), TrajectoryStateOnSurface::localError(), LogDebug, LocalTrajectoryError::matrix(), SiStripMatchedRecHit2D::monoHit(), TrajectorySeed::nHits(), PropagatorWithMaterial::propagate(), DetId::rawId(), TrajectorySeed::recHits(), RHBuilder, seed_plus, startingTSOS(), SiStripMatchedRecHit2D::stereoHit(), thePropagator, tracker, TrajectoryMeasurement::updatedState(), and useMatchedHits.

Referenced by run().

302  {
303 
304 
305  //The Hits with global y more than the seed are discarded
306  //The Hits correspondign to the seed are discarded
307  //At the end all the hits are sorted in y
308  vector<const TrackingRecHit*> allHits;
309 
311  TrajectorySeed::range hRange= seed.recHits();
313  float yref=0.;
314 
315  if (debug_info) cout << "SEED " << startingTSOS(seed) << endl;
316  if (debug_info) cout << "seed hits size " << seed.nHits() << endl;
317 
318 
319  //seed debugging:
320  // GeomDet *seedDet = tracker->idToDet(seed.startingState().detId());
321 
322 
323 // edm::LogInfo("CRackTrajectoryBuilder::SortHits") << "SEED " << seed.startingState();
324 
325  edm::LogInfo("CRackTrajectoryBuilder::SortHits") << "SEED " << startingTSOS(seed);
326 // edm::LogInfo("CRackTrajectoryBuilder::SortHits" << "seed hits size " << seed.nHits();
327 
328 
329 
330  // if (seed.nHits()<2)
331  // return allHits;
332 
333  float_t yMin=0.;
334  float_t yMax=0.;
335 
336  int seedHitSize= hRange.second - hRange.first;
337 
338  vector <int> detIDSeedMatched (seedHitSize);
339  vector <int> detIDSeedRphi (seedHitSize);
340  vector <int> detIDSeedStereo (seedHitSize);
341 
342  for (ihit = hRange.first;
343  ihit != hRange.second; ihit++) {
344 
345  // need to find track with lowest (seed_plus)/ highest y (seed_minus)
346  // split matched hits ...
347  const SiStripMatchedRecHit2D* matchedhit=dynamic_cast<const SiStripMatchedRecHit2D*>(&(*ihit));
348 
349  yref=RHBuilder->build(&(*ihit))->globalPosition().y();
350  if (ihit == hRange.first)
351  {
352  yMin = yref;
353  yMax = yref;
354  }
355 
356  if (matchedhit)
357  {
358  float_t yGlobRPhi = RHBuilder->build( (matchedhit->monoHit()))->globalPosition().y();
359  float_t yGlobStereo = RHBuilder->build( (matchedhit->stereoHit()))->globalPosition().y();
360 
361  if (debug_info) cout << "Rphi ..." << yGlobRPhi << endl;
362  if (debug_info) cout << "Stereo ..." << yGlobStereo << endl;
363 
364  if ( yGlobStereo < yMin ) yMin = yGlobStereo;
365  if ( yGlobRPhi < yMin ) yMin = yGlobRPhi;
366 
367  if ( yGlobStereo > yMax ) yMax = yGlobStereo;
368  if ( yGlobRPhi > yMax ) yMax = yGlobRPhi;
369 
370  detIDSeedMatched.push_back ( matchedhit->geographicalId().rawId() );
371  detIDSeedRphi.push_back ( matchedhit->monoHit()->geographicalId().rawId() );
372  detIDSeedStereo.push_back ( matchedhit->stereoHit()->geographicalId().rawId() );
373 
374  if (bAddSeedHits)
375  {
376  if (useMatchedHits)
377  {
378 
379  hits.push_back((RHBuilder->build(&(*ihit))));
380  }
381  else
382  {
383  if ( ( (yGlobRPhi > yGlobStereo ) && seed_plus ) || ((yGlobRPhi < yGlobStereo ) && !seed_plus ))
384  {
385  hits.push_back((RHBuilder->build(matchedhit->monoHit())));
386  hits.push_back((RHBuilder->build(matchedhit->stereoHit())));
387  }
388  else
389  {
390  hits.push_back((RHBuilder->build(matchedhit->stereoHit())));
391  hits.push_back((RHBuilder->build(matchedhit->monoHit())));
392 
393  }
394  }
395  }
396  }
397  else if (bAddSeedHits)
398  {
399  hits.push_back((RHBuilder->build(&(*ihit))));
400  detIDSeedRphi.push_back ( ihit->geographicalId().rawId() );
401  detIDSeedMatched.push_back ( -1 );
402  detIDSeedStereo.push_back ( -1 );
403 
404  }
405 
406  if ( yref < yMin ) yMin = yref;
407  if ( yref > yMax ) yMax = yref;
408 
409 // if (bAddSeedHits)
410 // hits.push_back((RHBuilder->build(&(*ihit))));
411 
412  LogDebug("CosmicTrackFinder")<<"SEED HITS"<<RHBuilder->build(&(*ihit))->globalPosition();
413  if (debug_info) cout <<"SEED HITS"<<RHBuilder->build(&(*ihit))->globalPosition() << endl;
414 // if (debug_info) cout <<"SEED HITS"<< seed.startingState().parameters() << endl;
415 
416 
417  }
418 
419  yref = (seed_plus) ? yMin : yMax;
420 
421  if ((&collpixel)!=0){
423  for(ipix=collpixel.data().begin();ipix!=collpixel.data().end();ipix++){
424  float ych= RHBuilder->build(&(*ipix))->globalPosition().y();
425  if ((seed_plus && (ych<yref)) || (!(seed_plus) && (ych>yref)))
426  allHits.push_back(&(*ipix));
427  }
428  }
429 
430 
431  if (useMatchedHits) // use matched
432  {
433  //add the matched hits ...
435 
436  if ((&collmatched)!=0){
437  for(istripm=collmatched.data().begin();istripm!=collmatched.data().end();istripm++){
438  float ych= RHBuilder->build(&(*istripm))->globalPosition().y();
439 
440  int cDetId=istripm->geographicalId().rawId();
441  bool noSeedDet = ( detIDSeedMatched.end() == find (detIDSeedMatched.begin(), detIDSeedMatched.end(), cDetId ) ) ;
442 
443  if ( noSeedDet )
444  if ((seed_plus && (ych<yref)) || (!(seed_plus) && (ych>yref)))
445  {
446  //if (debug_info) cout << "adding matched hit " << &(*istripm) << endl;
447  allHits.push_back(&(*istripm));
448  }
449  }
450  }
451 
452  //add the rpi hits, but only accept hits that are not matched hits
453  if ((&collrphi)!=0){
454  for(istrip=collrphi.data().begin();istrip!=collrphi.data().end();istrip++){
455  float ych= RHBuilder->build(&(*istrip))->globalPosition().y();
456  StripSubdetector monoDetId(istrip->geographicalId());
457  if (monoDetId.partnerDetId())
458  {
459  edm::LogInfo("CRackTrajectoryBuilder::SortHits") << "this det belongs to a glued det " << ych << endl;
460  continue;
461  }
462  int cDetId=istrip->geographicalId().rawId();
463  bool noSeedDet = ( detIDSeedRphi.end()== find (detIDSeedRphi.begin(), detIDSeedRphi.end(), cDetId ) ) ;
464  if (noSeedDet)
465  if ((seed_plus && (ych<yref)) || (!(seed_plus) && (ych>yref)))
466  {
467 
468  bool hitIsUnique = true;
469  //now
470  if ((&collmatched)!=0)
471  for(istripm=collmatched.data().begin();istripm!=collmatched.data().end();istripm++)
472  {
473  // if ( isDifferentStripReHit2D ( *istrip, * (istripm->stereoHit() ) ) == false)
474  if ( isDifferentStripReHit2D ( *istrip, * (istripm->monoHit() ) ) == false)
475  {
476  hitIsUnique = false;
477  edm::LogInfo("CRackTrajectoryBuilder::SortHits") << "rphi hit is in matched hits; y: " << ych << endl;
478  break;
479  }
480  } //end loop over all matched
481  if (hitIsUnique)
482  {
483  // if (debug_info) cout << "adding rphi hit " << &(*istrip) << endl;
484  allHits.push_back(&(*istrip));
485  }
486  }
487  }
488  }
489 
490 
491  //add the stereo hits except the hits that are in the matched collection
492  //update do not use unmatched rphi hist due to limitation of alignment framework
493  //if (!useMatchedHits)
494  //if ((&collstereo)!=0){
495  // for(istrip=collstereo.data().begin();istrip!=collstereo.data().end();istrip++){
496  // float ych= RHBuilder->build(&(*istrip))->globalPosition().y();
497  //
498  //
499  // int cDetId = istrip->geographicalId().rawId();
500  // bool noSeedDet = ( detIDSeedStereo.end()== find (detIDSeedStereo.begin(), detIDSeedStereo.end(), cDetId ) ) ;
501  //
502  // if (noSeedDet)
503  // if ((seed_plus && (ych<yref)) || (!(seed_plus) && (ych>yref)))
504  // {
505  //
506  // bool hitIsUnique = true;
507  // //now
508  // if ((&collmatched)!=0)
509  // for(istripm=collmatched.data().begin();istripm!=collmatched.data().end();istripm++)
510  // {
511  // if ( isDifferentStripReHit2D ( *istrip, * (istripm->stereoHit() ) ) == false)
512  // {
513  // hitIsUnique = false;
514  // edm::LogInfo("CRackTrajectoryBuilder::SortHits") << "stereo hit already in matched hits; y: " << ych << endl;
515  // break;
516  // }
517  // } //end loop over all stereo
518  // if (hitIsUnique)
519  // {
520  //
521  // // if (debug_info) cout << "now I am adding a stero hit, either noise or not in overlap ...!!!!" << endl;
522  // allHits.push_back(&(*istrip));
523  // }
524  // }
525  // }
526  //}
527  }
528  else // dont use matched ...
529  {
530 
531  if ((&collrphi)!=0){
532  for(istrip=collrphi.data().begin();istrip!=collrphi.data().end();istrip++){
533  float ych= RHBuilder->build(&(*istrip))->globalPosition().y();
534  if ((seed_plus && (ych<yref)) || (!(seed_plus) && (ych>yref)))
535  allHits.push_back(&(*istrip));
536  }
537  }
538 
539 
540  if ((&collstereo)!=0){
541  for(istrip=collstereo.data().begin();istrip!=collstereo.data().end();istrip++){
542  float ych= RHBuilder->build(&(*istrip))->globalPosition().y();
543  if ((seed_plus && (ych<yref)) || (!(seed_plus) && (ych>yref)))
544  allHits.push_back(&(*istrip));
545  }
546  }
547 
548  }
549 
550 
551 // if (seed_plus){
552 // stable_sort(allHits.begin(),allHits.end(),CompareDetY_plus(*tracker));
553 // }
554 // else {
555 // stable_sort(allHits.begin(),allHits.end(),CompareDetY_minus(*tracker));
556 // }
557 
558 
559  if (seed_plus){
560  stable_sort(allHits.begin(),allHits.end(),CompareHitY_plus(*tracker));
561  }
562  else {
563  stable_sort(allHits.begin(),allHits.end(),CompareHitY(*tracker));
564  }
565 
566 
567 
568  if (debug_info)
569  {
570  if (debug_info) cout << "all hits" << endl;
571 
572  //starting trajectory
573  Trajectory startingTraj = createStartingTrajectory(seed);
574 
575  if (debug_info) cout << "START " << startingTraj.lastMeasurement().updatedState() << endl;
576  if (debug_info) cout << "START Err" << startingTraj.lastMeasurement().updatedState().localError().matrix() << endl;
577 
578 
579  vector<const TrackingRecHit*>::iterator iHit;
580  for (iHit=allHits.begin(); iHit<allHits.end(); iHit++)
581  {
582  GlobalPoint gphit=RHBuilder->build(*iHit)->globalPosition();
583  if (debug_info) cout << "GH " << gphit << endl;
584 
585  // tracker->idToDet((*iHit)->geographicalId())->surface();
586 
587  TSOS prSt = thePropagator->propagate(startingTraj.lastMeasurement().updatedState(),
588  tracker->idToDet((*iHit)->geographicalId())->surface());
589 
590  if(prSt.isValid())
591  {
592  if (debug_info) cout << "PR " << prSt.globalPosition() << endl;
593  //if (debug_info) cout << "PR Err" << prSt.localError().matrix() << endl;
594 
595  }
596  else
597  {
598  if (debug_info) cout << "not valid" << endl;
599  }
600  }
601  if (debug_info) cout << "all hits end" << endl;
602  }
603 
604 
605  return allHits;
606 }
#define LogDebug(id)
Trajectory createStartingTrajectory(const TrajectorySeed &seed) const
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const SiStripRecHit2D * stereoHit() const
GlobalPoint globalPosition() const
const TransientTrackingRecHitBuilder * RHBuilder
TransientTrackingRecHit::RecHitContainer hits
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
recHitContainer::const_iterator const_iterator
edm::ESHandle< TrackerGeometry > tracker
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
TrajectoryMeasurement const & lastMeasurement() const
Definition: Trajectory.h:147
data_type const * data(size_t cell) const
std::pair< const_iterator, const_iterator > range
TrajectoryStateOnSurface updatedState() const
const AlgebraicSymMatrix55 & matrix() const
const LocalTrajectoryError & localError() const
TSOS startingTSOS(const TrajectorySeed &seed) const
virtual TrajectoryStateOnSurface propagate(const TrajectoryStateOnSurface &tsos, const Plane &plane) const
bool isDifferentStripReHit2D(const SiStripRecHit2D &hitA, const SiStripRecHit2D &hitB)
range recHits() const
PropagatorWithMaterial * thePropagator
unsigned int nHits() const
tuple cout
Definition: gather_cfg.py:41
DetId geographicalId() const
const SiStripRecHit2D * monoHit() const
TrajectoryStateOnSurface CRackTrajectoryBuilder::startingTSOS ( const TrajectorySeed seed) const
private

Definition at line 609 of file CRackTrajectoryBuilder.cc.

References magfield, TrajectorySeed::startingState(), TrajectoryStateTransform::transientState(), and tsTransform.

Referenced by seedMeasurements(), and SortHits().

610 {
611  PTrajectoryStateOnDet pState( seed.startingState());
612  const GeomDet* gdet = (&(*tracker))->idToDet(DetId(pState.detId()));
613  TSOS State= tsTransform.transientState( pState, &(gdet->surface()),
614  &(*magfield));
615  return State;
616 
617 }
edm::ESHandle< MagneticField > magfield
TrajectoryStateOnSurface transientState(const PTrajectoryStateOnDet &ts, const Surface *surface, const MagneticField *field) const
TrajectoryStateTransform tsTransform
Definition: DetId.h:20
PTrajectoryStateOnDet const & startingState() const
void CRackTrajectoryBuilder::updateTrajectory ( Trajectory traj,
const TM tm,
const TransientTrackingRecHit hit 
) const
private

Friends And Related Function Documentation

friend class CompareDetByTraj
friend

Definition at line 144 of file CRackTrajectoryBuilder.h.

Referenced by AddHit().

Member Data Documentation

double CRackTrajectoryBuilder::chi2cut
private

Definition at line 259 of file CRackTrajectoryBuilder.h.

Referenced by AddHit(), CRackTrajectoryBuilder(), and init().

edm::ParameterSet CRackTrajectoryBuilder::conf_
private

Definition at line 239 of file CRackTrajectoryBuilder.h.

Referenced by CRackTrajectoryBuilder(), and init().

bool CRackTrajectoryBuilder::debug_info
private

Definition at line 254 of file CRackTrajectoryBuilder.h.

Referenced by AddHit(), CRackTrajectoryBuilder(), run(), and SortHits().

bool CRackTrajectoryBuilder::fastPropagation
private

Definition at line 255 of file CRackTrajectoryBuilder.h.

Referenced by AddHit(), and CRackTrajectoryBuilder().

std::string CRackTrajectoryBuilder::geometry
private
TransientTrackingRecHit::RecHitContainer CRackTrajectoryBuilder::hits
private

Definition at line 262 of file CRackTrajectoryBuilder.h.

Referenced by AddHit(), qualityFilter(), run(), and SortHits().

edm::ESHandle<MagneticField> CRackTrajectoryBuilder::magfield
private

Definition at line 237 of file CRackTrajectoryBuilder.h.

Referenced by init(), and startingTSOS().

const TransientTrackingRecHitBuilder* CRackTrajectoryBuilder::RHBuilder
private

Definition at line 249 of file CRackTrajectoryBuilder.h.

Referenced by AddHit(), init(), run(), seedMeasurements(), and SortHits().

bool CRackTrajectoryBuilder::seed_plus
private

Definition at line 263 of file CRackTrajectoryBuilder.h.

Referenced by init(), run(), and SortHits().

Chi2MeasurementEstimator* CRackTrajectoryBuilder::theEstimator
private

Definition at line 248 of file CRackTrajectoryBuilder.h.

Referenced by AddHit(), init(), and run().

const KFTrajectoryFitter* CRackTrajectoryBuilder::theFitter
private

Definition at line 251 of file CRackTrajectoryBuilder.h.

Referenced by init(), and run().

int CRackTrajectoryBuilder::theMinHits
private

Definition at line 258 of file CRackTrajectoryBuilder.h.

Referenced by CRackTrajectoryBuilder(), and qualityFilter().

PropagatorWithMaterial* CRackTrajectoryBuilder::thePropagator
private

Definition at line 241 of file CRackTrajectoryBuilder.h.

Referenced by init(), innerState(), run(), and SortHits().

PropagatorWithMaterial* CRackTrajectoryBuilder::thePropagatorOp
private

Definition at line 242 of file CRackTrajectoryBuilder.h.

Referenced by init(), and run().

const KFTrajectorySmoother* CRackTrajectoryBuilder::theSmoother
private

Definition at line 250 of file CRackTrajectoryBuilder.h.

Referenced by init(), and run().

KFUpdator* CRackTrajectoryBuilder::theUpdator
private

Definition at line 247 of file CRackTrajectoryBuilder.h.

Referenced by AddHit(), init(), and run().

edm::ESHandle<TrackerGeometry> CRackTrajectoryBuilder::tracker
private

Definition at line 238 of file CRackTrajectoryBuilder.h.

Referenced by AddHit(), init(), run(), and SortHits().

std::vector<Trajectory> CRackTrajectoryBuilder::trajFit
private

Definition at line 260 of file CRackTrajectoryBuilder.h.

Referenced by run().

TrajectoryStateTransform CRackTrajectoryBuilder::tsTransform
private

Definition at line 240 of file CRackTrajectoryBuilder.h.

Referenced by startingTSOS().

bool CRackTrajectoryBuilder::useMatchedHits
private

Definition at line 256 of file CRackTrajectoryBuilder.h.

Referenced by CRackTrajectoryBuilder(), and SortHits().