Go to the documentation of this file.00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022 #include <vector>
00023 #include <boost/format.hpp>
00024
00025
00026 #include "GeneratorInterface/GenFilters/interface/ComphepSingletopFilter.h"
00027 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00028
00029
00030
00031
00032
00033
00034 ComphepSingletopFilter::ComphepSingletopFilter(const edm::ParameterSet& iConfig) {
00035 ptsep = iConfig.getParameter<double>("pTSep");
00036 }
00037
00038
00039 ComphepSingletopFilter::~ComphepSingletopFilter() {}
00040
00041
00042 void ComphepSingletopFilter::beginJob() {
00043 read22 = read23 = 0;
00044 pass22 = pass23 = 0;
00045 }
00046
00047
00048 void ComphepSingletopFilter::endJob() {
00049 using namespace std;
00050 cout << "Proc: 2-->2 2-->3 Total" << endl;
00051 cout << boost::format("Read: %9d %9d %9d") % read22 % read23 % (read22+read23)
00052 << endl;
00053 cout << boost::format("Pass: %9d %9d %9d") % pass22 % pass23 % (pass22+pass23)
00054 << endl;
00055 }
00056
00057
00058 bool ComphepSingletopFilter::filter(
00059 edm::Event& iEvent,
00060 const edm::EventSetup& iSetup) {
00061
00062 using namespace std;
00063 using namespace HepMC;
00064
00065 edm::Handle<edm::HepMCProduct> evt;
00066 iEvent.getByType(evt);
00067 const HepMC::GenEvent * myEvt = evt->GetEvent();
00068
00069
00070 const GenParticle * gp_clep = NULL;
00071
00072 for (GenEvent::particle_const_iterator it = myEvt->particles_begin();
00073 it != myEvt->particles_end(); ++it) {
00074 if ((*it)->status() == 3) {
00075 int abs_id = abs((*it)->pdg_id());
00076 if (abs_id==11 || abs_id==13 || abs_id==15) {
00077 gp_clep = *it;
00078 break;
00079 }
00080 }
00081 }
00082
00083 if (! gp_clep) {
00084 cerr << "ERROR: ComphepSingletopFilter: no charged lepton" << endl;
00085 return false;
00086 }
00087
00088 int id_bdec = 0;
00089 const GenParticle * gp_bdec = NULL;
00090 vector<const GenParticle *> vgp_bsec;
00091
00092 GenVertex * gv_hard = gp_clep->production_vertex();
00093
00094 for (GenVertex::particle_iterator it = gv_hard->particles_begin(children);
00095 it != gv_hard->particles_end(children); ++it) {
00096 int pdg_id = (*it)->pdg_id();
00097 if (abs(pdg_id) == 5) {
00098 if (pdg_id * (gp_clep->pdg_id()) < 0) {
00099 id_bdec = pdg_id;
00100 gp_bdec = *it;
00101 } else {
00102 vgp_bsec.push_back(*it);
00103 }
00104 }
00105 }
00106
00107 bool process22 = (vgp_bsec.size() == 0);
00108
00109 GenVertex * gv = NULL;
00110 if (process22) {
00111 for (GenVertex::particle_iterator it = gv_hard->particles_begin(parents);
00112 it != gv_hard->particles_end(parents); ++it) {
00113 if ((*it)->pdg_id() == id_bdec) {
00114 gv = (*it)->production_vertex();
00115 break;
00116 }
00117 }
00118 if (! gv) {
00119 cerr << "ERROR: ComphepSingletopFilter: HepMC inconsistency" << endl;
00120 myEvt->print();
00121 return false;
00122 }
00123 } else {
00124 gv = vgp_bsec.back()->end_vertex();
00125 }
00126 const GenParticle * gp;
00127 while (gv) {
00128 gp = NULL;
00129 for (GenVertex::particle_iterator it = gv->particles_begin(children);
00130 it != gv->particles_end(children); ++it) {
00131 if ((*it)->pdg_id() == -id_bdec) {
00132 if (!gp || (*it)->momentum().perp2() > gp->momentum().perp2()) {
00133 gp = *it;
00134 }
00135 }
00136 }
00137 if (gp) {
00138 gv = gp->end_vertex();
00139 vgp_bsec.push_back(gp);
00140 } else {
00141 gv = NULL;
00142 }
00143 }
00144
00145 if (vgp_bsec.size() == 0) {
00146 cerr << "ERROR: ComphepSingletopFilter: HepMC inconsistency" << endl;
00147 return false;
00148 }
00149
00150 double pt = vgp_bsec.back()->momentum().perp();
00151 bool pass;
00152 if (process22) {
00153 read22 += 1;
00154 pass = pt < ptsep;
00155 if (pass) pass22 += 1;
00156 } else {
00157 read23 += 1;
00158 pass = ptsep <= pt;
00159 if (pass) pass23 += 1;
00160 }
00161
00162 return pass;
00163 }
00164
00165