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
4 // Geometry
6 //
9 //
10 
14 
16  conf_(config),
17  theUpdator_()
18 {
19 
20  LogDebug("ConversionSeedFinder") << " CTOR " << "\n";
21 
22  theMeasurementTrackerName_=config.getParameter<std::string>("MeasurementTrackerName");
23 }
24 
25 
26 
28 
30 
31  //get the BeamSpot
32  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
33  evt.getByLabel("offlineBeamSpot",recoBeamSpotHandle);
34  theBeamSpot_ = *recoBeamSpotHandle;
35 
36  evt.getByLabel(edm::InputTag("MeasurementTrackerEvent"), theTrackerData_);
37 
38 }
39 
42  es.get<IdealMagneticFieldRecord>().get( theMF_ );
43 
44 
45  edm::ESHandle<MeasurementTracker> measurementTrackerHandle;
46  es.get<CkfComponentsRecord>().get(theMeasurementTrackerName_,measurementTrackerHandle);
47  theMeasurementTracker_ = measurementTrackerHandle.product();
48 
49  edm::ESHandle<Propagator> propagatorAlongMomHandle;
50  es.get<TrackingComponentsRecord>().get("alongMomElePropagator",propagatorAlongMomHandle);
51  thePropagatorAlongMomentum_ = &(*propagatorAlongMomHandle);
52 
53 
54  edm::ESHandle<Propagator> propagatorOppoToMomHandle;
55  es.get<TrackingComponentsRecord>().get("oppositeToMomElePropagator",propagatorOppoToMomHandle);
56  thePropagatorOppositeToMomentum_ = &(*propagatorOppoToMomHandle);
57 
58 
59 
60 
61 }
62 
63 
65 
66 
67  int charge;
68  //List the DetLayers crossed by a straight line from the centre of the
69  //detector to the supercluster position
70  // GlobalPoint vertex(0.,0.,0.);
72  charge=-1;
73  FreeTrajectoryState theStraightLineFTS = trackStateFromClusters(charge, vertex, alongMomentum, 1.);
74 
75  findLayers( theStraightLineFTS );
76 
77 
78 }
79 
81  PropagationDirection dir, float scaleFactor) const {
82 
83 
84 
85  double caloEnergy = theSCenergy_ * scaleFactor ;
86 
87  GlobalVector radiusCalo = theSCPosition_ - theOrigin ;
88 
89  GlobalVector momentumWithoutCurvature = radiusCalo.unit() * caloEnergy;
90 
91 
93  if(dir == alongMomentum) {
94  gtp = GlobalTrajectoryParameters(theOrigin, momentumWithoutCurvature, charge, &(*theMF_) ) ;
95  } else {
96  gtp = GlobalTrajectoryParameters(theSCPosition_, momentumWithoutCurvature, charge, &(*theMF_) ) ;
97  }
98 
99 
100 
101 
102  // now create error matrix
103  // dpos = 4mm/sqrt(E), dtheta = move vertex by 1sigma
104  float dpos = 0.4/sqrt(theSCenergy_);
105  dpos *= 2.;
106  float dphi = dpos/theSCPosition_.perp();
107  // float dp = 0.03 * sqrt(theCaloEnergy);
108  // float dp = theCaloEnergy / sqrt(12.); // for fun
109  float theta1 = theSCPosition_.theta();
110  float theta2 = atan2(double(theSCPosition_.perp()), theSCPosition_.z()-5.5);
111  float dtheta = theta1 - theta2;
113  m[0][0] = 1.; m[1][1] = dpos*dpos ; m[2][2] = dpos*dpos ;
114  m[3][3] = dphi*dphi ; m[4][4] = dtheta * dtheta ;
115 
117 
118  return fts ;
119 
120 
121 }
122 
124 
125 
126 
127  theLayerList_.clear();
128 
129 
131 
132  StartingLayerFinder starter(&prop, this->getMeasurementTracker() );
133 
134  LayerCollector collector(&prop, &starter, 5., 5.);
135 
136  theLayerList_ = collector.allLayers(traj);
137 
138 
139  for(unsigned int i = 0; i < theLayerList_.size(); ++i) {
140  printLayer(i);
141  }
142 
143 
144 }
145 
146 
148  const DetLayer * layer = theLayerList_[i];
149  if (layer->location() == GeomDetEnumerators::barrel ) {
150  // const BarrelDetLayer * barrelLayer = dynamic_cast<const BarrelDetLayer*>(layer);
151  //float r = barrelLayer->specificSurface().radius();
152  // std::cout << " barrel layer radius " << r << " " << barrelLayer->specificSurface().bounds().length()/2. << "\n";
153 
154  } else {
155  // const ForwardDetLayer * forwardLayer = dynamic_cast<const ForwardDetLayer*>(layer);
156  // float z = fabs(forwardLayer->surface().position().z());
157  // std::cout << " forward layer position " << z << " " << forwardLayer->specificSurface().innerRadius() << " " << forwardLayer->specificSurface().outerRadius() << "\n";
158  }
159 }
160 
161 
162 
163 
164 
165 
166 
167 
168 
169 
170 
#define LogDebug(id)
const Propagator * thePropagatorAlongMomentum_
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
const MeasurementTracker * getMeasurementTracker() const
const MeasurementTracker * theMeasurementTracker_
T perp() const
Definition: PV3DBase.h:72
virtual Location location() const =0
Which part of the detector (barrel, endcap)
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
PropagationDirection
double charge(const std::vector< uint8_t > &Ampls)
const Propagator * thePropagatorOppositeToMomentum_
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
std::string theMeasurementTrackerName_
void printLayer(int i) const
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
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:390
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.
const Point & position() const
position
Definition: BeamSpot.h:62
dbl *** dir
Definition: mlp_gen.cc:35
edm::Handle< MeasurementTrackerEvent > theTrackerData_