CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MCDecayingPionKaonFilter.cc
Go to the documentation of this file.
2 
4 #include <iostream>
5 
6 using namespace edm;
7 using namespace std;
8 
9 
11 token_(consumes<edm::HepMCProduct>(edm::InputTag(iConfig.getUntrackedParameter("moduleLabel",std::string("generator")),"unsmeared")))
12 {
13  //here do whatever other initialization is needed
14  vector<int> defpid ;
15  defpid.push_back(0) ;
16  particleID = iConfig.getUntrackedParameter< vector<int> >("ParticleID",defpid);
17  vector<double> defptmin ;
18  defptmin.push_back(0.);
19  ptMin = iConfig.getUntrackedParameter< vector<double> >("MinPt", defptmin);
20 
21  vector<double> defetamin ;
22  defetamin.push_back(-10.);
23  etaMin = iConfig.getUntrackedParameter< vector<double> >("MinEta", defetamin);
24  vector<double> defetamax ;
25  defetamax.push_back(10.);
26  etaMax = iConfig.getUntrackedParameter< vector<double> >("MaxEta", defetamax);
27 
28  vector<double> defDecayRadiusmin ;
29  defDecayRadiusmin.push_back(-10.);
30  decayRadiusMin = iConfig.getUntrackedParameter< vector<double> >("MinDecayRadius", defDecayRadiusmin);
31 
32  vector<double> defDecayRadiusmax ;
33  defDecayRadiusmax.push_back(1.e5);
34  decayRadiusMax = iConfig.getUntrackedParameter< vector<double> >("MaxDecayRadius", defDecayRadiusmax);
35 
36  vector<double> defDecayZmin ;
37  defDecayZmin.push_back(-1.e5);
38  decayZMin = iConfig.getUntrackedParameter< vector<double> >("MinDecayZ", defDecayZmin);
39 
40  vector<double> defDecayZmax ;
41  defDecayZmax.push_back(1.e5);
42  decayZMax = iConfig.getUntrackedParameter< vector<double> >("MaxDecayZ", defDecayZmax);
43 
44  ptMuMin = iConfig.getUntrackedParameter<double>("PtMuMin",0.);
45  // check for same size
46  if ( (ptMin.size() > 1 && particleID.size() != ptMin.size())
47  || (etaMin.size() > 1 && particleID.size() != etaMin.size())
48  || (etaMax.size() > 1 && particleID.size() != etaMax.size())
49  || (decayRadiusMin.size() > 1 && particleID.size() != decayRadiusMin.size())
50  || (decayRadiusMax.size() > 1 && particleID.size() != decayRadiusMax.size())
51  || (decayZMin.size() > 1 && particleID.size() != decayZMin.size())
52  || (decayZMax.size() > 1 && particleID.size() != decayZMax.size()) ) {
53  cout << "WARNING: MCPROCESSFILTER : size of MinPthat and/or MaxPthat not matching with ProcessID size!!" << endl;
54  }
55 
56  // if ptMin size smaller than particleID , fill up further with defaults
57  if (particleID.size() > ptMin.size() ){
58  vector<double> defptmin2 ;
59  for (unsigned int i = 0; i < particleID.size(); i++){ defptmin2.push_back(0.);}
60  ptMin = defptmin2;
61  }
62  // if etaMin size smaller than particleID , fill up further with defaults
63  if (particleID.size() > etaMin.size() ){
64  vector<double> defetamin2 ;
65  for (unsigned int i = 0; i < particleID.size(); i++){ defetamin2.push_back(-10.);}
66  etaMin = defetamin2;
67  }
68  // if etaMax size smaller than particleID , fill up further with defaults
69  if (particleID.size() > etaMax.size() ){
70  vector<double> defetamax2 ;
71  for (unsigned int i = 0; i < particleID.size(); i++){ defetamax2.push_back(10.);}
72  etaMax = defetamax2;
73  }
74 
75  // if decayRadiusMin size smaller than particleID , fill up further with defaults
76  if (particleID.size() > decayRadiusMin.size() ){
77  vector<double> decayRadiusmin2 ;
78  for (unsigned int i = 0; i < particleID.size(); i++){ decayRadiusmin2.push_back(-10.);}
79  decayRadiusMin = decayRadiusmin2;
80  }
81  // if decayRadiusMax size smaller than particleID , fill up further with defaults
82  if (particleID.size() > decayRadiusMax.size() ){
83  vector<double> decayRadiusmax2 ;
84  for (unsigned int i = 0; i < particleID.size(); i++){ decayRadiusmax2.push_back(1.e5);}
85  decayRadiusMax = decayRadiusmax2;
86  }
87 
88  // if decayZMin size smaller than particleID , fill up further with defaults
89  if (particleID.size() > decayZMin.size() ){
90  vector<double> decayZmin2 ;
91  for (unsigned int i = 0; i < particleID.size(); i++){ decayZmin2.push_back(-1.e5);}
92  decayZMin = decayZmin2;
93  }
94  // if decayZMax size smaller than particleID , fill up further with defaults
95  if (particleID.size() > decayZMax.size() ){
96  vector<double> decayZmax2 ;
97  for (unsigned int i = 0; i < particleID.size(); i++){ decayZmax2.push_back(1.e5);}
98  decayZMax = decayZmax2;
99  }
100 
101 }
102 
103 
105 {
106 
107  // do anything here that needs to be done at desctruction time
108  // (e.g. close files, deallocate resources etc.)
109 
110 }
111 
112 
113 // ------------ method called to skim the data ------------
115 {
116  using namespace edm;
117  bool accepted = false;
119  iEvent.getByToken(token_, evt);
120 
121  const HepMC::GenEvent * myGenEvent = evt->GetEvent();
122 
123 
124  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
125  p != myGenEvent->particles_end(); ++p ) {
126 
127 
128  for (unsigned int i = 0; i < particleID.size(); i++){
129  if (!((*p)->end_vertex())) continue;
130  if (particleID[i] != (*p)->pdg_id() && particleID[i] != 0) continue;
131 
132  if ( (*p)->momentum().perp() < ptMin[i] ) continue;
133 
134  if ( (*p)->momentum().eta() < etaMin[i] ) continue;
135  if ( (*p)->momentum().eta() > etaMax[i] ) continue;
136 
137  double decx = (*p)->end_vertex()->position().x();
138  double decy = (*p)->end_vertex()->position().y();
139  double decrad = sqrt(decx*decx+decy*decy);
140  if (decrad<decayRadiusMin[i] ) continue;
141  if (decrad>decayRadiusMax[i] ) continue;
142 
143  double decz = (*p)->end_vertex()->position().z();
144  if (decz<decayZMin[i] ) continue;
145  if (decz>decayZMax[i] ) continue;
146 
147 
148  if ((*p)->status()==2) {
149  for (HepMC::GenVertex::particle_iterator
150  vpdec= (*p)->end_vertex()->particles_begin(HepMC::children);
151  vpdec != (*p)->end_vertex()->particles_end(HepMC::children);
152  ++vpdec) {
153  if (abs((*vpdec)->pdg_id())!=13) continue;
154  if (fabs((*vpdec)->momentum().perp())>ptMuMin) {
155  accepted = true;
156  break;
157  }
158  }
159  } else if ((*p)->status()==1) {
160  accepted = true;
161  }
162 
163  }
164 
165 
166  }
167 
168  delete myGenEvent;
169 
170 
171  if (accepted){ return true; } else {return false;}
172 
173 }
174 
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
std::vector< double > decayZMax
edm::EDGetTokenT< edm::HepMCProduct > token_
std::vector< double > decayRadiusMin
MCDecayingPionKaonFilter(const edm::ParameterSet &)
std::vector< double > decayRadiusMax
int iEvent
Definition: GenABIO.cc:230
T sqrt(T t)
Definition: SSEVec.h:48
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< double > decayZMin
virtual bool filter(edm::Event &, const edm::EventSetup &)
tuple cout
Definition: gather_cfg.py:121