CMS 3D CMS Logo

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