test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFRecoTauDiscriminationAgainstElectronDeadECAL.cc
Go to the documentation of this file.
1 
30 
31 #include <TMath.h>
32 
33 using namespace reco;
34 
36 {
37  public:
40  moduleLabel_(cfg.getParameter<std::string>("@module_label")),
41  isFirstEvent_(true)
42  {
43  minStatus_ = cfg.getParameter<uint32_t>("minStatus");
44  dR_ = cfg.getParameter<double>("dR");
45 
46  verbosity_ = ( cfg.exists("verbosity") ) ?
47  cfg.getParameter<int>("verbosity") : 0;
48  }
50 
51  void beginEvent(const edm::Event& evt, const edm::EventSetup& es) override
52  {
53  updateBadTowers(es);
54  }
55 
56  double discriminate(const PFTauRef& pfTau) const override
57  {
58  if ( verbosity_ ) {
59  edm::LogPrint("PFTauAgainstEleDeadECAL") << "<PFRecoTauDiscriminationAgainstElectronDeadECAL::discriminate>:" ;
60  edm::LogPrint("PFTauAgainstEleDeadECAL") << " moduleLabel = " << moduleLabel_ ;
61  edm::LogPrint("PFTauAgainstEleDeadECAL") << "#badTowers = " << badTowers_.size() ;
62  edm::LogPrint("PFTauAgainstEleDeadECAL") << "tau: Pt = " << pfTau->pt() << ", eta = " << pfTau->eta() << ", phi = " << pfTau->phi() ;
63  }
64  double discriminator = 1.;
65  for ( std::vector<towerInfo>::const_iterator badTower = badTowers_.begin();
66  badTower != badTowers_.end(); ++badTower ) {
67  if ( deltaR(badTower->eta_, badTower->phi_, pfTau->eta(), pfTau->phi()) < dR_ ) {
68  if ( verbosity_ ) {
69  edm::LogPrint("PFTauAgainstEleDeadECAL") << " matches badTower: eta = " << badTower->eta_ << ", phi = " << badTower->phi_ ;
70  }
71  discriminator = 0.;
72  }
73  }
74  if ( verbosity_ ) {
75  edm::LogPrint("PFTauAgainstEleDeadECAL") << "--> discriminator = " << discriminator ;
76  }
77  return discriminator;
78  }
79 
80  private:
82  {
83  // NOTE: modified version of SUSY CAF code
84  // UserCode/SusyCAF/plugins/SusyCAF_EcalDeadChannels.cc
85  const uint32_t channelStatusId = es.get<EcalChannelStatusRcd>().cacheIdentifier();
86  const uint32_t caloGeometryId = es.get<CaloGeometryRecord>().cacheIdentifier();
87  const uint32_t idealGeometryId = es.get<IdealGeometryRecord>().cacheIdentifier();
88 
89  if ( !isFirstEvent_ && channelStatusId == channelStatusId_cache_ && caloGeometryId == caloGeometryId_cache_ && idealGeometryId == idealGeometryId_cache_ ) return;
90 
91  edm::ESHandle<EcalChannelStatus> channelStatus;
92  es.get<EcalChannelStatusRcd>().get(channelStatus);
93  channelStatusId_cache_ = channelStatusId;
94 
95  edm::ESHandle<CaloGeometry> caloGeometry;
96  es.get<CaloGeometryRecord>().get(caloGeometry);
97  caloGeometryId_cache_ = caloGeometryId;
98 
100  es.get<IdealGeometryRecord>().get(ttMap);
101  idealGeometryId_cache_ = idealGeometryId;
102 
103  std::map<uint32_t,unsigned> nBadCrystals, maxStatus;
104  std::map<uint32_t,double> sumEta, sumPhi;
105 
106  loopXtals<EBDetId>(nBadCrystals, maxStatus, sumEta, sumPhi, channelStatus.product(), caloGeometry.product(), ttMap.product());
107  loopXtals<EEDetId>(nBadCrystals, maxStatus, sumEta, sumPhi, channelStatus.product(), caloGeometry.product(), ttMap.product());
108 
109  badTowers_.clear();
110  for ( std::map<uint32_t, unsigned>::const_iterator it = nBadCrystals.begin();
111  it != nBadCrystals.end(); ++it ) {
112  uint32_t key = it->first;
113  badTowers_.push_back(towerInfo(key, it->second, maxStatus[key], sumEta[key]/it->second, sumPhi[key]/it->second));
114  }
115 
116  isFirstEvent_ = false;
117  }
118 
119  template <class Id>
120  void loopXtals(std::map<uint32_t, unsigned>& nBadCrystals,
121  std::map<uint32_t, unsigned>& maxStatus,
122  std::map<uint32_t, double>& sumEta,
123  std::map<uint32_t, double>& sumPhi ,
124  const EcalChannelStatus* channelStatus,
125  const CaloGeometry* caloGeometry,
126  const EcalTrigTowerConstituentsMap* ttMap) const
127  {
128  // NOTE: modified version of SUSY CAF code
129  // UserCode/SusyCAF/plugins/SusyCAF_EcalDeadChannels.cc
130  for ( int i = 0; i < Id::kSizeForDenseIndexing; ++i ) {
131  Id id = Id::unhashIndex(i);
132  if ( id == Id(0) ) continue;
133  EcalChannelStatusMap::const_iterator it = channelStatus->getMap().find(id.rawId());
134  unsigned status = ( it == channelStatus->end() ) ?
135  0 : (it->getStatusCode() & statusMask_);
136  if ( status >= minStatus_ ) {
137  const GlobalPoint& point = caloGeometry->getPosition(id);
138  uint32_t key = ttMap->towerOf(id);
139  maxStatus[key] = TMath::Max(status, maxStatus[key]);
140  ++nBadCrystals[key];
141  sumEta[key] += point.eta();
142  sumPhi[key] += point.phi();
143  }
144  }
145  }
146 
147  struct towerInfo
148  {
149  towerInfo(uint32_t id, unsigned nBad, unsigned maxStatus, double eta, double phi)
150  : id_(id),
151  nBad_(nBad),
152  maxStatus_(maxStatus),
153  eta_(eta),
154  phi_(phi)
155  {}
156  uint32_t id_;
157  unsigned nBad_;
158  unsigned maxStatus_;
159  double eta_;
160  double phi_;
161  };
162  typedef ROOT::Math::LorentzVector<ROOT::Math::PtEtaPhiE4D<double> > PolarLorentzVector;
163 
165  unsigned minStatus_;
166  double dR_;
167 
168  std::vector<towerInfo> badTowers_;
169  static const uint16_t statusMask_ = 0x1F;
170 
175 
177 };
178 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
tuple cfg
Definition: looper.py:293
const self & getMap() const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
ROOT::Math::LorentzVector< ROOT::Math::PtEtaPhiE4D< double > > PolarLorentzVector
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
bool exists(std::string const &parameterName) const
checks if a parameter exists
tuple discriminator
EcalTrigTowerDetId towerOf(const DetId &id) const
Get the tower id for this det id (or null if not known)
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
const GlobalPoint & getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:68
void loopXtals(std::map< uint32_t, unsigned > &nBadCrystals, std::map< uint32_t, unsigned > &maxStatus, std::map< uint32_t, double > &sumEta, std::map< uint32_t, double > &sumPhi, const EcalChannelStatus *channelStatus, const CaloGeometry *caloGeometry, const EcalTrigTowerConstituentsMap *ttMap) const
string key
FastSim: produces sample of signal events, overlayed with premixed minbias events.
T Max(T a, T b)
Definition: MathUtil.h:44
const T & get() const
Definition: EventSetup.h:56
std::vector< Item >::const_iterator const_iterator
T const * product() const
Definition: ESHandle.h:86
T eta() const
Definition: PV3DBase.h:76
const_iterator find(uint32_t rawId) const
const_iterator end() const
moduleLabel_(iConfig.getParameter< string >("@module_label"))
towerInfo(uint32_t id, unsigned nBad, unsigned maxStatus, double eta, double phi)
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
void beginEvent(const edm::Event &evt, const edm::EventSetup &es) override
tuple status
Definition: mps_update.py:57