CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/GeneratorInterface/GenFilters/src/BsJpsiPhiFilter.cc

Go to the documentation of this file.
00001 #include "GeneratorInterface/GenFilters/interface/BsJpsiPhiFilter.h"
00002 
00003 using namespace edm;
00004 using namespace std;
00005 using namespace HepMC;
00006 
00007 
00008 
00009 BsJpsiPhiFilter::BsJpsiPhiFilter(const edm::ParameterSet& iConfig)
00010 {
00011   label_ = iConfig.getUntrackedParameter("moduleLabel",std::string("generator"));
00012   hadronCuts.type = iConfig.getParameter< int >("hadronType");
00013   hadronCuts.etaMin = iConfig.getParameter<double>("hadronEtaMin");
00014   hadronCuts.etaMax = iConfig.getParameter<double>("hadronEtaMax");
00015   hadronCuts.ptMin = iConfig.getParameter<double>("hadronPtMin");
00016   leptonCuts.type = iConfig.getParameter< int >("leptonType");
00017   leptonCuts.etaMin = iConfig.getParameter<double>("leptonEtaMin");
00018   leptonCuts.etaMax = iConfig.getParameter<double>("leptonEtaMax");
00019   leptonCuts.ptMin = iConfig.getParameter<double>("leptonPtMin");
00020 
00021   noAccepted = 0;
00022 }
00023 
00024 
00025 BsJpsiPhiFilter::~BsJpsiPhiFilter()
00026 {  
00027   std::cout << "Total number of accepted events = " << noAccepted << std::endl;
00028 }
00029 
00030 /*
00031 HepMC::GenParticle * BsJpsiPhiFilter::findParticle(const GenPartVect genPartVect,
00032         const int requested_id)
00033 {
00034   for (GenPartVectIt p = genPartVect.begin(); p != genPartVect.end(); p++)
00035     {
00036       int event_particle_id = abs( (*p)->pdg_id() );
00037   cout << "isC "<<event_particle_id<<"\n";
00038       if (requested_id == event_particle_id) return *p;
00039     }
00040   return 0;
00041 }
00042 */
00043 
00044 HepMC::GenParticle * BsJpsiPhiFilter::findParticle(HepMC::GenVertex* vertex, 
00045                                                    const int requested_id)
00046 {
00047   for(std::vector<GenParticle*>::const_iterator p = vertex->particles_out_const_begin();
00048       p != vertex->particles_out_const_end(); p++)
00049     {
00050       int event_particle_id = abs( (*p)->pdg_id() );
00051       cout << "isC "<<event_particle_id<<"\n";
00052       if (requested_id == event_particle_id) return *p;
00053     }
00054   return 0;
00055 }
00056 
00057 HepMC::GenEvent::particle_const_iterator 
00058 BsJpsiPhiFilter::getNextBs(const HepMC::GenEvent::particle_const_iterator start, 
00059                            const HepMC::GenEvent::particle_const_iterator end)
00060 {
00061   HepMC::GenEvent::particle_const_iterator p;
00062   for (p = start; p != end; p++) 
00063     {
00064       int event_particle_id = abs( (*p)->pdg_id() );
00065 //   cout << "search "<<event_particle_id<<"\n";
00066       if (event_particle_id == 531) return p;
00067     }
00068   return p;  
00069 }
00070 
00071 
00072 bool BsJpsiPhiFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)
00073 {
00074   edm::Handle<HepMCProduct> evt;
00075   iEvent.getByLabel(label_, evt);
00076   
00077   const HepMC::GenEvent * generated_event = evt->GetEvent();
00078   //cout << "Start\n";
00079 
00080   bool event_passed = false;
00081   HepMC::GenEvent::particle_const_iterator bs = getNextBs(generated_event->particles_begin(), 
00082                                                     generated_event->particles_end());
00083   while (bs!=  generated_event->particles_end() ) {
00084 
00085     // vector< GenParticle * > bsChild = (*bs)->listChildren();
00086 
00087     //***
00088     HepMC::GenVertex* outVertex = (*bs)->end_vertex();
00089     //***
00090     
00091     GenParticle * jpsi = 0;
00092     GenParticle * phi = 0;
00093     // cout << "bs size "<<bsChild.size()<<endl;
00094     //***
00095     int numChildren = outVertex->particles_out_size();
00096     cout<< "bs size "<<numChildren<<endl;
00097     //***
00098     
00099     /*    if ((bsChild.size()==2) && ((jpsi = findParticle(bsChild, 443))!=0) && 
00100           ((phi = findParticle(bsChild, 333))!=0)) {
00101           cout << bsChild[0]->momentum()<<" "<<bsChild[0]->momentum().eta()
00102           <<" "<<bsChild[1]->momentum()<<" "<<bsChild[1]->momentum().eta()<<endl;
00103     */
00104     
00105     //***
00106     if( (numChildren==2) && ((jpsi = findParticle(outVertex, 443))!=0) && 
00107         ((phi = findParticle(outVertex, 333))!=0)) {
00108       
00109       cout << jpsi->momentum().rho()<<" "<<jpsi->momentum().eta()
00110            <<" "<<phi->momentum().rho() <<" "<<phi->momentum().eta()<<endl;
00111       cout <<"bs dec trouve"<<endl;
00112       if (cuts(phi, hadronCuts) && cuts(jpsi, leptonCuts)) {
00113         cout <<"OK trouve"<<endl;
00114         event_passed = true;
00115         break;
00116       }
00117     }
00118     bs = getNextBs(++bs, generated_event->particles_end());
00119   }
00120   
00121   if (event_passed) noAccepted++;
00122   cout << "End filter\n";
00123   
00124   delete generated_event; 
00125   
00126   return event_passed;
00127 }
00128 
00129 /*
00130 bool BsJpsiPhiFilter::cuts(const GenParticle * jpsi, const CutStruct cut)
00131 {
00132   GenPartVect psiChild = jpsi->listChildren();
00133     cout << psiChild[0]->pdg_id()<<" "<<psiChild[1]->pdg_id()<<endl;
00134   if (psiChild.size()==2 && (abs(psiChild[0]->pdg_id()) == cut.type) &&
00135     (abs(psiChild[1]->pdg_id()) == cut.type))
00136   {
00137     cout << psiChild[0]->momentum()<<" "<<psiChild[0]->momentum().eta()
00138     <<" "<<psiChild[1]->momentum()<<" "<<psiChild[1]->momentum().eta()<<endl;
00139     return ( (etaInRange(psiChild[0]->momentum().eta(), cut.etaMin, cut.etaMax)) &&
00140              (etaInRange(psiChild[1]->momentum().eta(), cut.etaMin, cut.etaMax)) &&
00141              (psiChild[0]->momentum().perp()> cut.ptMin) &&
00142              (psiChild[1]->momentum().perp()> cut.ptMin));
00143   }
00144   return false;
00145 }
00146 */
00147 
00148 bool BsJpsiPhiFilter::cuts(const GenParticle * jpsi, const CutStruct cut)
00149 {
00150   HepMC::GenVertex* myVertex = jpsi->end_vertex();
00151   int numChildren = myVertex->particles_out_size();
00152   std::vector<HepMC::GenParticle*> psiChild;
00153   for(std::vector<GenParticle*>::const_iterator p = myVertex->particles_out_const_begin();
00154       p != myVertex->particles_out_const_end(); p++) 
00155     psiChild.push_back((*p));
00156   
00157   if(numChildren>1) {
00158     cout << psiChild[0]->pdg_id()<<" "<<psiChild[1]->pdg_id()<<endl;
00159   
00160     if (psiChild.size()==2 && (abs(psiChild[0]->pdg_id()) == cut.type) &&
00161         (abs(psiChild[1]->pdg_id()) == cut.type))
00162       {
00163         cout << psiChild[0]->momentum().rho()<<" "<<psiChild[0]->momentum().eta()
00164              <<" "<<psiChild[1]->momentum().rho()<<" "<<psiChild[1]->momentum().eta()<<endl;
00165         return ( (etaInRange(psiChild[0]->momentum().eta(), cut.etaMin, cut.etaMax)) &&
00166                  (etaInRange(psiChild[1]->momentum().eta(), cut.etaMin, cut.etaMax)) &&
00167                  (psiChild[0]->momentum().perp()> cut.ptMin) &&
00168                  (psiChild[1]->momentum().perp()> cut.ptMin));
00169       }
00170     return false;
00171   }
00172   return false;
00173 }
00174 
00175 bool BsJpsiPhiFilter::etaInRange(float eta, float etamin, float etamax)
00176 {
00177   return ( (etamin < eta) && (eta < etamax) );
00178 }