Go to the documentation of this file.00001 #include "RecoTracker/SpecialSeedGenerators/interface/CosmicTrackingRegion.h"
00002 #include "TrackingTools/KalmanUpdators/interface/EtaPhiMeasurementEstimator.h"
00003 #include "RecoTracker/Record/interface/CkfComponentsRecord.h"
00004 #include "RecoTracker/MeasurementDet/interface/MeasurementTracker.h"
00005
00006 #include "TrackingTools/MeasurementDet/interface/LayerMeasurements.h"
00007 #include "TrackingTools/PatternTools/interface/TrajectoryMeasurement.h"
00008 #include "TrackingTools/TrajectoryParametrization/interface/GlobalTrajectoryParameters.h"
00009 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
00010 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00011 #include "TrackingTools/GeomPropagators/interface/AnalyticalPropagator.h"
00012 #include "TrackingTools/DetLayers/interface/DetLayer.h"
00013
00014 #include "DataFormats/GeometrySurface/interface/BoundPlane.h"
00015
00016 #include "FWCore/Framework/interface/ESHandle.h"
00017
00018 namespace {
00019 template <class T> T sqr( T t) {return t*t;}
00020 }
00021
00022 using namespace std;
00023 using namespace ctfseeding;
00024
00025
00026 TrackingRegion::Hits CosmicTrackingRegion::hits(const edm::Event& ev,
00027 const edm::EventSetup& es,
00028 const ctfseeding::SeedingLayer* layer) const
00029 {
00030
00031
00032
00033
00034
00035 TrackingRegion::Hits result;
00036
00037
00038 const DetLayer * detLayer = layer->detLayer();
00039 LogDebug("CosmicTrackingRegion") << "Looking at hits on subdet/layer " << layer->name();
00040 EtaPhiMeasurementEstimator est(0.3,0.3);
00041
00042
00043 edm::ESHandle<MagneticField> field;
00044 es.get<IdealMagneticFieldRecord>().get(field);
00045 const MagneticField * magField = field.product();
00046
00047
00048 const GlobalPoint vtx = origin();
00049 GlobalVector dir = direction();
00050 LogDebug("CosmicTrackingRegion") <<"The initial region characteristics are:" << "\n"
00051 <<" Origin = " << origin() << "\n"
00052 <<" Direction = " << direction() << "\n"
00053 <<" Eta = " << origin().eta() << "\n"
00054 <<" Phi = " << origin().phi();
00055
00056
00057 float phi = dir.phi();
00058 Surface::RotationType rot( sin(phi), -cos(phi), 0,
00059 0, 0, -1,
00060 cos(phi), sin(phi), 0);
00061
00062 Plane::PlanePointer surface = Plane::build(vtx, rot);
00063 FreeTrajectoryState fts( GlobalTrajectoryParameters(vtx, dir, 1, magField) );
00064 TrajectoryStateOnSurface tsos(fts, *surface);
00065 LogDebug("CosmicTrackingRegion")
00066 << "The state used to find measurement with the measurement tracker is:\n" << tsos;
00067
00068
00069 AnalyticalPropagator prop( magField, alongMomentum);
00070
00071
00072
00073
00074
00075 TrajectoryStateOnSurface stateOnLayer = prop.propagate( *tsos.freeState(),
00076 detLayer->surface());
00077
00078
00079 if (stateOnLayer.isValid()){
00080 LogDebug("CosmicTrackingRegion") << "The initial state propagates to the layer surface: \n" << stateOnLayer
00081 << "R = " << stateOnLayer.globalPosition().perp() << "\n"
00082 << "Eta = " << stateOnLayer.globalPosition().eta() << "\n"
00083 << "Phi = " << stateOnLayer.globalPosition().phi();
00084
00085 }
00086 else{
00087 LogDebug("CosmicTrackingRegion") << "The initial state does not propagate to the layer surface.";
00088 }
00089
00090
00091 typedef DetLayer::DetWithState DetWithState;
00092 vector<DetWithState> compatDets = detLayer->compatibleDets(tsos, prop, est);
00093 LogDebug("CosmicTrackingRegion") << "Compatible dets = " << compatDets.size();
00094
00095
00096
00097
00098
00099
00100 edm::ESHandle<MeasurementTracker> measurementTrackerESH;
00101 es.get<CkfComponentsRecord>().get(measurementTrackerName_,measurementTrackerESH);
00102 const MeasurementTracker * measurementTracker = measurementTrackerESH.product();
00103 measurementTracker->update(ev);
00104 LayerMeasurements lm(measurementTracker);
00105 vector<TrajectoryMeasurement> meas = lm.measurements(*detLayer, tsos, prop, est);
00106 LogDebug("CosmicTrackingRegion") << "Number of Trajectory measurements = " << meas.size()
00107 <<" but the last one is always an invalid hit, by construction.";
00108
00109
00110 typedef vector<TrajectoryMeasurement>::const_iterator IM;
00111
00112 for (IM im = meas.begin(); im != meas.end(); im++) {
00113 TrajectoryMeasurement::ConstRecHitPointer ptrHit = im->recHit();
00114
00115 if (ptrHit->isValid()) {
00116 LogDebug("CosmicTrackingRegion") << "Hit found in the region at position: "<<ptrHit->globalPosition();
00117 result.push_back( ptrHit );
00118 }
00119
00120 else LogDebug("CosmicTrackingRegion") << "No valid hit";
00121 }
00122
00123
00124
00125
00126
00127 return result;
00128 }
00129
00130