CMS 3D CMS Logo

GenLeadTrackFilter.cc
Go to the documentation of this file.
1 /***********************************************************
2 * GenLeadTrackFilter *
3 * ------------------ *
4 * *
5 * Original Author: Souvik Das, Cornell University *
6 * Created : 7 August 2009 *
7 * *
8 * Allows events which have at least one generator level *
9 * charged particle with pT greater than X GeV within *
10 * |eta| less than Y, where X and Y are specified in the *
11 * cfi configuration file. *
12 ***********************************************************/
13 
16 
17 using namespace edm;
18 using namespace reco;
19 using namespace math;
20 
22 {
23  hepMCProduct_label_ = iConfig.getParameter<InputTag>("HepMCProduct");
24  genLeadTrackPt_ = iConfig.getParameter<double>("GenLeadTrackPt");
25  genEta_ = iConfig.getParameter<double>("GenEta");
26 }
27 
29 {
30 }
31 
33 {
34  bool allow=false;
35 
37  iSetup.getData(pdt);
38 
40  iEvent.getByLabel(hepMCProduct_label_, genEvent);
41  if (genEvent.isValid())
42  {
43  float genLeadTrackPt=-100;
44  for (HepMC::GenEvent::particle_const_iterator iter=(*(genEvent->GetEvent())).particles_begin();
45  iter!=(*(genEvent->GetEvent())).particles_end();
46  ++iter)
47  {
48  HepMC::GenParticle* theParticle=*iter;
49  double pt=pow(pow(theParticle->momentum().px(),2)+pow(theParticle->momentum().py(),2), 0.5);
50  double charge=pdt->particle(theParticle->pdg_id())->charge();
51  if (theParticle->status()==1 &&
52  charge!=0 &&
53  fabs(theParticle->momentum().eta())<genEta_ &&
54  pt>genLeadTrackPt)
55  {
56  genLeadTrackPt=pt;
57  }
58  }
59  if (genLeadTrackPt>genLeadTrackPt_) allow=true; else allow=false;
60  }
61  else
62  {
63  std::cout<<"genEvent in not valid!"<<std::endl;
64  }
65  return allow;
66 }
67 
68 // ------------ method called once each job just before starting event loop ------------
69 
70 // ------------ method called once each job just after ending the event loop ------------
72 {
73 }
74 
75 //define this as a plug-in
T getParameter(std::string const &) const
void endJob() override
GenLeadTrackFilter(const edm::ParameterSet &)
bool getData(T &iHolder) const
Definition: EventSetup.h:111
~GenLeadTrackFilter() override
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool isValid() const
Definition: HandleBase.h:74
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:480
Definition: Error.h:16
bool filter(edm::Event &, const edm::EventSetup &) override
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
fixed size matrix
HLT enums.
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40