CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ConversionSeedFinder.cc
Go to the documentation of this file.
3 // Field
5 // Geometry
8 //
12 //
13 
17 
19  conf_(config),
20  theUpdator_()
21 {
22 
23  LogDebug("ConversionSeedFinder") << " CTOR " << "\n";
24 
25 
26 }
27 
28 
29 
31 
34 
35  //get the BeamSpot
36  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
37  evt.getByLabel("offlineBeamSpot",recoBeamSpotHandle);
38  theBeamSpot_ = *recoBeamSpotHandle;
39 
40 
41 
42 }
43 
46  es.get<IdealMagneticFieldRecord>().get( theMF_ );
47 
48 
49  edm::ESHandle<MeasurementTracker> measurementTrackerHandle;
50  es.get<CkfComponentsRecord>().get(measurementTrackerHandle);
51  theMeasurementTracker_ = measurementTrackerHandle.product();
52 
53  edm::ESHandle<Propagator> propagatorAlongMomHandle;
54  es.get<TrackingComponentsRecord>().get("alongMomElePropagator",propagatorAlongMomHandle);
55  thePropagatorAlongMomentum_ = &(*propagatorAlongMomHandle);
56 
57 
58  edm::ESHandle<Propagator> propagatorOppoToMomHandle;
59  es.get<TrackingComponentsRecord>().get("oppositeToMomElePropagator",propagatorOppoToMomHandle);
60  thePropagatorOppositeToMomentum_ = &(*propagatorOppoToMomHandle);
61 
62 
63 
64 
65 }
66 
67 
69 
70 
71  int charge;
72  //List the DetLayers crossed by a straight line from the centre of the
73  //detector to the supercluster position
74  // GlobalPoint vertex(0.,0.,0.);
76  charge=-1;
77  FreeTrajectoryState theStraightLineFTS = trackStateFromClusters(charge, vertex, alongMomentum, 1.);
78 
79  findLayers( theStraightLineFTS );
80 
81 
82 }
83 
85  PropagationDirection dir, float scaleFactor) const {
86 
87 
88 
89  double caloEnergy = theSCenergy_ * scaleFactor ;
90 
91  GlobalVector radiusCalo = theSCPosition_ - theOrigin ;
92 
93  GlobalVector momentumWithoutCurvature = radiusCalo.unit() * caloEnergy;
94 
95 
97  if(dir == alongMomentum) {
98  gtp = GlobalTrajectoryParameters(theOrigin, momentumWithoutCurvature, charge, &(*theMF_) ) ;
99  } else {
100  gtp = GlobalTrajectoryParameters(theSCPosition_, momentumWithoutCurvature, charge, &(*theMF_) ) ;
101  }
102 
103 
104 
105 
106  // now create error matrix
107  // dpos = 4mm/sqrt(E), dtheta = move vertex by 1sigma
108  float dpos = 0.4/sqrt(theSCenergy_);
109  dpos *= 2.;
110  float dphi = dpos/theSCPosition_.perp();
111  // float dp = 0.03 * sqrt(theCaloEnergy);
112  // float dp = theCaloEnergy / sqrt(12.); // for fun
113  float theta1 = theSCPosition_.theta();
114  float theta2 = atan2(double(theSCPosition_.perp()), theSCPosition_.z()-5.5);
115  float dtheta = theta1 - theta2;
116  AlgebraicSymMatrix m(5,1) ;
117  m[0][0] = 1.; m[1][1] = dpos*dpos ; m[2][2] = dpos*dpos ;
118  m[3][3] = dphi*dphi ; m[4][4] = dtheta * dtheta ;
119 
121 
122  return fts ;
123 
124 
125 }
126 
128 
129 
130 
131  theLayerList_.clear();
132 
133 
135 
136  StartingLayerFinder starter(&prop, this->getMeasurementTracker() );
137 
138  LayerCollector collector(&prop, &starter, 5., 5.);
139 
140  theLayerList_ = collector.allLayers(traj);
141 
142 
143  for(unsigned int i = 0; i < theLayerList_.size(); ++i) {
144  printLayer(i);
145  }
146 
147 
148 }
149 
150 
152  const DetLayer * layer = theLayerList_[i];
153  if (layer->location() == GeomDetEnumerators::barrel ) {
154  // const BarrelDetLayer * barrelLayer = dynamic_cast<const BarrelDetLayer*>(layer);
155  //float r = barrelLayer->specificSurface().radius();
156  // std::cout << " barrel layer radius " << r << " " << barrelLayer->specificSurface().bounds().length()/2. << "\n";
157 
158  } else {
159  // const ForwardDetLayer * forwardLayer = dynamic_cast<const ForwardDetLayer*>(layer);
160  // float z = fabs(forwardLayer->surface().position().z());
161  // std::cout << " forward layer position " << z << " " << forwardLayer->specificSurface().innerRadius() << " " << forwardLayer->specificSurface().outerRadius() << "\n";
162  }
163 }
164 
165 
166 
167 
168 
169 
170 
171 
172 
173 
174 
#define LogDebug(id)
const Propagator * thePropagatorAlongMomentum_
int i
Definition: DBlmapReader.cc:9
const MeasurementTracker * getMeasurementTracker() const
const MeasurementTracker * theMeasurementTracker_
T perp() const
Definition: PV3DBase.h:66
virtual Location location() const =0
Which part of the detector (barrel, endcap)
virtual void update(const edm::Event &) const
PropagationDirection
double charge(const std::vector< uint8_t > &Ampls)
const Propagator * thePropagatorOppositeToMomentum_
Geom::Theta< T > theta() const
Definition: PV3DBase.h:69
void printLayer(int i) const
T sqrt(T t)
Definition: SSEVec.h:28
T z() const
Definition: PV3DBase.h:58
std::vector< const DetLayer * > theLayerList_
const TrackingGeometry * geomTracker() const
const TrackingGeometry * theTrackerGeom_
edm::ESHandle< GeometricSearchTracker > theGeomSearchTracker_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
std::vector< const DetLayer * > allLayers(const FTS &aFts) const
Vector3DBase unit() const
Definition: Vector3DBase.h:57
edm::ESHandle< MagneticField > theMF_
void setEvent(const edm::Event &e)
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
FreeTrajectoryState trackStateFromClusters(int aCharge, const GlobalPoint &gpOrigine, PropagationDirection dir, float scaleFactor) const
void setEventSetup(const edm::EventSetup &es)
Initialize EventSetup objects at each event.
tuple config
Definition: cmsDriver.py:17
CLHEP::HepSymMatrix AlgebraicSymMatrix
const Point & position() const
position
Definition: BeamSpot.h:63
dbl *** dir
Definition: mlp_gen.cc:35