Go to the documentation of this file.00001 #include "TrackingTools/TrackRefitter/interface/TrackTransformer.h"
00002
00003 #include "FWCore/Framework/interface/EventSetup.h"
00004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006
00007 #include "Geometry/Records/interface/GlobalTrackingGeometryRecord.h"
00008 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00009
00010 #include "TrackingTools/Records/interface/TrackingComponentsRecord.h"
00011 #include "TrackingTools/TrackFitters/interface/TrajectoryFitter.h"
00012 #include "TrackingTools/PatternTools/interface/TrajectorySmoother.h"
00013
00014 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00015 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00016 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHitBuilder.h"
00017 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
00018 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00019 #include "TrackingTools/GeomPropagators/interface/Propagator.h"
00020
00021 #include "DataFormats/TrackReco/interface/Track.h"
00022 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00023 #include "DataFormats/DetId/interface/DetId.h"
00024
00025 using namespace std;
00026 using namespace edm;
00027
00029 TrackTransformer::TrackTransformer(const ParameterSet& parameterSet){
00030
00031
00032 string refitDirectionName = parameterSet.getParameter<string>("RefitDirection");
00033 theRefitDirection = RefitDirection(refitDirectionName);
00034
00035 theFitterName = parameterSet.getParameter<string>("Fitter");
00036 theSmootherName = parameterSet.getParameter<string>("Smoother");
00037 thePropagatorName = parameterSet.getParameter<string>("Propagator");
00038
00039 theTrackerRecHitBuilderName = parameterSet.getParameter<string>("TrackerRecHitBuilder");
00040 theMuonRecHitBuilderName = parameterSet.getParameter<string>("MuonRecHitBuilder");
00041
00042 theRPCInTheFit = parameterSet.getParameter<bool>("RefitRPCHits");
00043 theDoPredictionsOnly = parameterSet.getParameter<bool>("DoPredictionsOnly");
00044
00045 theCacheId_TC = theCacheId_GTG = theCacheId_MG = theCacheId_TRH = 0;
00046 }
00047
00049 TrackTransformer::~TrackTransformer(){}
00050
00051
00052 void TrackTransformer::setServices(const EventSetup& setup){
00053
00054 const std::string metname = "Reco|TrackingTools|TrackTransformer";
00055
00056 setup.get<TrajectoryFitter::Record>().get(theFitterName,theFitter);
00057 setup.get<TrajectoryFitter::Record>().get(theSmootherName,theSmoother);
00058
00059
00060 unsigned long long newCacheId_TC = setup.get<TrackingComponentsRecord>().cacheIdentifier();
00061
00062 if ( newCacheId_TC != theCacheId_TC ){
00063 LogTrace(metname) << "Tracking Component changed!";
00064 theCacheId_TC = newCacheId_TC;
00065 setup.get<TrackingComponentsRecord>().get(thePropagatorName,thePropagator);
00066 }
00067
00068
00069 unsigned long long newCacheId_GTG = setup.get<GlobalTrackingGeometryRecord>().cacheIdentifier();
00070 if ( newCacheId_GTG != theCacheId_GTG ) {
00071 LogTrace(metname) << "GlobalTrackingGeometry changed!";
00072 theCacheId_GTG = newCacheId_GTG;
00073 setup.get<GlobalTrackingGeometryRecord>().get(theTrackingGeometry);
00074 }
00075
00076
00077 unsigned long long newCacheId_MG = setup.get<IdealMagneticFieldRecord>().cacheIdentifier();
00078 if ( newCacheId_MG != theCacheId_MG ) {
00079 LogTrace(metname) << "Magnetic Field changed!";
00080 theCacheId_MG = newCacheId_MG;
00081 setup.get<IdealMagneticFieldRecord>().get(theMGField);
00082 }
00083
00084
00085 unsigned long long newCacheId_TRH = setup.get<TransientRecHitRecord>().cacheIdentifier();
00086 if ( newCacheId_TRH != theCacheId_TRH ) {
00087 theCacheId_TRH = newCacheId_TRH;
00088 LogTrace(metname) << "TransientRecHitRecord changed!";
00089 setup.get<TransientRecHitRecord>().get(theTrackerRecHitBuilderName,theTrackerRecHitBuilder);
00090 setup.get<TransientRecHitRecord>().get(theMuonRecHitBuilderName,theMuonRecHitBuilder);
00091 }
00092 }
00093
00094
00095 vector<Trajectory> TrackTransformer::transform(const reco::TrackRef& track) const {
00096 return transform(*track);
00097 }
00098
00099
00100 TransientTrackingRecHit::ConstRecHitContainer
00101 TrackTransformer::getTransientRecHits(const reco::TransientTrack& track) const {
00102
00103 TransientTrackingRecHit::ConstRecHitContainer result;
00104
00105 for (trackingRecHit_iterator hit = track.recHitsBegin(); hit != track.recHitsEnd(); ++hit) {
00106 if((*hit)->isValid()) {
00107 if ( (*hit)->geographicalId().det() == DetId::Tracker ) {
00108 result.push_back(theTrackerRecHitBuilder->build(&**hit));
00109 } else if ( (*hit)->geographicalId().det() == DetId::Muon ){
00110 if( (*hit)->geographicalId().subdetId() == 3 && !theRPCInTheFit){
00111 LogTrace("Reco|TrackingTools|TrackTransformer") << "RPC Rec Hit discarged";
00112 continue;
00113 }
00114 result.push_back(theMuonRecHitBuilder->build(&**hit));
00115 }
00116 }
00117 }
00118
00119 return result;
00120 }
00121
00122
00123 RefitDirection::GeometricalDirection
00124 TrackTransformer::checkRecHitsOrdering(TransientTrackingRecHit::ConstRecHitContainer& recHits) const {
00125
00126 if (!recHits.empty()){
00127 GlobalPoint first = trackingGeometry()->idToDet(recHits.front()->geographicalId())->position();
00128 GlobalPoint last = trackingGeometry()->idToDet(recHits.back()->geographicalId())->position();
00129
00130 double rFirst = first.mag();
00131 double rLast = last.mag();
00132 if(rFirst < rLast) return RefitDirection::insideOut;
00133 else if(rFirst > rLast) return RefitDirection::outsideIn;
00134 else{
00135 LogDebug("Reco|TrackingTools|TrackTransformer") << "Impossible to determine the rechits order" <<endl;
00136 return RefitDirection::undetermined;
00137 }
00138 }
00139 else{
00140 LogDebug("Reco|TrackingTools|TrackTransformer") << "Impossible to determine the rechits order" <<endl;
00141 return RefitDirection::undetermined;
00142 }
00143 }
00144
00145
00146
00147
00148
00149
00150
00151
00152
00153
00154
00156 vector<Trajectory> TrackTransformer::transform(const reco::Track& newTrack) const {
00157
00158 const std::string metname = "Reco|TrackingTools|TrackTransformer";
00159
00160 reco::TransientTrack track(newTrack,magneticField(),trackingGeometry());
00161
00162
00163 TransientTrackingRecHit::ConstRecHitContainer recHitsForReFit = getTransientRecHits(track);
00164
00165 return transform(track, recHitsForReFit);
00166 }
00167
00168
00170 vector<Trajectory> TrackTransformer::transform(const reco::TransientTrack track,
00171 TransientTrackingRecHit::ConstRecHitContainer recHitsForReFit) const {
00172
00173 const std::string metname = "Reco|TrackingTools|TrackTransformer";
00174
00175 if(recHitsForReFit.size() < 2) return vector<Trajectory>();
00176
00177
00178
00179
00180
00181
00182
00183
00184
00185
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208 RefitDirection::GeometricalDirection recHitsOrder = checkRecHitsOrdering(recHitsForReFit);
00209 LogTrace(metname) << "RH order (0-insideOut, 1-outsideIn): " << recHitsOrder;
00210
00211 PropagationDirection propagationDirection = theRefitDirection.propagationDirection();
00212
00213
00214 if(propagationDirection == anyDirection){
00215 GlobalVector momentum = track.innermostMeasurementState().globalMomentum();
00216 GlobalVector position = track.innermostMeasurementState().globalPosition() - GlobalPoint(0,0,0);
00217 RefitDirection::GeometricalDirection p = (momentum.x()*position.x() > 0 || momentum.y()*position.y() > 0) ? RefitDirection::insideOut : RefitDirection::outsideIn;
00218
00219 propagationDirection = p == theRefitDirection.geometricalDirection() ? alongMomentum : oppositeToMomentum;
00220 LogTrace(metname) << "P (0-insideOut, 1-outsideIn): " << p;
00221 LogTrace(metname) << "FD (0-OM, 1-AM, 2-ANY): " << propagationDirection;
00222 }
00223
00224
00225
00226 if(theRefitDirection.propagationDirection() != anyDirection){
00227 if((recHitsOrder == RefitDirection::insideOut && propagationDirection == oppositeToMomentum) ||
00228 (recHitsOrder == RefitDirection::outsideIn && propagationDirection == alongMomentum) )
00229 reverse(recHitsForReFit.begin(),recHitsForReFit.end());}
00230
00231
00232 else{
00233
00234 if(theRefitDirection.geometricalDirection() != recHitsOrder) reverse(recHitsForReFit.begin(),recHitsForReFit.end());
00235 }
00236
00237
00238
00239 TrajectoryStateOnSurface firstTSOS = track.innermostMeasurementState();
00240 unsigned int innerId = track.track().innerDetId();
00241 if(theRefitDirection.propagationDirection() != anyDirection){
00242 if(propagationDirection == oppositeToMomentum){
00243 innerId = track.track().outerDetId();
00244 firstTSOS = track.outermostMeasurementState();
00245 }
00246 }
00247 else {
00248
00249 if(theRefitDirection.geometricalDirection() == RefitDirection::outsideIn){
00250 innerId = track.track().outerDetId();
00251 firstTSOS = track.outermostMeasurementState();
00252 }
00253
00254 }
00255
00256
00257 if(!firstTSOS.isValid()){
00258 LogTrace(metname)<<"Error wrong initial state!"<<endl;
00259 return vector<Trajectory>();
00260 }
00261
00262 TrajectorySeed seed(PTrajectoryStateOnDet(),TrajectorySeed::recHitContainer(),propagationDirection);
00263
00264 if(recHitsForReFit.front()->geographicalId() != DetId(innerId)){
00265 LogTrace(metname)<<"Propagation occured"<<endl;
00266 firstTSOS = propagator()->propagate(firstTSOS, recHitsForReFit.front()->det()->surface());
00267 if(!firstTSOS.isValid()){
00268 LogTrace(metname)<<"Propagation error!"<<endl;
00269 return vector<Trajectory>();
00270 }
00271 }
00272
00273 if(theDoPredictionsOnly){
00274 Trajectory aTraj(seed,propagationDirection);
00275 TrajectoryStateOnSurface predTSOS = firstTSOS;
00276 for(TransientTrackingRecHit::ConstRecHitContainer::const_iterator ihit = recHitsForReFit.begin();
00277 ihit != recHitsForReFit.end(); ++ihit ) {
00278 predTSOS = propagator()->propagate(predTSOS, (*ihit)->det()->surface());
00279 if (predTSOS.isValid()) aTraj.push(TrajectoryMeasurement(predTSOS, *ihit));
00280 }
00281 return vector<Trajectory>(1, aTraj);
00282 }
00283
00284
00285 vector<Trajectory> trajectories = theFitter->fit(seed,recHitsForReFit,firstTSOS);
00286
00287 if(trajectories.empty()){
00288 LogTrace(metname)<<"No Track refitted!"<<endl;
00289 return vector<Trajectory>();
00290 }
00291
00292 Trajectory trajectoryBW = trajectories.front();
00293
00294 vector<Trajectory> trajectoriesSM = theSmoother->trajectories(trajectoryBW);
00295
00296 if(trajectoriesSM.empty()){
00297 LogTrace(metname)<<"No Track smoothed!"<<endl;
00298 return vector<Trajectory>();
00299 }
00300
00301 return trajectoriesSM;
00302
00303 }
00304
00305