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
PFMETFilter Class Reference

#include <PFMETFilter.h>

Inheritance diagram for PFMETFilter:
edm::EDFilter edm::ProducerBase edm::ProductRegistryHelper

Public Member Functions

virtual void beginJob ()
 
bool checkInput ()
 
virtual void endJob ()
 
virtual bool filter (edm::Event &, const edm::EventSetup &)
 
 PFMETFilter (const edm::ParameterSet &)
 
virtual ~PFMETFilter ()
 
- Public Member Functions inherited from edm::EDFilter
 EDFilter ()
 
virtual ~EDFilter ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
boost::function< void(const
BranchDescription &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 

Private Attributes

std::vector< std::string > collections_
 
double DeltaMEXsigma_
 
std::vector< int > doMax_
 
std::vector< int > doMin_
 
std::vector< double > max_
 
std::vector< double > min_
 
double sigma_a_
 
double sigma_b_
 
double sigma_c_
 
std::string TrueMET_
 
std::vector< std::string > variables_
 
bool verbose_
 

Additional Inherited Members

- Public Types inherited from edm::EDFilter
typedef EDFilter ModuleType
 
typedef WorkerT< EDFilterWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDFilter
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDFilter
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
template<class TProducer , class TMethod >
void callWhenNewProductsRegistered (TProducer *iProd, TMethod iMethod)
 

Detailed Description

Definition at line 11 of file PFMETFilter.h.

Constructor & Destructor Documentation

PFMETFilter::PFMETFilter ( const edm::ParameterSet iConfig)
explicit

Definition at line 12 of file PFMETFilter.cc.

References collections_, DeltaMEXsigma_, doMax_, doMin_, edm::ParameterSet::getParameter(), max_, min_, sigma_a_, sigma_b_, sigma_c_, TrueMET_, variables_, and verbose_.

13 {
14  collections_= iConfig.getParameter< std::vector<std::string> >("Collections");
15  variables_= iConfig.getParameter< std::vector<std::string> >("Variables");
16  min_ = iConfig.getParameter< std::vector<double> >("Mins");
17  max_ = iConfig.getParameter< std::vector<double> >("Maxs");
18  doMin_ = iConfig.getParameter< std::vector<int> >("DoMin");
19  doMax_ = iConfig.getParameter< std::vector<int> >("DoMax");
20  TrueMET_ = iConfig.getParameter<std::string>("TrueMET");
21  DeltaMEXsigma_ = iConfig.getParameter<double>("DeltaMEXsigma");
22  sigma_a_ = iConfig.getParameter<double>("sigma_a");
23  sigma_b_ = iConfig.getParameter<double>("sigma_b");
24  sigma_c_ = iConfig.getParameter<double>("sigma_c");
25  verbose_ = iConfig.getParameter<bool>("verbose");
26 }
T getParameter(std::string const &) const
std::vector< double > max_
Definition: PFMETFilter.h:26
double sigma_c_
Definition: PFMETFilter.h:36
std::vector< double > min_
Definition: PFMETFilter.h:25
std::vector< int > doMax_
Definition: PFMETFilter.h:28
double sigma_a_
Definition: PFMETFilter.h:34
std::string TrueMET_
Definition: PFMETFilter.h:32
double DeltaMEXsigma_
Definition: PFMETFilter.h:33
std::vector< std::string > collections_
Definition: PFMETFilter.h:23
std::vector< int > doMin_
Definition: PFMETFilter.h:27
double sigma_b_
Definition: PFMETFilter.h:35
std::vector< std::string > variables_
Definition: PFMETFilter.h:24
bool verbose_
Definition: PFMETFilter.h:37
PFMETFilter::~PFMETFilter ( )
virtual

Definition at line 28 of file PFMETFilter.cc.

28 {}

Member Function Documentation

void PFMETFilter::beginJob ( void  )
virtual

Reimplemented from edm::EDFilter.

Definition at line 70 of file PFMETFilter.cc.

71 {
72  //std::cout << "FL: beginJob" << std::endl;
73 }
bool PFMETFilter::checkInput ( )

Definition at line 30 of file PFMETFilter.cc.

References collections_, gather_cfg::cout, doMax_, doMin_, max_, min_, and variables_.

Referenced by filter().

31 {
32  if (collections_.size()!=min_.size())
33  {
34  std::cout << "Error: in PFMETFilter: collections_.size()!=min_.size()" << std::endl;
35  std::cout << "collections_.size() = " << collections_.size() << std::endl;
36  std::cout << "min_.size() = " << min_.size() << std::endl;
37  return false;
38  }
39  if (collections_.size()!=max_.size())
40  {
41  std::cout << "Error: in PFMETFilter: collections_.size()!=max_.size()" << std::endl;
42  std::cout << "collections_.size() = " << collections_.size() << std::endl;
43  std::cout << "max_.size() = " << max_.size() << std::endl;
44  return false;
45  }
46  if (collections_.size()!=doMin_.size())
47  {
48  std::cout << "Error: in PFMETFilter: collections_.size()!=min_.size()" << std::endl;
49  std::cout << "collections_.size() = " << collections_.size() << std::endl;
50  std::cout << "doMin_.size() = " << doMin_.size() << std::endl;
51  return false;
52  }
53  if (collections_.size()!=doMax_.size())
54  {
55  std::cout << "Error: in PFMETFilter: collections_.size()!=min_.size()" << std::endl;
56  std::cout << "collections_.size() = " << collections_.size() << std::endl;
57  std::cout << "doMax_.size() = " << doMax_.size() << std::endl;
58  return false;
59  }
60  if (collections_.size()!=variables_.size())
61  {
62  std::cout << "Error: in PFMETFilter: collections_.size()!=variables_.size()" << std::endl;
63  std::cout << "collections_.size() = " << collections_.size() << std::endl;
64  std::cout << "variables_.size() = " << variables_.size() << std::endl;
65  return false;
66  }
67  return true;
68 }
std::vector< double > max_
Definition: PFMETFilter.h:26
std::vector< double > min_
Definition: PFMETFilter.h:25
std::vector< int > doMax_
Definition: PFMETFilter.h:28
std::vector< std::string > collections_
Definition: PFMETFilter.h:23
std::vector< int > doMin_
Definition: PFMETFilter.h:27
tuple cout
Definition: gather_cfg.py:121
std::vector< std::string > variables_
Definition: PFMETFilter.h:24
void PFMETFilter::endJob ( void  )
virtual

Reimplemented from edm::EDFilter.

Definition at line 282 of file PFMETFilter.cc.

283 {
284  //std::cout << "FL: endJob" << std::endl;
285 }
bool PFMETFilter::filter ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
virtual

Implements edm::EDFilter.

Definition at line 75 of file PFMETFilter.cc.

References checkInput(), collections_, gather_cfg::cout, delta, DeltaMEXsigma_, doMax_, doMin_, reco::Candidate::et(), reco::Candidate::eta(), edm::Event::getByLabel(), M_PI, max_, CaloMET_cfi::met, min_, reco::Candidate::phi(), edm::Handle< T >::product(), reco::Candidate::px(), reco::Candidate::py(), sigma_a_, sigma_b_, sigma_c_, findQualityFiles::size, createPayload::skip, mathSSE::sqrt(), reco::MET::sumEt(), TrueMET_, variables_, and verbose_.

77 {
78  //std::cout << "FL: filter" << std::endl;
79  //std::cout << "FL: Mins = " << min_ << std::endl;
80 
81  if (!checkInput()) return true; // no filtering !
82 
83  bool skip=false;
84 
85  for (unsigned int varc=0;varc<collections_.size();++varc)
86  {
87  //std::cout << "FL: var[" << varc << "] = " << collections_[varc] << std::endl;
88  //std::cout << "FL: var[0] = " << collections_[0] << std::endl;
89 
90  // if the collection is collection1-collection2:
91  const unsigned int minuspos=collections_[varc].find("-");
92  if (minuspos<collections_[varc].size())
93  {
94  std::string collection1;
95  collection1.assign(collections_[varc],0,minuspos);
96  //std::cout << "collection1 = " << collection1 << std::endl;
97  std::string collection2;
98  collection2.assign(collections_[varc],minuspos+1,collections_[varc].size());
99  //std::cout << "collection2 = " << collection2 << std::endl;
100 
101  const edm::View<reco::Candidate> *var1;
103  bool isVar1 = iEvent.getByLabel(collection1, var1_hnd);
104  if ( !isVar1 ) {
105  std::cout << "Warning : no " << collection1 << " in input !" << std::endl;
106  return false;
107  }
108  var1 = var1_hnd.product();
109  const reco::Candidate *var10 = &(*var1)[0];
110  //std::cout << "FL: var10.pt = " << var10->et() << std::endl;
111  double coll_var1;
112  if (variables_[varc]=="et") coll_var1=var10->et();
113  else if (variables_[varc]=="phi") coll_var1=var10->phi();
114  else if (variables_[varc]=="eta") coll_var1=var10->eta();
115  else
116  {
117  std::cout << "Error: PFMETFilter: variable unknown: " << variables_[varc] << std::endl;
118  return true;
119  }
120  //std::cout << "FL: coll_var1[" << variables_[varc] << "] = " << coll_var1 << std::endl;
121 
122  const edm::View<reco::Candidate> *var2;
124  bool isVar2 = iEvent.getByLabel(collection2, var2_hnd);
125  if ( !isVar2 ) {
126  std::cout << "Warning : no " << collection2 << " in input !" << std::endl;
127  return false;
128  }
129  var2 = var2_hnd.product();
130  const reco::Candidate *var20 = &(*var2)[0];
131  //std::cout << "FL: var20.pt = " << var20->et() << std::endl;
132  double coll_var2;
133  if (variables_[varc]=="et") coll_var2=var20->et();
134  else if (variables_[varc]=="phi") coll_var2=var20->phi();
135  else if (variables_[varc]=="eta") coll_var2=var20->eta();
136  else
137  {
138  std::cout << "Error: PFMETFilter: variable unknown: " << variables_[varc] << std::endl;
139  return true;
140  }
141  //std::cout << "FL: coll_var2[" << variables_[varc] << "] = " << coll_var2 << std::endl;
142  //std::cout << "FL: max_[varc] = " << max_[varc] << std::endl;
143  //std::cout << "FL: min_[varc] = " << min_[varc] << std::endl;
144 
145  // Delta computation
146  double delta=coll_var1-coll_var2;
147  if (variables_[varc]=="phi")
148  {
149  if (coll_var1 > M_PI) coll_var1 -= ceil((coll_var1 - M_PI) / (2 * M_PI)) * 2 * M_PI;
150  if (coll_var1 <= - M_PI) coll_var1 += ceil((coll_var1 + M_PI) / (-2. * M_PI)) * 2. * M_PI;
151  if (coll_var2 > M_PI) coll_var2 -= ceil((coll_var2 - M_PI) / (2 * M_PI)) * 2 * M_PI;
152  if (coll_var2 <= - M_PI) coll_var2 += ceil((coll_var2 + M_PI) / (-2. * M_PI)) * 2 * M_PI;
153 
154  double deltaphi=-999.0;
155  if (fabs(coll_var1 - coll_var2)<M_PI)
156  {
157  deltaphi=(coll_var1-coll_var2);
158  }
159  else
160  {
161  if ((coll_var1-coll_var2)>0.0)
162  {
163  deltaphi=(2*M_PI - fabs(coll_var1 - coll_var2));
164  }
165  else
166  {
167  deltaphi=-(2*M_PI - fabs(coll_var1 - coll_var2));
168  }
169  }
170  delta=deltaphi;
171  }
172 
173  // cuts
174  if (doMin_[varc] && doMax_[varc] && max_[varc]<min_[varc])
175  {
176  if (delta>max_[varc] && delta<min_[varc]) skip=true;
177  }
178  else
179  {
180  if (doMin_[varc] && delta<min_[varc]) skip=true;
181  if (doMax_[varc] && delta>max_[varc]) skip=true;
182  }
183  //std::cout << "skip = " << skip << std::endl;
184  }
185  else
186  {
187  // get the variable:
188  const edm::View<reco::Candidate> *var0;
190  bool isVar0 = iEvent.getByLabel(collections_[varc], var0_hnd);
191  if ( !isVar0 ) {
192  std::cout << "Warning : no " << collections_[varc] << " in input !" << std::endl;
193  return false;
194  }
195  var0 = var0_hnd.product();
196  const reco::Candidate *var00 = &(*var0)[0];
197  //std::cout << "FL: var00.pt = " << var00->et() << std::endl;
198  double coll_var;
199  if (variables_[varc]=="et") coll_var=var00->et();
200  else if (variables_[varc]=="phi") coll_var=var00->phi();
201  else if (variables_[varc]=="eta") coll_var=var00->eta();
202  else if (variables_[varc]=="DeltaMEXcut")
203  {
204 
205  const edm::View<reco::Candidate> *truevar0;
207  bool istrueVar0 = iEvent.getByLabel(TrueMET_, truevar0_hnd);
208  if ( !istrueVar0 ) {
209  std::cout << "Warning : no " << TrueMET_ << " in input !" << std::endl;
210  return false;
211  }
212  truevar0 = truevar0_hnd.product();
213  const reco::Candidate *truevar00 = &(*truevar0)[0];
214 
215  const double DeltaMEX=var00->px()-truevar00->px();
216  const double DeltaMEY=var00->py()-truevar00->py();
217  const double cutvalc=sqrt(DeltaMEX*DeltaMEX+DeltaMEY*DeltaMEY);
218  const reco::MET* met=static_cast<const reco::MET*>(truevar00);
219  const double SETc=met->sumEt();
220  //std::cout << "FL: SETc = " << SETc << std::endl;
221  const double sigmac=sigma_a_+sigma_b_*sqrt(SETc)+sigma_c_*SETc;
222  if (cutvalc>DeltaMEXsigma_*sigmac)
223  {
224  if (verbose_)
225  {
226  std::cout << "DeltaMET = " << var00->et()-truevar00->et() << std::endl;
227  std::cout << "trueSET = " << SETc << std::endl;
228  std::cout << "pfMET = " << var00->et() << std::endl;
229  std::cout << "trueMET = " << truevar00->et() << std::endl;
230  std::cout << "DeltaMEX = " << DeltaMEX << std::endl;
231  std::cout << "DeltaMEY = " << DeltaMEY << std::endl;
232  std::cout << "cutvalc = " << cutvalc << std::endl;
233  std::cout << "sigmac = " << sigmac << std::endl;
234  std::cout << "cutvalc/sigmac = " << cutvalc/sigmac << std::endl;
235  }
236  return true;
237  }
238  else
239  {
240  if (verbose_ && (var00->et()-truevar00->et())>300.0)
241  {
242  std::cout << "EVENT NOT KEPT:" << std::endl;
243  std::cout << "DeltaMET = " << var00->et()-truevar00->et() << std::endl;
244  std::cout << "SETc = " << SETc << std::endl;
245  std::cout << "pfMET = " << var00->et() << std::endl;
246  std::cout << "trueMET = " << truevar00->et() << std::endl;
247  std::cout << "DeltaMEX = " << DeltaMEX << std::endl;
248  std::cout << "DeltaMEY = " << DeltaMEY << std::endl;
249  std::cout << "cutvalc = " << cutvalc << std::endl;
250  std::cout << "sigmac = " << sigmac << std::endl;
251  std::cout << "cutvalc/sigmac = " << cutvalc/sigmac << std::endl;
252  }
253  return false;
254  }
255  }
256  else
257  {
258  std::cout << "Error: PFMETFilter: variable unknown: " << variables_[varc] << std::endl;
259  return true;
260  }
261  //std::cout << "FL: coll_var[" << variables_[varc] << "] = " << coll_var << std::endl;
262  //std::cout << "FL: max_[varc] = " << max_[varc] << std::endl;
263  //std::cout << "FL: min_[varc] = " << min_[varc] << std::endl;
264 
265  // cuts
266  if (doMin_[varc] && doMax_[varc] && max_[varc]<min_[varc])
267  {
268  if (coll_var>max_[varc] && coll_var<min_[varc]) skip=true;
269  }
270  else
271  {
272  if (doMin_[varc] && coll_var<min_[varc]) skip=true;
273  if (doMax_[varc] && coll_var>max_[varc]) skip=true;
274  }
275  //std::cout << "skip = " << skip << std::endl;
276  }
277  }
278  //std::cout << "final skip = " << skip << std::endl;
279  return !skip;
280 }
dbl * delta
Definition: mlp_gen.cc:36
virtual double et() const =0
transverse energy
std::vector< double > max_
Definition: PFMETFilter.h:26
double sigma_c_
Definition: PFMETFilter.h:36
std::vector< double > min_
Definition: PFMETFilter.h:25
std::vector< int > doMax_
Definition: PFMETFilter.h:28
double sumEt() const
Definition: MET.h:48
double sigma_a_
Definition: PFMETFilter.h:34
Definition: MET.h:32
T sqrt(T t)
Definition: SSEVec.h:46
bool checkInput()
Definition: PFMETFilter.cc:30
std::string TrueMET_
Definition: PFMETFilter.h:32
double DeltaMEXsigma_
Definition: PFMETFilter.h:33
virtual double py() const =0
y coordinate of momentum vector
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
virtual double px() const =0
x coordinate of momentum vector
#define M_PI
Definition: BFit3D.cc:3
std::vector< std::string > collections_
Definition: PFMETFilter.h:23
std::vector< int > doMin_
Definition: PFMETFilter.h:27
T const * product() const
Definition: Handle.h:74
double sigma_b_
Definition: PFMETFilter.h:35
tuple cout
Definition: gather_cfg.py:121
std::vector< std::string > variables_
Definition: PFMETFilter.h:24
tuple size
Write out results.
bool verbose_
Definition: PFMETFilter.h:37
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity

Member Data Documentation

std::vector<std::string> PFMETFilter::collections_
private

Definition at line 23 of file PFMETFilter.h.

Referenced by checkInput(), filter(), and PFMETFilter().

double PFMETFilter::DeltaMEXsigma_
private

Definition at line 33 of file PFMETFilter.h.

Referenced by filter(), and PFMETFilter().

std::vector<int> PFMETFilter::doMax_
private

Definition at line 28 of file PFMETFilter.h.

Referenced by checkInput(), filter(), and PFMETFilter().

std::vector<int> PFMETFilter::doMin_
private

Definition at line 27 of file PFMETFilter.h.

Referenced by checkInput(), filter(), and PFMETFilter().

std::vector<double> PFMETFilter::max_
private

Definition at line 26 of file PFMETFilter.h.

Referenced by checkInput(), filter(), and PFMETFilter().

std::vector<double> PFMETFilter::min_
private

Definition at line 25 of file PFMETFilter.h.

Referenced by checkInput(), filter(), and PFMETFilter().

double PFMETFilter::sigma_a_
private

Definition at line 34 of file PFMETFilter.h.

Referenced by filter(), and PFMETFilter().

double PFMETFilter::sigma_b_
private

Definition at line 35 of file PFMETFilter.h.

Referenced by filter(), and PFMETFilter().

double PFMETFilter::sigma_c_
private

Definition at line 36 of file PFMETFilter.h.

Referenced by filter(), and PFMETFilter().

std::string PFMETFilter::TrueMET_
private

Definition at line 32 of file PFMETFilter.h.

Referenced by filter(), and PFMETFilter().

std::vector<std::string> PFMETFilter::variables_
private

Definition at line 24 of file PFMETFilter.h.

Referenced by checkInput(), filter(), and PFMETFilter().

bool PFMETFilter::verbose_
private

Definition at line 37 of file PFMETFilter.h.

Referenced by filter(), and PFMETFilter().