test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTHcalMETNoiseCleaner.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Class: HLTHcalMETNoiseCleaner
4 //
12 //
13 // Original Author: Alexander Mott
14 // Created: Mon Nov 21 11:32:00 CEST 2011
15 //
16 //
17 //
18 
20 
26 
28 
35 
39 
40 #include <iostream>
41 #include <string>
42 #include <fstream>
43 #include <TVector3.h>
44 #include <TLorentzVector.h>
45 //#include <Point.h>
46 
48  : HcalNoiseRBXCollectionTag_(iConfig.getParameter<edm::InputTag>("HcalNoiseRBXCollection")),
49  CaloMetCollectionTag_(iConfig.getParameter<edm::InputTag>("CaloMetCollection")),
50  CaloMetCut_(iConfig.getParameter<double>("CaloMetCut")),
51  severity_(iConfig.getParameter<int> ("severity")),
52  maxNumRBXs_(iConfig.getParameter<int>("maxNumRBXs")),
53  numRBXsToConsider_(iConfig.getParameter<int>("numRBXsToConsider")),
54  accept2NoiseRBXEvents_(iConfig.getParameter<bool>("accept2NoiseRBXEvents")),
55  needEMFCoincidence_(iConfig.getParameter<bool>("needEMFCoincidence")),
56  minRBXEnergy_(iConfig.getParameter<double>("minRBXEnergy")),
57  minRatio_(iConfig.getParameter<double>("minRatio")),
58  maxRatio_(iConfig.getParameter<double>("maxRatio")),
59  minHPDHits_(iConfig.getParameter<int>("minHPDHits")),
60  minRBXHits_(iConfig.getParameter<int>("minRBXHits")),
61  minHPDNoOtherHits_(iConfig.getParameter<int>("minHPDNoOtherHits")),
62  minZeros_(iConfig.getParameter<int>("minZeros")),
63  minHighEHitTime_(iConfig.getParameter<double>("minHighEHitTime")),
64  maxHighEHitTime_(iConfig.getParameter<double>("maxHighEHitTime")),
65  maxRBXEMF_(iConfig.getParameter<double>("maxRBXEMF")),
66  minRecHitE_(iConfig.getParameter<double>("minRecHitE")),
67  minLowHitE_(iConfig.getParameter<double>("minLowHitE")),
68  minHighHitE_(iConfig.getParameter<double>("minHighHitE")),
69  minR45HitE_(5.0),
70  TS4TS5EnergyThreshold_(iConfig.getParameter<double>("TS4TS5EnergyThreshold"))
71 {
72 
73  std::vector<double> TS4TS5UpperThresholdTemp = iConfig.getParameter<std::vector<double> >("TS4TS5UpperThreshold");
74  std::vector<double> TS4TS5UpperCutTemp = iConfig.getParameter<std::vector<double> >("TS4TS5UpperCut");
75  std::vector<double> TS4TS5LowerThresholdTemp = iConfig.getParameter<std::vector<double> >("TS4TS5LowerThreshold");
76  std::vector<double> TS4TS5LowerCutTemp = iConfig.getParameter<std::vector<double> >("TS4TS5LowerCut");
77 
78  for(int i = 0; i < (int)TS4TS5UpperThresholdTemp.size() && i < (int)TS4TS5UpperCutTemp.size(); i++)
79  TS4TS5UpperCut_.push_back(std::pair<double, double>(TS4TS5UpperThresholdTemp[i], TS4TS5UpperCutTemp[i]));
80  sort(TS4TS5UpperCut_.begin(), TS4TS5UpperCut_.end());
81 
82  for(int i = 0; i < (int)TS4TS5LowerThresholdTemp.size() && i < (int)TS4TS5LowerCutTemp.size(); i++)
83  TS4TS5LowerCut_.push_back(std::pair<double, double>(TS4TS5LowerThresholdTemp[i], TS4TS5LowerCutTemp[i]));
84  sort(TS4TS5LowerCut_.begin(), TS4TS5LowerCut_.end());
85 
86  m_theCaloMetToken = consumes<reco::CaloMETCollection>(CaloMetCollectionTag_);
87  m_theHcalNoiseToken = consumes<reco::HcalNoiseRBXCollection>(HcalNoiseRBXCollectionTag_);
88 
89  if(iConfig.existsAs<double>("minR45HitE"))
90  minR45HitE_ = iConfig.getParameter<double>("minR45HitE");
91 
92  produces<reco::CaloMETCollection>();
93 }
94 
95 
97 
98 void
101  desc.add<edm::InputTag>("HcalNoiseRBXCollection",edm::InputTag("hltHcalNoiseInfoProducer"));
102  desc.add<edm::InputTag>("CaloMetCollection",edm::InputTag("hltMet"));
103  desc.add<double>("CaloMetCut",0.0);
104  desc.add<int>("severity",1);
105  desc.add<int>("maxNumRBXs",2);
106  desc.add<int>("numRBXsToConsider",2);
107  desc.add<bool>("accept2NoiseRBXEvents",true);
108  desc.add<bool>("needEMFCoincidence",true);
109  desc.add<double>("minRBXEnergy",50.0);
110  desc.add<double>("minRatio",-999.);
111  desc.add<double>("maxRatio",999.);
112  desc.add<int>("minHPDHits",17);
113  desc.add<int>("minRBXHits",999);
114  desc.add<int>("minHPDNoOtherHits",10);
115  desc.add<int>("minZeros",10);
116  desc.add<double>("minHighEHitTime",-9999.0);
117  desc.add<double>("maxHighEHitTime",9999.0);
118  desc.add<double>("maxRBXEMF",0.02);
119  desc.add<double>("minRecHitE",1.5);
120  desc.add<double>("minLowHitE",10.0);
121  desc.add<double>("minHighHitE",25.0);
122  desc.add<double>("minR45HitE",5.0);
123  desc.add<double>("TS4TS5EnergyThreshold",50.0);
124 
125  double TS4TS5UpperThresholdArray[5] = {70, 90, 100, 400, 4000 };
126  double TS4TS5UpperCutArray[5] = {1, 0.8, 0.75, 0.72, 0.72};
127  double TS4TS5LowerThresholdArray[7] = {100, 120, 150, 200, 300, 400, 500};
128  double TS4TS5LowerCutArray[7] = {-1, -0.7, -0.4, -0.2, -0.08, 0, 0.1};
129  std::vector<double> TS4TS5UpperThreshold(TS4TS5UpperThresholdArray, TS4TS5UpperThresholdArray+5);
130  std::vector<double> TS4TS5UpperCut(TS4TS5UpperCutArray, TS4TS5UpperCutArray+5);
131  std::vector<double> TS4TS5LowerThreshold(TS4TS5LowerThresholdArray, TS4TS5LowerThresholdArray+7);
132  std::vector<double> TS4TS5LowerCut(TS4TS5LowerCutArray, TS4TS5LowerCutArray+7);
133 
134  desc.add<std::vector<double> >("TS4TS5UpperThreshold", TS4TS5UpperThreshold);
135  desc.add<std::vector<double> >("TS4TS5UpperCut", TS4TS5UpperCut);
136  desc.add<std::vector<double> >("TS4TS5LowerThreshold", TS4TS5LowerThreshold);
137  desc.add<std::vector<double> >("TS4TS5LowerCut", TS4TS5LowerCut);
138  descriptions.add("hltHcalMETNoiseCleaner",desc);
139 }
140 
141 //
142 // member functions
143 //
144 
146 {
147  using namespace reco;
148 
149  //output collection
150  std::auto_ptr<CaloMETCollection> CleanedMET(new CaloMETCollection);
151 
152  //get the calo MET / MHT
154  iEvent.getByToken(m_theCaloMetToken,met_h);
155 
156  if(not met_h.isValid() or met_h->size()==0 or met_h->front().pt()<0){ //No Valid MET, don't do anything and accept the event
157  return true; // we shouldn't get here, but lets not crash
158  }
159 
160  reco::CaloMET inCaloMet = met_h->front();
161 
162 
163  // in this case, do not filter anything
164  if(severity_==0){
165  CleanedMET->push_back(inCaloMet);
166  iEvent.put(CleanedMET);
167  return true;
168  }
169 
170  // get the RBXs produced by RecoMET/METProducers/HcalNoiseInfoProducer
172  iEvent.getByToken(m_theHcalNoiseToken,rbxs_h);
173  if(!rbxs_h.isValid()) {
174  edm::LogError("DataNotFound") << "HLTHcalMETNoiseCleaner: Could not find HcalNoiseRBXCollection product named "
175  << HcalNoiseRBXCollectionTag_ << "." << std::endl;
176  CleanedMET->push_back(inCaloMet);
177  iEvent.put(CleanedMET);
178  return true; // no valid RBXs
179  }
180 
181  // create a sorted set of the RBXs, ordered by energy
183  for(HcalNoiseRBXCollection::const_iterator it=rbxs_h->begin(); it!=rbxs_h->end(); ++it) {
184  const HcalNoiseRBX &rbx=(*it);
187  data.insert(d);
188  }
189  //if 0 RBXs are in the list, just accept
190  if(data.size()<1){
191  CleanedMET->push_back(inCaloMet);
192  iEvent.put(CleanedMET);
193  return true;
194  }
195  // data is now sorted by RBX energy
196  // only consider top N=numRBXsToConsider_ energy RBXs
197  int cntr=0;
198  int nNoise=0;
199 
200  TVector3 metVec;
201  metVec.SetPtEtaPhi(met_h->front().pt(), 0, met_h->front().phi() );
202 
203  TVector3 noiseHPDVector(0,0,0);
204  TVector3 secondHPDVector(0,0,0);
205  for(noisedataset_t::const_iterator it=data.begin();
206  it!=data.end() && cntr<numRBXsToConsider_;
207  it++, cntr++) {
208  bool isNoise=false;
209  bool passFilter=true;
210  bool passEMF=true;
211  if(it->energy()>minRBXEnergy_) {
212  if(it->validRatio() && it->ratio()<minRatio_) passFilter=false;
213  else if(it->validRatio() && it->ratio()>maxRatio_) passFilter=false;
214  else if(it->numHPDHits()>=minHPDHits_) passFilter=false;
215  else if(it->numRBXHits()>=minRBXHits_) passFilter=false;
216  else if(it->numHPDNoOtherHits()>=minHPDNoOtherHits_) passFilter=false;
217  else if(it->numZeros()>=minZeros_) passFilter=false;
218  else if(it->minHighEHitTime()<minHighEHitTime_) passFilter=false;
219  else if(it->maxHighEHitTime()>maxHighEHitTime_) passFilter=false;
220  else if(!it->PassTS4TS5()) passFilter=false;
221 
222  if(it->RBXEMF()<maxRBXEMF_){
223  passEMF=false;
224  }
225  }
226 
227  if((needEMFCoincidence_ && !passEMF && !passFilter) ||
228  (!needEMFCoincidence_ && !passFilter)) { // check for noise
229  LogDebug("") << "HLTHcalMETNoiseCleaner debug: Found a noisy RBX: "
230  << "energy=" << it->energy() << "; "
231  << "ratio=" << it->ratio() << "; "
232  << "# RBX hits=" << it->numRBXHits() << "; "
233  << "# HPD hits=" << it->numHPDHits() << "; "
234  << "# Zeros=" << it->numZeros() << "; "
235  << "min time=" << it->minHighEHitTime() << "; "
236  << "max time=" << it->maxHighEHitTime() << "; "
237  << "passTS4TS5=" << it->PassTS4TS5() << "; "
238  << "RBX EMF=" << it->RBXEMF()
239  << std::endl;
240  nNoise++;
241  isNoise=true;
242  }// OK, checked for noise
243 
244  //------------First Noisy RBX-----------------------
245  if(isNoise && nNoise==1){
246  edm::RefVector<CaloTowerCollection> noiseTowers = it->rbxTowers();
248  // get the energy vector for this RBX from the calotowers
249  for( noiseTowersIt = noiseTowers.begin(); noiseTowersIt != noiseTowers.end(); noiseTowersIt++){
250  TVector3 towerVec;
251  towerVec.SetPtEtaPhi((*noiseTowersIt)->pt(),(*noiseTowersIt)->eta(),(*noiseTowersIt)->phi());
252  noiseHPDVector+=towerVec; // add this tower to the vector for the RBX
253  }
254  if(noiseHPDVector.Mag()>0) noiseHPDVector.SetPtEtaPhi(noiseHPDVector.Pt(),0,noiseHPDVector.Phi()); // make the noise transverse
255  else noiseHPDVector.SetPtEtaPhi(0,0,0);
256  }
257  //-----------FOUND a SECOND NOISY RBX-------------------
258  if(isNoise && cntr > 0){
259  CleanedMET->push_back(inCaloMet);
260  iEvent.put(CleanedMET);
261  return accept2NoiseRBXEvents_; // don't try to clean these for the moment, just keep or throw away
262  }
263  //----------LEADING RBX is NOT NOISY--------------------
264  if(!isNoise && cntr == 0){
265  CleanedMET->push_back(inCaloMet);
266  iEvent.put(CleanedMET);
267  return true; // don't reject the event if the leading RBX isn't noise
268  }
269  //-----------SUBLEADING RBX is NOT NOISY: STORE INFO----
270  if(!isNoise && nNoise>0){ //second RBX isn't noisy (and first one was), so clean
271  edm::RefVector<CaloTowerCollection> noiseTowers = it->rbxTowers();
273  for( noiseTowersIt = noiseTowers.begin(); noiseTowersIt != noiseTowers.end(); noiseTowersIt++){
274  TVector3 towerVec;
275  towerVec.SetPtEtaPhi((*noiseTowersIt)->pt(),(*noiseTowersIt)->eta(),(*noiseTowersIt)->phi());
276  secondHPDVector+=towerVec; // add this tower to the vector for the RBX
277  }
278  if(secondHPDVector.Mag()>0) secondHPDVector.SetPtEtaPhi(secondHPDVector.Pt(),0,secondHPDVector.Phi()); // make the second transverse
279  else secondHPDVector.SetPtEtaPhi(0,0,0);
280  break;
281  }
282  } // end RBX loop
283 
284  if(noiseHPDVector.Mag()==0){
285  CleanedMET->push_back(inCaloMet);
286  iEvent.put(CleanedMET);
287  return true; // don't reject the event if the leading RBX isn't noise
288  }
289 
290  //********************************************************************************
291  //The Event gets here only if it had exactly 1 noisy RBX in the lead position
292  //********************************************************************************
293 
294  float METsumet = met_h->front().energy();
295 
296  metVec+=noiseHPDVector;
297 
298  float ZMETsumet = METsumet-noiseHPDVector.Mag();
299  float ZMETpt = metVec.Pt();
300  float ZMETphi = metVec.Phi();
301 
302  //put the second RBX vector in the eta phi position of the leading RBX vector
303 
304  float SMETsumet = 0;
305  float SMETpt = 0;
306  float SMETphi = 0;
307  if(secondHPDVector.Mag()>0.){
308  secondHPDVector.SetPtEtaPhi(secondHPDVector.Pt(),noiseHPDVector.Eta(),noiseHPDVector.Phi());
309  metVec-= secondHPDVector;
310  SMETsumet = METsumet-noiseHPDVector.Mag();
311  SMETpt = metVec.Pt();
312  SMETphi = metVec.Phi();
313  }
314  //Get the maximum MET:
315  float CorMetSumEt,CorMetPt,CorMetPhi;
316  if(ZMETpt>SMETpt){
317  CorMetSumEt = ZMETsumet;
318  CorMetPt = ZMETpt;
319  CorMetPhi = ZMETphi;
320  }else{
321  CorMetSumEt = SMETsumet;
322  CorMetPt = SMETpt;
323  CorMetPhi = SMETphi;
324  }
325 
326  reco::CaloMET corMet = BuildCaloMet(CorMetSumEt,CorMetPt,CorMetPhi);
327  CleanedMET->push_back(corMet);
328  iEvent.put(CleanedMET);
329 
330  return (corMet.pt() > CaloMetCut_);
331 }
332 
334  // Instantiate the container to hold the calorimeter specific information
335 
336  typedef math::XYZPoint Point;
338 
340  // Initialise the container
341  specific.MaxEtInEmTowers = 0.0; // Maximum energy in EM towers
342  specific.MaxEtInHadTowers = 0.0; // Maximum energy in HCAL towers
343  specific.HadEtInHO = 0.0; // Hadronic energy fraction in HO
344  specific.HadEtInHB = 0.0; // Hadronic energy in HB
345  specific.HadEtInHF = 0.0; // Hadronic energy in HF
346  specific.HadEtInHE = 0.0; // Hadronic energy in HE
347  specific.EmEtInEB = 0.0; // Em energy in EB
348  specific.EmEtInEE = 0.0; // Em energy in EE
349  specific.EmEtInHF = 0.0; // Em energy in HF
350  specific.EtFractionHadronic = 0.0; // Hadronic energy fraction
351  specific.EtFractionEm = 0.0; // Em energy fraction
352  specific.CaloSETInpHF = 0.0; // CaloSET in HF+
353  specific.CaloSETInmHF = 0.0; // CaloSET in HF-
354  specific.CaloMETInpHF = 0.0; // CaloMET in HF+
355  specific.CaloMETInmHF = 0.0; // CaloMET in HF-
356  specific.CaloMETPhiInpHF = 0.0; // CaloMET-phi in HF+
357  specific.CaloMETPhiInmHF = 0.0; // CaloMET-phi in HF-
358  specific.METSignificance = 0.0;
359 
360  TLorentzVector p4TL;
361  p4TL.SetPtEtaPhiM(pt,0.,phi,0.);
362  const LorentzVector p4(p4TL.X(),p4TL.Y(),0,p4TL.T());
363  const Point vtx( 0.0, 0.0, 0.0 );
364  reco::CaloMET specificmet( specific, sumet, p4, vtx );
365  return specificmet;
366  }
367 
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
dictionary specific
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventIDconst &, edm::Timestampconst & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
std::vector< std::pair< double, double > > TS4TS5UpperCut_
edm::EDGetTokenT< reco::CaloMETCollection > m_theCaloMetToken
std::set< CommonHcalNoiseRBXData, noisedatacomp > noisedataset_t
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:253
std::pair< double, double > Point
Definition: CaloEllipse.h:18
reco::CaloMET BuildCaloMet(float sumet, float pt, float phi)
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:248
virtual bool filter(edm::Event &, const edm::EventSetup &)
Collection of Calo MET.
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
tuple d
Definition: ztail.py:151
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
double p4[4]
Definition: TauolaWrapper.h:92
HLTHcalMETNoiseCleaner(const edm::ParameterSet &)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::InputTag HcalNoiseRBXCollectionTag_
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
edm::EDGetTokenT< reco::HcalNoiseRBXCollection > m_theHcalNoiseToken
void add(std::string const &label, ParameterSetDescription const &psetDescription)
math::PtEtaPhiELorentzVectorF LorentzVector
std::vector< std::pair< double, double > > TS4TS5LowerCut_