CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/ElectroWeakAnalysis/ZMuMu/plugins/ZMuMuAnalyzer_cynematics.cc

Go to the documentation of this file.
00001 /* \class ZMuMuAnalyzer_cynematics
00002  *
00003  * Z->mu+m- standard analysis for cross section 
00004  * measurements. Take as input the output of the
00005  * standard EWK skim: zToMuMu
00006  * 
00007  * Produces mass spectra and other histograms for
00008  * the samples in input:
00009  *
00010  *  + Z -> mu+mu-, both muons are "global" muons
00011  *  + Z -> mu+mu-, one muons is "global" muons, one unmatched tracks
00012  *  + Z -> mu+mu-, one muons is "global" muons, one unmatched stand-alone muon
00013  * 
00014  *
00015  * \author Michele de Gruttola, 
00016  * \modified by Davide Piccolo, INFN Naples to include gerarchyc selection of Z and histos as a finction of eta pt phi
00017  * 
00018  *
00019  * \id $Id: ZMuMuAnalyzer_cynematics.cc,v 1.6 2010/02/19 02:46:29 wmtan Exp $
00020  *
00021  */
00022 #include "FWCore/Framework/interface/EDAnalyzer.h"
00023 #include "FWCore/Utilities/interface/InputTag.h"
00024 #include "FWCore/Framework/interface/Event.h"
00025 #include "DataFormats/Common/interface/Handle.h"
00026 #include "DataFormats/Candidate/interface/OverlapChecker.h"
00027 #include "DataFormats/Candidate/interface/Particle.h"
00028 #include "DataFormats/Candidate/interface/Candidate.h"
00029 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
00030 #include "DataFormats/Candidate/interface/CandMatchMap.h"
00031 #include "FWCore/ServiceRegistry/interface/Service.h"
00032 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00034 #include "DataFormats/Common/interface/AssociationVector.h"
00035 #include "DataFormats/Candidate/interface/OverlapChecker.h"
00036 #include "DataFormats/Math/interface/LorentzVector.h"
00037 #include "DataFormats/TrackReco/interface/Track.h"
00038 #include "TFile.h"
00039 #include "TH1.h"
00040 #include "TH2.h"
00041 #include <iostream>
00042 #include <iterator>
00043 #include <sstream>
00044 using namespace edm;
00045 using namespace std;
00046 using namespace reco;
00047 
00048 typedef edm::AssociationVector<reco::CandidateRefProd, std::vector<double> > IsolationCollection;
00049 
00050 class ZMuMuAnalyzer_cynematics : public edm::EDAnalyzer {
00051 public:
00052   ZMuMuAnalyzer_cynematics(const edm::ParameterSet& pset);
00053 private:
00054   virtual void analyze(const edm::Event& event, const edm::EventSetup& setup);
00055   bool isContained(const Candidate &, const Candidate &);
00056   virtual void endJob();
00057   
00058   OverlapChecker overlap_;
00059   InputTag zMuMu_;
00060   InputTag zMuTrack_;
00061   InputTag zMuStandAlone_;
00062   InputTag  muIso_, trackIso_, standAloneIso_;
00063   InputTag  zMuMuMap_ ,zMuTrackMap_, zMuStandAloneMap_;
00064   double isocut_, etacut_, ptcut_,ptSTAcut_,  minZmass_, maxZmass_;
00065   TH1D * h_zMuMu_numberOfCand, * h_zMuMu_numberOfCand_passed, * h_zMuMu_numberOfCand_ptpassed,* h_zMuMu_numberOfCand_etapassed,
00066     * h_zMuMu_numberOfCand_masspassed, * h_zMuMu_numberOfCand_isopassed, * h_zMuMu_numberOfCand_ptetapassed, 
00067     * h_zMuMu_numberOfCand_ptetamasspassed, * h_zMuMu_mass_, * h_zMuSingleTrack_mass_, * h_zMuSingleStandAlone_mass_,
00068     * h_zMuSingleStandAloneOverlap_mass_, * h_zMuMuMatched_mass_, * h_zMuSingleTrackMatched_mass_,
00069     * h_zMuSingleStandAloneMatched_mass_, 
00070     * h_zMuSingleStandAloneOverlapMatched_mass_;
00071 
00072   TH1D * h_zMuSta_numberOfCand,* h_zMuSta_numberOfCand_passed,* h_zMuSta_MCmatched_numberOfCand_passed, 
00073     * h_zMuSta_numberOfCand_notcontained, 
00074     * h_zMuTrack_numberOfCand, * h_zMuTrack_numberOfCand_notcontained, * h_zMuTrack_numberOfCand_passed,
00075     * h_zMuTrack_MCmatched_numberOfCand_passed; 
00076   TH2D * h_OneSta_mass;
00077 
00078   double etamin, etamax, phimin, phimax, ptmin, ptmax;
00079   int numberOfIntervals;        // number of intervals in which to divide cynematic variables
00080   double binEta,binPhi, binPt; 
00081   vector<TH1D *>  hmumu_eta, hmusta_eta, hmutrack_eta;
00082   vector<TH1D *>  hmumu_phi, hmusta_phi, hmutrack_phi;
00083   vector<TH1D *>  hmumu_pt, hmusta_pt, hmutrack_pt;
00084 
00085 };
00086 
00087 ZMuMuAnalyzer_cynematics::ZMuMuAnalyzer_cynematics(const edm::ParameterSet& pset) : 
00088   zMuMu_( pset.getParameter<InputTag>( "zMuMu" ) ),
00089   zMuTrack_( pset.getParameter<InputTag>( "zMuTrack" ) ),
00090   zMuStandAlone_( pset.getParameter<InputTag>( "zMuStandAlone" ) ),
00091   muIso_( pset.getParameter<InputTag>( "muIso" ) ),
00092   trackIso_( pset.getParameter<InputTag>( "trackIso" ) ),
00093   standAloneIso_( pset.getParameter<InputTag>( "standAloneIso" ) ),
00094   zMuMuMap_( pset.getParameter<InputTag>( "zMuMuMap" ) ),
00095   zMuTrackMap_( pset.getParameter<InputTag>( "zMuTrackMap" ) ),
00096   zMuStandAloneMap_( pset.getParameter<InputTag>( "zMuStandAloneMap" ) ),
00097   isocut_( pset.getParameter<double>( "isocut" ) ),
00098   etacut_( pset.getParameter<double>( "etacut" ) ),
00099   ptcut_( pset.getParameter<double>( "ptcut" ) ),
00100   ptSTAcut_( pset.getParameter<double>( "ptSTAcut" ) ),  
00101   minZmass_( pset.getParameter<double>( "minZmass" )),
00102   maxZmass_( pset.getParameter<double>( "maxZmass" )) {
00103   
00104   Service<TFileService> fs;
00105   h_zMuMu_numberOfCand = fs->make<TH1D>("ZMuMunumberOfCand","number of ZMuMu cand",10, -.5, 9.5);
00106   h_zMuMu_numberOfCand_passed = fs->make<TH1D>("ZMuMunumberOfCandpassed","number of ZMuMu cand selected",10, -.5, 9.5);
00107   h_zMuMu_numberOfCand_ptpassed = fs->make<TH1D>("ZMuMunumberOfCandptpassed","number of ZMuMu cand after pt cut selected",10, -.5, 9.5);
00108   h_zMuMu_numberOfCand_etapassed = fs->make<TH1D>("ZMuMunumberOfCandetapassed","number of ZMuMu cand after eta cut selected",10, -.5, 9.5);
00109   h_zMuMu_numberOfCand_masspassed = fs->make<TH1D>("ZMuMunumberOfCandmasspassed","number of ZMuMu cand after mass cut selected",10, -.5, 9.5);
00110   h_zMuMu_numberOfCand_isopassed = fs->make<TH1D>("ZMuMunumberOfCandisopassed","number of ZMuMu cand after iso cut selected",10, -.5, 9.5);
00111   h_zMuMu_numberOfCand_ptetapassed = fs->make<TH1D>("ZMuMunumberOfCandptetapassed","number of ZMuMu cand after pt & eta cut selected",10, -.5, 9.5);
00112   h_zMuMu_numberOfCand_ptetamasspassed = fs->make<TH1D>("ZMuMunumberOfCandptetamaspassed","number of ZMuMu cand after pt & eta & mass cut selected",10, -.5, 9.5);
00113 
00114 
00115   h_zMuMu_mass_ = fs->make<TH1D>( "ZMuMumass", "ZMuMu mass(GeV)", 200,  0., 200. );
00116   h_zMuSingleTrack_mass_ = fs->make<TH1D>( "ZMuSingleTrackmass", "ZMuSingleTrack mass(GeV)", 100,  0., 200. );
00117   h_zMuSingleStandAlone_mass_ = fs->make<TH1D>( "ZMuSingleStandAlonemass", "ZMuSingleStandAlone mass(GeV)", 50,  0., 200. );
00118   h_zMuSingleStandAloneOverlap_mass_ = fs->make<TH1D>( "ZMuSingleStandAloneOverlapmass", "ZMuSingleStandAloneOverlap  mass(GeV)", 50,  0., 200. );
00119   
00120   
00121   h_zMuMuMatched_mass_ = fs->make<TH1D>( "ZMuMuMatchedmass", "ZMuMu Matched  mass(GeV)", 200,  0., 200. );
00122   h_zMuSingleTrackMatched_mass_ = fs->make<TH1D>( "ZMuSingleTrackmassMatched", "ZMuSingleTrackMatched mass(GeV)", 100,  0., 200. );
00123   h_zMuSingleStandAloneMatched_mass_ = fs->make<TH1D>( "ZMuSingleStandAlonemassMatched", "ZMuSingleStandAloneMatched mass(GeV)", 50,  0., 200. );
00124   h_zMuSingleStandAloneOverlapMatched_mass_ = fs->make<TH1D>( "ZMuSingleStandAloneOverlapmassMatched", "ZMuSingleStandAloneMatched Overlap  mass(GeV)", 50,  0., 200. );
00125   
00126   h_zMuSta_numberOfCand = fs->make<TH1D>("ZMuStanumberOfCand","number of ZMuSta cand (if ZMuMu not selected)",10, -.5, 9.5);
00127   h_OneSta_mass = fs->make<TH2D>("ZOneMuStaMass","inv. mass of ZMuSta1 vs ZMuSta2 when one ZMuSta has been found (if ZMuMu not selected)",100, 0., 400, 100, 0., 400.);
00128   h_zMuSta_numberOfCand_notcontained = fs->make<TH1D>("ZMuStanumberOfCandnotcontained","number of independent ZMuSta cand (if ZMuMu not selected)",10, -.5, 9.5);
00129   h_zMuSta_numberOfCand_passed = fs->make<TH1D>("ZMuStanumberOfCandpassed","number of ZMuSta cand selected (if ZMuMu not selected)",10, -.5, 9.5);
00130   h_zMuSta_MCmatched_numberOfCand_passed = fs->make<TH1D>("ZMuStaMCmatchedNumberOfCandpassed","number of ZMuSta MC matched cand selected (if ZMuMu not selected)",10, -.5, 9.5);
00131   h_zMuTrack_numberOfCand = fs->make<TH1D>("ZMuTranumberOfCand","number of ZMuTrack cand (if ZMuMu and ZMuSTa not selected)",10, -.5, 9.5);
00132   h_zMuTrack_numberOfCand_notcontained = fs->make<TH1D>("ZMuTranumberOfCandnotcontaind","number of indeendent ZMuTrack cand (if ZMuMu and ZMuSTa not selected)",10, -.5, 9.5);
00133   h_zMuTrack_numberOfCand_passed = fs->make<TH1D>("ZMuTranumberOfCandpassed","number of ZMuTrack cand selected (if ZMuMu and ZMuSta not selected)",10, -.5, 9.5);
00134   h_zMuTrack_MCmatched_numberOfCand_passed = fs->make<TH1D>("ZMuTraMCmacthedNumberOfCandpassed","number of ZMuTrack MC matched cand selected (if ZMuMu and ZMuSta not selected)",10, -.5, 9.5);
00135 
00136 
00137   // creating histograms for each Pt, eta, phi interval
00138 
00139   etamin = -etacut_;
00140   etamax = etacut_;
00141   phimin = -3.1415;
00142   phimax = 3.1415;
00143   ptmin = ptcut_;
00144   ptmax = 100;
00145   numberOfIntervals = 8;        // number of intervals in which to divide cynematic variables
00146   binEta = (etamax - etamin)/numberOfIntervals; 
00147   binPhi = (phimax - phimin)/numberOfIntervals; 
00148   binPt = (ptmax - ptmin)/numberOfIntervals; 
00149   TFileDirectory etaDirectory = fs->mkdir("etaIntervals");   // in this directory will be saved all the histos of different eta intervals
00150   TFileDirectory phiDirectory = fs->mkdir("phiIntervals");   // in this directory will be saved all the histos of different phi intervals
00151   TFileDirectory ptDirectory = fs->mkdir("ptIntervals");   // in this directory will be saved all the histos of different pt intervals
00152 
00153   // eta histograms creation
00154 
00155   for (int i=0;i<numberOfIntervals;i++) {
00156     double range0 = etamin + i*binEta;
00157     double range1= range0 + binEta;
00158     char a[30], b[50];
00159     sprintf(a,"zmumu_etaRange%d",i);
00160     sprintf(b,"zmumu mass eta Range %f to %f",range0,range1);
00161     hmumu_eta.push_back(etaDirectory.make<TH1D>(a,b,200,0.,200.));
00162     char asta[30], bsta[50];
00163     sprintf(asta,"zmusta_etaRange%d",i);
00164     sprintf(bsta,"zmusta mass eta Range %f to %f",range0,range1);
00165     hmusta_eta.push_back(etaDirectory.make<TH1D>(asta,bsta,50,0.,200.));
00166     char atk[30], btk[50];
00167     sprintf(atk,"zmutrack_etaRange%d",i);
00168     sprintf(btk,"zmutrack mass eta Range %f to %f",range0,range1);
00169     hmutrack_eta.push_back(etaDirectory.make<TH1D>(atk,btk,100,0.,200.));
00170   } 
00171 
00172   // phi histograms creation
00173 
00174   for (int i=0;i<numberOfIntervals;i++) {
00175     double range0 = phimin + i*binPhi;
00176     double range1= range0 + binPhi;
00177     char a[30], b[50];
00178     sprintf(a,"zmumu_phiRange%d",i);
00179     sprintf(b,"zmumu mass phi Range %f to %f",range0,range1);
00180     hmumu_phi.push_back(phiDirectory.make<TH1D>(a,b,200,0.,200.));
00181     char asta[30], bsta[50];
00182     sprintf(asta,"zmusta_phiRange%d",i);
00183     sprintf(bsta,"zmusta mass phi Range %f to %f",range0,range1);
00184     hmusta_phi.push_back(phiDirectory.make<TH1D>(asta,bsta,50,0.,200.));
00185     char atk[30], btk[50];
00186     sprintf(atk,"zmutrack_phiRange%d",i);
00187     sprintf(btk,"zmutrack mass phi Range %f to %f",range0,range1);
00188     hmutrack_phi.push_back(phiDirectory.make<TH1D>(atk,btk,100,0.,200.));
00189   } 
00190 
00191   // pt histograms creation
00192 
00193   for (int i=0;i<numberOfIntervals;i++) {
00194     double range0 = ptmin + i*binPt;
00195     double range1= range0 + binPt;
00196     char a[30], b[50];
00197     sprintf(a,"zmumu_ptRange%d",i);
00198     sprintf(b,"zmumu mass pt Range %f to %f",range0,range1);
00199     hmumu_pt.push_back(ptDirectory.make<TH1D>(a,b,200,0.,200.));
00200     char asta[30], bsta[50];
00201     sprintf(asta,"zmusta_ptRange%d",i);
00202     sprintf(bsta,"zmusta mass pt Range %f to %f",range0,range1);
00203     hmusta_pt.push_back(ptDirectory.make<TH1D>(asta,bsta,50,0.,200.));
00204     char atk[30], btk[50];
00205     sprintf(atk,"zmutrack_ptRange%d",i);
00206     sprintf(btk,"zmutrack mass pt Range %f to %f",range0,range1);
00207     hmutrack_pt.push_back(ptDirectory.make<TH1D>(atk,btk,100,0.,200.));
00208   } 
00209  }
00210 
00211 void ZMuMuAnalyzer_cynematics::analyze(const edm::Event& event, const edm::EventSetup& setup) {
00212   Handle<CandidateCollection> zMuMu;
00213   event.getByLabel(zMuMu_, zMuMu);
00214   Handle<CandidateCollection> zMuTrack;
00215   event.getByLabel( zMuTrack_, zMuTrack );
00216   Handle<CandidateCollection> zMuStandAlone;
00217   event.getByLabel( zMuStandAlone_, zMuStandAlone );  
00218 
00219   unsigned int nZMuMu = zMuMu->size();
00220   unsigned int nZTrackMu = zMuTrack->size();
00221   unsigned int nZStandAloneMu = zMuStandAlone->size();
00222   //  static const double zMass = 91.1876; // PDG Z mass
00223 
00224   cout << "++++++++++++++++++++++++++" << endl;
00225   cout << "nZMuMu = " << nZMuMu << endl;
00226   cout << "nZTrackMu = " << nZTrackMu << endl;
00227   cout << "nZStandAloneMu = " << nZStandAloneMu << endl;
00228   cout << "++++++++++++++++++++++++++" << endl;
00229 
00230   // ZMuMu counters
00231 
00232   int ZMuMu_passed = 0;
00233   int ZMuMu_ptcut_counter = 0;
00234   int ZMuMu_etacut_counter = 0;
00235   int ZMuMu_masscut_counter = 0;
00236   int ZMuMu_isocut_counter = 0;
00237   int ZMuMu_ptetacut_counter = 0;
00238   int ZMuMu_ptetamasscut_counter = 0;
00239   int ZMuMu_allcut_counter = 0;
00240 
00241   // ZMuTrack counters
00242 
00243   int ZMuTrack_passed = 0;
00244   int ZMuTrack_notcontained = 0;
00245   int ZMuTrack_MCmatched_passed = 0;
00246 
00247   // ZMuStandalone counters
00248   int ZMuStandalone_notcontained = 0;
00249   int ZMuStandalone_passed = 0;
00250   int ZMuStandalone_MCmatched_passed = 0;
00251 
00252   Handle<CandMatchMap> zMuMuMap;
00253   if( nZMuMu > 0 ) {
00254     event.getByLabel(zMuMuMap_, zMuMuMap);
00255   }
00256 
00257   Handle<CandMatchMap> zMuTrackMap;
00258   if( nZTrackMu > 0 ) {
00259     event.getByLabel( zMuTrackMap_, zMuTrackMap );
00260   }
00261 
00262   Handle<CandMatchMap> zMuStandAloneMap;
00263   if( nZStandAloneMu > 0 ) {
00264     event.getByLabel( zMuStandAloneMap_, zMuStandAloneMap );  
00265   }    
00266 
00267   Handle<IsolationCollection> muIso;
00268   event.getByLabel(muIso_, muIso);
00269   ProductID muIsoId = muIso->keyProduct().id();
00270   Handle<IsolationCollection> trackIso;
00271   event.getByLabel(trackIso_, trackIso);
00272   ProductID trackIsoId = trackIso->keyProduct().id();
00273   
00274   Handle<IsolationCollection> standAloneIso;
00275   event.getByLabel(standAloneIso_, standAloneIso);
00276   ProductID standAloneIsoId = standAloneIso->keyProduct().id();
00277 
00278   if (nZMuMu > 0) {
00279     // double mass = 1000000.;
00280     for( unsigned int i = 0; i < nZMuMu; i++ ) {
00281       bool ptcutAccept = false;
00282       bool etacutAccept = false;
00283       bool masscutAccept = false;
00284       bool isocutAccept = false;
00285       const Candidate & zmmCand = (*zMuMu)[ i ];
00286       CandidateRef CandRef(zMuMu,i);
00287       CandidateRef lep1 = zmmCand.daughter( 0 )->masterClone().castTo<CandidateRef>();
00288       CandidateRef lep2 = zmmCand.daughter( 1 )->masterClone().castTo<CandidateRef>();
00289       
00290       const  double iso1 = muIso->value( lep1.key() );  
00291       const  double iso2 = muIso->value( lep2.key() );  
00292       
00293       double m = zmmCand.mass();
00294       // check single cuts
00295 
00296       if (lep1->pt()>ptcut_ && lep2->pt()>ptcut_) ptcutAccept = true;
00297       if (fabs(lep1->eta())<etacut_ && fabs(lep2->eta())<etacut_) etacutAccept = true;
00298       if (m>minZmass_ && m<maxZmass_) masscutAccept = true;
00299       if (iso1 < isocut_ && iso2 <isocut_) isocutAccept = true;
00300 
00301 
00302       if (ptcutAccept) ZMuMu_ptcut_counter++;
00303       if (etacutAccept) ZMuMu_etacut_counter++;
00304       if (masscutAccept) ZMuMu_masscut_counter++;
00305       if (isocutAccept) ZMuMu_isocut_counter++;
00306 
00307       // check sequencial cuts
00308 
00309       if (ptcutAccept && etacutAccept) {
00310         ZMuMu_ptetacut_counter++;
00311         if (masscutAccept) {
00312           ZMuMu_ptetamasscut_counter++;
00313           if (isocutAccept) {
00314             ZMuMu_passed++;}
00315         }
00316       }
00317 
00318       if (ptcutAccept && etacutAccept && masscutAccept && isocutAccept)  {
00319         ZMuMu_allcut_counter++;
00320         h_zMuMu_mass_->Fill( m );         
00321 
00322         // check the cynematics to fill correct histograms
00323         for (int j=0;j<numberOfIntervals;j++) {
00324           bool statusBinEta = false;
00325           bool statusBinPhi = false;
00326           bool statusBinPt  = false;
00327           double range0 = etamin + j*binEta;
00328           double range1= range0 + binEta;
00329           double range0phi = phimin + j*binPhi;
00330           double range1phi= range0phi + binPhi;
00331           double range0pt = ptmin + j*binPt;
00332           double range1pt = range0pt + binPt;
00333           // eta histograms
00334           if (lep1->eta()>=range0 && lep1->eta()<range1)
00335             {
00336               hmumu_eta[j]->Fill(m);
00337               statusBinEta = true;
00338             }
00339           if (lep2->eta()>=range0 && lep2->eta()<range1 && !statusBinEta){
00340             hmumu_eta[j]->Fill(m);                               // If eta1 is in the same bin of eta2 fill just once
00341           }
00342           // phi histograms
00343           if (lep1->phi()>=range0phi && lep1->phi()<range1phi)
00344             {
00345               hmumu_phi[j]->Fill(m);
00346               statusBinPhi = true;
00347             }
00348           if (lep2->phi()>=range0phi && lep2->phi()<range1phi && !statusBinPhi){
00349             hmumu_phi[j]->Fill(m);                               // If phi1 is in the same bin of phi2 fill just once
00350           }
00351           // pt histograms
00352           if (lep1->pt()>=range0pt && lep1->pt()<range1pt)
00353             {
00354               hmumu_pt[j]->Fill(m);
00355               statusBinPt = true;
00356             }
00357           if (lep2->pt()>=range0pt && lep2->pt()<range1pt && !statusBinPt){
00358             hmumu_pt[j]->Fill(m);                               // If pt1 is in the same bin of pt2 fill just once
00359           }
00360         }
00361         
00362         CandMatchMap::const_iterator m0 = zMuMuMap->find(CandRef);
00363         if( m0 != zMuMuMap->end()) {                                            // the Z is matched to MC thruth
00364             h_zMuMuMatched_mass_->Fill( m );
00365         }
00366       }
00367     }
00368   }
00369 
00370   h_zMuMu_numberOfCand->Fill(nZMuMu);                                             // number of Z cand found per event
00371   h_zMuMu_numberOfCand_passed->Fill(ZMuMu_allcut_counter);                        // number of Z cand after all cuts found per event
00372   h_zMuMu_numberOfCand_ptpassed->Fill(ZMuMu_ptcut_counter);                       // number of Z cand afer pt cut found per event
00373   h_zMuMu_numberOfCand_etapassed->Fill(ZMuMu_etacut_counter);                     // number of Z cand afer eta cut found per event
00374   h_zMuMu_numberOfCand_masspassed->Fill(ZMuMu_masscut_counter);                   // number of Z cand afer mass cut found per event
00375   h_zMuMu_numberOfCand_isopassed->Fill(ZMuMu_isocut_counter);                     // number of Z cand afer iso cut found per event
00376   h_zMuMu_numberOfCand_ptetapassed->Fill(ZMuMu_ptetacut_counter);                 // number of Z cand afer pt&eta cut found per event
00377   h_zMuMu_numberOfCand_ptetamasspassed->Fill(ZMuMu_ptetamasscut_counter);         // number of Z cand afer pt&eta&mass cut found per event
00378 
00379   
00380   //ZmuSingleStandAlone (check MuStandalone if MuMu has not been selected by cuts)
00381   //  cout << "ZMuMuanalyzer : n of zMuMu " << nZMuMu << " passed " << ZMuMu_passed << "     n. of zStaMu " << nZStandAloneMu << endl;
00382 
00383   if (ZMuMu_passed == 0 && nZStandAloneMu>0 ) {
00384     //      unsigned int index = 1000;
00385     for( unsigned int j = 0; j < nZStandAloneMu; j++ ) {
00386       const Candidate & zsmCand = (*zMuStandAlone)[ j ];
00387       bool skipZ = false;
00388       for( unsigned int i = 0; i < nZMuMu; i++ ) {              // chek if the ZMuSTandalone is contained in a ZMuMu 
00389         const Candidate & zmmCand = (*zMuMu)[ i ];        // if yes .. the event has to be skipped
00390         if (isContained(zmmCand,zsmCand)) skipZ=true;
00391       }
00392       if (!skipZ) {                                       // ZSMuSTandalone not contained in a ZMuMu
00393         ZMuStandalone_notcontained++;
00394         CandidateRef CandRef(zMuStandAlone,j);
00395         CandidateRef lep1 = zsmCand.daughter( 0 )->masterClone().castTo<CandidateRef>();
00396         CandidateRef lep2 = zsmCand.daughter( 1 )->masterClone().castTo<CandidateRef>();
00397         
00398         ProductID id1 = lep1.id();
00399         ProductID id2 = lep2.id();
00400         double iso1 = -1;
00401         double iso2 = -1;
00402         
00403         if( id1 == muIsoId ) 
00404           iso1 = muIso->value( lep1.key() );    
00405         else if ( id1 == standAloneIsoId )
00406           iso1 = standAloneIso->value( lep1.key() );    
00407         
00408         if( id2 == muIsoId ) 
00409           iso2 = muIso->value( lep2.key() );    
00410         else if ( id2 == standAloneIsoId )
00411           iso2 = standAloneIso->value( lep2.key() );    
00412         
00413         double ms = zsmCand.mass();
00414         if (lep1->pt()>ptSTAcut_ && lep2->pt()>ptSTAcut_ &&    
00415             fabs(lep1->eta())<etacut_ && fabs(lep2->eta())<etacut_ &&
00416             ms>minZmass_ && ms<maxZmass_ && iso1<isocut_ && iso2 <isocut_) {
00417           h_zMuSingleStandAlone_mass_->Fill( ms );
00418           ZMuStandalone_passed++;
00419           // check the cynematics to fill correct histograms
00420           for (int j=0;j<numberOfIntervals;j++) {
00421             double range0 = etamin + j*binEta;
00422             double range1= range0 + binEta;
00423             double range0phi = phimin + j*binPhi;
00424             double range1phi= range0phi + binPhi;
00425             double range0pt = ptmin + j*binPt;
00426             double range1pt = range0pt + binPt;
00427 
00428             // check which muon is a standalone (standalone means that there is a reference missing.)
00429             if ((lep1->get<TrackRef,reco::StandAloneMuonTag>()).isNull()) 
00430               {
00431                 if (lep1->eta()>=range0 && lep1->eta()<range1)          hmusta_eta[j]->Fill(ms);
00432                 if (lep1->phi()>=range0phi && lep1->phi()<range1phi)    hmusta_phi[j]->Fill(ms);
00433                 if (lep1->pt()>=range0pt && lep1->pt()<range1pt)        hmusta_pt[j]->Fill(ms);
00434               }
00435             if ((lep2->get<TrackRef,reco::StandAloneMuonTag>()).isNull()) 
00436               {
00437                 if (lep2->eta()>=range0 && lep2->eta()<range1)          hmusta_eta[j]->Fill(ms);
00438                 if (lep2->phi()>=range0phi && lep2->phi()<range1phi)    hmusta_phi[j]->Fill(ms);
00439                 if (lep2->pt()>=range0pt && lep2->pt()<range1pt)        hmusta_pt[j]->Fill(ms);
00440               }
00441 
00442           }
00443           CandMatchMap::const_iterator m0 = zMuStandAloneMap->find(CandRef);
00444           if( m0 != zMuStandAloneMap->end()) {
00445            ZMuStandalone_MCmatched_passed++;
00446             h_zMuSingleStandAloneMatched_mass_->Fill( ms );
00447           }
00448         }
00449       }
00450     }
00451     h_zMuSta_numberOfCand->Fill(nZStandAloneMu);                    // number of ZMuStandalone cand found per event (no higher priority Z selected)
00452     h_zMuSta_numberOfCand_notcontained->Fill(ZMuStandalone_notcontained);
00453     h_zMuSta_numberOfCand_passed->Fill(ZMuStandalone_passed);        // number of ZMuSTa cand after all cuts found per event (no higher prioriy Z selected)
00454     h_zMuSta_MCmatched_numberOfCand_passed->Fill(ZMuStandalone_MCmatched_passed);   // number of ZMuSTa MC matched cand after all cuts found per event (no higher prioriy Z selected)
00455 
00456   }
00457 
00458   //ZmuSingleTRack  (check MuTrack if MuMu has not been selected)
00459   if (ZMuMu_passed == 0 && ZMuStandalone_passed == 0 && nZTrackMu>0) {
00460     for( unsigned int j = 0; j < nZTrackMu; j++ ) {
00461       const Candidate & ztmCand = (*zMuTrack)[ j ];
00462       bool skipZ = false;
00463       for( unsigned int i = 0; i < nZMuMu; i++ ) {              // chek if the ZMuTrack is contained in a ZMuMu 
00464         const Candidate & zmmCand = (*zMuMu)[ i ];        // if yes .. the event has to be skipped
00465         if (isContained(zmmCand,ztmCand)) skipZ=true;
00466       }
00467       if (!skipZ) {
00468         ZMuTrack_notcontained++;
00469         CandidateRef CandRef(zMuTrack,j);
00470         CandidateRef lep1 = ztmCand.daughter( 0 )->masterClone().castTo<CandidateRef>();
00471         CandidateRef lep2 = ztmCand.daughter( 1 )->masterClone().castTo<CandidateRef>();
00472         
00473         ProductID id1 = lep1.id();
00474         ProductID id2 = lep2.id();
00475         double iso1 = -1;
00476         double iso2 = -1;
00477         
00478         if( id1 == muIsoId ) 
00479           iso1 = muIso->value( lep1.key() );    
00480         else if ( id1 == trackIsoId )
00481           iso1 = trackIso->value( lep1.key() ); 
00482         
00483         if( id2 == muIsoId ) 
00484           iso2 = muIso->value( lep2.key() );    
00485         else if ( id2 == trackIsoId )
00486           iso2 = trackIso->value( lep2.key() ); 
00487         
00488         double mt = ztmCand.mass();
00489         if (lep1->pt()>ptcut_ && lep2->pt()>ptcut_ &&    
00490             fabs(lep1->eta())<etacut_ && fabs(lep2->eta())<etacut_ &&
00491             mt>minZmass_ && mt<maxZmass_ && iso1<isocut_ && iso2 <isocut_) {
00492           h_zMuSingleTrack_mass_->Fill( mt );
00493           ZMuTrack_passed++;            
00494 
00495           // check the cynematics to fill correct histograms
00496           for (int j=0;j<numberOfIntervals;j++) {
00497             double range0 = etamin + j*binEta;
00498             double range1= range0 + binEta;
00499             double range0phi = phimin + j*binPhi;
00500             double range1phi= range0phi + binPhi;
00501             double range0pt = ptmin + j*binPt;
00502             double range1pt = range0pt + binPt;
00503 
00504             // check which muon is a track only (track only means that there is a reference missing.)
00505             if ((lep1->get<TrackRef,reco::StandAloneMuonTag>()).isNull()) 
00506               {
00507                 if (lep1->eta()>=range0 && lep1->eta()<range1)          hmutrack_eta[j]->Fill(mt);
00508                 if (lep1->phi()>=range0phi && lep1->phi()<range1phi)    hmutrack_phi[j]->Fill(mt);
00509                 if (lep1->pt()>=range0pt && lep1->pt()<range1pt)        hmutrack_pt[j]->Fill(mt);
00510               }
00511             if ((lep2->get<TrackRef,reco::StandAloneMuonTag>()).isNull()) 
00512               {
00513                 if (lep2->eta()>=range0 && lep2->eta()<range1)          hmutrack_eta[j]->Fill(mt);
00514                 if (lep2->phi()>=range0phi && lep2->phi()<range1phi)    hmutrack_phi[j]->Fill(mt);
00515                 if (lep2->pt()>=range0pt && lep2->pt()<range1pt)        hmutrack_pt[j]->Fill(mt);
00516               }
00517           }
00518           CandMatchMap::const_iterator m0 = zMuTrackMap->find(CandRef);
00519           if( m0 != zMuTrackMap->end()) {
00520             ZMuTrack_MCmatched_passed++;
00521             h_zMuSingleTrackMatched_mass_->Fill( mt );
00522           }
00523         }
00524       }
00525     }
00526     h_zMuTrack_numberOfCand->Fill(nZTrackMu);                     // number of ZMuTrack cand found per event (no higher priority Z selected)
00527     h_zMuTrack_numberOfCand_notcontained->Fill(ZMuTrack_notcontained); // number of ZMuTrack cand not cntained in ZMuMu (no higher priority Z selected)
00528 
00529     h_zMuTrack_numberOfCand_passed->Fill(ZMuTrack_passed);        // number of ZMuTrack cand after all cuts found per event (no higher priority Z selected)
00530 
00531     h_zMuTrack_MCmatched_numberOfCand_passed->Fill(ZMuTrack_MCmatched_passed);
00532  
00533   }   
00534 }
00535 
00536 bool ZMuMuAnalyzer_cynematics::isContained(const Candidate & obj1, const Candidate & obj2)
00537 {  
00538   // check if a candidate obj2 is different from obj1  (assume that obj1 is a ZMuMu and obj2 is any other type)
00539   // (for example a Z can be done with two global muons, or with a global muon plus a standalone muon.
00540   // if the standalone muon is part of the second global muon in fact this is the same Z)
00541 
00542   const int maxd = 10;
00543   const Candidate * daughters1[maxd];
00544   const Candidate * daughters2[maxd];
00545   TrackRef trackerTrack1[maxd];
00546   TrackRef stAloneTrack1[maxd];
00547   TrackRef globalTrack1[maxd];
00548   TrackRef trackerTrack2[maxd];
00549   TrackRef stAloneTrack2[maxd];
00550   TrackRef globalTrack2[maxd];
00551   bool flag;
00552   unsigned int nd1 = obj1.numberOfDaughters();
00553   unsigned int nd2 = obj2.numberOfDaughters();
00554   unsigned int matched=0;
00555 
00556   for( unsigned int i = 0; i < nd1; ++ i ) {
00557     daughters1[i] = obj1.daughter( i );
00558     trackerTrack1[i] = daughters1[i]->get<TrackRef>();
00559     stAloneTrack1[i] = daughters1[i]->get<TrackRef,reco::StandAloneMuonTag>();
00560     globalTrack1[i]  = daughters1[i]->get<TrackRef,reco::CombinedMuonTag>();
00561 
00562     /*********************************************** just used for debug ********************
00563     if (trackerTrack1[i].isNull()) 
00564       cout << "in ZMuMu daughter " << i << " tracker ref non found " << endl;
00565     else
00566       cout << "in ZMuMu daughter " << i << " tracker ref FOUND" 
00567            << " id: " << trackerTrack1[i].id() << ", index: " << trackerTrack1[i].key() 
00568            << endl;
00569     if (stAloneTrack1[i].isNull()) 
00570       cout << "in ZMuMu daughter " << i << " stalone ref non found " << endl;
00571     else
00572       cout << "in ZMuMu daughter " << i << " stalone ref FOUND" 
00573            << " id: " << stAloneTrack1[i].id() << ", index: " << stAloneTrack1[i].key() 
00574            << endl;
00575     
00576     if (globalTrack1[i].isNull()) 
00577       cout << "in ZMuMu daughter " << i << " global ref non found " << endl;
00578     else
00579       cout << "in ZMuMu daughter " << i << " global ref FOUND"  
00580            << " id: " << globalTrack1[i].id() << ", index: " << globalTrack1[i].key() 
00581            << endl;
00582     */
00583   }
00584   for( unsigned int i = 0; i < nd2; ++ i ) {
00585     daughters2[i] = obj2.daughter( i );
00586     trackerTrack2[i] = daughters2[i]->get<TrackRef>();
00587     stAloneTrack2[i] = daughters2[i]->get<TrackRef,reco::StandAloneMuonTag>();
00588     globalTrack2[i]  = daughters2[i]->get<TrackRef,reco::CombinedMuonTag>();
00589 
00590     /******************************************** just used for debug ************
00591     if (trackerTrack2[i].isNull()) 
00592       cout << "in ZMuSta daughter " << i << " tracker ref non found " << endl;
00593     else
00594       cout << "in ZMuSta daughter " << i << " tracker ref FOUND"  
00595            << " id: " << trackerTrack2[i].id() << ", index: " << trackerTrack2[i].key() 
00596            << endl;
00597     if (stAloneTrack2[i].isNull()) 
00598       cout << "in ZMuSta daughter " << i << " standalone ref non found " << endl;
00599     else
00600       cout << "in ZMuSta daughter " << i << " standalone ref FOUND" 
00601            << " id: " << stAloneTrack2[i].id() << ", index: " << stAloneTrack2[i].key() 
00602            << endl;
00603     
00604     if (globalTrack2[i].isNull()) 
00605       cout << "in ZMuSta daughter " << i << " global ref non found " << endl;
00606     else
00607       cout << "in ZMuSta daughter " << i << " global ref FOUND" 
00608            << " id: " << globalTrack2[i].id() << ", index: " << globalTrack2[i].key() 
00609            << endl;
00610            
00611    */  
00612   }
00613   if (nd1 != nd2) 
00614     {
00615       cout << "ZMuMuAnalyzer::isContained WARNING n.of daughters different " << nd1 << "  " << nd2 << endl;
00616     }
00617   else
00618     {
00619       for (unsigned int i = 0; i < nd1; i++) {
00620         flag = false;
00621         for (unsigned int j = 0; j < nd2; j++) {           // if the obj2 is a standalone the trackref is alwais in the trackerTRack position
00622           if ( ((trackerTrack2[i].id()==trackerTrack1[j].id()) && (trackerTrack2[i].key()==trackerTrack1[j].key())) ||
00623                ((trackerTrack2[i].id()==stAloneTrack1[j].id()) && (trackerTrack2[i].key()==stAloneTrack1[j].key())) ) {
00624             flag = true;
00625           }
00626         }
00627         if (flag) matched++;
00628       }
00629     }
00630   if (matched==nd1) // return true if all the childrens of the ZMuMu have a children matched in ZMuXX
00631     return true;
00632   else 
00633     return false;
00634 }
00635    
00636 void ZMuMuAnalyzer_cynematics::endJob() {
00637 
00638   // candidate analysis
00639   // ZMuMu
00640   double Nzmmc = h_zMuMu_numberOfCand->GetEntries();
00641   double Nzmmc_0Z = h_zMuMu_numberOfCand->GetBinContent(1);
00642   double Nzmmc_1Z = h_zMuMu_numberOfCand->GetBinContent(2);
00643   double Nzmmc_moreZ = Nzmmc-Nzmmc_0Z-Nzmmc_1Z;
00644   double Nzmmc_passed_0Z = h_zMuMu_numberOfCand_passed->GetBinContent(1);  
00645   double Nzmmc_passed_1Z = h_zMuMu_numberOfCand_passed->GetBinContent(2);  
00646   double Nzmmc_passed_moreZ = Nzmmc-Nzmmc_passed_0Z-Nzmmc_passed_1Z;
00647   double Nzmmc_ptpassed_0Z = h_zMuMu_numberOfCand_ptpassed->GetBinContent(1);
00648   double Nzmmc_ptpassed_1Z = h_zMuMu_numberOfCand_ptpassed->GetBinContent(2);
00649   double Nzmmc_etapassed_0Z = h_zMuMu_numberOfCand_etapassed->GetBinContent(1);
00650   double Nzmmc_etapassed_1Z = h_zMuMu_numberOfCand_etapassed->GetBinContent(2);
00651   double Nzmmc_masspassed_0Z = h_zMuMu_numberOfCand_masspassed->GetBinContent(1);
00652   double Nzmmc_masspassed_1Z = h_zMuMu_numberOfCand_masspassed->GetBinContent(2);
00653   double Nzmmc_isopassed_0Z = h_zMuMu_numberOfCand_isopassed->GetBinContent(1);
00654   double Nzmmc_isopassed_1Z = h_zMuMu_numberOfCand_isopassed->GetBinContent(2);
00655   double Nzmmc_ptetapassed_0Z = h_zMuMu_numberOfCand_ptetapassed->GetBinContent(1);
00656   double Nzmmc_ptetapassed_1Z = h_zMuMu_numberOfCand_ptetapassed->GetBinContent(2);
00657   double Nzmmc_ptetamasspassed_0Z = h_zMuMu_numberOfCand_ptetamasspassed->GetBinContent(1);
00658   double Nzmmc_ptetamasspassed_1Z = h_zMuMu_numberOfCand_ptetamasspassed->GetBinContent(2);
00659   double Nzmmc_ptpassed_moreZ = Nzmmc-Nzmmc_ptpassed_0Z-Nzmmc_ptpassed_1Z;
00660   double Nzmmc_etapassed_moreZ = Nzmmc-Nzmmc_etapassed_0Z-Nzmmc_etapassed_1Z;
00661   double Nzmmc_masspassed_moreZ = Nzmmc-Nzmmc_masspassed_0Z-Nzmmc_masspassed_1Z;
00662   double Nzmmc_isopassed_moreZ = Nzmmc-Nzmmc_isopassed_0Z-Nzmmc_isopassed_1Z;
00663   double Nzmmc_ptetapassed_moreZ = Nzmmc-Nzmmc_ptetapassed_0Z-Nzmmc_ptetapassed_1Z;
00664   double Nzmmc_ptetamasspassed_moreZ = Nzmmc-Nzmmc_ptetamasspassed_0Z-Nzmmc_ptetamasspassed_1Z;
00665   double Nzmsc = h_zMuSta_numberOfCand->GetEntries();
00666   double Nzmsc_0Z = h_zMuSta_numberOfCand->GetBinContent(1);
00667   double Nzmsc_1Z = h_zMuSta_numberOfCand->GetBinContent(2);
00668   double Nzmsc_moreZ = Nzmsc - Nzmsc_0Z - Nzmsc_1Z;
00669   double Nzmsc_notcontained_0Z = h_zMuSta_numberOfCand_notcontained->GetBinContent(1);
00670   double Nzmsc_notcontained_1Z = h_zMuSta_numberOfCand_notcontained->GetBinContent(2);
00671   double Nzmsc_notcontained_moreZ = Nzmsc-Nzmsc_notcontained_0Z-Nzmsc_notcontained_1Z;
00672   double Nzmsc_passed_0Z = h_zMuSta_numberOfCand_passed->GetBinContent(1);
00673   double Nzmsc_passed_1Z = h_zMuSta_numberOfCand_passed->GetBinContent(2);
00674   double Nzmsc_passed_moreZ = Nzmsc - Nzmsc_passed_0Z - Nzmsc_passed_1Z;
00675   double Nzmsc_MCmatched_passed_0Z = h_zMuSta_MCmatched_numberOfCand_passed->GetBinContent(1);
00676   double Nzmsc_MCmatched_passed_1Z = h_zMuSta_MCmatched_numberOfCand_passed->GetBinContent(2);
00677   double Nzmsc_MCmatched_passed_moreZ = Nzmsc - Nzmsc_MCmatched_passed_0Z - Nzmsc_MCmatched_passed_1Z;
00678   double Nzmtc = h_zMuTrack_numberOfCand->GetEntries();
00679   double Nzmtc_0Z = h_zMuTrack_numberOfCand->GetBinContent(1);
00680   double Nzmtc_1Z = h_zMuTrack_numberOfCand->GetBinContent(2);
00681   double Nzmtc_moreZ = Nzmtc - Nzmtc_0Z - Nzmtc_1Z;
00682   double Nzmtc_notcontained_0Z = h_zMuTrack_numberOfCand_notcontained->GetBinContent(1);
00683   double Nzmtc_notcontained_1Z = h_zMuTrack_numberOfCand_notcontained->GetBinContent(2);
00684   double Nzmtc_notcontained_moreZ = Nzmtc-Nzmtc_notcontained_0Z-Nzmtc_notcontained_1Z;
00685   double Nzmtc_passed_0Z = h_zMuTrack_numberOfCand_passed->GetBinContent(1);
00686   double Nzmtc_passed_1Z = h_zMuTrack_numberOfCand_passed->GetBinContent(2);
00687   double Nzmtc_passed_moreZ = Nzmtc - Nzmtc_passed_0Z - Nzmtc_passed_1Z;
00688   double Nzmtc_MCmatched_passed_0Z = h_zMuTrack_MCmatched_numberOfCand_passed->GetBinContent(1);
00689   double Nzmtc_MCmatched_passed_1Z = h_zMuTrack_MCmatched_numberOfCand_passed->GetBinContent(2);
00690   double Nzmtc_MCmatched_passed_moreZ = Nzmtc - Nzmtc_MCmatched_passed_0Z - Nzmtc_MCmatched_passed_1Z;
00691 
00692   cout << "--------------- Statistics ----------------------------------------------------------" << endl;
00693   cout << "n of ZMuMu entries   ...................................................... " << Nzmmc << endl;
00694   cout << "n of ZMuMu events with 0 cand ............................................. " << Nzmmc_0Z << endl;
00695   cout << "n of ZMuMu events with 1 cand ............................................. " << Nzmmc_1Z << endl;
00696   cout << "n of ZMuMu events with 2 or more cand ..................................... " << Nzmmc_moreZ << endl << endl ;
00697 
00698   cout << "n of ZMuMu events not selected by cuts .................................... " << Nzmmc_passed_0Z << endl;
00699   cout << "n of ZMuMu events with 1 cand selected by cuts ............................ " << Nzmmc_passed_1Z << endl;
00700   cout << "n of ZMuMu events with 2 or more cand elected by cuts ..................... " << Nzmmc_passed_moreZ << endl<< endl ;
00701 
00702   cout << "n of ZMuMu events not selected by pt cut .................................. " << Nzmmc_ptpassed_0Z << endl;
00703   cout << "n of ZMuMu events with 1 cand selected by pt cut .......................... " << Nzmmc_ptpassed_1Z << endl;
00704   cout << "n of ZMuMu events with 2 or more cand elected by pt cut ................... " << Nzmmc_ptpassed_moreZ << endl<< endl ;
00705 
00706   cout << "n of ZMuMu events not selected by eta cut ................................. " << Nzmmc_etapassed_0Z << endl;
00707   cout << "n of ZMuMu events with 1 cand selected by eta cut ......................... " << Nzmmc_etapassed_1Z << endl;
00708   cout << "n of ZMuMu events with 2 or more cand elected by eta cut .................. " << Nzmmc_etapassed_moreZ << endl<< endl ;
00709 
00710   cout << "n of ZMuMu events not selected by mass cut ................................ " << Nzmmc_masspassed_0Z << endl;
00711   cout << "n of ZMuMu events with 1 cand selected by mass cut ........................ " << Nzmmc_masspassed_1Z << endl;
00712   cout << "n of ZMuMu events with 2 or more cand elected by mass cut ................. " << Nzmmc_masspassed_moreZ << endl<< endl ;
00713 
00714   cout << "n of ZMuMu events not selected by iso cut ................................. " << Nzmmc_isopassed_0Z << endl;
00715   cout << "n of ZMuMu events with 1 cand selected iso cut ............................ " << Nzmmc_isopassed_1Z << endl;
00716   cout << "n of ZMuMu events with 2 or more cand elected iso cut ..................... " << Nzmmc_isopassed_moreZ << endl<< endl ;
00717 
00718   cout << "n of ZMuMu events not selected by pt and eta cut .......................... " << Nzmmc_ptetapassed_0Z << endl;
00719   cout << "n of ZMuMu events with 1 cand selected by pt and eta cut .................. " << Nzmmc_ptetapassed_1Z << endl;
00720   cout << "n of ZMuMu events with 2 or more cand elected by pt and eta cut ........... " << Nzmmc_ptetapassed_moreZ << endl<< endl ;
00721 
00722   cout << "n of ZMuMu events not selected by pt and eta and mass cut ................. " << Nzmmc_ptetamasspassed_0Z << endl;
00723   cout << "n of ZMuMu events with 1 cand selected by pt and eta and mass cut ......... " << Nzmmc_ptetamasspassed_1Z << endl;
00724   cout << "n of ZMuMu events with 2 or more cand elected by pt and eta and mass cut .. " << Nzmmc_ptetamasspassed_moreZ << endl<< endl ;
00725 
00726   cout << "................When No ZMuMu are selected.................................." << endl;
00727   cout << "n of ZMuSta entries ....................................................... " << Nzmsc << endl;
00728   cout << "n of ZMuSta events with 0 cand ............................................ " << Nzmsc_0Z << endl;
00729   cout << "n of ZMuSta events with 1 cand ............................................ " << Nzmsc_1Z << endl;
00730   cout << "n of ZMuSta events with 2 or more cand .................................... " << Nzmsc_moreZ << endl<< endl ;
00731 
00732   cout << "n of ZMuSta not contained events with 0 cand .............................. " << Nzmsc_notcontained_0Z << endl;
00733   cout << "n of ZMuSta events not contained with 1 cand .............................. " << Nzmsc_notcontained_1Z << endl;
00734   cout << "n of ZMuSta events no contained with 2 or more cand ....................... " << Nzmsc_notcontained_moreZ << endl<< endl ;
00735 
00736   cout << "n of ZMuSta cand not selectd by cuts ...................................... " << Nzmsc_passed_0Z << endl;
00737   cout << "n of ZMuSta events with 1 cand selected by cuts ........................... " << Nzmsc_passed_1Z << endl;
00738   cout << "n of ZMuSta events with 2 or more cand selected by cuts ................... " << Nzmsc_passed_moreZ << endl<< endl ;
00739 
00740   cout << "n of ZMuSta MCmatched cand not selectd by cuts ............................ " << Nzmsc_MCmatched_passed_0Z << endl;
00741   cout << "n of ZMuSta MCmatched events with 1 cand selected by cuts ................. " << Nzmsc_MCmatched_passed_1Z << endl;
00742   cout << "n of ZMuSta MCmatched events with 2 or more cand selected by cuts ......... " << Nzmsc_MCmatched_passed_moreZ << endl<< endl ;
00743 
00744   cout << "...............When no ZMuMu and ZMuSta are selcted........................." << endl;
00745   cout << "n of ZMuTrack entries   ................................................... " << Nzmtc << endl;
00746   cout << "n of ZMuTrack events with 0 cand ..........................................." << Nzmtc_0Z << endl;
00747   cout << "n of ZMuTrack events with 1 cand ..........................................." << Nzmtc_1Z << endl;
00748   cout << "n of ZMuTrack events with 2 or more cand ..................................." << Nzmtc_moreZ << endl<< endl ;
00749 
00750   cout << "n of ZMuTrack not contained events with 0 cand ............................ " << Nzmtc_notcontained_0Z << endl;
00751   cout << "n of ZMuTrack events not contained with 1 cand ............................ " << Nzmtc_notcontained_1Z << endl;
00752   cout << "n of ZMuTrack events no contained with 2 or more cand ..................... " << Nzmtc_notcontained_moreZ << endl<< endl ;
00753 
00754   cout << "n of ZMuTrack cand not selectd by cuts ....................................." << Nzmtc_passed_0Z << endl;
00755   cout << "n of ZMuTrack events with 1 cand selected by cuts .........................." << Nzmtc_passed_1Z << endl;
00756   cout << "n of ZMuTrack events with 2 or more cand selected by cuts .................." << Nzmtc_passed_moreZ << endl<< endl ;
00757 
00758   cout << "n of ZMuTrack MCmatched cand not selectd by cuts .......................... " << Nzmtc_MCmatched_passed_0Z << endl;
00759   cout << "n of ZMuTrcak MCmatched events with 1 cand selected by cuts ............... " << Nzmtc_MCmatched_passed_1Z << endl;
00760   cout << "n of ZMuTrack MCmatched events with 2 or more cand selected by cuts ....... " << Nzmtc_MCmatched_passed_moreZ << endl;
00761 
00762   cout << "------------------------------------------------------------------------------------------" << endl;
00763 
00764   double Nzmm = h_zMuMu_mass_->GetEntries() ;
00765   double Nzsm = h_zMuSingleStandAlone_mass_->GetEntries()  ;
00766   double Nzsnom = h_zMuSingleStandAloneOverlap_mass_->GetEntries()  ;
00767   double Nztm = h_zMuSingleTrack_mass_->GetEntries();
00768   
00769   double NzmmMatch = h_zMuMuMatched_mass_->GetEntries() ;
00770   double NzsmMatch = h_zMuSingleStandAloneMatched_mass_->GetEntries()  ;
00771   double NzsnomMatch = h_zMuSingleStandAloneOverlapMatched_mass_->GetEntries()  ;
00772   double NztmMatch = h_zMuSingleTrackMatched_mass_->GetEntries();
00773 
00774   cout<<"-- N SingleTrackMu = "<<Nztm<<endl;
00775   cout<<"-----N SinglStandAloneMu = "<<Nzsm<<endl;
00776   cout<<"-----N SingleStandAloneOverlapMu = "<<Nzsnom<<endl;
00777   cout<<"------- N MuMu = "<<Nzmm<<endl;
00778   
00779   cout<<"-- N SingleTrackMuMatched = "<<NztmMatch<<endl;
00780   cout<<"-----N SinglStandAloneMuMatched = "<<NzsmMatch<<endl;
00781   cout<<"-----N SingleStandAloneOverlapMuMatched  = "<<NzsnomMatch<<endl;
00782   cout<<"------- N MuMu Matched  = "<<NzmmMatch<<endl;
00783 }
00784   
00785 #include "FWCore/Framework/interface/MakerMacros.h"
00786 
00787 DEFINE_FWK_MODULE(ZMuMuAnalyzer_cynematics);
00788