Go to the documentation of this file.00001
00013 #include "RecoMuon/GlobalTrackingTools/interface/MuonTrackingRegionBuilder.h"
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024 #include "FWCore/Framework/interface/Event.h"
00025 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00026 #include "DataFormats/Common/interface/Handle.h"
00027 #include "RecoMuon/TrackingTools/interface/MuonServiceProxy.h"
00028 #include "RecoTracker/TkTrackingRegions/interface/RectangularEtaPhiTrackingRegion.h"
00029 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00030 #include "TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h"
00031 #include "DataFormats/VertexReco/interface/Vertex.h"
00032 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00033 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
00034 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00035
00036 using namespace std;
00037
00038
00039
00040
00041
00042 void MuonTrackingRegionBuilder::init(const MuonServiceProxy* service) { theService= service;}
00043 MuonTrackingRegionBuilder::MuonTrackingRegionBuilder(const edm::ParameterSet& par)
00044 {
00045 build(par);
00046 }
00047 void MuonTrackingRegionBuilder::build(const edm::ParameterSet& par){
00048
00049 theBeamSpotTag = par.getParameter<edm::InputTag>("beamSpot");
00050 theVertexCollTag = par.getParameter<edm::InputTag>("vertexCollection");
00051
00052
00053 theNsigmaEta = par.getParameter<double>("Rescale_eta");
00054 theNsigmaPhi = par.getParameter<double>("Rescale_phi");
00055 theNsigmaDz = par.getParameter<double>("Rescale_Dz");
00056 theTkEscapePt = par.getParameter<double>("EscapePt");
00057
00058
00059 theEtaRegionPar1 = par.getParameter<double>("EtaR_UpperLimit_Par1");
00060 theEtaRegionPar2 = par.getParameter<double>("EtaR_UpperLimit_Par2");
00061 thePhiRegionPar1 = par.getParameter<double>("PhiR_UpperLimit_Par1");
00062 thePhiRegionPar2 = par.getParameter<double>("PhiR_UpperLimit_Par2");
00063
00064 useVertex = par.getParameter<bool>("UseVertex");
00065 useFixedRegion = par.getParameter<bool>("UseFixedRegion");
00066
00067
00068 thePhiFixed = par.getParameter<double>("Phi_fixed");
00069 theEtaFixed = par.getParameter<double>("Eta_fixed");
00070
00071 thePhiMin = par.getParameter<double>("Phi_min");
00072 theEtaMin = par.getParameter<double>("Eta_min");
00073 theHalfZ = par.getParameter<double>("DeltaZ_Region");
00074 theDeltaR = par.getParameter<double>("DeltaR");
00075
00076
00077
00078 theOnDemand = par.getParameter<double>("OnDemand");
00079 theMeasurementTrackerName = par.getParameter<std::string>("MeasurementTrackerName");
00080 }
00081
00082
00083
00084
00085
00086 RectangularEtaPhiTrackingRegion*
00087 MuonTrackingRegionBuilder::region(const reco::TrackRef& track) const {
00088
00089 return region(*track);
00090
00091 }
00092
00093
00094
00095
00096
00097 void MuonTrackingRegionBuilder::setEvent(const edm::Event& event) {
00098
00099 theEvent = &event;
00100
00101 }
00102
00103
00104
00105
00106
00107 RectangularEtaPhiTrackingRegion*
00108 MuonTrackingRegionBuilder::region(const reco::Track& staTrack) const {
00109
00110
00111 TSCPBuilderNoMaterial tscpBuilder;
00112 TrajectoryStateTransform tsTransform;
00113 FreeTrajectoryState muFTS = tsTransform.initialFreeState(staTrack,&*theService->magneticField());
00114
00115 LogDebug("MuonTrackingRegionBuilder")<<"from state: "<<muFTS;
00116
00117
00118 GlobalVector dirVector(muFTS.momentum());
00119
00120
00121 const math::XYZVector& mo = staTrack.innerMomentum();
00122 GlobalVector mom(mo.x(),mo.y(),mo.z());
00123 if ( staTrack.p() > 1.5 ) {
00124 mom = dirVector;
00125 }
00126
00127
00128 GlobalPoint vertexPos(0.0,0.0,0.0);
00129 double deltaZatVTX = 0.0;
00130
00131
00132 edm::Handle<reco::BeamSpot> bsHandle;
00133 bool bsHandleFlag = theEvent->getByLabel(theBeamSpotTag, bsHandle);
00134
00135
00136
00137 if ( bsHandleFlag && !useVertex ) {
00138 const reco::BeamSpot& bs = *bsHandle;
00139 vertexPos = GlobalPoint(bs.x0(), bs.y0(), bs.z0());
00140 } else {
00141
00142 edm::Handle<reco::VertexCollection> vertexCollection;
00143 bool vtxHandleFlag = theEvent->getByLabel(theVertexCollTag,vertexCollection);
00144
00145 if ( vtxHandleFlag && !vertexCollection->empty() ) {
00146
00147 reco::VertexCollection::const_iterator vtx = vertexCollection->begin();
00148 vertexPos = GlobalPoint(vtx->x(),vtx->y(),vtx->z());
00149
00150 deltaZatVTX = vtx->zError() * theNsigmaDz;
00151 }
00152 }
00153
00154
00155
00156
00157
00158 double deta = 0.4;
00159 double dphi = 0.6;
00160
00161
00162 if ( bsHandleFlag ) {
00163
00164 const reco::BeamSpot& bs = *bsHandle;
00165
00166 TSCBLBuilderNoMaterial tscblBuilder;
00167 TrajectoryStateClosestToBeamLine tscbl = tscblBuilder(muFTS,bs);
00168
00169
00170
00171 if(tscbl.isValid()){
00172
00173 PerigeeConversions tspConverter;
00174 PerigeeTrajectoryError trackPerigeeErrors = tspConverter.ftsToPerigeeError(tscbl.trackStateAtPCA());
00175 GlobalVector pTrack = tscbl.trackStateAtPCA().momentum();
00176
00177
00178 double deltaTheta = trackPerigeeErrors.thetaError();
00179 double theta = pTrack.theta();
00180 double sin_theta = sin(theta);
00181
00182
00183 deta = theNsigmaEta*(1/fabs(sin_theta))*deltaTheta;
00184 dphi = theNsigmaPhi*(trackPerigeeErrors.phiError());
00185
00186 }
00187 }
00188
00189
00190
00191
00192
00193 double region_dEta = 0;
00194 double region_dPhi = 0;
00195 double eta=0; double phi=0;
00196
00197
00198 float pt = fabs(mom.perp());
00199
00200 if ( pt <= 10. ) {
00201
00202 float acoeff_Phi = (thePhiRegionPar2 - thePhiRegionPar1)/5;
00203 float acoeff_Eta = (theEtaRegionPar2 - theEtaRegionPar1)/5;
00204
00205 eta = theEtaRegionPar1 + (acoeff_Eta)*(mom.perp()-5.);
00206 phi = thePhiRegionPar1 + (acoeff_Phi)*(mom.perp()-5.) ;
00207 }
00208
00209 if ( pt > 10. && pt < 100. ) {
00210 eta = theEtaRegionPar2;
00211 phi = thePhiRegionPar2;
00212 }
00213
00214 if ( pt >= 100. ) {
00215
00216 float acoeff_Phi = (thePhiRegionPar1 - thePhiRegionPar2)/900;
00217 float acoeff_Eta = (theEtaRegionPar1 - theEtaRegionPar2)/900;
00218
00219 eta = theEtaRegionPar2 + (acoeff_Eta)*(mom.perp()-100.);
00220 phi = thePhiRegionPar2 + (acoeff_Phi)*(mom.perp()-100.);
00221 }
00222
00223
00224 double region_dPhi1 = min(phi,dphi);
00225 double region_dEta1 = min(eta,deta);
00226
00227
00228 region_dPhi = max(thePhiMin,region_dPhi1);
00229 region_dEta = max(theEtaMin,region_dEta1);
00230
00231 float deltaZ = 0.0;
00232
00233 if ( useVertex ) {
00234 deltaZ = max(theHalfZ,deltaZatVTX);
00235 } else {
00236 deltaZ = theHalfZ;
00237 }
00238
00239 float deltaR = theDeltaR;
00240 double minPt = max(theTkEscapePt,mom.perp()*0.6);
00241
00242 RectangularEtaPhiTrackingRegion* region = 0;
00243
00244 if (useFixedRegion) {
00245 region_dEta = theEtaFixed;
00246 region_dPhi = thePhiFixed;
00247 }
00248
00249 region = new RectangularEtaPhiTrackingRegion(dirVector, vertexPos,
00250 minPt, deltaR,
00251 deltaZ, region_dEta, region_dPhi,
00252 theOnDemand,
00253 true,
00254 theMeasurementTrackerName);
00255
00256 LogDebug("MuonTrackingRegionBuilder")<<"the region parameters are:\n"
00257 <<"\n dirVector: "<<dirVector
00258 <<"\n vertexPos: "<<vertexPos
00259 <<"\n minPt: "<<minPt
00260 <<"\n deltaR:"<<deltaR
00261 <<"\n deltaZ:"<<deltaZ
00262 <<"\n region_dEta:"<<region_dEta
00263 <<"\n region_dPhi:"<<region_dPhi
00264 <<"\n on demand parameter: "<<theOnDemand;
00265
00266
00267 return region;
00268
00269 }