CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/RecoMuon/MuonIdentification/src/MuonTimingFiller.cc

Go to the documentation of this file.
00001 //
00002 // Package:    MuonTimingFiller
00003 // Class:      MuonTimingFiller
00004 // 
00012 //
00013 // Original Author:  Piotr Traczyk, CERN
00014 //         Created:  Mon Mar 16 12:27:22 CET 2009
00015 // $Id: MuonTimingFiller.cc,v 1.13 2011/02/24 15:41:53 farrell3 Exp $
00016 //
00017 //
00018 
00019 
00020 // system include files
00021 #include <memory>
00022 
00023 // user include files
00024 #include "FWCore/Framework/interface/Frameworkfwd.h"
00025 #include "FWCore/Framework/interface/EDProducer.h"
00026 
00027 #include "FWCore/Framework/interface/Event.h"
00028 #include "FWCore/Framework/interface/MakerMacros.h"
00029 
00030 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00031 
00032 #include "DataFormats/MuonReco/interface/Muon.h" 
00033 #include "DataFormats/MuonReco/interface/MuonFwd.h" 
00034 #include "DataFormats/MuonReco/interface/MuonTimeExtra.h"
00035 #include "DataFormats/MuonReco/interface/MuonTimeExtraMap.h"
00036 
00037 #include "RecoMuon/MuonIdentification/interface/MuonTimingFiller.h"
00038 #include "RecoMuon/MuonIdentification/interface/TimeMeasurementSequence.h"
00039 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
00040 
00041 //
00042 // constructors and destructor
00043 //
00044 MuonTimingFiller::MuonTimingFiller(const edm::ParameterSet& iConfig)
00045 {
00046    // Load parameters for the DTTimingExtractor
00047    edm::ParameterSet dtTimingParameters = iConfig.getParameter<edm::ParameterSet>("DTTimingParameters");
00048    theDTTimingExtractor_ = new DTTimingExtractor(dtTimingParameters);
00049 
00050    // Load parameters for the CSCTimingExtractor
00051    edm::ParameterSet cscTimingParameters = iConfig.getParameter<edm::ParameterSet>("CSCTimingParameters");
00052    theCSCTimingExtractor_ = new CSCTimingExtractor(cscTimingParameters);
00053    
00054    errorEB_ = iConfig.getParameter<double>("ErrorEB");
00055    errorEE_ = iConfig.getParameter<double>("ErrorEE");
00056    ecalEcut_ = iConfig.getParameter<double>("EcalEnergyCut");
00057    
00058    useDT_ = iConfig.getParameter<bool>("UseDT");
00059    useCSC_ = iConfig.getParameter<bool>("UseCSC");
00060    useECAL_ = iConfig.getParameter<bool>("UseECAL");
00061    
00062 }
00063 
00064 
00065 MuonTimingFiller::~MuonTimingFiller()
00066 {
00067   if (theDTTimingExtractor_) delete theDTTimingExtractor_;
00068   if (theCSCTimingExtractor_) delete theCSCTimingExtractor_;
00069 }
00070 
00071 
00072 //
00073 // member functions
00074 //
00075 
00076 void 
00077 MuonTimingFiller::fillTiming( const reco::Muon& muon, reco::MuonTimeExtra& dtTime, reco::MuonTimeExtra& cscTime, reco::MuonTimeExtra& combinedTime, edm::Event& iEvent, const edm::EventSetup& iSetup )
00078 {
00079   TimeMeasurementSequence dtTmSeq,cscTmSeq;
00080      
00081   if ( !(muon.combinedMuon().isNull()) ) {
00082     theDTTimingExtractor_->fillTiming(dtTmSeq, muon.combinedMuon(), iEvent, iSetup);
00083     theCSCTimingExtractor_->fillTiming(cscTmSeq, muon.combinedMuon(), iEvent, iSetup);
00084   } else
00085     if ( !(muon.standAloneMuon().isNull()) ) {
00086       theDTTimingExtractor_->fillTiming(dtTmSeq, muon.standAloneMuon(), iEvent, iSetup);
00087       theCSCTimingExtractor_->fillTiming(cscTmSeq, muon.standAloneMuon(), iEvent, iSetup);
00088     }
00089   
00090   // Fill DT-specific timing information block     
00091   fillTimeFromMeasurements(dtTmSeq, dtTime);
00092 
00093   // Fill CSC-specific timing information block     
00094   fillTimeFromMeasurements(cscTmSeq, cscTime);
00095        
00096   // Combine the TimeMeasurementSequences from all subdetectors
00097   TimeMeasurementSequence combinedTmSeq;
00098   combineTMSequences(muon,dtTmSeq,cscTmSeq,combinedTmSeq);
00099   // add ECAL info
00100   if (useECAL_) addEcalTime(muon,combinedTmSeq);
00101 
00102   // Fill the master timing block
00103   fillTimeFromMeasurements(combinedTmSeq, combinedTime);
00104     
00105   LogTrace("MuonTime") << "Global 1/beta: " << combinedTime.inverseBeta() << " +/- " << combinedTime.inverseBetaErr()<<std::endl;
00106   LogTrace("MuonTime") << "  Free 1/beta: " << combinedTime.freeInverseBeta() << " +/- " << combinedTime.freeInverseBetaErr()<<std::endl;
00107   LogTrace("MuonTime") << "  Vertex time (in-out): " << combinedTime.timeAtIpInOut() << " +/- " << combinedTime.timeAtIpInOutErr()
00108                        << "  # of points: " << combinedTime.nDof() <<std::endl;
00109   LogTrace("MuonTime") << "  Vertex time (out-in): " << combinedTime.timeAtIpOutIn() << " +/- " << combinedTime.timeAtIpOutInErr()<<std::endl;
00110   LogTrace("MuonTime") << "  direction: "   << combinedTime.direction() << std::endl;
00111      
00112 }
00113 
00114 
00115 void 
00116 MuonTimingFiller::fillTimeFromMeasurements( TimeMeasurementSequence tmSeq, reco::MuonTimeExtra &muTime ) {
00117 
00118   std::vector <double> x,y;
00119   double invbeta=0, invbetaerr=0;
00120   double vertexTime=0, vertexTimeErr=0, vertexTimeR=0, vertexTimeRErr=0;    
00121   double freeBeta, freeBetaErr, freeTime, freeTimeErr;
00122 
00123   if (tmSeq.dstnc.size()<=1) return;
00124 
00125   for (unsigned int i=0;i<tmSeq.dstnc.size();i++) {
00126     invbeta+=(1.+tmSeq.local_t0.at(i)/tmSeq.dstnc.at(i)*30.)*tmSeq.weightInvbeta.at(i)/tmSeq.totalWeightInvbeta;
00127     x.push_back(tmSeq.dstnc.at(i)/30.);
00128     y.push_back(tmSeq.local_t0.at(i)+tmSeq.dstnc.at(i)/30.);
00129     vertexTime+=tmSeq.local_t0.at(i)*tmSeq.weightVertex.at(i)/tmSeq.totalWeightVertex;
00130     vertexTimeR+=(tmSeq.local_t0.at(i)+2*tmSeq.dstnc.at(i)/30.)*tmSeq.weightVertex.at(i)/tmSeq.totalWeightVertex;
00131   }
00132 
00133   double diff;
00134   for (unsigned int i=0;i<tmSeq.dstnc.size();i++) {
00135     diff=(1.+tmSeq.local_t0.at(i)/tmSeq.dstnc.at(i)*30.)-invbeta;
00136     invbetaerr+=diff*diff*tmSeq.weightInvbeta.at(i);
00137     diff=tmSeq.local_t0.at(i)-vertexTime;
00138     vertexTimeErr+=diff*diff*tmSeq.weightVertex.at(i);
00139     diff=tmSeq.local_t0.at(i)+2*tmSeq.dstnc.at(i)/30.-vertexTimeR;
00140     vertexTimeRErr+=diff*diff*tmSeq.weightVertex.at(i);
00141   }
00142   
00143   double cf = 1./(tmSeq.dstnc.size()-1);
00144   invbetaerr=sqrt(invbetaerr/tmSeq.totalWeightInvbeta*cf);
00145   vertexTimeErr=sqrt(vertexTimeErr/tmSeq.totalWeightVertex*cf);
00146   vertexTimeRErr=sqrt(vertexTimeRErr/tmSeq.totalWeightVertex*cf);
00147 
00148   muTime.setInverseBeta(invbeta);
00149   muTime.setInverseBetaErr(invbetaerr);
00150   muTime.setTimeAtIpInOut(vertexTime);
00151   muTime.setTimeAtIpInOutErr(vertexTimeErr);
00152   muTime.setTimeAtIpOutIn(vertexTimeR);
00153   muTime.setTimeAtIpOutInErr(vertexTimeRErr);
00154       
00155   rawFit(freeBeta, freeBetaErr, freeTime, freeTimeErr, x, y);
00156 
00157   muTime.setFreeInverseBeta(freeBeta);
00158   muTime.setFreeInverseBetaErr(freeBetaErr);
00159     
00160   muTime.setNDof(tmSeq.dstnc.size());
00161 }
00162 
00163 void 
00164 MuonTimingFiller::combineTMSequences( const reco::Muon& muon, 
00165                                       TimeMeasurementSequence dtSeq, 
00166                                       TimeMeasurementSequence cscSeq, 
00167                                       TimeMeasurementSequence &cmbSeq ) {
00168                                         
00169   if (useDT_) for (unsigned int i=0;i<dtSeq.dstnc.size();i++) {
00170     cmbSeq.dstnc.push_back(dtSeq.dstnc.at(i));
00171     cmbSeq.local_t0.push_back(dtSeq.local_t0.at(i));
00172     cmbSeq.weightVertex.push_back(dtSeq.weightVertex.at(i));
00173     cmbSeq.weightInvbeta.push_back(dtSeq.weightInvbeta.at(i));
00174 
00175     cmbSeq.totalWeightVertex+=dtSeq.weightVertex.at(i);
00176     cmbSeq.totalWeightInvbeta+=dtSeq.weightInvbeta.at(i);
00177   }
00178 
00179   if (useCSC_) for (unsigned int i=0;i<cscSeq.dstnc.size();i++) {
00180     cmbSeq.dstnc.push_back(cscSeq.dstnc.at(i));
00181     cmbSeq.local_t0.push_back(cscSeq.local_t0.at(i));
00182     cmbSeq.weightVertex.push_back(cscSeq.weightVertex.at(i));
00183     cmbSeq.weightInvbeta.push_back(cscSeq.weightInvbeta.at(i));
00184 
00185     cmbSeq.totalWeightVertex+=cscSeq.weightVertex.at(i);
00186     cmbSeq.totalWeightInvbeta+=cscSeq.weightInvbeta.at(i);
00187   }
00188 }
00189 
00190 
00191 void 
00192 MuonTimingFiller::addEcalTime( const reco::Muon& muon, 
00193                                TimeMeasurementSequence &cmbSeq ) {
00194 
00195   reco::MuonEnergy muonE;
00196   if (muon.isEnergyValid())  
00197     muonE = muon.calEnergy();
00198   
00199   // Cut on the crystal energy and restrict to the ECAL barrel for now
00200 //  if (muonE.emMax<ecalEcut_ || fabs(muon.eta())>1.5) return;    
00201   if (muonE.emMax<ecalEcut_) return;    
00202   
00203   // A simple parametrization of the error on the ECAL time measurement
00204   double emErr;
00205   if (muonE.ecal_id.subdetId()==EcalBarrel) emErr= errorEB_/muonE.emMax; else
00206     emErr=errorEE_/muonE.emMax;
00207   double hitWeight = 1/(emErr*emErr);
00208   double hitDist=muonE.ecal_position.r();
00209         
00210   cmbSeq.local_t0.push_back(muonE.ecal_time);
00211   cmbSeq.weightVertex.push_back(hitWeight);
00212   cmbSeq.weightInvbeta.push_back(hitDist*hitDist*hitWeight/(30.*30.));
00213 
00214   cmbSeq.dstnc.push_back(hitDist);
00215   
00216   cmbSeq.totalWeightVertex+=hitWeight;
00217   cmbSeq.totalWeightInvbeta+=hitDist*hitDist*hitWeight/(30.*30.);
00218                                       
00219 }
00220 
00221 
00222 
00223 void 
00224 MuonTimingFiller::rawFit(double &a, double &da, double &b, double &db, const std::vector<double> hitsx, const std::vector<double> hitsy) {
00225 
00226   double s=0,sx=0,sy=0,x,y;
00227   double sxx=0,sxy=0;
00228 
00229   a=b=0;
00230   if (hitsx.size()==0) return;
00231   if (hitsx.size()==1) {
00232     b=hitsy[0];
00233   } else {
00234     for (unsigned int i = 0; i != hitsx.size(); i++) {
00235       x=hitsx[i];
00236       y=hitsy[i];
00237       sy += y;
00238       sxy+= x*y;
00239       s += 1.;
00240       sx += x;
00241       sxx += x*x;
00242     }
00243 
00244     double d = s*sxx - sx*sx;
00245     b = (sxx*sy- sx*sxy)/ d;
00246     a = (s*sxy - sx*sy) / d;
00247     da = sqrt(sxx/d);
00248     db = sqrt(s/d);
00249   }
00250 }
00251