CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CosmicTrackFinder.h
Go to the documentation of this file.
1 #ifndef CosmicTrackFinder_h
2 #define CosmicTrackFinder_h
3 
4 // Package: RecoTracker/SingleTrackPattern
5 // Class: CosmicTrackFinder
6 // Original Author: Michele Pioppi-INFN perugia
7 
8 
19 
20 namespace cms
21 {
23  public:
25  Trajectory *t2){
26  AnalHits(t1->recHits());
27  unsigned int alay=nlay;
28  AnalHits(t2->recHits());
29  unsigned int blay=nlay;
30  if (alay!=blay) return alay > blay;
31  if (t1->foundHits() != t2->foundHits())
32  return t1->foundHits()> t2->foundHits();
33  return t1->chiSquared()< t2->chiSquared();
34  // std::cout<<"chi "<<t1.chiSquared()<<" "<<t2.chiSquared()<<std::endl;
35  // return false;
36  }
38  ltob1=false; ltob2=false; ltib1=false; ltib2=false;
39  std::vector< ConstReferenceCountingPointer< TransientTrackingRecHit> >::const_iterator hit;
40  // ConstRecHitIterator hit;
41  for(hit=hits.begin();hit!=hits.end();hit++){
42  unsigned int iid=(*hit)->hit()->geographicalId().rawId();
43 
44  int sub=(iid>>25)&0x7 ;
45  int lay=(iid>>16) & 0xF;
46  if ((lay==1)&&(sub==3)) ltib1=true;
47  if ((lay==2)&&(sub==3)) ltib2=true;
48  if ((lay==1)&&(sub==5)) ltob1=true;
49  if ((lay==2)&&(sub==5)) ltob2=true;
50  }
52 
53  }
54 
55  private:
57  unsigned int nlay;
58 
59  };
61  public:
63  Trajectory *t2){
64  if (t1->foundHits() != t2->foundHits())
65  return t1->foundHits()> t2->foundHits();
66  return t1->chiSquared()< t2->chiSquared();
67  }
68  };
70  {
71 
73  public:
74 
75  explicit CosmicTrackFinder(const edm::ParameterSet& conf);
76 
77  virtual ~CosmicTrackFinder();
78 
79  virtual void produce(edm::Event& e, const edm::EventSetup& c);
80 
81  private:
85  std::string geometry;
86  bool trinevents;
87  //bool useHitsSplitting_;
88  };
89 }
90 
91 #endif
int foundHits() const
Definition: Trajectory.h:190
bool operator()(Trajectory *t1, Trajectory *t2)
TrajectoryStateOnSurface TSOS
CosmicTrajectoryBuilder cosmicTrajectoryBuilder_
void AnalHits(std::vector< ConstReferenceCountingPointer< TransientTrackingRecHit > > hits)
ConstRecHitContainer recHits(bool splitting=false) const
Definition: Trajectory.cc:67
CRackTrajectoryBuilder crackTrajectoryBuilder_
bool operator()(Trajectory *t1, Trajectory *t2)
CosmicTrackFinder(const edm::ParameterSet &conf)
tuple conf
Definition: dbtoconf.py:185
virtual void produce(edm::Event &e, const edm::EventSetup &c)
edm::ParameterSet conf_
double chiSquared() const
Definition: Trajectory.h:208