CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_2_9_HLT1_bphpatch4/src/DQMOffline/JetMET/src/METAnalyzer.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2011/08/12 15:29:52 $
00005  *  $Revision: 1.38.6.4 $
00006  *  \author A.Apresyan - Caltech
00007  *          K.Hatakeyama - Baylor
00008  */
00009 
00010 #include "DQMOffline/JetMET/interface/METAnalyzer.h"
00011 #include "DataFormats/Common/interface/Handle.h"
00012 #include "FWCore/Common/interface/TriggerNames.h"
00013 
00014 #include "DataFormats/Math/interface/LorentzVector.h"
00015 
00016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00017 
00018 #include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
00019 #include "DataFormats/CaloTowers/interface/CaloTowerDetId.h"
00020 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
00021 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00022 
00023 #include "DataFormats/Math/interface/LorentzVector.h"
00024 
00025 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutSetupFwd.h"
00026 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutSetup.h"
00027 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
00028 
00029 #include <string>
00030 using namespace edm;
00031 using namespace reco;
00032 using namespace math;
00033 
00034 // ***********************************************************
00035 METAnalyzer::METAnalyzer(const edm::ParameterSet& pSet) {
00036 
00037   parameters = pSet;
00038 
00039   edm::ParameterSet highptjetparms = parameters.getParameter<edm::ParameterSet>("highPtJetTrigger");
00040   edm::ParameterSet lowptjetparms  = parameters.getParameter<edm::ParameterSet>("lowPtJetTrigger" );
00041   edm::ParameterSet minbiasparms   = parameters.getParameter<edm::ParameterSet>("minBiasTrigger"  );
00042   edm::ParameterSet highmetparms   = parameters.getParameter<edm::ParameterSet>("highMETTrigger"  );
00043   edm::ParameterSet lowmetparms    = parameters.getParameter<edm::ParameterSet>("lowMETTrigger"   );
00044   edm::ParameterSet eleparms       = parameters.getParameter<edm::ParameterSet>("eleTrigger"      );
00045   edm::ParameterSet muonparms      = parameters.getParameter<edm::ParameterSet>("muonTrigger"     );
00046 
00047   _HighPtJetEventFlag = new GenericTriggerEventFlag( highptjetparms );
00048   _LowPtJetEventFlag  = new GenericTriggerEventFlag( lowptjetparms  );
00049   _MinBiasEventFlag   = new GenericTriggerEventFlag( minbiasparms   );
00050   _HighMETEventFlag   = new GenericTriggerEventFlag( highmetparms   );
00051   _LowMETEventFlag    = new GenericTriggerEventFlag( lowmetparms    );
00052   _EleEventFlag       = new GenericTriggerEventFlag( eleparms       );
00053   _MuonEventFlag      = new GenericTriggerEventFlag( muonparms      );
00054 
00055 }
00056 
00057 // ***********************************************************
00058 METAnalyzer::~METAnalyzer() { 
00059 
00060   delete _HighPtJetEventFlag;
00061   delete _LowPtJetEventFlag;
00062   delete _MinBiasEventFlag;
00063   delete _HighMETEventFlag;
00064   delete _LowMETEventFlag;
00065   delete _EleEventFlag;
00066   delete _MuonEventFlag;
00067 
00068 }
00069 
00070 void METAnalyzer::beginJob(DQMStore * dbe) {
00071 
00072   evtCounter = 0;
00073   metname = "METAnalyzer";
00074 
00075   // trigger information
00076   HLTPathsJetMBByName_ = parameters.getParameter<std::vector<std::string > >("HLTPathsJetMB");
00077 
00078   theCleaningParameters = parameters.getParameter<ParameterSet>("CleaningParameters"),
00079 
00080   //Trigger parameters
00081   gtTag          = theCleaningParameters.getParameter<edm::InputTag>("gtLabel");
00082   _techTrigsAND  = theCleaningParameters.getParameter<std::vector<unsigned > >("techTrigsAND");
00083   _techTrigsOR   = theCleaningParameters.getParameter<std::vector<unsigned > >("techTrigsOR");
00084   _techTrigsNOT  = theCleaningParameters.getParameter<std::vector<unsigned > >("techTrigsNOT");
00085 
00086   _doHLTPhysicsOn = theCleaningParameters.getParameter<bool>("doHLTPhysicsOn");
00087   _hlt_PhysDec    = theCleaningParameters.getParameter<std::string>("HLT_PhysDec");
00088 
00089   _tightBHFiltering     = theCleaningParameters.getParameter<bool>("tightBHFiltering");
00090   _tightJetIDFiltering  = theCleaningParameters.getParameter<int>("tightJetIDFiltering");
00091   _tightHcalFiltering   = theCleaningParameters.getParameter<bool>("tightHcalFiltering");
00092 
00093   // ==========================================================
00094   //DCS information
00095   // ==========================================================
00096   DCSFilter = new JetMETDQMDCSFilter(parameters.getParameter<ParameterSet>("DCSFilter"));
00097 
00098   //Vertex requirements
00099   _doPVCheck          = theCleaningParameters.getParameter<bool>("doPrimaryVertexCheck");
00100   vertexTag  = theCleaningParameters.getParameter<edm::InputTag>("vertexLabel");
00101 
00102   if (_doPVCheck) {
00103     _nvtx_min        = theCleaningParameters.getParameter<int>("nvtx_min");
00104     _nvtxtrks_min    = theCleaningParameters.getParameter<int>("nvtxtrks_min");
00105     _vtxndof_min     = theCleaningParameters.getParameter<int>("vtxndof_min");
00106     _vtxchi2_max     = theCleaningParameters.getParameter<double>("vtxchi2_max");
00107     _vtxz_max        = theCleaningParameters.getParameter<double>("vtxz_max");
00108   }
00109 
00110   // MET information
00111   theMETCollectionLabel       = parameters.getParameter<edm::InputTag>("METCollectionLabel");
00112   _source                     = parameters.getParameter<std::string>("Source");
00113 
00114   if (theMETCollectionLabel.label() == "tcMet" ) {
00115     inputTrackLabel         = parameters.getParameter<edm::InputTag>("InputTrackLabel");    
00116     inputMuonLabel          = parameters.getParameter<edm::InputTag>("InputMuonLabel");
00117     inputElectronLabel      = parameters.getParameter<edm::InputTag>("InputElectronLabel");
00118     inputBeamSpotLabel      = parameters.getParameter<edm::InputTag>("InputBeamSpotLabel");
00119   }
00120 
00121   // Other data collections
00122   theJetCollectionLabel       = parameters.getParameter<edm::InputTag>("JetCollectionLabel");
00123   HcalNoiseRBXCollectionTag   = parameters.getParameter<edm::InputTag>("HcalNoiseRBXCollection");
00124   HcalNoiseSummaryTag         = parameters.getParameter<edm::InputTag>("HcalNoiseSummary");
00125   BeamHaloSummaryTag          = parameters.getParameter<edm::InputTag>("BeamHaloSummaryLabel");
00126 
00127   // misc
00128   _verbose      = parameters.getParameter<int>("verbose");
00129   _etThreshold  = parameters.getParameter<double>("etThreshold"); // MET threshold
00130   _allhist      = parameters.getParameter<bool>("allHist");       // Full set of monitoring histograms
00131   _allSelection = parameters.getParameter<bool>("allSelection");  // Plot with all sets of event selection
00132   _cleanupSelection = parameters.getParameter<bool>("cleanupSelection");  // Plot with all sets of event selection
00133 
00134   _FolderName              = parameters.getUntrackedParameter<std::string>("FolderName");
00135 
00136   _highPtJetThreshold = parameters.getParameter<double>("HighPtJetThreshold"); // High Pt Jet threshold
00137   _lowPtJetThreshold  = parameters.getParameter<double>("LowPtJetThreshold");   // Low Pt Jet threshold
00138   _highMETThreshold   = parameters.getParameter<double>("HighMETThreshold");     // High MET threshold
00139   _lowMETThreshold    = parameters.getParameter<double>("LowMETThreshold");       // Low MET threshold
00140 
00141   //
00142   jetID = new reco::helper::JetIDHelper(parameters.getParameter<ParameterSet>("JetIDParams"));
00143 
00144   // DQStore stuff
00145   LogTrace(metname)<<"[METAnalyzer] Parameters initialization";
00146   std::string DirName = _FolderName+_source;
00147   dbe->setCurrentFolder(DirName);
00148 
00149   hmetME = dbe->book1D("metReco", "metReco", 4, 1, 5);
00150   hmetME->setBinLabel(2,"MET",1);
00151 
00152   _dbe = dbe;
00153 
00154   _FolderNames.push_back("All");
00155   _FolderNames.push_back("BasicCleanup");
00156   _FolderNames.push_back("ExtraCleanup");
00157   _FolderNames.push_back("HcalNoiseFilter");
00158   _FolderNames.push_back("HcalNoiseFilterTight");
00159   _FolderNames.push_back("JetIDMinimal");
00160   _FolderNames.push_back("JetIDLoose");
00161   _FolderNames.push_back("JetIDTight");
00162   _FolderNames.push_back("BeamHaloIDTightPass");
00163   _FolderNames.push_back("BeamHaloIDLoosePass");
00164   _FolderNames.push_back("Triggers");
00165   _FolderNames.push_back("PV");
00166 
00167   for (std::vector<std::string>::const_iterator ic = _FolderNames.begin(); 
00168        ic != _FolderNames.end(); ic++){
00169     if (*ic=="All")                  bookMESet(DirName+"/"+*ic);
00170     if (_cleanupSelection){
00171     if (*ic=="BasicCleanup")         bookMESet(DirName+"/"+*ic);
00172     if (*ic=="ExtraCleanup")         bookMESet(DirName+"/"+*ic);
00173     }
00174     if (_allSelection){
00175       if (*ic=="HcalNoiseFilter")      bookMESet(DirName+"/"+*ic);
00176       if (*ic=="HcalNoiseFilterTight") bookMESet(DirName+"/"+*ic);
00177       if (*ic=="JetIDMinimal")         bookMESet(DirName+"/"+*ic);
00178       if (*ic=="JetIDLoose")           bookMESet(DirName+"/"+*ic);
00179       if (*ic=="JetIDTight")           bookMESet(DirName+"/"+*ic);
00180       if (*ic=="BeamHaloIDTightPass")  bookMESet(DirName+"/"+*ic);
00181       if (*ic=="BeamHaloIDLoosePass")  bookMESet(DirName+"/"+*ic);
00182       if (*ic=="Triggers")             bookMESet(DirName+"/"+*ic);
00183       if (*ic=="PV")                   bookMESet(DirName+"/"+*ic);
00184     }
00185   }
00186 }
00187 
00188 // ***********************************************************
00189 void METAnalyzer::endJob() {
00190 
00191   delete jetID;
00192   delete DCSFilter;
00193 
00194 }
00195 
00196 // ***********************************************************
00197 void METAnalyzer::bookMESet(std::string DirName)
00198 {
00199 
00200   bool bLumiSecPlot=false;
00201   if (DirName.find("All")!=std::string::npos) bLumiSecPlot=true;
00202 
00203   bookMonitorElement(DirName,bLumiSecPlot);
00204 
00205   if ( _HighPtJetEventFlag->on() ) {
00206     bookMonitorElement(DirName+"/"+"HighPtJet",false);
00207     hTriggerName_HighPtJet = _dbe->bookString("triggerName_HighPtJet", _hlt_HighPtJet);
00208   }  
00209 
00210   if ( _LowPtJetEventFlag->on() ) {
00211     bookMonitorElement(DirName+"/"+"LowPtJet",false);
00212     hTriggerName_LowPtJet = _dbe->bookString("triggerName_LowPtJet", _hlt_LowPtJet);
00213   }
00214 
00215   if ( _MinBiasEventFlag->on() ) {
00216     bookMonitorElement(DirName+"/"+"MinBias",false);
00217     hTriggerName_MinBias = _dbe->bookString("triggerName_MinBias", _hlt_MinBias);
00218   }
00219 
00220   if ( _HighMETEventFlag->on() ) {
00221     bookMonitorElement(DirName+"/"+"HighMET",false);
00222     hTriggerName_HighMET = _dbe->bookString("triggerName_HighMET", _hlt_HighMET);
00223   }
00224 
00225   if ( _LowMETEventFlag->on() ) {
00226     bookMonitorElement(DirName+"/"+"LowMET",false);
00227     hTriggerName_LowMET = _dbe->bookString("triggerName_LowMET", _hlt_LowMET);
00228   }
00229 
00230   if ( _EleEventFlag->on() ) {
00231     bookMonitorElement(DirName+"/"+"Ele",false);
00232     hTriggerName_Ele = _dbe->bookString("triggerName_Ele", _hlt_Ele);
00233   }
00234 
00235   if ( _MuonEventFlag->on() ) {
00236     bookMonitorElement(DirName+"/"+"Muon",false);
00237     hTriggerName_Muon = _dbe->bookString("triggerName_Muon", _hlt_Muon);
00238   }
00239 }
00240 
00241 // ***********************************************************
00242 void METAnalyzer::bookMonitorElement(std::string DirName, bool bLumiSecPlot=false)
00243 {
00244 
00245   if (_verbose) std::cout << "bookMonitorElement " << DirName << std::endl;
00246   _dbe->setCurrentFolder(DirName);
00247  
00248   //removed for optimizations//hNevents            = _dbe->book1D("METTask_Nevents", "METTask_Nevents"   ,1,0,1);
00249   hMEx                = _dbe->book1D("METTask_MEx",   "METTask_MEx"   ,200,-500,500);
00250   hMEx->setAxisTitle("MEx [GeV]",1);
00251   hMEy                = _dbe->book1D("METTask_MEy",   "METTask_MEy"   ,200,-500,500);
00252   hMEy->setAxisTitle("MEy [GeV]",1);
00253   //removed for optimizations//hEz                 = _dbe->book1D("METTask_Ez",    "METTask_Ez"    ,500,-500,500);
00254   //removed for optimizations//hEz->setAxisTitle("MEz [GeV]",1);
00255   hMETSig             = _dbe->book1D("METTask_METSig","METTask_METSig",51,0,51);
00256   hMETSig->setAxisTitle("CaloMETSig",1);
00257   hMET                = _dbe->book1D("METTask_MET",   "METTask_MET"   ,200,0,1000);
00258   hMET->setAxisTitle("MET [GeV]",1);
00259   hMETPhi             = _dbe->book1D("METTask_METPhi","METTask_METPhi",60,-TMath::Pi(),TMath::Pi());
00260   hMETPhi->setAxisTitle("METPhi [rad]",1);
00261   hSumET              = _dbe->book1D("METTask_SumET", "METTask_SumET" ,400,0,2000);
00262   hSumET->setAxisTitle("SumET [GeV]",1);
00263 
00264   hMET_logx           = _dbe->book1D("METTask_MET_logx",   "METTask_MET_logx"   ,40,-1.,7.);
00265   hMET_logx->setAxisTitle("log(MET) [GeV]",1);
00266   hSumET_logx         = _dbe->book1D("METTask_SumET_logx", "METTask_SumET_logx" ,40,-1.,7.);
00267   hSumET_logx->setAxisTitle("log(SumET) [GeV]",1);
00268 
00269   //removed for optimizations//hMETIonFeedbck      = _dbe->book1D("METTask_METIonFeedbck", "METTask_METIonFeedbck" ,500,0,1000);
00270   //removed for optimizations//hMETIonFeedbck->setAxisTitle("MET [GeV]",1);
00271   //removed for optimizations//hMETHPDNoise        = _dbe->book1D("METTask_METHPDNoise",   "METTask_METHPDNoise"   ,500,0,1000);
00272   //removed for optimizations//hMETHPDNoise->setAxisTitle("MET [GeV]",1);
00273   //removed for optimizations//hMETRBXNoise        = _dbe->book1D("METTask_METRBXNoise",   "METTask_METRBXNoise"   ,500,0,1000);
00274   //removed for optimizations//hMETRBXNoise->setAxisTitle("MET [GeV]",1);
00275 
00276   if (_allhist){
00277     if (bLumiSecPlot){
00278       hMExLS = _dbe->book2D("METTask_MEx_LS","METTask_MEx_LS",200,-200,200,50,0.,500.);
00279       hMExLS->setAxisTitle("MEx [GeV]",1);
00280       hMExLS->setAxisTitle("Lumi Section",2);
00281       hMEyLS = _dbe->book2D("METTask_MEy_LS","METTask_MEy_LS",200,-200,200,50,0.,500.);
00282       hMEyLS->setAxisTitle("MEy [GeV]",1);
00283       hMEyLS->setAxisTitle("Lumi Section",2);
00284     }
00285   }
00286 
00287   if (theMETCollectionLabel.label() == "tcMet" ) {
00288     htrkPt    = _dbe->book1D("METTask_trackPt", "METTask_trackPt", 50, 0, 500);
00289     htrkEta   = _dbe->book1D("METTask_trackEta", "METTask_trackEta", 60, -3.0, 3.0);
00290     htrkNhits = _dbe->book1D("METTask_trackNhits", "METTask_trackNhits", 50, 0, 50);
00291     htrkChi2  = _dbe->book1D("METTask_trackNormalizedChi2", "METTask_trackNormalizedChi2", 20, 0, 20);
00292     htrkD0    = _dbe->book1D("METTask_trackD0", "METTask_trackd0", 50, -1, 1);
00293     helePt    = _dbe->book1D("METTask_electronPt", "METTask_electronPt", 50, 0, 500);
00294     heleEta   = _dbe->book1D("METTask_electronEta", "METTask_electronEta", 60, -3.0, 3.0);
00295     heleHoE   = _dbe->book1D("METTask_electronHoverE", "METTask_electronHoverE", 25, 0, 0.5);
00296     hmuPt     = _dbe->book1D("METTask_muonPt", "METTask_muonPt", 50, 0, 500);
00297     hmuEta    = _dbe->book1D("METTask_muonEta", "METTask_muonEta", 60, -3.0, 3.0);
00298     hmuNhits  = _dbe->book1D("METTask_muonNhits", "METTask_muonNhits", 50, 0, 50);
00299     hmuChi2   = _dbe->book1D("METTask_muonNormalizedChi2", "METTask_muonNormalizedChi2", 20, 0, 20);
00300     hmuD0     = _dbe->book1D("METTask_muonD0", "METTask_muonD0", 50, -1, 1);
00301   }
00302 
00303   hMExCorrection       = _dbe->book1D("METTask_MExCorrection", "METTask_MExCorrection", 100, -500.0,500.0);
00304   hMEyCorrection       = _dbe->book1D("METTask_MEyCorrection", "METTask_MEyCorrection", 100, -500.0,500.0);
00305   hMuonCorrectionFlag  = _dbe->book1D("METTask_CorrectionFlag","METTask_CorrectionFlag", 5, -0.5, 4.5);
00306 
00307 }
00308 
00309 // ***********************************************************
00310 void METAnalyzer::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup)
00311 {
00312   if ( _HighPtJetEventFlag->on() ) _HighPtJetEventFlag->initRun( iRun, iSetup );
00313   if ( _LowPtJetEventFlag ->on() ) _LowPtJetEventFlag ->initRun( iRun, iSetup );
00314   if ( _MinBiasEventFlag  ->on() ) _MinBiasEventFlag  ->initRun( iRun, iSetup );
00315   if ( _HighMETEventFlag  ->on() ) _HighMETEventFlag  ->initRun( iRun, iSetup );
00316   if ( _LowMETEventFlag   ->on() ) _LowMETEventFlag   ->initRun( iRun, iSetup );
00317   if ( _EleEventFlag      ->on() ) _EleEventFlag      ->initRun( iRun, iSetup );
00318   if ( _MuonEventFlag     ->on() ) _MuonEventFlag     ->initRun( iRun, iSetup );
00319 
00320 }
00321 
00322 // ***********************************************************
00323 void METAnalyzer::endRun(const edm::Run& iRun, const edm::EventSetup& iSetup, DQMStore * dbe)
00324 {
00325   
00326   //
00327   //--- Check the time length of the Run from the lumi section plots
00328 
00329   std::string dirName = _FolderName+_source+"/";
00330   _dbe->setCurrentFolder(dirName);
00331 
00332   TH1F* tlumisec;
00333 
00334   MonitorElement *meLumiSec = _dbe->get("aaa");
00335   meLumiSec = _dbe->get("JetMET/lumisec");
00336 
00337   int totlsec=0;
00338   double totltime=0.;
00339   if ( meLumiSec->getRootObject() ) {
00340     tlumisec = meLumiSec->getTH1F();
00341     for (int i=0; i<500; i++){
00342       if (tlumisec->GetBinContent(i+1)) totlsec++;
00343     }
00344     totltime = double(totlsec*90); // one lumi sec ~ 90 (sec)
00345   }
00346 
00347   if (totltime==0.) totltime=1.; 
00348 
00349   //
00350   //--- Make the integrated plots with rate (Hz)
00351 
00352   for (std::vector<std::string>::const_iterator ic = _FolderNames.begin(); ic != _FolderNames.end(); ic++)
00353     {
00354 
00355       std::string DirName;
00356       DirName = dirName+*ic;
00357 
00358       makeRatePlot(DirName,totltime);
00359       if ( _HighPtJetEventFlag->on() ) 
00360         makeRatePlot(DirName+"/"+"triggerName_HighJetPt",totltime);
00361       if ( _LowPtJetEventFlag->on() ) 
00362         makeRatePlot(DirName+"/"+"triggerName_LowJetPt",totltime);
00363       if ( _MinBiasEventFlag->on() ) 
00364         makeRatePlot(DirName+"/"+"triggerName_MinBias",totltime);
00365       if ( _HighMETEventFlag->on() ) 
00366         makeRatePlot(DirName+"/"+"triggerName_HighMET",totltime);
00367       if ( _LowMETEventFlag->on() ) 
00368         makeRatePlot(DirName+"/"+"triggerName_LowMET",totltime);
00369       if ( _EleEventFlag->on() ) 
00370         makeRatePlot(DirName+"/"+"triggerName_Ele",totltime);
00371       if ( _MuonEventFlag->on() ) 
00372         makeRatePlot(DirName+"/"+"triggerName_Muon",totltime);
00373     }
00374 }
00375 
00376 // ***********************************************************
00377 void METAnalyzer::makeRatePlot(std::string DirName, double totltime)
00378 {
00379 
00380   _dbe->setCurrentFolder(DirName);
00381   MonitorElement *meMET = _dbe->get(DirName+"/"+"METTask_MET");
00382 
00383   TH1F* tMET;
00384   TH1F* tMETRate;
00385 
00386   if ( meMET )
00387     if ( meMET->getRootObject() ) {
00388       tMET     = meMET->getTH1F();
00389       
00390       // Integral plot & convert number of events to rate (hz)
00391       tMETRate = (TH1F*) tMET->Clone("METTask_METRate");
00392       for (int i = tMETRate->GetNbinsX()-1; i>=0; i--){
00393         tMETRate->SetBinContent(i+1,tMETRate->GetBinContent(i+2)+tMET->GetBinContent(i+1));
00394       }
00395       for (int i = 0; i<tMETRate->GetNbinsX(); i++){
00396         tMETRate->SetBinContent(i+1,tMETRate->GetBinContent(i+1)/double(totltime));
00397       }      
00398 
00399       tMETRate->SetName("METTask_METRate");
00400       tMETRate->SetTitle("METTask_METRate");
00401       hMETRate      = _dbe->book1D("METTask_METRate",tMETRate);
00402     }
00403 }
00404 
00405 // ***********************************************************
00406 void METAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup, 
00407                             const edm::TriggerResults& triggerResults) {
00408 
00409   if (_verbose) std::cout << "METAnalyzer analyze" << std::endl;
00410 
00411   std::string DirName = _FolderName+_source;
00412 
00413   LogTrace(metname)<<"[METAnalyzer] Analyze MET";
00414 
00415   hmetME->Fill(2);
00416 
00417   // ==========================================================  
00418   // Trigger information 
00419   //
00420   _trig_JetMB=0;
00421   _trig_HighPtJet=0;
00422   _trig_LowPtJet=0;
00423   _trig_MinBias=0;
00424   _trig_HighMET=0;
00425   _trig_LowMET=0;
00426   _trig_Ele=0;
00427   _trig_Muon=0;
00428   _trig_PhysDec=0;
00429   if(&triggerResults) {   
00430 
00432 
00433     //
00434     //
00435     // Check how many HLT triggers are in triggerResults 
00436     int ntrigs = triggerResults.size();
00437     if (_verbose) std::cout << "ntrigs=" << ntrigs << std::endl;
00438     
00439     //
00440     //
00441     // If index=ntrigs, this HLT trigger doesn't exist in the HLT table for this data.
00442     const edm::TriggerNames & triggerNames = iEvent.triggerNames(triggerResults);
00443     
00444     //
00445     //
00446     // count number of requested Jet or MB HLT paths which have fired
00447     for (unsigned int i=0; i!=HLTPathsJetMBByName_.size(); i++) {
00448       unsigned int triggerIndex = triggerNames.triggerIndex(HLTPathsJetMBByName_[i]);
00449       if (triggerIndex<triggerResults.size()) {
00450         if (triggerResults.accept(triggerIndex)) {
00451           _trig_JetMB++;
00452         }
00453       }
00454     }
00455     // for empty input vectors (n==0), take all HLT triggers!
00456     if (HLTPathsJetMBByName_.size()==0) _trig_JetMB=triggerResults.size()-1;
00457 
00458     if ( _HighPtJetEventFlag->on() && _HighPtJetEventFlag->accept( iEvent, iSetup ) ) 
00459       _trig_HighPtJet=1;
00460 
00461     if ( _LowPtJetEventFlag->on() && _LowPtJetEventFlag->accept( iEvent, iSetup ) ) 
00462       _trig_LowPtJet=1;
00463     
00464     if ( _MinBiasEventFlag->on() && _MinBiasEventFlag->accept( iEvent, iSetup ) ) 
00465       _trig_MinBias=1;
00466     
00467     if ( _HighMETEventFlag->on() && _HighMETEventFlag->accept( iEvent, iSetup ) ) 
00468       _trig_HighMET=1;
00469     
00470     if ( _LowMETEventFlag->on() && _LowMETEventFlag->accept( iEvent, iSetup ) ) 
00471       _trig_LowMET=1;
00472     
00473     if ( _EleEventFlag->on() && _EleEventFlag->accept( iEvent, iSetup ) ) 
00474       _trig_Ele=1;
00475     
00476     if ( _MuonEventFlag->on() && _MuonEventFlag->accept( iEvent, iSetup ) ) 
00477       _trig_Muon=1;
00478     
00479     if (triggerNames.triggerIndex(_hlt_PhysDec)   != triggerNames.size() &&
00480         triggerResults.accept(triggerNames.triggerIndex(_hlt_PhysDec)))   _trig_PhysDec=1;
00481     
00482   } else {
00483 
00484     edm::LogInfo("MetAnalyzer") << "TriggerResults::HLT not found, "
00485       "automatically select events"; 
00486 
00487     // TriggerResults object not found. Look at all events.    
00488     _trig_JetMB=1;
00489   }
00490 
00491   // ==========================================================
00492   // MET information
00493   
00494   // **** Get the MET container  
00495   edm::Handle<reco::METCollection> metcoll;
00496   iEvent.getByLabel(theMETCollectionLabel, metcoll);
00497   
00498   if(!metcoll.isValid()) {
00499     std::cout<<"Unable to find MET results for MET collection "<<theMETCollectionLabel<<std::endl;
00500     return;
00501   }
00502 
00503   const METCollection *metcol = metcoll.product();
00504   const MET *met;
00505   met = &(metcol->front());
00506     
00507   LogTrace(metname)<<"[METAnalyzer] Call to the MET analyzer";
00508 
00509   // ==========================================================
00510   // TCMET 
00511 
00512   if (theMETCollectionLabel.label() == "tcMet" ) {
00513 
00514     iEvent.getByLabel(inputMuonLabel, muon_h);
00515     iEvent.getByLabel(inputTrackLabel, track_h);
00516     iEvent.getByLabel(inputElectronLabel, electron_h);
00517     iEvent.getByLabel(inputBeamSpotLabel, beamSpot_h);
00518     iEvent.getByLabel("muonTCMETValueMapProducer" , "muCorrData", tcMet_ValueMap_Handle);
00519 
00520     if(!muon_h.isValid())     edm::LogInfo("OutputInfo") << "falied to retrieve muon data require by MET Task";
00521     if(!track_h.isValid())    edm::LogInfo("OutputInfo") << "falied to retrieve track data require by MET Task";
00522     if(!electron_h.isValid()) edm::LogInfo("OutputInfo") << "falied to retrieve electron data require by MET Task";
00523     if(!beamSpot_h.isValid()) edm::LogInfo("OutputInfo") << "falied to retrieve beam spot data require by MET Task";
00524 
00525     bspot = ( beamSpot_h.isValid() ) ? beamSpot_h->position() : math::XYZPoint(0, 0, 0);
00526     
00527   }
00528 
00529   // ==========================================================
00530   // HcalNoiseSummary
00531   //
00532 
00533   edm::Handle<HcalNoiseRBXCollection> HRBXCollection;
00534   iEvent.getByLabel(HcalNoiseRBXCollectionTag,HRBXCollection);
00535   if (!HRBXCollection.isValid()) {
00536     LogDebug("") << "METAnalyzer: Could not find HcalNoiseRBX Collection" << std::endl;
00537     if (_verbose) std::cout << "METAnalyzer: Could not find HcalNoiseRBX Collection" << std::endl;
00538   }
00539   
00540   edm::Handle<HcalNoiseSummary> HNoiseSummary;
00541   iEvent.getByLabel(HcalNoiseSummaryTag,HNoiseSummary);
00542   if (!HNoiseSummary.isValid()) {
00543     LogDebug("") << "METAnalyzer: Could not find Hcal NoiseSummary product" << std::endl;
00544     if (_verbose) std::cout << "METAnalyzer: Could not find Hcal NoiseSummary product" << std::endl;
00545   }
00546 
00547   edm::Handle<reco::CaloJetCollection> caloJets;
00548   iEvent.getByLabel(theJetCollectionLabel, caloJets);
00549   if (!caloJets.isValid()) {
00550     LogDebug("") << "METAnalyzer: Could not find jet product" << std::endl;
00551     if (_verbose) std::cout << "METAnalyzer: Could not find jet product" << std::endl;
00552   }
00553 
00554   // ==========================================================
00555   // MET sanity check
00556 
00557   //   if (_source=="MET") validateMET(*met, tcCandidates);
00558   
00559   // ==========================================================
00560   // JetID 
00561 
00562   if (_verbose) std::cout << "JetID starts" << std::endl;
00563   
00564   //
00565   // --- Minimal cuts
00566   //
00567   bool bJetIDMinimal=true;
00568   for (reco::CaloJetCollection::const_iterator cal = caloJets->begin(); 
00569        cal!=caloJets->end(); ++cal){
00570     jetID->calculate(iEvent, *cal);
00571     if (cal->pt()>10.){
00572       if (fabs(cal->eta())<=2.6 && 
00573           cal->emEnergyFraction()<=0.01) bJetIDMinimal=false;
00574     }
00575   }
00576 
00577   //
00578   // --- Loose cuts, not  specific for now!
00579   //
00580   bool bJetIDLoose=true;
00581   for (reco::CaloJetCollection::const_iterator cal = caloJets->begin(); 
00582        cal!=caloJets->end(); ++cal){ 
00583     jetID->calculate(iEvent, *cal);
00584     if (_verbose) std::cout << jetID->n90Hits() << " " 
00585                             << jetID->restrictedEMF() << " "
00586                             << cal->pt() << std::endl;
00587     if (cal->pt()>10.){
00588       //
00589       // for all regions
00590       if (jetID->n90Hits()<2)  bJetIDLoose=false; 
00591       if (jetID->fHPD()>=0.98) bJetIDLoose=false; 
00592       //if (jetID->restrictedEMF()<0.01) bJetIDLoose=false; 
00593       //
00594       // for non-forward
00595       if (fabs(cal->eta())<2.55){
00596         if (cal->emEnergyFraction()<=0.01) bJetIDLoose=false; 
00597       }
00598       // for forward
00599       else {
00600         if (cal->emEnergyFraction()<=-0.9) bJetIDLoose=false; 
00601         if (cal->pt()>80.){
00602         if (cal->emEnergyFraction()>= 1.0) bJetIDLoose=false; 
00603         }
00604       } // forward vs non-forward
00605     }   // pt>10 GeV/c
00606   }     // calor-jets loop
00607  
00608   //
00609   // --- Tight cuts
00610   //
00611   bool bJetIDTight=true;
00612   bJetIDTight=bJetIDLoose;
00613   for (reco::CaloJetCollection::const_iterator cal = caloJets->begin(); 
00614        cal!=caloJets->end(); ++cal){
00615     jetID->calculate(iEvent, *cal);
00616     if (cal->pt()>25.){
00617       //
00618       // for all regions
00619       if (jetID->fHPD()>=0.95) bJetIDTight=false; 
00620       //
00621       // for 1.0<|eta|<1.75
00622       if (fabs(cal->eta())>=1.00 && fabs(cal->eta())<1.75){
00623         if (cal->pt()>80. && cal->emEnergyFraction()>=1.) bJetIDTight=false; 
00624       }
00625       //
00626       // for 1.75<|eta|<2.55
00627       else if (fabs(cal->eta())>=1.75 && fabs(cal->eta())<2.55){
00628         if (cal->pt()>80. && cal->emEnergyFraction()>=1.) bJetIDTight=false; 
00629       }
00630       //
00631       // for 2.55<|eta|<3.25
00632       else if (fabs(cal->eta())>=2.55 && fabs(cal->eta())<3.25){
00633         if (cal->pt()< 50.                   && cal->emEnergyFraction()<=-0.3) bJetIDTight=false; 
00634         if (cal->pt()>=50. && cal->pt()< 80. && cal->emEnergyFraction()<=-0.2) bJetIDTight=false; 
00635         if (cal->pt()>=80. && cal->pt()<340. && cal->emEnergyFraction()<=-0.1) bJetIDTight=false; 
00636         if (cal->pt()>=340.                  && cal->emEnergyFraction()<=-0.1 
00637                                              && cal->emEnergyFraction()>=0.95) bJetIDTight=false; 
00638       }
00639       //
00640       // for 3.25<|eta|
00641       else if (fabs(cal->eta())>=3.25){
00642         if (cal->pt()< 50.                   && cal->emEnergyFraction()<=-0.3
00643                                              && cal->emEnergyFraction()>=0.90) bJetIDTight=false; 
00644         if (cal->pt()>=50. && cal->pt()<130. && cal->emEnergyFraction()<=-0.2
00645                                              && cal->emEnergyFraction()>=0.80) bJetIDTight=false; 
00646         if (cal->pt()>=130.                  && cal->emEnergyFraction()<=-0.1 
00647                                              && cal->emEnergyFraction()>=0.70) bJetIDTight=false; 
00648       }
00649     }   // pt>10 GeV/c
00650   }     // calor-jets loop
00651   
00652   if (_verbose) std::cout << "JetID ends" << std::endl;
00653 
00654 
00655   // ==========================================================
00656   // HCAL Noise filter
00657   
00658   bool bHcalNoiseFilter      = HNoiseSummary->passLooseNoiseFilter();
00659   bool bHcalNoiseFilterTight = HNoiseSummary->passTightNoiseFilter();
00660 
00661   // ==========================================================
00662   // Get BeamHaloSummary
00663   edm::Handle<BeamHaloSummary> TheBeamHaloSummary ;
00664   iEvent.getByLabel(BeamHaloSummaryTag, TheBeamHaloSummary) ;
00665 
00666   if (!TheBeamHaloSummary.isValid()) {
00667     std::cout << "BeamHaloSummary doesn't exist" << std::endl;
00668   }
00669 
00670   bool bBeamHaloIDTightPass = true;
00671   bool bBeamHaloIDLoosePass = true;
00672 
00673   if(!TheBeamHaloSummary.isValid()) {
00674 
00675   const BeamHaloSummary TheSummary = (*TheBeamHaloSummary.product() );
00676 
00677   if( !TheSummary.EcalLooseHaloId()  && !TheSummary.HcalLooseHaloId() && 
00678       !TheSummary.CSCLooseHaloId()   && !TheSummary.GlobalLooseHaloId() )
00679     bBeamHaloIDLoosePass = false;
00680 
00681   if( !TheSummary.EcalTightHaloId()  && !TheSummary.HcalTightHaloId() && 
00682       !TheSummary.CSCTightHaloId()   && !TheSummary.GlobalTightHaloId() )
00683     bBeamHaloIDTightPass = false;
00684 
00685   }
00686 
00687   // ==========================================================
00688   //Vertex information
00689 
00690   bool bPrimaryVertex = true;
00691   if(_doPVCheck){
00692     bPrimaryVertex = false;
00693     Handle<VertexCollection> vertexHandle;
00694 
00695     iEvent.getByLabel(vertexTag, vertexHandle);
00696 
00697     if (!vertexHandle.isValid()) {
00698       LogDebug("") << "CaloMETAnalyzer: Could not find vertex collection" << std::endl;
00699       if (_verbose) std::cout << "CaloMETAnalyzer: Could not find vertex collection" << std::endl;
00700     }
00701     
00702     if ( vertexHandle.isValid() ){
00703       VertexCollection vertexCollection = *(vertexHandle.product());
00704       int vertex_number     = vertexCollection.size();
00705       VertexCollection::const_iterator v = vertexCollection.begin();
00706       double vertex_chi2    = v->normalizedChi2();
00707       double vertex_ndof    = v->ndof();
00708       bool   fakeVtx        = v->isFake();
00709       double vertex_Z       = v->z();
00710       
00711       if (  !fakeVtx
00712             && vertex_number>=_nvtx_min
00713             && vertex_ndof   >_vtxndof_min
00714             && vertex_chi2   <_vtxchi2_max
00715             && fabs(vertex_Z)<_vtxz_max ) bPrimaryVertex = true;
00716     }
00717   }
00718   // ==========================================================
00719 
00720   edm::Handle< L1GlobalTriggerReadoutRecord > gtReadoutRecord;
00721   iEvent.getByLabel( gtTag, gtReadoutRecord);
00722 
00723   if (!gtReadoutRecord.isValid()) {
00724     LogDebug("") << "CaloMETAnalyzer: Could not find GT readout record" << std::endl;
00725     if (_verbose) std::cout << "CaloMETAnalyzer: Could not find GT readout record product" << std::endl;
00726   }
00727   
00728   bool bTechTriggers    = true;
00729   bool bTechTriggersAND = true;
00730   bool bTechTriggersOR  = false;
00731   bool bTechTriggersNOT = false;
00732 
00733   if (gtReadoutRecord.isValid()) {
00734     const TechnicalTriggerWord&  technicalTriggerWordBeforeMask = gtReadoutRecord->technicalTriggerWord();
00735     
00736     if (_techTrigsAND.size() == 0)
00737       bTechTriggersAND = true;
00738     else
00739       for (unsigned ttr = 0; ttr != _techTrigsAND.size(); ttr++) {
00740         bTechTriggersAND = bTechTriggersAND && technicalTriggerWordBeforeMask.at(_techTrigsAND.at(ttr));
00741       }
00742     
00743     if (_techTrigsAND.size() == 0)
00744       bTechTriggersOR = true;
00745     else
00746       for (unsigned ttr = 0; ttr != _techTrigsOR.size(); ttr++) {
00747         bTechTriggersOR = bTechTriggersOR || technicalTriggerWordBeforeMask.at(_techTrigsOR.at(ttr));
00748       }
00749     if (_techTrigsNOT.size() == 0)
00750       bTechTriggersNOT = false;
00751     else
00752       for (unsigned ttr = 0; ttr != _techTrigsNOT.size(); ttr++) {
00753         bTechTriggersNOT = bTechTriggersNOT || technicalTriggerWordBeforeMask.at(_techTrigsNOT.at(ttr));
00754       }
00755   }
00756   else
00757     {
00758       bTechTriggersAND = true;
00759       bTechTriggersOR  = true;
00760       bTechTriggersNOT = false;
00761     }
00762   
00763   if (_techTrigsAND.size()==0)
00764     bTechTriggersAND = true;
00765   if (_techTrigsOR.size()==0)
00766     bTechTriggersOR = true;
00767   if (_techTrigsNOT.size()==0)
00768     bTechTriggersNOT = false;
00769   
00770   bTechTriggers = bTechTriggersAND && bTechTriggersOR && !bTechTriggersNOT;
00771   
00772   // ==========================================================
00773   // Reconstructed MET Information - fill MonitorElements
00774   
00775   bool bHcalNoise   = bHcalNoiseFilter;
00776   bool bBeamHaloID  = bBeamHaloIDLoosePass;
00777   bool bJetID       = true;
00778 
00779   bool bPhysicsDeclared = true;
00780   if(_doHLTPhysicsOn) bPhysicsDeclared =_trig_PhysDec;
00781 
00782 
00783   if      (_tightHcalFiltering)     bHcalNoise  = bHcalNoiseFilterTight;
00784   if      (_tightBHFiltering)       bBeamHaloID = bBeamHaloIDTightPass;
00785 
00786   if      (_tightJetIDFiltering==1)  bJetID      = bJetIDMinimal;
00787   else if (_tightJetIDFiltering==2)  bJetID      = bJetIDLoose;
00788   else if (_tightJetIDFiltering==3)  bJetID      = bJetIDTight;
00789   else if (_tightJetIDFiltering==-1) bJetID      = true;
00790 
00791   bool bBasicCleanup = bTechTriggers && bPrimaryVertex && bPhysicsDeclared;
00792   bool bExtraCleanup = bBasicCleanup && bHcalNoise && bJetID && bBeamHaloID;
00793 
00794   //std::string DirName = _FolderName+_source;
00795   
00796   for (std::vector<std::string>::const_iterator ic = _FolderNames.begin(); 
00797        ic != _FolderNames.end(); ic++){
00798     if (*ic=="All")                                             fillMESet(iEvent, DirName+"/"+*ic, *met);
00799     if (DCSFilter->filter(iEvent, iSetup)) {
00800     if (_cleanupSelection){
00801     if (*ic=="BasicCleanup" && bBasicCleanup)                   fillMESet(iEvent, DirName+"/"+*ic, *met);
00802     if (*ic=="ExtraCleanup" && bExtraCleanup)                   fillMESet(iEvent, DirName+"/"+*ic, *met);
00803     }
00804     if (_allSelection) {
00805       if (*ic=="HcalNoiseFilter"      && bHcalNoiseFilter )       fillMESet(iEvent, DirName+"/"+*ic, *met);
00806       if (*ic=="HcalNoiseFilterTight" && bHcalNoiseFilterTight )  fillMESet(iEvent, DirName+"/"+*ic, *met);
00807       if (*ic=="JetIDMinimal"         && bJetIDMinimal)           fillMESet(iEvent, DirName+"/"+*ic, *met);
00808       if (*ic=="JetIDLoose"           && bJetIDLoose)             fillMESet(iEvent, DirName+"/"+*ic, *met);
00809       if (*ic=="JetIDTight"           && bJetIDTight)             fillMESet(iEvent, DirName+"/"+*ic, *met);
00810       if (*ic=="BeamHaloIDTightPass"  && bBeamHaloIDTightPass)    fillMESet(iEvent, DirName+"/"+*ic, *met);
00811       if (*ic=="BeamHaloIDLoosePass"  && bBeamHaloIDLoosePass)    fillMESet(iEvent, DirName+"/"+*ic, *met);
00812       if (*ic=="Triggers"             && bTechTriggers)           fillMESet(iEvent, DirName+"/"+*ic, *met);
00813       if (*ic=="PV"                   && bPrimaryVertex)          fillMESet(iEvent, DirName+"/"+*ic, *met);
00814     }
00815     } // DCS
00816   }
00817 }
00818 
00819 // ***********************************************************
00820 void METAnalyzer::fillMESet(const edm::Event& iEvent, std::string DirName, 
00821                               const reco::MET& met)
00822 {
00823 
00824   _dbe->setCurrentFolder(DirName);
00825 
00826   bool bLumiSecPlot=false;
00827   if (DirName.find("All")) bLumiSecPlot=true;
00828 
00829   if (_trig_JetMB) fillMonitorElement(iEvent,DirName,"",met, bLumiSecPlot);
00830   if (_trig_HighPtJet)
00831     fillMonitorElement(iEvent,DirName,"HighPtJet",met,false);
00832   if (_trig_LowPtJet)
00833     fillMonitorElement(iEvent,DirName,"LowPtJet",met,false);
00834   if (_trig_MinBias)
00835     fillMonitorElement(iEvent,DirName,"MinBias",met,false);
00836   if (_trig_HighMET)
00837     fillMonitorElement(iEvent,DirName,"HighMET",met,false);
00838   if (_trig_LowMET)
00839     fillMonitorElement(iEvent,DirName,"LowMET",met,false);
00840   if (_trig_Ele)
00841     fillMonitorElement(iEvent,DirName,"Ele",met,false);
00842   if (_trig_Muon)
00843     fillMonitorElement(iEvent,DirName,"Muon",met,false);
00844 
00845 }
00846 
00847 // ***********************************************************
00848 void METAnalyzer::fillMonitorElement(const edm::Event& iEvent, std::string DirName, 
00849                                          std::string TriggerTypeName, 
00850                                          const reco::MET& met, bool bLumiSecPlot)
00851 {
00852 
00853   if (TriggerTypeName=="HighPtJet") {
00854     if (!selectHighPtJetEvent(iEvent)) return;
00855   }
00856   else if (TriggerTypeName=="LowPtJet") {
00857     if (!selectLowPtJetEvent(iEvent)) return;
00858   }
00859   else if (TriggerTypeName=="HighMET") {
00860     if (met.pt()<_highMETThreshold) return;
00861   }
00862   else if (TriggerTypeName=="LowMET") {
00863     if (met.pt()<_lowMETThreshold) return;
00864   }
00865   else if (TriggerTypeName=="Ele") {
00866     if (!selectWElectronEvent(iEvent)) return;
00867   }
00868   else if (TriggerTypeName=="Muon") {
00869     if (!selectWMuonEvent(iEvent)) return;
00870   }
00871   
00872 // Reconstructed MET Information
00873   double SumET  = met.sumEt();
00874   double METSig = met.mEtSig();
00875   //removed for optimizations//double Ez     = met.e_longitudinal();
00876   double MET    = met.pt();
00877   double MEx    = met.px();
00878   double MEy    = met.py();
00879   double METPhi = met.phi();
00880 
00881   //
00882   int myLuminosityBlock;
00883   //  myLuminosityBlock = (evtCounter++)/1000;
00884   myLuminosityBlock = iEvent.luminosityBlock();
00885   //
00886 
00887   if (TriggerTypeName!="") DirName = DirName +"/"+TriggerTypeName;
00888 
00889   if (_verbose) std::cout << "_etThreshold = " << _etThreshold << std::endl;
00890   if (SumET>_etThreshold){
00891     
00892     hMEx    = _dbe->get(DirName+"/"+"METTask_MEx");     if (hMEx           && hMEx->getRootObject())     hMEx          ->Fill(MEx);
00893     hMEy    = _dbe->get(DirName+"/"+"METTask_MEy");     if (hMEy           && hMEy->getRootObject())     hMEy          ->Fill(MEy);
00894     hMET    = _dbe->get(DirName+"/"+"METTask_MET");     if (hMET           && hMET->getRootObject())     hMET          ->Fill(MET);
00895     hMETPhi = _dbe->get(DirName+"/"+"METTask_METPhi");  if (hMETPhi        && hMETPhi->getRootObject())  hMETPhi       ->Fill(METPhi);
00896     hSumET  = _dbe->get(DirName+"/"+"METTask_SumET");   if (hSumET         && hSumET->getRootObject())   hSumET        ->Fill(SumET);
00897     hMETSig = _dbe->get(DirName+"/"+"METTask_METSig");  if (hMETSig        && hMETSig->getRootObject())  hMETSig       ->Fill(METSig);
00898     //removed for optimizations//hEz     = _dbe->get(DirName+"/"+"METTask_Ez");      if (hEz            && hEz->getRootObject())      hEz           ->Fill(Ez);
00899 
00900     hMET_logx   = _dbe->get(DirName+"/"+"METTask_MET_logx");    if (hMET_logx      && hMET_logx->getRootObject())    hMET_logx->Fill(log10(MET));
00901     hSumET_logx = _dbe->get(DirName+"/"+"METTask_SumET_logx");  if (hSumET_logx    && hSumET_logx->getRootObject())  hSumET_logx->Fill(log10(SumET));
00902 
00903     //removed for optimizations//hMETIonFeedbck = _dbe->get(DirName+"/"+"METTask_METIonFeedbck");  if (hMETIonFeedbck && hMETIonFeedbck->getRootObject())  hMETIonFeedbck->Fill(MET);
00904     //removed for optimizations//hMETHPDNoise   = _dbe->get(DirName+"/"+"METTask_METHPDNoise");    if (hMETHPDNoise   && hMETHPDNoise->getRootObject())    hMETHPDNoise->Fill(MET);
00905        
00906     if (_allhist){
00907       if (bLumiSecPlot){
00908         hMExLS = _dbe->get(DirName+"/"+"METTask_MExLS"); if (hMExLS  &&  hMExLS->getRootObject())   hMExLS->Fill(MEx,myLuminosityBlock);
00909         hMEyLS = _dbe->get(DirName+"/"+"METTask_MEyLS"); if (hMEyLS  &&  hMEyLS->getRootObject())   hMEyLS->Fill(MEy,myLuminosityBlock);
00910       }
00911     } // _allhist
00912 
00914     if (theMETCollectionLabel.label() == "tcMet" ) {
00915     
00916       if(track_h.isValid()) {
00917         for( edm::View<reco::Track>::const_iterator trkit = track_h->begin(); trkit != track_h->end(); trkit++ ) {
00918           htrkPt    = _dbe->get(DirName+"/"+"METTask_trackPt");     if (htrkPt    && htrkPt->getRootObject())     htrkPt->Fill( trkit->pt() );
00919           htrkEta   = _dbe->get(DirName+"/"+"METTask_trackEta");    if (htrkEta   && htrkEta->getRootObject())    htrkEta->Fill( trkit->eta() );
00920           htrkNhits = _dbe->get(DirName+"/"+"METTask_trackNhits");  if (htrkNhits && htrkNhits->getRootObject())  htrkNhits->Fill( trkit->numberOfValidHits() );
00921           htrkChi2  = _dbe->get(DirName+"/"+"METTask_trackChi2");   if (htrkChi2  && htrkChi2->getRootObject())   htrkChi2->Fill( trkit->chi2() / trkit->ndof() );
00922           double d0 = -1 * trkit->dxy( bspot );
00923           htrkD0    = _dbe->get(DirName+"/"+"METTask_trackD0");     if (htrkD0 && htrkD0->getRootObject())        htrkD0->Fill( d0 );
00924         }
00925       }
00926       
00927       if(electron_h.isValid()) {
00928         for( edm::View<reco::GsfElectron>::const_iterator eleit = electron_h->begin(); eleit != electron_h->end(); eleit++ ) {
00929           helePt  = _dbe->get(DirName+"/"+"METTask_electronPt");   if (helePt  && helePt->getRootObject())   helePt->Fill( eleit->p4().pt() );  
00930           heleEta = _dbe->get(DirName+"/"+"METTask_electronEta");  if (heleEta && heleEta->getRootObject())  heleEta->Fill( eleit->p4().eta() );
00931           heleHoE = _dbe->get(DirName+"/"+"METTask_electronHoE");  if (heleHoE && heleHoE->getRootObject())  heleHoE->Fill( eleit->hadronicOverEm() );
00932         }
00933       }
00934       
00935       if(muon_h.isValid()) {
00936         for( reco::MuonCollection::const_iterator muonit = muon_h->begin(); muonit != muon_h->end(); muonit++ ) {      
00937           const reco::TrackRef siTrack = muonit->innerTrack();      
00938           hmuPt    = _dbe->get(DirName+"/"+"METTask_muonPt");     if (hmuPt    && hmuPt->getRootObject())  hmuPt   ->Fill( muonit->p4().pt() );
00939           hmuEta   = _dbe->get(DirName+"/"+"METTask_muonEta");    if (hmuEta   && hmuEta->getRootObject())  hmuEta  ->Fill( muonit->p4().eta() );
00940           hmuNhits = _dbe->get(DirName+"/"+"METTask_muonNhits");  if (hmuNhits && hmuNhits->getRootObject())  hmuNhits->Fill( siTrack.isNonnull() ? siTrack->numberOfValidHits() : -999 );
00941           hmuChi2  = _dbe->get(DirName+"/"+"METTask_muonChi2");   if (hmuChi2  && hmuChi2->getRootObject())  hmuChi2 ->Fill( siTrack.isNonnull() ? siTrack->chi2()/siTrack->ndof() : -999 );
00942           double d0 = siTrack.isNonnull() ? -1 * siTrack->dxy( bspot) : -999;
00943           hmuD0    = _dbe->get(DirName+"/"+"METTask_muonD0");     if (hmuD0    && hmuD0->getRootObject())  hmuD0->Fill( d0 );
00944         }
00945         
00946         const unsigned int nMuons = muon_h->size();      
00947         for( unsigned int mus = 0; mus < nMuons; mus++ ) {
00948           reco::MuonRef muref( muon_h, mus);
00949           reco::MuonMETCorrectionData muCorrData = (*tcMet_ValueMap_Handle)[muref];
00950           hMExCorrection      = _dbe->get(DirName+"/"+"METTask_MExCorrection");       if (hMExCorrection      && hMExCorrection->getRootObject())       hMExCorrection-> Fill(muCorrData.corrY());
00951           hMEyCorrection      = _dbe->get(DirName+"/"+"METTask_MEyCorrection");       if (hMEyCorrection      && hMEyCorrection->getRootObject())       hMEyCorrection-> Fill(muCorrData.corrX());
00952           hMuonCorrectionFlag = _dbe->get(DirName+"/"+"METTask_CorrectionFlag");      if (hMuonCorrectionFlag && hMuonCorrectionFlag->getRootObject())  hMuonCorrectionFlag-> Fill(muCorrData.type());
00953         }
00954       }
00955     }
00956 
00957   } // et threshold cut
00958 
00959 }
00960 
00961 // ***********************************************************
00962 bool METAnalyzer::selectHighPtJetEvent(const edm::Event& iEvent){
00963 
00964   bool return_value=false;
00965 
00966   edm::Handle<reco::CaloJetCollection> caloJets;
00967   iEvent.getByLabel(theJetCollectionLabel, caloJets);
00968   if (!caloJets.isValid()) {
00969     LogDebug("") << "METAnalyzer: Could not find jet product" << std::endl;
00970     if (_verbose) std::cout << "METAnalyzer: Could not find jet product" << std::endl;
00971   }
00972 
00973   for (reco::CaloJetCollection::const_iterator cal = caloJets->begin(); 
00974        cal!=caloJets->end(); ++cal){
00975     if (cal->pt()>_highPtJetThreshold){
00976       return_value=true;
00977     }
00978   }
00979   
00980   return return_value;
00981 }
00982 
00983 // // ***********************************************************
00984 bool METAnalyzer::selectLowPtJetEvent(const edm::Event& iEvent){
00985 
00986   bool return_value=false;
00987 
00988   edm::Handle<reco::CaloJetCollection> caloJets;
00989   iEvent.getByLabel(theJetCollectionLabel, caloJets);
00990   if (!caloJets.isValid()) {
00991     LogDebug("") << "METAnalyzer: Could not find jet product" << std::endl;
00992     if (_verbose) std::cout << "METAnalyzer: Could not find jet product" << std::endl;
00993   }
00994 
00995   for (reco::CaloJetCollection::const_iterator cal = caloJets->begin(); 
00996        cal!=caloJets->end(); ++cal){
00997     if (cal->pt()>_lowPtJetThreshold){
00998       return_value=true;
00999     }
01000   }
01001 
01002   return return_value;
01003 
01004 }
01005 
01006 
01007 // ***********************************************************
01008 bool METAnalyzer::selectWElectronEvent(const edm::Event& iEvent){
01009 
01010   bool return_value=true;
01011 
01012   /*
01013     W-electron event selection comes here
01014    */
01015 
01016   return return_value;
01017 
01018 }
01019 
01020 // ***********************************************************
01021 bool METAnalyzer::selectWMuonEvent(const edm::Event& iEvent){
01022 
01023   bool return_value=true;
01024 
01025   /*
01026     W-muon event selection comes here
01027    */
01028 
01029   return return_value;
01030 
01031 }
01032