CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/DQM/Physics/src/EwkDQM.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2012/10/25 11:54:49 $
00005  *  $Revision: 1.22 $
00006  *  \author Valentina Gori, University of Firenze
00007  */
00008 
00009 #include "DQM/Physics/src/EwkDQM.h"
00010 
00011 #include <vector>
00012 #include <string>
00013 #include <cmath>
00014 
00015 #include "DataFormats/Candidate/interface/Candidate.h"
00016 #include "FWCore/Framework/interface/Event.h"
00017 #include "FWCore/Framework/interface/EventSetup.h"
00018 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00019 
00020 #include "DQMServices/Core/interface/DQMStore.h"
00021 #include "DQMServices/Core/interface/MonitorElement.h"
00022 #include "FWCore/ServiceRegistry/interface/Service.h"
00023 
00024 #include "DataFormats/Common/interface/Handle.h"
00025 
00026 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00027 #include "FWCore/Common/interface/TriggerNames.h"
00028 
00029 // Physics Objects
00030 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00031 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00032 #include "DataFormats/MuonReco/interface/Muon.h"
00033 #include "DataFormats/MuonReco/interface/MuonFwd.h"
00034 #include "DataFormats/JetReco/interface/Jet.h"
00035 #include "DataFormats/METReco/interface/MET.h"
00036 #include "DataFormats/TrackReco/interface/Track.h"
00037 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00038 #include "DataFormats/VertexReco/interface/Vertex.h"
00039 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00040 #include "DataFormats/Math/interface/LorentzVector.h"
00041 
00042 #include "TLorentzVector.h"
00043 
00044 using namespace std;
00045 using namespace edm;
00046 using namespace reco;
00047 
00048 
00049 
00050 // EwkDQM::EwkDQM(const ParameterSet& parameters) {
00051 EwkDQM::EwkDQM(const ParameterSet& parameters) {
00052   eJetMin_     = parameters.getUntrackedParameter<double>("EJetMin", 999999.);
00053 
00054   // riguardare questa sintassi
00055   // Get parameters from configuration file
00056   theElecTriggerPathToPass_ =
00057     parameters.getParameter<std::vector<string> >("elecTriggerPathToPass");
00058   theMuonTriggerPathToPass_ =
00059     parameters.getParameter<std::vector<string> >("muonTriggerPathToPass");
00060   //   eleTrigPathNames_ =
00061   //     parameters.getUntrackedParameter<std::vector<std::string> >("eleTrigPathNames");
00062   //   muTrigPathNames_ =
00063   //     parameters.getUntrackedParameter< std::vector<std::string> >("muTrigPathNames");
00064   theTriggerResultsCollection_ =
00065     parameters.getParameter<InputTag>("triggerResultsCollection");
00066   theMuonCollectionLabel_ = parameters.getParameter<InputTag>("muonCollection");
00067   theElectronCollectionLabel_ =
00068     parameters.getParameter<InputTag>("electronCollection");
00069   //  theCaloJetCollectionLabel_   = parameters.getParameter<InputTag>("caloJetCollection");
00070   thePFJetCollectionLabel_ =
00071     parameters.getParameter<InputTag>("PFJetCollection");
00072   theCaloMETCollectionLabel_ = parameters.getParameter<InputTag>("caloMETCollection");
00073 
00074   // just to initialize
00075   isValidHltConfig_ = false;
00076 
00077   h_vertex_number = 0;
00078   h_vertex_chi2 = 0;
00079   h_vertex_numTrks = 0;
00080   h_vertex_sumTrks = 0;
00081   h_vertex_d0 = 0;
00082 
00083   h_jet_count = 0;
00084   h_jet_et = 0;
00085   h_jet_pt = 0;
00086   h_jet_eta = 0;
00087   h_jet_phi = 0;
00088   h_jet2_et = 0;
00089   // h_jet2_pt = 0;
00090   h_jet2_eta = 0;
00091   h_jet2_phi = 0;
00092 
00093   h_e1_et = 0;
00094   h_e2_et = 0;
00095   h_e1_eta = 0;
00096   h_e2_eta = 0;
00097   h_e1_phi = 0;
00098   h_e2_phi = 0;
00099 
00100   h_m1_pt = 0;
00101   h_m2_pt = 0;
00102   h_m1_eta = 0;
00103   h_m2_eta = 0;
00104   h_m1_phi = 0;
00105   h_m2_phi = 0;
00106 
00107   // h_t1_et = 0;
00108   // h_t1_eta = 0;
00109   // h_t1_phi = 0;
00110 
00111   h_met = 0;
00112   h_met_phi = 0;
00113 
00114   h_e_invWMass = 0;
00115   h_m_invWMass = 0;
00116   h_mumu_invMass = 0;
00117   h_ee_invMass = 0;
00118 
00119   theDbe = Service<DQMStore>().operator->();
00120 }
00121 
00122 EwkDQM::~EwkDQM() {
00123 }
00124 
00125 
00126 void EwkDQM::beginJob() {
00127   char chtitle[256] = "";
00128   const size_t title_s = sizeof(chtitle);
00129 
00130   logTraceName = "EwkAnalyzer";
00131 
00132   LogTrace(logTraceName) << "Parameters initialization";
00133   theDbe->setCurrentFolder("Physics/EwkDQM");  // Use folder with name of PAG
00134 
00135   const float pi = 4*atan(1);
00136 
00137   // Keep the number of plots and number of bins to a minimum!
00138   h_vertex_number = theDbe->book1D("vertex_number",
00139                                    "Number of event vertices in collection",
00140                                    10, -0.5, 9.5);
00141   h_vertex_chi2  = theDbe->book1D("vertex_chi2",
00142                                   "Event Vertex #chi^{2}/n.d.o.f.",
00143                                   20, 0.0, 2.0);
00144   h_vertex_numTrks = theDbe->book1D("vertex_numTrks",
00145                                     "Event Vertex, number of tracks",
00146                                     20, -0.5, 59.5);
00147   h_vertex_sumTrks = theDbe->book1D("vertex_sumTrks",
00148                                     "Event Vertex, sum of track pt",
00149                                     20, 0.0, 100.0);
00150   h_vertex_d0 = theDbe->book1D("vertex_d0",
00151                                "Event Vertex d0", 20, 0.0, 0.05);
00152 
00153   snprintf(chtitle, title_s, "Number of %s (E_{T} > 15 GeV);Number of Jets",
00154            thePFJetCollectionLabel_.label().data());
00155   h_jet_count = theDbe->book1D("jet_count", chtitle, 8, -0.5, 7.5);
00156 
00157   snprintf(chtitle, title_s, "Leading jet E_{T} (from %s);E_{T}(1^{st} jet) (GeV)",
00158            thePFJetCollectionLabel_.label().data());
00159   h_jet_et = theDbe->book1D("jet_et", chtitle,    20, 0., 200.0);
00160 
00161   snprintf(chtitle, title_s, "Leading jet p_{T} (from %s);p_{T}(1^{st} jet) (GeV/c)",
00162            thePFJetCollectionLabel_.label().data());
00163   h_jet_pt = theDbe->book1D("jet_pt", chtitle,  20, 0., 200.0);
00164 
00165   snprintf(chtitle, title_s,  "Leading jet #eta (from %s); #eta (1^{st} jet)",
00166            thePFJetCollectionLabel_.label().data());
00167   h_jet_eta = theDbe->book1D("jet_eta", chtitle, 20, -10., 10.0);
00168   snprintf(chtitle, title_s, "Leading jet #phi (from %s); #phi(1^{st} jet)",
00169            thePFJetCollectionLabel_.label().data());
00170   h_jet_phi = theDbe->book1D("jet_phi", chtitle, 22, -1.1*pi, 1.1*pi);
00171 
00172   snprintf(chtitle, title_s, "2^{nd} leading jet E_{T} (from %s);E_{T}(2^{nd} jet) (GeV)",
00173            thePFJetCollectionLabel_.label().data());
00174   h_jet2_et = theDbe->book1D("jet2_et", chtitle,  20, 0., 200.0);
00175   // snprintf(chtitle, title_s, "2^{nd} leading jet p_{T} (from %s);p_{T}(2^{nd} jet) (GeV/c)",
00176   // thePFJetCollectionLabel_.label().data());
00177   // h_jet2_pt = theDbe->book1D("jet2_pt", chtitle,  20, 0., 200.0);
00178 
00179   snprintf(chtitle, title_s,  "2^{nd} leading jet #eta (from %s); #eta (2^{nd} jet)",
00180            thePFJetCollectionLabel_.label().data());
00181   h_jet2_eta = theDbe->book1D("jet2_eta", chtitle,  20, -10., 10.0);
00182 
00183   snprintf(chtitle, title_s, "2^{nd} leading jet #phi (from %s); #phi(2^{nd} jet)",
00184            thePFJetCollectionLabel_.label().data());
00185   h_jet2_phi = theDbe->book1D("jet2_phi", chtitle,  22, -1.1*pi, 1.1*pi);
00186 
00187   h_e1_et  = theDbe->book1D("e1_et",  "E_{T} of Leading Electron;E_{T} (GeV)",
00188                             20, 0.0, 100.0);
00189   h_e2_et  = theDbe->book1D("e2_et",  "E_{T} of Second Electron;E_{T} (GeV)",
00190                             20, 0.0, 100.0);
00191   h_e1_eta = theDbe->book1D("e1_eta", "#eta of Leading Electron;#eta",
00192                             20, -4.0, 4.0);
00193   h_e2_eta = theDbe->book1D("e2_eta", "#eta of Second Electron;#eta",
00194                             20, -4.0, 4.0);
00195   h_e1_phi = theDbe->book1D("e1_phi", "#phi of Leading Electron;#phi",
00196                             22, -1.1*pi, 1.1*pi);
00197   h_e2_phi = theDbe->book1D("e2_phi", "#phi of Second Electron;#phi",
00198                             22, -1.1*pi, 1.1*pi);
00199   h_m1_pt  = theDbe->book1D("m1_pt",  "p_{T} of Leading Muon;p_{T}(1^{st} #mu) (GeV)",
00200                             20, 0.0, 100.0);
00201   h_m2_pt  = theDbe->book1D("m2_pt",  "p_{T} of Second Muon;p_{T}(2^{nd} #mu) (GeV)",
00202                             20, 0.0, 100.0);
00203   h_m1_eta = theDbe->book1D("m1_eta", "#eta of Leading Muon;#eta(1^{st} #mu)",
00204                             20, -4.0, 4.0);
00205   h_m2_eta = theDbe->book1D("m2_eta", "#eta of Second Muon;#eta(2^{nd} #mu)",
00206                             20, -4.0, 4.0);
00207   h_m1_phi = theDbe->book1D("m1_phi", "#phi of Leading Muon;#phi(1^{st} #mu)",
00208                             20, (-1. - 1./10.)*pi, (1. + 1./10.)*pi);
00209   h_m2_phi = theDbe->book1D("m2_phi", "#phi of Second Muon;#phi(2^{nd} #mu)",
00210                             20, (-1. - 1./10.)*pi, (1. + 1./10.)*pi);
00211   //  h_t1_et = theDbe->book1D("t1_et", "E_{T} of Leading Tau;E_{T} (GeV)",
00212   //  20, 0.0 , 100.0);
00213   //  h_t1_eta = theDbe->book1D("t1_eta", "#eta of Leading Tau;#eta",
00214   //  20, -4.0, 4.0);
00215   //  h_t1_phi = theDbe->book1D("t1_phi", "#phi of Leading Tau;#phi",
00216   //  20, -4.0, 4.0);
00217   snprintf(chtitle, title_s, "Missing E_{T} (%s); GeV",
00218            theCaloMETCollectionLabel_.label().data());
00219   h_met          = theDbe->book1D("met", chtitle, 20, 0.0, 100);
00220   h_met_phi      = theDbe->book1D("met_phi", "Missing E_{T} #phi;#phi(MET)",
00221                                   22, (-1. - 1./10.)*pi, (1. + 1./10.)*pi);
00222 
00223   h_e_invWMass   = theDbe->book1D("we_invWMass", "W-> e #nu Transverse Mass;M_{T} (GeV)",
00224                                   20, 0.0, 140.0);
00225   h_m_invWMass   = theDbe->book1D("wm_invWMass", "W-> #mu #nu Transverse Mass;M_{T} (GeV)",
00226                                   20, 0.0, 140.0);
00227   h_mumu_invMass = theDbe->book1D("z_mm_invMass", "#mu#mu Invariant Mass;InvMass (GeV)",
00228                                   20, 40.0, 140.0);
00229   h_ee_invMass   = theDbe->book1D("z_ee_invMass", "ee Invariant Mass;InvMass (Gev)",
00230                                   20, 40.0, 140.0);
00231 }
00232 
00236 void EwkDQM::beginRun(const edm::Run& theRun, const edm::EventSetup& theSetup) {
00237   // passed as parameter to HLTConfigProvider::init(), not yet used
00238   bool isConfigChanged = false;
00239 
00240   // isValidHltConfig_ used to short-circuit analyze() in case of problems
00241   const std::string hltProcessName(theTriggerResultsCollection_.process());
00242   isValidHltConfig_ = hltConfigProvider_.init(theRun, theSetup,
00243                                               hltProcessName, isConfigChanged);
00244 }
00245 
00246 
00247 void EwkDQM::analyze(const Event& iEvent, const EventSetup& iSetup) {
00248   // short-circuit if hlt problems
00249   if (!isValidHltConfig_)
00250     return;
00251 
00252   LogTrace(logTraceName) << "Analysis of event # ";
00253   // Did it pass certain HLT path?
00254   Handle<TriggerResults> HLTresults;
00255   iEvent.getByLabel(theTriggerResultsCollection_, HLTresults);
00256   if (!HLTresults.isValid())
00257     return;
00258 
00259   const edm::TriggerNames & trigNames = iEvent.triggerNames(*HLTresults);
00260 
00261   // a temporary, until we have a list of triggers of interest
00262   std::vector<std::string> eleTrigPathNames;
00263   std::vector<std::string> muTrigPathNames;
00264 
00265   // eleTrigPathNames.push_back(theElecTriggerPathToPass_);
00266   // muTrigPathNames.push_back(theMuonTriggerPathToPass_);
00267   // end of temporary
00268 
00269   bool passed_electron_HLT = false;
00270   bool passed_muon_HLT     = false;
00271   for (unsigned int i = 0; i < HLTresults->size(); i++) {
00272     const std::string trigName = trigNames.triggerName(i);
00273     // check if triggerName matches electronPath
00274     for (unsigned int index = 0;
00275          index < theElecTriggerPathToPass_.size() && !passed_electron_HLT;
00276          index++) {
00277       // 0 if found, pos if not
00278       size_t trigPath = trigName.find(theElecTriggerPathToPass_[index]);
00279       if (trigPath == 0) {
00280         //      cout << "MuonTrigger passed (=trigName): " << trigName <<endl;
00281         passed_electron_HLT = HLTresults->accept(i);
00282       }
00283     }
00284     // check if triggerName matches muonPath
00285     for (unsigned int index = 0;
00286          index < theMuonTriggerPathToPass_.size() && !passed_muon_HLT;
00287          index++) {
00288       // 0 if found, pos if not
00289       size_t trigPath = trigName.find(theMuonTriggerPathToPass_[index]);
00290       if (trigPath == 0) {
00291         //      cout << "MuonTrigger passed (=trigName): " << trigName <<endl;
00292         passed_muon_HLT = HLTresults->accept(i);
00293       }
00294     }
00295   }
00296 
00297   // we are interested in events with a valid electron or muon
00298   if (!(passed_electron_HLT || passed_muon_HLT))
00299     return;
00300 
00302   // Vertex information
00303   Handle<VertexCollection> vertexHandle;
00304   iEvent.getByLabel("offlinePrimaryVertices", vertexHandle);
00305   if (!vertexHandle.isValid())
00306     return;
00307   VertexCollection vertexCollection = *(vertexHandle.product());
00308   VertexCollection::const_iterator v = vertexCollection.begin();
00309   int vertex_number     = vertexCollection.size();
00310   double vertex_chi2    = v->normalizedChi2();  // v->chi2();
00311   double vertex_d0      = sqrt(v->x()*v->x()+v->y()*v->y());
00312   double vertex_numTrks = v->tracksSize();
00313   double vertex_sumTrks = 0.0;
00314   // std::cout << "vertex_d0=" << vertex_d0 << "\n";
00315   // double vertex_ndof    = v->ndof();cout << "ndof="<<vertex_ndof<<endl;
00316   for (Vertex::trackRef_iterator vertex_curTrack = v->tracks_begin();
00317        vertex_curTrack != v->tracks_end(); vertex_curTrack++)
00318     vertex_sumTrks += (*vertex_curTrack)->pt();
00319 
00321   // Missing ET
00322   Handle< View<MET> > caloMETCollection;
00323   iEvent.getByLabel(theCaloMETCollectionLabel_, caloMETCollection);
00324   if (!caloMETCollection.isValid())
00325     return;
00326   float missing_et = caloMETCollection->begin()->et();
00327   float met_phi    = caloMETCollection->begin()->phi();
00328 
00329 
00331   // grab "gaussian sum fitting" electrons
00332   Handle<GsfElectronCollection> electronCollection;
00333   iEvent.getByLabel(theElectronCollectionLabel_, electronCollection);
00334   if (!electronCollection.isValid())
00335     return;
00336 
00337   // Find the highest and 2nd highest electron
00338   float electron_et   = -8.0;
00339   float electron_eta  = -8.0;
00340   float electron_phi  = -8.0;
00341   float electron2_et  = -9.0;
00342   float electron2_eta = -9.0;
00343   float electron2_phi = -9.0;
00344   float ee_invMass    = -9.0;
00345   TLorentzVector e1, e2;
00346 
00347   // If it passed electron HLT and the collection was found, find electrons near Z mass
00348   if (passed_electron_HLT) {
00349     for (reco::GsfElectronCollection::const_iterator recoElectron = electronCollection->begin();
00350          recoElectron != electronCollection->end(); recoElectron++) {
00351       // Require electron to pass some basic cuts
00352       if (recoElectron->et() < 20 || fabs(recoElectron->eta()) > 2.5)
00353         continue;
00354 
00355       // Tighter electron cuts
00356       if (recoElectron->deltaPhiSuperClusterTrackAtVtx() > 0.58 ||
00357           recoElectron->deltaEtaSuperClusterTrackAtVtx() > 0.01 ||
00358           recoElectron->sigmaIetaIeta() > 0.027)
00359         continue;
00360 
00361       if (recoElectron->et() > electron_et) {
00362         electron2_et  = electron_et;  // 2nd highest gets values from current highest
00363         electron2_eta = electron_eta;
00364         electron2_phi = electron_phi;
00365         electron_et   = recoElectron->et();  // 1st highest gets values from new highest
00366         electron_eta  = recoElectron->eta();
00367         electron_phi  = recoElectron->phi();
00368         e1 = TLorentzVector(recoElectron->momentum().x(),
00369                             recoElectron->momentum().y(),
00370                             recoElectron->momentum().z(),
00371                             recoElectron->p());
00372       } else if (recoElectron->et() > electron2_et) {
00373         electron2_et  = recoElectron->et();
00374         electron2_eta = recoElectron->eta();
00375         electron2_phi = recoElectron->phi();
00376         e2 = TLorentzVector(recoElectron->momentum().x(),
00377                             recoElectron->momentum().y(),
00378                             recoElectron->momentum().z(),
00379                             recoElectron->p());
00380       }
00381     }  // end of loop over electrons
00382     if (electron2_et > 0.0) {
00383       TLorentzVector pair = e1+e2;
00384       ee_invMass = pair.M();
00385     }
00386   }  // end of "are electrons valid"
00388 
00389 
00390 
00392   // Take the STA muon container
00393   Handle<MuonCollection> muonCollection;
00394   iEvent.getByLabel(theMuonCollectionLabel_, muonCollection);
00395   if (!muonCollection.isValid())
00396     return;
00397 
00398   // Find the highest pt muons
00399   float mm_invMass = -9.0;
00400   float muon_pt    = -9.0;
00401   float muon_eta   = -9.0;
00402   float muon_phi   = -9.0;
00403   float muon2_pt   = -9.0;
00404   float muon2_eta  = -9.0;
00405   float muon2_phi  = -9.0;
00406   TLorentzVector m1, m2;
00407 
00408   if (passed_muon_HLT) {
00409     for (reco::MuonCollection::const_iterator recoMuon = muonCollection->begin();
00410          recoMuon != muonCollection->end(); recoMuon++) {
00411       // Require muon to pass some basic cuts
00412       if (recoMuon->pt() < 20 || !recoMuon->isGlobalMuon())
00413         continue;
00414       // Some tighter muon cuts
00415       if (recoMuon->globalTrack()->normalizedChi2() > 10)
00416         continue;
00417 
00418       if (recoMuon->pt() > muon_pt) {
00419         muon2_pt  = muon_pt;  // 2nd highest gets values from current highest
00420         muon2_eta = muon_eta;
00421         muon2_phi = muon_phi;
00422         muon_pt   = recoMuon->pt();  // 1st highest gets values from new highest
00423         muon_eta  = recoMuon->eta();
00424         muon_phi  = recoMuon->phi();
00425         m1 = TLorentzVector(recoMuon->momentum().x(),
00426                             recoMuon->momentum().y(),
00427                             recoMuon->momentum().z(),
00428                             recoMuon->p());
00429       } else if (recoMuon->pt() > muon2_pt) {
00430         muon2_pt  = recoMuon->pt();
00431         muon2_eta = recoMuon->eta();
00432         muon2_phi = recoMuon->phi();
00433         m2 = TLorentzVector(recoMuon->momentum().x(),
00434                             recoMuon->momentum().y(),
00435                             recoMuon->momentum().z(),
00436                             recoMuon->p());
00437       }
00438     }
00439   }
00440   if (muon2_pt > 0.0) {
00441     TLorentzVector pair = m1+m2;
00442     mm_invMass = pair.M();
00443   }
00445 
00446 
00448   // Find the highest et jet
00449 
00450   //  Handle<CaloJetCollection> caloJetCollection;
00451   Handle<View<Jet> > PFJetCollection;
00452   //  iEvent.getByLabel (theCaloJetCollectionLabel,caloJetCollection);
00453   iEvent.getByLabel(thePFJetCollectionLabel_, PFJetCollection);
00454   //  if ( !caloJetCollection.isValid() ) return;
00455   if (!PFJetCollection.isValid())
00456     return;
00457 
00458   unsigned int muonCollectionSize = muonCollection->size();
00459   // unsigned int jetCollectionSize = jetCollection->size();
00460   unsigned int PFJetCollectionSize = PFJetCollection->size();
00461   int jet_count = 0;
00462   // int LEADJET=-1;  double max_pt=0;
00463 
00464 
00465   float jet_et    = -80.0;
00466   float jet_pt    = -80.0;  // prova
00467   float jet_eta   = -80.0;  // now USED
00468   float jet_phi   = -80.0;  // now USED
00469   float jet2_et   = -90.0;
00470   float jet2_eta  = -90.0;  // now USED
00471   float jet2_phi  = -90.0;  // now USED
00472   //  for (CaloJetCollection::const_iterator i_calojet = caloJetCollection->begin();
00473   //       i_calojet != caloJetCollection->end(); i_calojet++) {
00474   //  for (PFJetCollection::const_iterator i_pfjet = PFJetCollection->begin();
00475   //       i_pfjet != PFJetCollection->end(); i_pfjet++) {
00476   //  float jet_current_et = i_calojet->et();
00477   //  float jet_current_et = i_pfjet->et();            // e` identico a jet.et()
00478   //    jet_count++;
00479 
00480   // cleaning: va messo prima del riempimento dell'istogramma // This is in order to use PFJets
00481   for (unsigned int i = 0; i < PFJetCollectionSize; i++) {
00482     const Jet& jet = PFJetCollection->at(i);
00483     // la classe "jet" viene definita qui!!!
00484     double minDistance = 99999;
00485     for (unsigned int j = 0; j < muonCollectionSize; j++) {
00486       const Muon& mu = muonCollection->at(j);
00487       double distance = sqrt((mu.eta() - jet.eta()) * (mu.eta() - jet.eta()) +
00488                              (mu.phi() - jet.phi()) * (mu.phi() - jet.phi()));
00489       if (minDistance > distance)
00490         minDistance = distance;
00491     }
00492     if (minDistance < 0.3)
00493       continue;  // 0.3 is the isolation cone around the muon
00494     // se la distanza muone-cono del jet e` minore di 0.3, passo avanti e non conteggio il mio jet
00495 
00496     // If it overlaps with ELECTRON, it is not a jet
00497     if (electron_et > 0.0 &&
00498         fabs(jet.eta() - electron_eta) < 0.2 &&
00499         calcDeltaPhi(jet.phi(), electron_phi) < 0.2)
00500       continue;
00501     if (electron2_et > 0.0 &&
00502         fabs(jet.eta() - electron2_eta) < 0.2 &&
00503         calcDeltaPhi(jet.phi(), electron2_phi) < 0.2)
00504       continue;
00505 
00506     // provo a cambiare la parte degli elettroni in modo simmetrico alla parte per i muoni
00507 
00508     // ...
00509     // ...
00510 
00511 
00512     // if it has too low Et, throw away
00513     if (jet.et() < eJetMin_)
00514       continue;
00515     jet_count++;
00516 
00517     // ovvero: incrementa jet_count se:
00518     //   - non c'e un muone entro 0.3 di distanza dal cono del jet;
00519     //   - se il jet non si sovrappone ad un elettrone;
00520     //   - se l'energia trasversa e` maggiore della soglia impostata (15?)
00521 
00522     // if(jet.et()>max_pt) { LEADJET=i; max_pt=jet.et();}
00523     // se l'energia del jet e` maggiore di max_pt, diventa "i"
00524     // l'indice del jet piu` energetico e max_pt la sua energia
00525 
00526     // riguardare questo!!!
00527     // fino ad ora, jet_et era inizializzato a -8.0
00528     if (jet.et() > jet_et) {
00529       jet2_et  = jet_et;  // 2nd highest jet gets et from current highest
00530       // perche` prende l'energia del primo jet??
00531       jet2_eta = jet_eta;  // now USED
00532       jet2_phi = jet_phi;  // now USED
00533       // jet_et = i_calojet->et(); // current highest jet gets
00534       // et from the new highest
00535       jet_et   = jet.et();  // current highest jet gets et from the new highest
00536       // ah, ok! lo riaggiorna solo dopo!
00537       jet_pt   = jet.pt();   // e` il pT del leading jet
00538       jet_eta  = jet.eta();  // now USED
00539       jet_phi  = jet.phi() * (Geom::pi() / 180.);  // now USED
00540     } else if (jet.et() > jet2_et) {
00541       //      jet2_et  = i_calojet->et();
00542       jet2_et  = jet.et();
00543       //      jet2_eta = i_calojet->eta();  // UNUSED
00544       //      jet2_phi = i_calojet->phi();  // UNUSED
00545       jet2_eta = jet.eta();  // now USED
00546       jet2_phi = jet.phi();  // now USED
00547     }
00548     // questo elseif funziona
00549   }
00551 
00552 
00553 
00555   //                 Fill Histograms                                            //
00557 
00558   bool fill_e1  = false;
00559   bool fill_e2  = false;
00560   bool fill_m1  = false;
00561   bool fill_m2  = false;
00562   bool fill_met = false;
00563 
00564   // Was Z->ee found?
00565   if (ee_invMass > 0.0) {
00566     h_ee_invMass->Fill(ee_invMass);
00567     fill_e1 = true;
00568     fill_e2 = true;
00569   }
00570 
00571   // Was Z->mu mu found?
00572   if (mm_invMass > 0.0) {
00573     h_mumu_invMass->Fill(mm_invMass);
00574     fill_m1 = true;
00575     fill_m2 = true;
00576     h_jet2_et ->Fill(jet2_et);
00577   }
00578 
00579   // Was W->e nu found?
00580   if (electron_et > 0.0 && missing_et > 20.0) {
00581     float dphiW  = fabs(met_phi-electron_phi);
00582     float W_mt_e = sqrt(2 * missing_et * electron_et * (1 - cos(dphiW)));
00583     h_e_invWMass->Fill(W_mt_e);
00584     fill_e1  = true;
00585     fill_met = true;
00586   }
00587 
00588   // Was W->mu nu found?
00589   if (muon_pt > 0.0 && missing_et > 20.0) {
00590     float dphiW = fabs(met_phi - muon_phi);
00591     float W_mt_m = sqrt(2 * missing_et * muon_pt * (1 - cos(dphiW)));
00592     h_m_invWMass->Fill(W_mt_m);
00593     fill_m1  = true;
00594     fill_met = true;
00595   }
00596 
00597   if (jet_et > -10.0) {
00598     h_jet_et->Fill(jet_et);
00599     h_jet_count->Fill(jet_count);
00600   }
00601 
00602   if (jet_pt > 0.) {
00603     h_jet_pt->Fill(jet_pt);
00604   }
00605 
00606   if (jet_eta > -50.) {
00607     h_jet_eta->Fill(jet_eta);
00608   }
00609 
00610   if (jet_phi > -10.) {
00611     h_jet_phi->Fill(jet_phi);
00612   }
00613 
00614   if (jet2_et > -10.0) {
00615     h_jet2_et->Fill(jet2_et);
00616   }
00617 
00618   // if (jet2_pt>0.) {
00619   //  h_jet2_pt   ->Fill(jet2_pt);
00620   // }
00621 
00622   if (jet2_eta > -50.) {
00623     h_jet2_eta->Fill(jet2_eta);
00624   }
00625 
00626   if (jet2_phi > -10.) {
00627     h_jet2_phi->Fill(jet2_phi);
00628   }
00629 
00630 
00631 
00632   if (fill_e1 || fill_m1) {
00633     h_vertex_number->Fill(vertex_number);
00634     h_vertex_chi2->Fill(vertex_chi2);
00635     h_vertex_d0->Fill(vertex_d0);
00636     h_vertex_numTrks->Fill(vertex_numTrks);
00637     h_vertex_sumTrks->Fill(vertex_sumTrks);
00638   }
00639 
00640   if (fill_e1) {
00641     h_e1_et->Fill(electron_et);
00642     h_e1_eta->Fill(electron_eta);
00643     h_e1_phi->Fill(electron_phi);
00644   }
00645   if (fill_e2) {
00646     h_e2_et->Fill(electron2_et);
00647     h_e2_eta->Fill(electron2_eta);
00648     h_e2_phi->Fill(electron2_phi);
00649   }
00650   if (fill_m1) {
00651     h_m1_pt->Fill(muon_pt);
00652     h_m1_eta->Fill(muon_eta);
00653     h_m1_phi->Fill(muon_phi);
00654   }
00655   if (fill_m2) {
00656     h_m2_pt->Fill(muon2_pt);
00657     h_m2_eta->Fill(muon2_eta);
00658     h_m2_phi->Fill(muon2_phi);
00659   }
00660   if (fill_met) {
00661     h_met->Fill(missing_et);
00662     h_met_phi->Fill(met_phi);
00663   }
00665 }
00666 
00667 
00668 void EwkDQM::endJob(void) {}
00669 
00670 
00671 // This always returns only a positive deltaPhi
00672 double EwkDQM::calcDeltaPhi(double phi1, double phi2) {
00673   double deltaPhi = phi1 - phi2;
00674 
00675   if (deltaPhi < 0)
00676     deltaPhi = -deltaPhi;
00677 
00678   if (deltaPhi > 3.1415926)
00679     deltaPhi = 2 * 3.1415926 - deltaPhi;
00680 
00681   return deltaPhi;
00682 }