CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
GenericDauHepMCFilter Class Reference

#include <GenericDauHepMCFilter.h>

Inheritance diagram for GenericDauHepMCFilter:
BaseHepMCFilter

Public Member Functions

virtual bool filter (const HepMC::GenEvent *evt)
 
 GenericDauHepMCFilter (const edm::ParameterSet &)
 
 ~GenericDauHepMCFilter ()
 
- Public Member Functions inherited from BaseHepMCFilter
 BaseHepMCFilter ()
 
virtual ~BaseHepMCFilter ()
 

Private Attributes

bool chargeconju
 
std::vector< int > dauIDs
 
double maxetacut
 
double maxptcut
 
double minetacut
 
double minptcut
 
int ndaughters
 
int particleID
 

Detailed Description

Description: Filter events using MotherId and ChildrenIds infos

Implementation: <Notes on="" implementation>="">

Definition at line 35 of file GenericDauHepMCFilter.h.

Constructor & Destructor Documentation

GenericDauHepMCFilter::GenericDauHepMCFilter ( const edm::ParameterSet iConfig)

Definition at line 12 of file GenericDauHepMCFilter.cc.

12  :
13  particleID(iConfig.getParameter<int>("ParticleID")),
14  chargeconju(iConfig.getParameter<bool>("ChargeConjugation")),
15  ndaughters(iConfig.getParameter<int>("NumberDaughters")),
16  dauIDs(iConfig.getParameter<std::vector<int> >("DaughterIDs")),
17  minptcut(iConfig.getParameter<double>("MinPt")),
18  maxptcut(iConfig.getParameter<double>("MaxPt")),
19  minetacut(iConfig.getParameter<double>("MinEta")),
20  maxetacut(iConfig.getParameter<double>("MaxEta"))
21 {
22 
23 }
T getParameter(std::string const &) const
std::vector< int > dauIDs
GenericDauHepMCFilter::~GenericDauHepMCFilter ( )

Definition at line 26 of file GenericDauHepMCFilter.cc.

27 {
28 
29 }

Member Function Documentation

bool GenericDauHepMCFilter::filter ( const HepMC::GenEvent *  evt)
virtual

Implements BaseHepMCFilter.

Definition at line 37 of file GenericDauHepMCFilter.cc.

References chargeconju, dauIDs, i, maxetacut, maxptcut, minetacut, minptcut, ndaughters, AlCaHLTBitMon_ParallelJobs::p, and particleID.

38 {
39 
40  for ( HepMC::GenEvent::particle_const_iterator p = evt->particles_begin();
41  p != evt->particles_end(); ++p ) {
42 
43  if( (*p)->pdg_id() != particleID ) continue ;
44 
45  int ndauac = 0;
46  int ndau = 0;
47  if ( (*p)->end_vertex() ) {
48  for ( HepMC::GenVertex::particle_iterator
49  des=(*p)->end_vertex()->particles_begin(HepMC::children);
50  des != (*p)->end_vertex()->particles_end(HepMC::children);
51  ++des ) {
52  ++ndau;
53  // cout << " -> Daughter = " << (*des)->pdg_id() << endl;
54  for( unsigned int i=0; i<dauIDs.size(); ++i) {
55  if( (*des)->pdg_id() != dauIDs[i] ) continue ;
56  if( (*des)->momentum().perp() > minptcut &&
57  (*des)->momentum().perp() < maxptcut &&
58  (*des)->momentum().eta() > minetacut &&
59  (*des)->momentum().eta() < maxetacut ) {
60  ++ndauac;
61  // cout << " -> pT = " << (*des)->momentum().perp() << endl;
62  break;
63  }
64  }
65  }
66  }
67  if( ndau == ndaughters && ndauac == ndaughters ) {
68  return true;
69  }
70 
71  }
72 
73 
74  if (chargeconju) {
75 
76  for ( HepMC::GenEvent::particle_const_iterator p = evt->particles_begin();
77  p != evt->particles_end(); ++p ) {
78 
79  if( (*p)->pdg_id() != -particleID ) continue ;
80  int ndauac = 0;
81  int ndau = 0;
82  if ( (*p)->end_vertex() ) {
83  for ( HepMC::GenVertex::particle_iterator
84  des=(*p)->end_vertex()->particles_begin(HepMC::children);
85  des != (*p)->end_vertex()->particles_end(HepMC::children);
86  ++des ) {
87  ++ndau;
88  for( unsigned int i=0; i<dauIDs.size(); ++i) {
89  bool has_antipart = !(dauIDs[i]==22 || dauIDs[i]==23);
90  int IDanti = has_antipart ? -dauIDs[i] : dauIDs[i];
91  if( (*des)->pdg_id() != IDanti ) continue ;
92  if( (*des)->momentum().perp() > minptcut &&
93  (*des)->momentum().perp() < maxptcut &&
94  (*des)->momentum().eta() > minetacut &&
95  (*des)->momentum().eta() < maxetacut ) {
96  ++ndauac;
97  break;
98  }
99  }
100  }
101  }
102  if( ndau == ndaughters && ndauac == ndaughters ) {
103  return true;
104  }
105  }
106 
107  }
108 
109  return false;
110 
111 }
int i
Definition: DBlmapReader.cc:9
std::vector< int > dauIDs

Member Data Documentation

bool GenericDauHepMCFilter::chargeconju
private

Definition at line 48 of file GenericDauHepMCFilter.h.

Referenced by filter().

std::vector<int> GenericDauHepMCFilter::dauIDs
private

Definition at line 50 of file GenericDauHepMCFilter.h.

Referenced by filter().

double GenericDauHepMCFilter::maxetacut
private

Definition at line 54 of file GenericDauHepMCFilter.h.

Referenced by filter().

double GenericDauHepMCFilter::maxptcut
private

Definition at line 52 of file GenericDauHepMCFilter.h.

Referenced by filter().

double GenericDauHepMCFilter::minetacut
private

Definition at line 53 of file GenericDauHepMCFilter.h.

Referenced by filter().

double GenericDauHepMCFilter::minptcut
private

Definition at line 51 of file GenericDauHepMCFilter.h.

Referenced by filter().

int GenericDauHepMCFilter::ndaughters
private

Definition at line 49 of file GenericDauHepMCFilter.h.

Referenced by filter().

int GenericDauHepMCFilter::particleID
private

Definition at line 47 of file GenericDauHepMCFilter.h.

Referenced by filter().