CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/Validation/EventGenerator/plugins/HiggsValidation.cc

Go to the documentation of this file.
00001 /*class HiggsValidation
00002  *  
00003  *  Class to fill dqm monitor elements from existing EDM file
00004  *
00005  *  $Date: 2012/08/12 16:13:29 $
00006  *  $Revision: 1.1 $
00007  */
00008  
00009 #include "Validation/EventGenerator/interface/HiggsValidation.h"
00010 
00011 #include "FWCore/Framework/interface/MakerMacros.h"  
00012 
00013 #include "CLHEP/Units/defs.h"
00014 #include "CLHEP/Units/PhysicalConstants.h"
00015 
00016 #include "DataFormats/Math/interface/LorentzVector.h"
00017 
00018 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00019 #include "Validation/EventGenerator/interface/PdtPdgMini.h"
00020 
00021 using namespace edm;
00022 
00023 HiggsValidation::HiggsValidation(const edm::ParameterSet& iPSet): 
00024   _wmanager(iPSet),
00025   hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
00026   particle_id(iPSet.getParameter<int>("pdg_id")),
00027   particle_name(iPSet.getParameter<std::string>("particleName"))
00028 {    
00029   dbe = 0;
00030   dbe = edm::Service<DQMStore>().operator->();
00031   
00032   monitoredDecays = new MonitoredDecays(iPSet);
00033 
00034 }
00035 
00036 HiggsValidation::~HiggsValidation() {}
00037 
00038 void HiggsValidation::beginJob()
00039 {
00040   if(dbe){
00042     TString dir="Generator/";
00043     dir+=particle_name;
00044     dbe->setCurrentFolder(dir.Data());
00045     
00046     // Number of analyzed events
00047     nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
00048     
00049     //decay type
00050     
00051     std::string channel = particle_name+"_DecayChannels";
00052     HiggsDecayChannels = dbe->book1D(channel.c_str(),(particle_name+" decay channels").c_str(),monitoredDecays->size(),0,monitoredDecays->size());
00053     
00054     for(size_t i = 0; i < monitoredDecays->size(); ++i){
00055       HiggsDecayChannels->setBinLabel(1+i,monitoredDecays->channel(i));
00056     }
00057   }
00058 
00059   //Kinematics 
00060   Higgs_pt = dbe->book1D((particle_name+"_pt"),(particle_name+" p_{t}"),50,0,250);
00061   Higgs_eta  = dbe->book1D((particle_name+"_eta"),(particle_name+" #eta"),50,-5,5); 
00062   Higgs_mass = dbe->book1D((particle_name+"_m"),(particle_name+" M"),500,0,500);
00063 
00064   int idx=0;
00065   for(unsigned int i=0;i<monitoredDecays->NDecayParticles();i++){
00066     HiggsDecayProd_pt.push_back(dbe->book1D((monitoredDecays->ConvertIndex(idx)+"_pt"),(monitoredDecays->ConvertIndex(idx)+" p_{t}"),50,0,250));
00067     HiggsDecayProd_eta.push_back(dbe->book1D((monitoredDecays->ConvertIndex(idx)+"_eta"),(monitoredDecays->ConvertIndex(idx)+" #eta"),50,-5,5));
00068     idx++;
00069   }
00070 
00071   return;
00072 }
00073 
00074 void HiggsValidation::endJob(){
00075   return;
00076 }
00077 
00078 void HiggsValidation::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup)
00079 {
00081   iSetup.getData( fPDGTable );
00082   return;
00083 }
00084 void HiggsValidation::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
00085 void HiggsValidation::analyze(const edm::Event& iEvent,const edm::EventSetup& iSetup)
00086 { 
00087   double weight = _wmanager.weight(iEvent);  
00088   nEvt->Fill(0.5,weight);
00089   
00090   //Gathering the HepMCProduct information
00091   edm::Handle<HepMCProduct> evt;
00092   iEvent.getByLabel(hepmcCollection_, evt);
00093   
00094   //Get EVENT
00095   HepMC::GenEvent *myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
00096 
00097   // loop over all particles
00098   bool filled = false;
00099   for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); 
00100       iter!= myGenEvent->particles_end() && !filled; ++iter) {
00101     if(particle_id == fabs((*iter)->pdg_id())){
00102       std::vector<HepMC::GenParticle*> decayprod;
00103       int channel = findHiggsDecayChannel(*iter,decayprod);
00104       HiggsDecayChannels->Fill(channel,weight);
00105       Higgs_pt->Fill((*iter)->momentum().perp(),weight); 
00106       Higgs_eta->Fill((*iter)->momentum().eta(),weight);
00107       Higgs_mass->Fill((*iter)->momentum().m(),weight);
00108       for(unsigned int i=0;i<decayprod.size();i++){
00109         int idx=monitoredDecays->isDecayParticle(decayprod.at(i)->pdg_id());
00110         if(0<=idx && idx<=(int)HiggsDecayProd_pt.size()){
00111           HiggsDecayProd_pt.at(idx)->Fill(decayprod.at(i)->momentum().perp(),weight);
00112           HiggsDecayProd_eta.at(idx)->Fill(decayprod.at(i)->momentum().eta(),weight);
00113         }
00114       }
00115       filled = true;
00116     }
00117   }
00118   
00119   delete myGenEvent;
00120   
00121 }//analyze
00122 
00123 int HiggsValidation::findHiggsDecayChannel(const HepMC::GenParticle* genParticle,std::vector<HepMC::GenParticle*> &decayprod){
00124   
00125   if(genParticle->status() == 1) return monitoredDecays->stable();
00126   
00127   std::vector<int> children;
00128   if ( genParticle->end_vertex() ) {
00129     HepMC::GenVertex::particle_iterator des;
00130     for(des = genParticle->end_vertex()->particles_begin(HepMC::descendants);
00131         des!= genParticle->end_vertex()->particles_end(HepMC::descendants);++des ) {
00132       
00133       if((*des)->pdg_id() == genParticle->pdg_id()) continue;
00134       
00135       HepMC::GenVertex::particle_iterator mother = (*des)->production_vertex()->particles_begin(HepMC::parents);
00136       if((*mother)->pdg_id() == genParticle->pdg_id()){
00137         children.push_back((*des)->pdg_id());
00138         decayprod.push_back((*des));
00139       }
00140     }
00141   }
00142   
00143   if(children.size() == 2 && children.at(0) != 0 && children.at(1) != 0) return monitoredDecays->position(children.at(0),children.at(1));
00144   return monitoredDecays->undetermined();
00145 }
00146 
00147 //define this as a plug-in
00148 DEFINE_FWK_MODULE(HiggsValidation);