CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SeedForPhotonConversionFromQuadruplets.h
Go to the documentation of this file.
1 #ifndef SeedForPhotonConversionFromQuadruplets_H
2 #define SeedForPhotonConversionFromQuadruplets_H
3 
9 
13 
15 public:
16  static const int cotTheta_Max=99999;
17 
19  thePropagatorLabel(cfg.getParameter<std::string>("propagator")),
20  theBOFFMomentum(cfg.existsAs<double>("SeedMomentumForBOFF") ? cfg.getParameter<double>("SeedMomentumForBOFF") : 5.0)
21  {}
22 
24  const std::string & propagator = "PropagatorWithMaterial", double seedMomentumForBOFF = -5.0)
25  : thePropagatorLabel(propagator), theBOFFMomentum(seedMomentumForBOFF) { }
26 
27  //dtor
29 
30  const TrajectorySeed * trajectorySeed( TrajectorySeedCollection & seedCollection,
31  const SeedingHitSet & phits,
32  const SeedingHitSet & mhits,
33  const TrackingRegion & region,
34  const edm::EventSetup& es,
35  std::stringstream& ss, std::vector<Quad> & quadV,
36  edm::ParameterSet& SeedComparitorPSet,
38 
39 
41  double verySimpleFit(int size, double* ax, double* ay, double* e2y, double& p0, double& e2p0, double& p1);
43 
44  //
45  // Some utility methods added by sguazz
46  void stupidPrint(std::string s,float* d);
47  void stupidPrint(std::string s,double* d);
48  void stupidPrint(const char* s,GlobalPoint* d);
49  void stupidPrint(const char* s,GlobalPoint* d, int n);
50  void bubbleSortVsPhi(GlobalPoint arr[], int n, GlobalPoint vtx);
51  void bubbleReverseSortVsPhi(GlobalPoint arr[], int n, GlobalPoint vtx);
52  //
53  //
54 
55 
56 
57 
58  protected:
59 
60  bool checkHit(
63  const edm::EventSetup& es) const { return true; }
64 
66  const SeedingHitSet & hits,
67  const GlobalPoint & vertexPos,
68  const edm::EventSetup& es,
69  const float cotTheta) const;
70 
72  const GlobalVector& vertexBounds,
73  float ptMin,
74  float sinTheta) const;
75 
76  const TrajectorySeed * buildSeed(
77  TrajectorySeedCollection & seedCollection,
78  const SeedingHitSet & hits,
79  const FreeTrajectoryState & fts,
80  const edm::EventSetup& es,
81  bool apply_dzCut,
82  const TrackingRegion &region) const;
83 
84  bool buildSeedBool(
85  TrajectorySeedCollection & seedCollection,
86  const SeedingHitSet & hits,
87  const FreeTrajectoryState & fts,
88  const edm::EventSetup& es,
89  bool apply_dzCut,
90  const TrackingRegion & region,
91  double dzcut) const;
92 
95  const TrajectoryStateOnSurface &state) const;
96 
97  bool similarQuadExist(Quad & thisQuad, std::vector<Quad>& quadV);
98 
100 
101 
102 protected:
105  double kPI_;
106 
107  std::stringstream * pss;
109 };
110 #endif
CurvilinearTrajectoryError initialError(const GlobalVector &vertexBounds, float ptMin, float sinTheta) const
bool buildSeedBool(TrajectorySeedCollection &seedCollection, const SeedingHitSet &hits, const FreeTrajectoryState &fts, const edm::EventSetup &es, bool apply_dzCut, const TrackingRegion &region, double dzcut) const
void bubbleReverseSortVsPhi(GlobalPoint arr[], int n, GlobalPoint vtx)
bool checkHit(const TrajectoryStateOnSurface &, const TransientTrackingRecHit::ConstRecHitPointer &hit, const edm::EventSetup &es) const
void bubbleSortVsPhi(GlobalPoint arr[], int n, GlobalPoint vtx)
SeedForPhotonConversionFromQuadruplets(const std::string &propagator="PropagatorWithMaterial", double seedMomentumForBOFF=-5.0)
double getSqrEffectiveErrorOnZ(const TransientTrackingRecHit::ConstRecHitPointer &hit, const TrackingRegion &region)
GlobalTrajectoryParameters initialKinematic(const SeedingHitSet &hits, const GlobalPoint &vertexPos, const edm::EventSetup &es, const float cotTheta) const
double verySimpleFit(int size, double *ax, double *ay, double *e2y, double &p0, double &e2p0, double &p1)
std::vector< TrajectorySeed > TrajectorySeedCollection
TransientTrackingRecHit::RecHitPointer refitHit(const TransientTrackingRecHit::ConstRecHitPointer &hit, const TrajectoryStateOnSurface &state) const
const TrajectorySeed * buildSeed(TrajectorySeedCollection &seedCollection, const SeedingHitSet &hits, const FreeTrajectoryState &fts, const edm::EventSetup &es, bool apply_dzCut, const TrackingRegion &region) const
double simpleGetSlope(const TransientTrackingRecHit::ConstRecHitPointer &ohit, const TransientTrackingRecHit::ConstRecHitPointer &nohit, const TransientTrackingRecHit::ConstRecHitPointer &ihit, const TransientTrackingRecHit::ConstRecHitPointer &nihit, const TrackingRegion &region, double &cotTheta, double &z0)
double DeltaPhiManual(math::XYZVector v1, math::XYZVector v2)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
const TrajectorySeed * trajectorySeed(TrajectorySeedCollection &seedCollection, const SeedingHitSet &phits, const SeedingHitSet &mhits, const TrackingRegion &region, const edm::EventSetup &es, std::stringstream &ss, std::vector< Quad > &quadV, edm::ParameterSet &SeedComparitorPSet, edm::ParameterSet &QuadCutPSet)
char state
Definition: procUtils.cc:75
double p1[4]
Definition: TauolaWrapper.h:89
bool similarQuadExist(Quad &thisQuad, std::vector< Quad > &quadV)
Definition: Quad.h:4
tuple size
Write out results.