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