CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/DQM/TrackerMonitorTrack/src/MonitorTrackResiduals.cc

Go to the documentation of this file.
00001 #include "FWCore/Framework/interface/ESHandle.h"
00002 #include "FWCore/Framework/interface/Event.h"
00003 #include "FWCore/ServiceRegistry/interface/Service.h"
00004 #include "CalibTracker/Records/interface/SiStripDetCablingRcd.h"
00005 #include "CalibFormats/SiStripObjects/interface/SiStripDetCabling.h"
00006 #include "DataFormats/SiStripDetId/interface/SiStripSubStructure.h"
00007 #include "DQM/SiStripCommon/interface/SiStripFolderOrganizer.h"
00008 #include "DQM/SiStripCommon/interface/SiStripHistoId.h"
00009 #include "DQM/TrackerMonitorTrack/interface/MonitorTrackResiduals.h"
00010 #include "Geometry/CommonTopologies/interface/StripTopology.h"
00011 #include "RecoTracker/TransientTrackingRecHit/interface/TkTransientTrackingRecHitBuilder.h"
00012 #include "TrackingTools/TrackFitters/interface/TrajectoryFitter.h"
00013 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00014 #include "DataFormats/GeometryCommonDetAlgo/interface/MeasurementVector.h"
00015 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00016 #include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
00017 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00018 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00019 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00020 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00021 #include "TrackingTools/TrackFitters/interface/TrajectoryStateCombiner.h"
00022 #include "Alignment/OfflineValidation/interface/TrackerValidationVariables.h"
00023 #include "Alignment/TrackerAlignment/interface/TrackerAlignableId.h"
00024 #include "DQMServices/Core/interface/DQMStore.h"
00025 #include "DQMServices/Core/interface/MonitorElement.h"
00026 #include "CommonTools/TriggerUtils/interface/GenericTriggerEventFlag.h"
00027 
00028 MonitorTrackResiduals::MonitorTrackResiduals(const edm::ParameterSet& iConfig) 
00029    : dqmStore_( edm::Service<DQMStore>().operator->() )
00030    , conf_(iConfig), m_cacheID_(0) 
00031    , genTriggerEventFlag_(new GenericTriggerEventFlag(iConfig)) {
00032 }
00033 
00034 MonitorTrackResiduals::~MonitorTrackResiduals() {
00035   if (genTriggerEventFlag_) delete genTriggerEventFlag_; 
00036 }
00037 
00038 
00039 void MonitorTrackResiduals::beginJob(void) {
00040   ModOn = conf_.getParameter<bool>("Mod_On");
00041   reset_me_after_each_run = conf_.getParameter<bool>("ResetAfterRun");
00042 }
00043 
00044 void MonitorTrackResiduals::beginRun(edm::Run const& run, edm::EventSetup const& iSetup) {
00045   
00046   
00047   unsigned long long cacheID = iSetup.get<SiStripDetCablingRcd>().cacheIdentifier();
00048   if (m_cacheID_ != cacheID) {
00049     m_cacheID_ = cacheID;     
00050     this->createMEs(iSetup);
00051   }
00052   if(reset_me_after_each_run) {
00053     if(ModOn) {
00054       for(std::map<int32_t, MonitorElement*>::const_iterator it = HitResidual.begin(), 
00055             itEnd = HitResidual.end(); it!= itEnd;++it) {
00056         this->resetModuleMEs(it->first);
00057         this->resetLayerMEs(folder_organizer.GetSubDetAndLayer(it->first));
00058       }
00059     } else {
00060       for(std::map< std::pair<std::string,int32_t>, MonitorElement*>::const_iterator it = m_SubdetLayerResiduals.begin(), 
00061             itEnd = m_SubdetLayerResiduals.end(); it!= itEnd;++it) {
00062         this->resetLayerMEs(it->first);
00063       }
00064     } // end if-else Module level on
00065   } // end reset after run
00066 
00067   // Initialize the GenericTriggerEventFlag
00068   if ( genTriggerEventFlag_->on() ) genTriggerEventFlag_->initRun( run, iSetup );
00069 }
00070 
00071 void MonitorTrackResiduals::createMEs(const edm::EventSetup& iSetup){
00072   Parameters = conf_.getParameter<edm::ParameterSet>("TH1ResModules");
00073   int32_t i_residuals_Nbins =  Parameters.getParameter<int32_t>("Nbinx");
00074   double d_residual_xmin = Parameters.getParameter<double>("xmin");
00075   double d_residual_xmax = Parameters.getParameter<double>("xmax");
00076   Parameters = conf_.getParameter<edm::ParameterSet>("TH1NormResModules");
00077   int32_t i_normres_Nbins =  Parameters.getParameter<int32_t>("Nbinx");
00078   double d_normres_xmin = Parameters.getParameter<double>("xmin");
00079   double d_normres_xmax = Parameters.getParameter<double>("xmax");
00080 
00081   
00082   // use SistripHistoId for producing histogram id (and title)
00083   SiStripHistoId hidmanager;
00084   folder_organizer.setSiStripFolder(); // top SiStrip folder
00085   
00086   // take from eventSetup the SiStripDetCabling object
00087   edm::ESHandle<SiStripDetCabling> tkmechstruct;
00088   iSetup.get<SiStripDetCablingRcd>().get(tkmechstruct);
00089   
00090   // get list of active detectors from SiStripDetCabling
00091   std::vector<uint32_t> activeDets;
00092   activeDets.clear(); // just in case
00093   tkmechstruct->addActiveDetectorsRawIds(activeDets);
00094   
00095   // use SiStripSubStructure for selecting certain regions
00096   SiStripSubStructure substructure;
00097   std::vector<uint32_t> DetIds = activeDets;
00098   
00099   // book histo per each detector module
00100   for (std::vector<uint32_t>::const_iterator DetItr=activeDets.begin(),  
00101          DetItrEnd = activeDets.end(); DetItr!=DetItrEnd; ++DetItr)
00102     {
00103       uint ModuleID = (*DetItr);
00104       
00105       // is this a StripModule?
00106       if( SiStripDetId(ModuleID).subDetector() != 0 ) {
00107         
00108         folder_organizer.setDetectorFolder(ModuleID); //  detid sets appropriate detector folder                
00109         // Book module histogramms? 
00110         if (ModOn) { 
00111           std::string hid = hidmanager.createHistoId("HitResiduals","det",ModuleID);
00112           std::string normhid = hidmanager.createHistoId("NormalizedHitResiduals","det",ModuleID);      
00113           HitResidual[ModuleID] = dqmStore_->book1D(hid, hid, 
00114                                                     i_residuals_Nbins,d_residual_xmin,d_residual_xmax);
00115           HitResidual[ModuleID]->setAxisTitle("(x_{pred} - x_{rec})' [cm]");
00116           NormedHitResiduals[ModuleID] = dqmStore_->book1D(normhid, normhid, 
00117                                                            i_normres_Nbins,d_normres_xmin,d_normres_xmax);
00118           NormedHitResiduals[ModuleID]->setAxisTitle("(x_{pred} - x_{rec})'/#sigma");
00119         }
00120         // book layer level histogramms
00121         std::pair<std::string,int32_t> subdetandlayer = folder_organizer.GetSubDetAndLayer(ModuleID);
00122         folder_organizer.setLayerFolder(ModuleID,subdetandlayer.second);
00123         if(! m_SubdetLayerResiduals[subdetandlayer ] ) {
00124           // book histogramms on layer level, check for barrel for correct labeling
00125           std::string histoname = Form(subdetandlayer.first.find("B") != std::string::npos ? 
00126                                        "HitResiduals_%s__Layer__%d" : "HitResiduals_%s__wheel__%d" ,
00127                                        subdetandlayer.first.c_str(),std::abs(subdetandlayer.second));
00128           std::string normhistoname = 
00129             Form(subdetandlayer.first.find("B") != std::string::npos ? 
00130                  "NormalizedHitResidual_%s__Layer__%d" : "NormalizedHitResidual_%s__wheel__%d" ,
00131                  subdetandlayer.first.c_str(),std::abs(subdetandlayer.second));
00132           m_SubdetLayerResiduals[subdetandlayer] = 
00133             dqmStore_->book1D(histoname.c_str(),histoname.c_str(),
00134                               i_residuals_Nbins,d_residual_xmin,d_residual_xmax);
00135           m_SubdetLayerResiduals[subdetandlayer]->setAxisTitle("(x_{pred} - x_{rec})' [cm]");
00136           m_SubdetLayerNormedResiduals[subdetandlayer] = 
00137             dqmStore_->book1D(normhistoname.c_str(),normhistoname.c_str(),
00138                               i_normres_Nbins,d_normres_xmin,d_normres_xmax);
00139           m_SubdetLayerNormedResiduals[subdetandlayer]->setAxisTitle("(x_{pred} - x_{rec})'/#sigma");
00140         } 
00141       } // end 'is strip module'
00142     } // end loop over activeDets
00143 }
00144 
00145 
00146 void MonitorTrackResiduals::resetModuleMEs(int32_t modid) {
00147   HitResidual[modid]->Reset();
00148   NormedHitResiduals[modid]->Reset();
00149 }
00150 
00151 void MonitorTrackResiduals::resetLayerMEs(const std::pair<std::string, int32_t> &subdetandlayer) {
00152   m_SubdetLayerResiduals      [subdetandlayer]->Reset();
00153   m_SubdetLayerNormedResiduals[subdetandlayer]->Reset();
00154 }
00155 
00156 
00157 
00158 void MonitorTrackResiduals::endRun(const edm::Run&, const edm::EventSetup&){
00159 }
00160 
00161 
00162 void MonitorTrackResiduals::endJob(void){
00163 
00164   //dqmStore_->showDirStructure();
00165   bool outputMEsInRootFile = conf_.getParameter<bool>("OutputMEsInRootFile");
00166   std::string outputFileName = conf_.getParameter<std::string>("OutputFileName");
00167   if(outputMEsInRootFile){
00168     dqmStore_->save(outputFileName);
00169   }
00170 }
00171 
00172 
00173 void MonitorTrackResiduals::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
00174 
00175   // Filter out events if Trigger Filtering is requested                                                                                           
00176   if (genTriggerEventFlag_->on()&& ! genTriggerEventFlag_->accept( iEvent, iSetup) ) return;
00177 
00178   TrackerValidationVariables avalidator_(iSetup,conf_);
00179   std::vector<TrackerValidationVariables::AVHitStruct> v_hitstruct;
00180   avalidator_.fillHitQuantities(iEvent,v_hitstruct);
00181   for (std::vector<TrackerValidationVariables::AVHitStruct>::const_iterator it = v_hitstruct.begin(),
00182        itEnd = v_hitstruct.end(); it != itEnd; ++it) {
00183     uint RawId = it->rawDetId;
00184     
00185     // fill if hit belongs to StripDetector and its error is not zero
00186     if( it->resXprimeErr != 0 && SiStripDetId(RawId).subDetector()  != 0 )  {
00187       if (ModOn && HitResidual[RawId]) { 
00188         HitResidual[RawId]->Fill(it->resXprime);
00189         NormedHitResiduals[RawId]->Fill(it->resXprime/it->resXprimeErr);
00190       }
00191       std::pair<std::string, int32_t> subdetandlayer = folder_organizer.GetSubDetAndLayer(RawId);
00192       if(m_SubdetLayerResiduals[subdetandlayer]) {
00193         m_SubdetLayerResiduals[subdetandlayer]->Fill(it->resXprime);
00194         m_SubdetLayerNormedResiduals[subdetandlayer]->Fill(it->resXprime/it->resXprimeErr);
00195       }
00196     }
00197   }
00198 
00199 }
00200 
00201 
00202 
00203 DEFINE_FWK_MODULE(MonitorTrackResiduals);
00204