CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/Alignment/OfflineValidation/src/TrackerValidationVariables.cc

Go to the documentation of this file.
00001 
00002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00003 #include "FWCore/Framework/interface/Frameworkfwd.h"
00004 #include "FWCore/Framework/interface/EDAnalyzer.h"
00005 #include "FWCore/Framework/interface/MakerMacros.h"
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007 #include "FWCore/Framework/interface/Event.h"
00008 #include "FWCore/Framework/interface/EventSetup.h"
00009 
00010 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00011 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00012 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00013 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
00014 #include "Geometry/CommonTopologies/interface/RadialStripTopology.h"
00015 
00016 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00017 #include "MagneticField/Engine/interface/MagneticField.h"
00018 
00019 #include "TrackingTools/TrackFitters/interface/TrajectoryStateCombiner.h"
00020 #include "TrackingTools/PatternTools/interface/TrajectoryMeasurement.h"
00021 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00022 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
00023 #include "TrackingTools/GeomPropagators/interface/AnalyticalPropagator.h"
00024 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00025 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
00026 
00027 #include "CondFormats/Alignment/interface/Definitions.h"
00028 
00029 #include "DataFormats/Math/interface/deltaPhi.h"
00030 #include "DataFormats/TrackReco/interface/Track.h"
00031 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00032 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00033 #include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
00034 #include "DataFormats/GeometryCommonDetAlgo/interface/MeasurementPoint.h"
00035 #include "DataFormats/GeometryCommonDetAlgo/interface/MeasurementError.h"
00036 #include "DataFormats/GeometrySurface/interface/TrapezoidalPlaneBounds.h"
00037 #include "DataFormats/GeometrySurface/interface/RectangularPlaneBounds.h"
00038 
00039 #include "Alignment/TrackerAlignment/interface/TrackerAlignableId.h"
00040 
00041 #include "Alignment/OfflineValidation/interface/TrackerValidationVariables.h"
00042 
00043 #include "TMath.h"
00044 
00045 TrackerValidationVariables::TrackerValidationVariables()
00046 {
00047 
00048 }
00049 
00050 TrackerValidationVariables::TrackerValidationVariables(const edm::EventSetup& es, const edm::ParameterSet& iSetup) 
00051   : conf_(iSetup)
00052 {
00053   es.get<TrackerDigiGeometryRecord>().get(tkGeom_);
00054   es.get<IdealMagneticFieldRecord>().get(magneticField_);
00055 }
00056 
00057 TrackerValidationVariables::~TrackerValidationVariables()
00058 {
00059 
00060 }
00061 
00062 void
00063 TrackerValidationVariables::fillHitQuantities(const Trajectory* trajectory, std::vector<AVHitStruct> & v_avhitout)
00064 {
00065   TrajectoryStateCombiner tsoscomb;
00066 
00067   const std::vector<TrajectoryMeasurement> &tmColl = trajectory->measurements();
00068   for(std::vector<TrajectoryMeasurement>::const_iterator itTraj = tmColl.begin();
00069       itTraj != tmColl.end();
00070       ++itTraj) {
00071     
00072     if (!itTraj->updatedState().isValid()) continue;
00073     
00074     TrajectoryStateOnSurface tsos = tsoscomb( itTraj->forwardPredictedState(), itTraj->backwardPredictedState() );
00075     TransientTrackingRecHit::ConstRecHitPointer hit = itTraj->recHit();
00076     
00077     if(!hit->isValid() || hit->geographicalId().det() != DetId::Tracker) continue;
00078     
00079     AVHitStruct hitStruct;
00080     const DetId& hit_detId = hit->geographicalId();
00081     unsigned int IntRawDetID = (hit_detId.rawId());     
00082     unsigned int IntSubDetID = (hit_detId.subdetId());
00083     
00084     if(IntSubDetID == 0) continue;
00085     
00086     //first calculate residuals in cartesian coordinates in the local module coordinate system
00087     
00088     LocalPoint lPHit = hit->localPosition();
00089     LocalPoint lPTrk = tsos.localPosition();
00090     LocalVector lVTrk = tsos.localDirection();
00091 
00092     hitStruct.localAlpha = atan2(lVTrk.x(), lVTrk.z()); // wrt. normal tg(alpha)=x/z
00093     hitStruct.localBeta  = atan2(lVTrk.y(), lVTrk.z()); // wrt. normal tg(beta)= y/z
00094 
00095     //LocalError errHit = hit->localPositionError();
00096     // adding APE to hitError
00097     AlgebraicROOTObject<2>::SymMatrix mat = asSMatrix<2>(hit->parametersError());
00098     LocalError errHit = LocalError( mat(0,0),mat(0,1),mat(1,1) );
00099     LocalError errTrk = tsos.localError().positionError();
00100     
00101     //check for negative error values: track error can have negative value, if matrix inversion fails (very rare case)
00102     //hit error should always give positive values
00103     if(errHit.xx()<0. || errHit.yy()<0. || errTrk.xx()<0. || errTrk.yy()<0.){
00104       edm::LogError("TrackerValidationVariables") << "@SUB=TrackerValidationVariables::fillHitQuantities"
00105                                                   << "One of the squared error methods gives negative result"
00106                                                   << "\n\terrHit.xx()\terrHit.yy()\terrTrk.xx()\terrTrk.yy()"
00107                                                   << "\n\t" << errHit.xx()
00108                                                   << "\t" << errHit.yy()
00109                                                   << "\t" << errTrk.xx()
00110                                                   << "\t" << errTrk.yy();
00111       continue;
00112     }
00113     
00114     align::LocalVector res = lPTrk - lPHit;
00115     
00116     float resXErr = std::sqrt( errHit.xx() + errTrk.xx() );
00117     float resYErr = std::sqrt( errHit.yy() + errTrk.yy() );
00118     
00119     hitStruct.resX = res.x();
00120     hitStruct.resY = res.y();
00121     hitStruct.resErrX = resXErr;
00122     hitStruct.resErrY = resYErr;
00123 
00124     // hitStruct.localX = lPhit.x();
00125     // hitStruct.localY = lPhit.y();
00126     // EM: use predictions for local coordinates
00127     hitStruct.localX = lPTrk.x();
00128     hitStruct.localY = lPTrk.y();
00129 
00130     // now calculate residuals taking global orientation of modules and radial topology in TID/TEC into account
00131     float resXprime(999.F), resYprime(999.F);
00132     float resXatTrkY(999.F);
00133     float resXprimeErr(999.F), resYprimeErr(999.F);
00134     
00135     if(hit->detUnit()){ // is it a single physical module?
00136       const GeomDetUnit& detUnit = *(hit->detUnit());
00137       float uOrientation(-999.F), vOrientation(-999.F);
00138       float resXTopol(999.F), resYTopol(999.F);
00139       float resXatTrkYTopol(999.F);       
00140 
00141       const Surface& surface = hit->detUnit()->surface();
00142       const BoundPlane& boundplane = hit->detUnit()->surface();
00143       const Bounds& bound = boundplane.bounds();
00144       
00145       float length = 0;
00146       float width = 0;
00147 
00148       LocalPoint lPModule(0.,0.,0.), lUDirection(1.,0.,0.), lVDirection(0.,1.,0.);
00149       GlobalPoint gPModule    = surface.toGlobal(lPModule),
00150                   gUDirection = surface.toGlobal(lUDirection),
00151                   gVDirection = surface.toGlobal(lVDirection);
00152       
00153       if (IntSubDetID == PixelSubdetector::PixelBarrel ||
00154           IntSubDetID == StripSubdetector::TIB || 
00155           IntSubDetID == StripSubdetector::TOB) {
00156         
00157         uOrientation = deltaPhi(gUDirection.phi(),gPModule.phi()) >= 0. ? +1.F : -1.F;
00158         vOrientation = gVDirection.z() - gPModule.z() >= 0 ? +1.F : -1.F;
00159         resXTopol = res.x();
00160         resXatTrkYTopol = res.x();
00161         resYTopol = res.y();
00162         resXprimeErr = resXErr;
00163         resYprimeErr = resYErr;
00164 
00165         const RectangularPlaneBounds *rectangularBound = dynamic_cast<const RectangularPlaneBounds*>(&bound);
00166         if (rectangularBound!=NULL) {
00167           hitStruct.inside = rectangularBound->inside(lPTrk);
00168           length = rectangularBound->length();
00169           width = rectangularBound->width();
00170           hitStruct.localXnorm = 2*hitStruct.localX/width;
00171           hitStruct.localYnorm = 2*hitStruct.localY/length;
00172         } else {
00173           throw cms::Exception("Geometry Error")
00174             << "[TrackerValidationVariables] Cannot cast bounds to RectangularPlaneBounds as expected for TPB, TIB and TOB";
00175         }
00176 
00177       } else if (IntSubDetID == PixelSubdetector::PixelEndcap) {
00178         
00179         uOrientation = gUDirection.perp() - gPModule.perp() >= 0 ? +1.F : -1.F;
00180         vOrientation = deltaPhi(gVDirection.phi(),gPModule.phi()) >= 0. ? +1.F : -1.F;
00181         resXTopol = res.x();
00182         resXatTrkYTopol = res.x();
00183         resYTopol = res.y();
00184         resXprimeErr = resXErr;
00185         resYprimeErr = resYErr;
00186         
00187         const RectangularPlaneBounds *rectangularBound = dynamic_cast<const RectangularPlaneBounds*>(&bound);
00188         if (rectangularBound!=NULL) {
00189           hitStruct.inside = rectangularBound->inside(lPTrk);
00190           length = rectangularBound->length();
00191           width = rectangularBound->width();
00192           hitStruct.localXnorm = 2*hitStruct.localX/width;
00193           hitStruct.localYnorm = 2*hitStruct.localY/length;
00194         } else {
00195           throw cms::Exception("Geometry Error")
00196             << "[TrackerValidationVariables] Cannot cast bounds to RectangularPlaneBounds as expected for TPE";
00197         }
00198 
00199       } else if (IntSubDetID == StripSubdetector::TID ||
00200                  IntSubDetID == StripSubdetector::TEC) {
00201         
00202         uOrientation = deltaPhi(gUDirection.phi(),gPModule.phi()) >= 0. ? +1.F : -1.F;
00203         vOrientation = gVDirection.perp() - gPModule.perp() >= 0. ? +1.F : -1.F;
00204         
00205         if (!dynamic_cast<const RadialStripTopology*>(&detUnit.type().topology()))continue;
00206         const RadialStripTopology& topol = dynamic_cast<const RadialStripTopology&>(detUnit.type().topology());
00207         
00208         MeasurementPoint measHitPos = topol.measurementPosition(lPHit);
00209         MeasurementPoint measTrkPos = topol.measurementPosition(lPTrk);
00210         
00211         MeasurementError measHitErr = topol.measurementError(lPHit,errHit);
00212         MeasurementError measTrkErr = topol.measurementError(lPTrk,errTrk);
00213         
00214         if (measHitErr.uu()<0. ||
00215             measHitErr.vv()<0. ||
00216             measTrkErr.uu()<0. ||
00217             measTrkErr.vv()<0.){
00218           edm::LogError("TrackerValidationVariables") << "@SUB=TrackerValidationVariables::fillHitQuantities"
00219                                                       << "One of the squared error methods gives negative result"
00220                                                       << "\n\tmeasHitErr.uu()\tmeasHitErr.vv()\tmeasTrkErr.uu()\tmeasTrkErr.vv()"
00221                                                       << "\n\t" << measHitErr.uu()
00222                                                       << "\t" << measHitErr.vv()
00223                                                       << "\t" << measTrkErr.uu()
00224                                                       << "\t" << measTrkErr.vv();
00225           continue;
00226         }
00227           
00228         float localStripLengthHit = topol.localStripLength(lPHit);
00229         float localStripLengthTrk = topol.localStripLength(lPTrk);
00230         float phiHit = topol.stripAngle(measHitPos.x());
00231         float phiTrk = topol.stripAngle(measTrkPos.x());
00232         float r_0 = topol.originToIntersection();
00233         
00234         resXTopol = (phiTrk-phiHit)*r_0;
00235 //      resXTopol = (tan(phiTrk)-tan(phiHit))*r_0;
00236         
00237     LocalPoint LocalHitPosCor= topol.localPosition(MeasurementPoint(measHitPos.x(), measTrkPos.y()));                       
00238        resXatTrkYTopol = lPTrk.x() - LocalHitPosCor.x();  
00239 
00240     
00241         //resYTopol = measTrkPos.y()*localStripLengthTrk - measHitPos.y()*localStripLengthHit;
00242         float cosPhiHit(cos(phiHit)), cosPhiTrk(cos(phiTrk)),
00243               sinPhiHit(sin(phiHit)), sinPhiTrk(sin(phiTrk));
00244         float l_0 = r_0 - topol.detHeight()/2;
00245         resYTopol = measTrkPos.y()*localStripLengthTrk - measHitPos.y()*localStripLengthHit + l_0*(1/cosPhiTrk - 1/cosPhiHit);
00246         
00247         resXprimeErr = std::sqrt(measHitErr.uu()+measTrkErr.uu())*topol.angularWidth()*r_0;
00248         //resYprimeErr = std::sqrt(measHitErr.vv()*localStripLengthHit*localStripLengthHit + measTrkErr.vv()*localStripLengthTrk*localStripLengthTrk);
00249         float helpSummand = l_0*l_0*topol.angularWidth()*topol.angularWidth()*(sinPhiHit*sinPhiHit/pow(cosPhiHit,4)*measHitErr.uu()
00250                                                                                + sinPhiTrk*sinPhiTrk/pow(cosPhiTrk,4)*measTrkErr.uu() );
00251         resYprimeErr = std::sqrt(measHitErr.vv()*localStripLengthHit*localStripLengthHit
00252                                  + measTrkErr.vv()*localStripLengthTrk*localStripLengthTrk + helpSummand );  
00253 
00254 
00255         const TrapezoidalPlaneBounds *trapezoidalBound = dynamic_cast < const TrapezoidalPlaneBounds * >(& bound);
00256         if (trapezoidalBound!=NULL) {
00257           hitStruct.inside = trapezoidalBound->inside(lPTrk);
00258           length = trapezoidalBound->length();
00259           width  = trapezoidalBound->width();
00260           //float widthAtHalfLength = trapezoidalBound->widthAtHalfLength();
00261 
00262           //    int yAxisOrientation=trapezoidalBound->yAxisOrientation(); 
00263           // for trapezoidal shape modules, scale with as function of local y coordinate 
00264           //    float widthAtlocalY=width-(1-yAxisOrientation*2*lPTrk.y()/length)*(width-widthAtHalfLength); 
00265           //    hitStruct.localXnorm = 2*hitStruct.localX/widthAtlocalY;  
00266           hitStruct.localXnorm = 2*hitStruct.localX/width; 
00267           hitStruct.localYnorm = 2*hitStruct.localY/length;
00268         } else {
00269           throw cms::Exception("Geometry Error")
00270             << "[TrackerValidationVariables] Cannot cast bounds to TrapezoidalPlaneBounds as expected for TID and TEC";
00271         }
00272 
00273       } else {
00274         edm::LogWarning("TrackerValidationVariables") << "@SUB=TrackerValidationVariables::fillHitQuantities" 
00275                                                       << "No valid tracker subdetector " << IntSubDetID;
00276         continue;
00277       }  
00278       
00279       resXprime = resXTopol*uOrientation;
00280       resXatTrkY = resXatTrkYTopol;
00281       resYprime = resYTopol*vOrientation;
00282       
00283     } else { // not a detUnit, so must be a virtual 2D-Module
00284       // FIXME: at present only for det units residuals are calculated and filled in the hitStruct
00285       // But in principle this method should also be useable for the gluedDets (2D modules in TIB, TID, TOB, TEC)
00286       // In this case, only orientation should be taken into account for primeResiduals, but not the radial topology
00287       // At present, default values (999.F) are given out
00288     }
00289     
00290     hitStruct.resXprime = resXprime;
00291     hitStruct.resXatTrkY = resXatTrkY;
00292     hitStruct.resYprime = resYprime;
00293     hitStruct.resXprimeErr = resXprimeErr;
00294     hitStruct.resYprimeErr = resYprimeErr;
00295     
00296     hitStruct.rawDetId = IntRawDetID;
00297     hitStruct.phi = tsos.globalDirection().phi();
00298     hitStruct.eta = tsos.globalDirection().eta();
00299     
00300     v_avhitout.push_back(hitStruct);
00301   } 
00302 }
00303 
00304 void
00305 TrackerValidationVariables::fillTrackQuantities(const edm::Event& event, std::vector<AVTrackStruct> & v_avtrackout)
00306 {
00307   edm::InputTag TrjTrackTag = conf_.getParameter<edm::InputTag>("Tracks");
00308   edm::Handle<TrajTrackAssociationCollection> TrajTracksMap;
00309   event.getByLabel(TrjTrackTag, TrajTracksMap);
00310   LogDebug("TrackerValidationVariables") << "TrajTrack collection size " << TrajTracksMap->size();
00311   
00312   const Trajectory* trajectory;
00313   const reco::Track* track;
00314   
00315   for ( TrajTrackAssociationCollection::const_iterator iPair = TrajTracksMap->begin();
00316         iPair != TrajTracksMap->end();
00317         ++iPair) {
00318     
00319     trajectory = &(*(*iPair).key);
00320     track = &(*(*iPair).val);
00321     
00322     AVTrackStruct trackStruct;
00323     
00324     trackStruct.p = track->p();
00325     trackStruct.pt = track->pt();
00326     trackStruct.ptError = track->ptError();
00327     trackStruct.px = track->px();
00328     trackStruct.py = track->py();
00329     trackStruct.pz = track->pz();
00330     trackStruct.eta = track->eta();
00331     trackStruct.phi = track->phi();
00332     trackStruct.chi2 = track->chi2();
00333     trackStruct.chi2Prob= TMath::Prob(track->chi2(),track->ndof());
00334     trackStruct.normchi2 = track->normalizedChi2();
00335     GlobalPoint gPoint(track->vx(), track->vy(), track->vz());
00336     double theLocalMagFieldInInverseGeV = magneticField_->inInverseGeV(gPoint).z();
00337     trackStruct.kappa = -track->charge()*theLocalMagFieldInInverseGeV/track->pt();
00338     trackStruct.charge = track->charge();
00339     trackStruct.d0 = track->d0();
00340     trackStruct.dz = track->dz();
00341     trackStruct.numberOfValidHits = track->numberOfValidHits();
00342     trackStruct.numberOfLostHits = track->numberOfLostHits();
00343     
00344     fillHitQuantities(trajectory, trackStruct.hits);
00345     
00346     v_avtrackout.push_back(trackStruct);
00347   }
00348 }
00349 
00350 void
00351 TrackerValidationVariables::fillHitQuantities(const edm::Event& event, std::vector<AVHitStruct> & v_avhitout)
00352 {
00353   edm::Handle<std::vector<Trajectory> > trajCollectionHandle;
00354   event.getByLabel(conf_.getParameter<std::string>("trajectoryInput"), trajCollectionHandle);
00355   
00356   LogDebug("TrackerValidationVariables") << "trajColl->size(): " << trajCollectionHandle->size() ;
00357 
00358   for (std::vector<Trajectory>::const_iterator it = trajCollectionHandle->begin(), itEnd = trajCollectionHandle->end(); 
00359        it!=itEnd;
00360        ++it) {
00361     
00362     fillHitQuantities(&(*it), v_avhitout);
00363   }
00364 }