CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoMET/METProducers/src/MuonTCMETValueMapProducer.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    MuonTCMETValueMapProducer
00004 // Class:      MuonTCMETValueMapProducer
00005 // 
00013 //
00014 // Original Author:  Frank Golf
00015 //         Created:  Sun Mar 15 11:33:20 CDT 2009
00016 // $Id: MuonTCMETValueMapProducer.cc,v 1.9 2011/02/21 10:09:21 benhoob Exp $
00017 //
00018 //
00019 
00020 // system include files
00021 #include <memory>
00022 
00023 // user include files
00024 #include "RecoMET/METProducers/interface/MuonTCMETValueMapProducer.h"
00025 
00026 #include "RecoMET/METAlgorithms/interface/TCMETAlgo.h"
00027 
00028 #include "DataFormats/MuonReco/interface/Muon.h"
00029 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00030 #include "DataFormats/TrackReco/interface/Track.h"
00031 #include "DataFormats/TrackReco/interface/TrackBase.h"
00032 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00033 #include "DataFormats/Common/interface/ValueMap.h" 
00034 #include "DataFormats/MuonReco/interface/MuonMETCorrectionData.h"
00035 #include "DataFormats/GeometrySurface/interface/Plane.h"
00036 #include "DataFormats/GeometrySurface/interface/Cylinder.h"
00037 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00038 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00039 #include "DataFormats/Math/interface/Point3D.h"
00040 
00041 #include "MagneticField/Engine/interface/MagneticField.h"
00042 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00043 
00044 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00045 #include "TrackingTools/GeomPropagators/interface/AnalyticalPropagator.h"
00046 
00047 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00048 #include "FWCore/Framework/interface/MakerMacros.h"
00049 
00050 #include "TH2D.h"
00051 #include "TVector3.h"
00052 #include "TMath.h"
00053 
00054 typedef math::XYZTLorentzVector LorentzVector;
00055 typedef math::XYZPoint Point;
00056 
00057 namespace cms {
00058   MuonTCMETValueMapProducer::MuonTCMETValueMapProducer(const edm::ParameterSet& iConfig) {
00059   
00060     produces<edm::ValueMap<reco::MuonMETCorrectionData> > ("muCorrData");
00061 
00062     // get input collections
00063     muonInputTag_     = iConfig.getParameter<edm::InputTag>("muonInputTag"    );
00064     beamSpotInputTag_ = iConfig.getParameter<edm::InputTag>("beamSpotInputTag");
00065     vertexInputTag_   = iConfig.getParameter<edm::InputTag>("vertexInputTag");
00066 
00067     rfType_     = iConfig.getParameter<int>("rf_type");
00068 
00069     nLayers_                = iConfig.getParameter<int>      ("nLayers");
00070     nLayersTight_           = iConfig.getParameter<int>      ("nLayersTight");
00071     vertexNdof_             = iConfig.getParameter<int>      ("vertexNdof");
00072     vertexZ_                = iConfig.getParameter<double>   ("vertexZ");
00073     vertexRho_              = iConfig.getParameter<double>   ("vertexRho");
00074     vertexMaxDZ_            = iConfig.getParameter<double>   ("vertexMaxDZ");
00075     maxpt_eta20_            = iConfig.getParameter<double>   ("maxpt_eta20");
00076     maxpt_eta25_            = iConfig.getParameter<double>   ("maxpt_eta25");
00077 
00078     // get configuration parameters
00079     maxTrackAlgo_    = iConfig.getParameter<int>("trackAlgo_max");
00080     maxd0cut_        = iConfig.getParameter<double>("d0_max"       );
00081     minpt_           = iConfig.getParameter<double>("pt_min"       );
00082     maxpt_           = iConfig.getParameter<double>("pt_max"       );
00083     maxeta_          = iConfig.getParameter<double>("eta_max"      );
00084     maxchi2_         = iConfig.getParameter<double>("chi2_max"     );
00085     minhits_         = iConfig.getParameter<double>("nhits_min"    );
00086     maxPtErr_        = iConfig.getParameter<double>("ptErr_max"    );
00087 
00088     trkQuality_      = iConfig.getParameter<std::vector<int> >("track_quality");
00089     trkAlgos_        = iConfig.getParameter<std::vector<int> >("track_algos"  );
00090     maxchi2_tight_   = iConfig.getParameter<double>("chi2_max_tight");
00091     minhits_tight_   = iConfig.getParameter<double>("nhits_min_tight");
00092     maxPtErr_tight_  = iConfig.getParameter<double>("ptErr_max_tight");
00093     usePvtxd0_       = iConfig.getParameter<bool>("usePvtxd0");
00094     d0cuta_          = iConfig.getParameter<double>("d0cuta");
00095     d0cutb_          = iConfig.getParameter<double>("d0cutb");
00096 
00097     muon_dptrel_  = iConfig.getParameter<double>("muon_dptrel");
00098     muond0_     = iConfig.getParameter<double>("d0_muon"    );
00099     muonpt_     = iConfig.getParameter<double>("pt_muon"    );
00100     muoneta_    = iConfig.getParameter<double>("eta_muon"   );
00101     muonchi2_   = iConfig.getParameter<double>("chi2_muon"  );
00102     muonhits_   = iConfig.getParameter<double>("nhits_muon" );
00103     muonGlobal_   = iConfig.getParameter<bool>("global_muon");
00104     muonTracker_  = iConfig.getParameter<bool>("tracker_muon");
00105     muonDeltaR_ = iConfig.getParameter<double>("deltaR_muon");
00106     useCaloMuons_ = iConfig.getParameter<bool>("useCaloMuons");
00107     muonMinValidStaHits_ = iConfig.getParameter<int>("muonMinValidStaHits");
00108 
00109     response_function = 0;
00110     tcmetAlgo_=new TCMETAlgo();
00111   }
00112 
00113   MuonTCMETValueMapProducer::~MuonTCMETValueMapProducer()
00114   {
00115  
00116     // do anything here that needs to be done at desctruction time
00117     // (e.g. close files, deallocate resources etc.)
00118     delete tcmetAlgo_;
00119   }
00120 
00121   //
00122   // member functions
00123   //
00124 
00125   // ------------ method called to produce the data  ------------
00126   void MuonTCMETValueMapProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
00127   
00128     //get input collections
00129     iEvent.getByLabel(muonInputTag_    , muon_h    );
00130     iEvent.getByLabel(beamSpotInputTag_, beamSpot_h);
00131 
00132     //get vertex collection
00133     hasValidVertex = false;
00134     if( usePvtxd0_ ){
00135       iEvent.getByLabel( vertexInputTag_  , VertexHandle  );
00136        
00137       if( VertexHandle.isValid() ) {
00138         vertexColl = VertexHandle.product();
00139         hasValidVertex = isValidVertex();
00140       } 
00141     }
00142 
00143     //get the Bfield
00144     edm::ESHandle<MagneticField> theMagField;
00145     bool haveBfield = true;
00146     if( !theMagField.isValid() ) haveBfield = false;
00147     iSetup.get<IdealMagneticFieldRecord>().get(theMagField);
00148     bField = theMagField.product();
00149 
00150     //make a ValueMap of ints => flags for 
00151     //met correction. The values and meanings of the flags are :
00152     // flag==0 --->    The muon is not used to correct the MET by default
00153     // flag==1 --->    The muon is used to correct the MET. The Global pt is used.
00154     // flag==2 --->    The muon is used to correct the MET. The tracker pt is used.
00155     // flag==3 --->    The muon is used to correct the MET. The standalone pt is used.
00156     // flag==4 --->    The muon is used to correct the MET as pion using the tcMET ZSP RF.
00157     // flag==5 --->    The muon is used to correct the MET.  The default fit is used; i.e. we get the pt from muon->pt().
00158     // In general, the flag should never be 3. You do not want to correct the MET using
00159     // the pt measurement from the standalone system (unless you really know what you're 
00160     // doing
00161 
00162     std::auto_ptr<edm::ValueMap<reco::MuonMETCorrectionData> > vm_muCorrData(new edm::ValueMap<reco::MuonMETCorrectionData>());
00163 
00164     std::vector<reco::MuonMETCorrectionData> v_muCorrData;
00165 
00166     unsigned int nMuons = muon_h->size();
00167 
00168     for (unsigned int iMu = 0; iMu < nMuons; iMu++) {
00169 
00170       const reco::Muon* mu = &(*muon_h)[iMu];
00171       double deltax = 0.0;
00172       double deltay = 0.0;
00173 
00174       reco::MuonMETCorrectionData muMETCorrData(reco::MuonMETCorrectionData::NotUsed, deltax, deltay);
00175     
00176       reco::TrackRef mu_track;
00177       if( mu->isGlobalMuon() || mu->isTrackerMuon() || mu->isCaloMuon() )
00178         mu_track = mu->innerTrack();
00179       else {
00180         v_muCorrData.push_back( muMETCorrData );
00181         continue;
00182       }
00183 
00184       // figure out depositions muons would make if they were treated as pions
00185       if( isGoodTrack( mu ) ) {
00186 
00187         if( mu_track->pt() < minpt_ ) 
00188           muMETCorrData = reco::MuonMETCorrectionData(reco::MuonMETCorrectionData::TreatedAsPion, deltax, deltay);
00189 
00190         else {
00191           int bin_index   = response_function->FindBin( mu_track->eta(), mu_track->pt() );
00192           double response = response_function->GetBinContent( bin_index );
00193 
00194           TVector3 outerTrkPosition = propagateTrack( mu );
00195 
00196           deltax = response * mu_track->p() * sin( outerTrkPosition.Theta() ) * cos( outerTrkPosition.Phi() );
00197           deltay = response * mu_track->p() * sin( outerTrkPosition.Theta() ) * sin( outerTrkPosition.Phi() );
00198 
00199           muMETCorrData = reco::MuonMETCorrectionData(reco::MuonMETCorrectionData::TreatedAsPion, deltax, deltay);
00200         }
00201       }
00202 
00203       // figure out muon flag
00204       if( isGoodMuon( mu ) )
00205         v_muCorrData.push_back( reco::MuonMETCorrectionData(reco::MuonMETCorrectionData::MuonCandidateValuesUsed, deltax, deltay) );
00206 
00207       else if( useCaloMuons_ && isGoodCaloMuon( mu, iMu ) )
00208         v_muCorrData.push_back( reco::MuonMETCorrectionData(reco::MuonMETCorrectionData::MuonCandidateValuesUsed, deltax, deltay) );
00209 
00210       else v_muCorrData.push_back( muMETCorrData );
00211     }
00212     
00213     edm::ValueMap<reco::MuonMETCorrectionData>::Filler dataFiller(*vm_muCorrData);
00214 
00215     dataFiller.insert( muon_h, v_muCorrData.begin(), v_muCorrData.end());
00216     dataFiller.fill();
00217     
00218     iEvent.put(vm_muCorrData, "muCorrData");    
00219   }
00220   
00221   // ------------ method called once each job just before starting event loop  ------------
00222   void MuonTCMETValueMapProducer::beginJob()
00223   {
00224 
00225     if( rfType_ == 1 )
00226                  response_function = tcmetAlgo_->getResponseFunction_fit();
00227     else if( rfType_ == 2 )
00228                  response_function = tcmetAlgo_->getResponseFunction_mode();
00229   }
00230 
00231   // ------------ method called once each job just after ending the event loop  ------------
00232   void MuonTCMETValueMapProducer::endJob() {
00233   }
00234 
00235   // ------------ check is muon is a good muon  ------------
00236   bool MuonTCMETValueMapProducer::isGoodMuon( const reco::Muon* muon ) {
00237     double d0    = -999;
00238     double nhits = 0;
00239     double chi2  = 999;  
00240 
00241     // get d0 corrected for beam spot
00242     bool haveBeamSpot = true;
00243     if( !beamSpot_h.isValid() ) haveBeamSpot = false;
00244 
00245     if( muonGlobal_  && !muon->isGlobalMuon()  ) return false;
00246     if( muonTracker_ && !muon->isTrackerMuon() ) return false;
00247 
00248     const reco::TrackRef siTrack     = muon->innerTrack();
00249     const reco::TrackRef globalTrack = muon->globalTrack();
00250 
00251     Point bspot = haveBeamSpot ? beamSpot_h->position() : Point(0,0,0);
00252     if( siTrack.isNonnull() ) nhits = siTrack->numberOfValidHits();
00253     if( globalTrack.isNonnull() ) {
00254       d0   = -1 * globalTrack->dxy( bspot );
00255       chi2 = globalTrack->normalizedChi2();
00256     }
00257 
00258     if( fabs( d0 ) > muond0_ )                          return false;
00259     if( muon->pt() < muonpt_ )                          return false;
00260     if( fabs( muon->eta() ) > muoneta_ )                return false;
00261     if( nhits < muonhits_ )                             return false;
00262     if( chi2 > muonchi2_ )                              return false;
00263     if( globalTrack->hitPattern().numberOfValidMuonHits() < muonMinValidStaHits_ ) return false;
00264 
00265     //reject muons with tracker dpt/pt > X
00266     if( !siTrack.isNonnull() )                                return false;
00267     if( siTrack->ptError() / siTrack->pt() > muon_dptrel_ )   return false;
00268 
00269     else return true;
00270   }
00271 
00272   // ------------ check if muon is a good calo muon  ------------
00273   bool MuonTCMETValueMapProducer::isGoodCaloMuon( const reco::Muon* muon, const unsigned int index ) {
00274 
00275     if( muon->pt() < 10 ) return false;
00276 
00277     if( !isGoodTrack( muon ) ) return false;
00278 
00279     const reco::TrackRef inputSiliconTrack = muon->innerTrack();
00280     if( !inputSiliconTrack.isNonnull() ) return false;
00281 
00282     //check if it is in the vicinity of a global or tracker muon
00283     unsigned int nMuons = muon_h->size();
00284     for (unsigned int iMu = 0; iMu < nMuons; iMu++) {
00285 
00286       if( iMu == index ) continue;
00287 
00288       const reco::Muon* mu = &(*muon_h)[iMu];
00289 
00290       const reco::TrackRef testSiliconTrack = mu->innerTrack();
00291       if( !testSiliconTrack.isNonnull() ) continue;
00292 
00293       double deltaEta = inputSiliconTrack.get()->eta() - testSiliconTrack.get()->eta();
00294       double deltaPhi = acos( cos( inputSiliconTrack.get()->phi() - testSiliconTrack.get()->phi() ) );
00295       double deltaR   = TMath::Sqrt( deltaEta * deltaEta + deltaPhi * deltaPhi );
00296 
00297       if( deltaR < muonDeltaR_ ) return false;
00298     }
00299 
00300     return true;
00301   }
00302 
00303   // ------------ check if track is good  ------------
00304   bool MuonTCMETValueMapProducer::isGoodTrack( const reco::Muon* muon ) {
00305     double d0    = -999;
00306 
00307     const reco::TrackRef siTrack = muon->innerTrack();
00308     if (!siTrack.isNonnull())
00309       return false;
00310 
00311     if( hasValidVertex ){
00312       //get d0 corrected for primary vertex
00313             
00314       const Point pvtx = Point(vertexColl->begin()->x(),
00315                                vertexColl->begin()->y(), 
00316                                vertexColl->begin()->z());
00317             
00318       d0 = -1 * siTrack->dxy( pvtx );
00319             
00320       double dz = siTrack->dz( pvtx );
00321             
00322       if( fabs( dz ) < vertexMaxDZ_ ){
00323               
00324         //get d0 corrected for pvtx
00325         d0 = -1 * siTrack->dxy( pvtx );
00326               
00327       }else{
00328               
00329         // get d0 corrected for beam spot
00330         bool haveBeamSpot = true;
00331         if( !beamSpot_h.isValid() ) haveBeamSpot = false;
00332               
00333         Point bspot = haveBeamSpot ? beamSpot_h->position() : Point(0,0,0);
00334         d0 = -1 * siTrack->dxy( bspot );
00335               
00336       }
00337     }else{
00338        
00339       // get d0 corrected for beam spot
00340       bool haveBeamSpot = true;
00341       if( !beamSpot_h.isValid() ) haveBeamSpot = false;
00342        
00343       Point bspot = haveBeamSpot ? beamSpot_h->position() : Point(0,0,0);
00344       d0 = -1 * siTrack->dxy( bspot );
00345     }
00346 
00347     if( siTrack->algo() < maxTrackAlgo_ ){
00348       //1st 4 tracking iterations (pT-dependent d0 cut)
00349        
00350       float d0cut = sqrt(std::pow(d0cuta_,2) + std::pow(d0cutb_/siTrack->pt(),2)); 
00351       if(d0cut > maxd0cut_) d0cut = maxd0cut_;
00352        
00353       if( fabs( d0 ) > d0cut )            return false;    
00354       if( nLayers( siTrack ) < nLayers_ ) return false;
00355     }
00356     else{
00357       //last 2 tracking iterations (tighten chi2, nhits, pt error cuts)
00358      
00359       if( siTrack->normalizedChi2() > maxchi2_tight_ )               return false;
00360       if( siTrack->numberOfValidHits() < minhits_tight_ )            return false;
00361       if( (siTrack->ptError() / siTrack->pt()) > maxPtErr_tight_ )   return false;
00362       if( nLayers( siTrack ) < nLayersTight_ )                       return false;
00363     }
00364 
00365     if( siTrack->numberOfValidHits() < minhits_ )                         return false;
00366     if( siTrack->normalizedChi2() > maxchi2_ )                            return false;
00367     if( fabs( siTrack->eta() ) > maxeta_ )                                return false;
00368     if( siTrack->pt() > maxpt_ )                                          return false;
00369     if( (siTrack->ptError() / siTrack->pt()) > maxPtErr_ )                return false;
00370     if( fabs( siTrack->eta() ) > 2.5 && siTrack->pt() > maxpt_eta25_ )    return false;
00371     if( fabs( siTrack->eta() ) > 2.0 && siTrack->pt() > maxpt_eta20_ )    return false;
00372 
00373     int cut = 0;          
00374     for( unsigned int i = 0; i < trkQuality_.size(); i++ ) {
00375 
00376       cut |= (1 << trkQuality_.at(i));
00377     }
00378 
00379     if( !( (siTrack->qualityMask() & cut) == cut ) ) return false;
00380           
00381     bool isGoodAlgo = false;    
00382     if( trkAlgos_.size() == 0 ) isGoodAlgo = true;
00383     for( unsigned int i = 0; i < trkAlgos_.size(); i++ ) {
00384 
00385       if( siTrack->algo() == trkAlgos_.at(i) ) isGoodAlgo = true;
00386     }
00387 
00388     if( !isGoodAlgo ) return false;
00389           
00390     return true;
00391   }
00392 
00393   // ------------ propagate track to calorimeter face  ------------
00394   TVector3 MuonTCMETValueMapProducer::propagateTrack( const reco::Muon* muon) {
00395 
00396     TVector3 outerTrkPosition;
00397 
00398     outerTrkPosition.SetPtEtaPhi( 999., -10., 2 * TMath::Pi() );
00399 
00400     const reco::TrackRef track = muon->innerTrack();
00401 
00402     if( !track.isNonnull() ) {
00403       return outerTrkPosition;
00404     }
00405 
00406     GlobalPoint  tpVertex ( track->vx(), track->vy(), track->vz() );
00407     GlobalVector tpMomentum ( track.get()->px(), track.get()->py(), track.get()->pz() );
00408     int tpCharge ( track->charge() );
00409 
00410     FreeTrajectoryState fts ( tpVertex, tpMomentum, tpCharge, bField);
00411 
00412     const double zdist = 314.;
00413 
00414     const double radius = 130.;
00415 
00416     const double corner = 1.479;
00417 
00418     Plane::PlanePointer lendcap = Plane::build( Plane::PositionType (0, 0, -zdist), Plane::RotationType () );
00419     Plane::PlanePointer rendcap = Plane::build( Plane::PositionType (0, 0, zdist),  Plane::RotationType () );
00420 
00421     Cylinder::CylinderPointer barrel = Cylinder::build( Cylinder::PositionType (0, 0, 0), Cylinder::RotationType (), radius);
00422 
00423     AnalyticalPropagator myAP (bField, alongMomentum, 2*M_PI);
00424 
00425     TrajectoryStateOnSurface tsos;
00426 
00427     if( track.get()->eta() < -corner ) {
00428       tsos = myAP.propagate( fts, *lendcap);
00429     }
00430     else if( fabs(track.get()->eta()) < corner ) {
00431       tsos = myAP.propagate( fts, *barrel);
00432     }
00433     else if( track.get()->eta() > corner ) {
00434       tsos = myAP.propagate( fts, *rendcap);
00435     }
00436 
00437     if( tsos.isValid() )
00438       outerTrkPosition.SetXYZ( tsos.globalPosition().x(), tsos.globalPosition().y(), tsos.globalPosition().z() );
00439 
00440     else 
00441       outerTrkPosition.SetPtEtaPhi( 999., -10., 2 * TMath::Pi() );
00442 
00443     return outerTrkPosition;
00444   }
00445 
00446   // ------------ single pion response function from fit  ------------
00447 
00448   int MuonTCMETValueMapProducer::nLayers(const reco::TrackRef track){
00449     const reco::HitPattern& p = track->hitPattern();
00450     return p.trackerLayersWithMeasurement();
00451   }
00452 
00453   //--------------------------------------------------------------------
00454 
00455   bool MuonTCMETValueMapProducer::isValidVertex(){
00456     
00457     if( vertexColl->begin()->isFake()                ) return false;
00458     if( vertexColl->begin()->ndof() < vertexNdof_    ) return false;
00459     if( fabs( vertexColl->begin()->z() ) > vertexZ_  ) return false;
00460     if( sqrt( std::pow( vertexColl->begin()->x() , 2 ) + std::pow( vertexColl->begin()->y() , 2 ) ) > vertexRho_ ) return false;
00461     
00462     return true;
00463     
00464   }
00465 }
00466