CMS 3D CMS Logo

InterestingDetIdCollectionProducer.cc
Go to the documentation of this file.
2 
6 
7 
10 
13 
17 
21 
22 
24 {
25 
26  recHitsToken_ =
27  consumes<EcalRecHitCollection>(iConfig.getParameter< edm::InputTag > ("recHitsLabel"));
29  consumes<reco::BasicClusterCollection>(iConfig.getParameter< edm::InputTag >("basicClustersLabel"));
30 
31  interestingDetIdCollection_ = iConfig.getParameter<std::string>("interestingDetIdCollection");
32 
33  minimalEtaSize_ = iConfig.getParameter<int> ("etaSize");
34  minimalPhiSize_ = iConfig.getParameter<int> ("phiSize");
35  if ( minimalPhiSize_ % 2 == 0 || minimalEtaSize_ % 2 == 0)
36  edm::LogError("InterestingDetIdCollectionProducerError") << "Size of eta/phi should be odd numbers";
37 
38  //register your products
39  produces< DetIdCollection > (interestingDetIdCollection_) ;
40 
41  severityLevel_ = iConfig.getParameter<int>("severityLevel");
42  keepNextToDead_ = iConfig.getParameter<bool>("keepNextToDead");
43  keepNextToBoundary_ = iConfig.getParameter<bool>("keepNextToBoundary");
44 }
45 
46 
48 {
49  edm::ESHandle<CaloTopology> theCaloTopology;
50  iSetup.get<CaloTopologyRecord>().get(theCaloTopology);
51  caloTopology_ = &(*theCaloTopology);
52 
54  iSetup.get<EcalSeverityLevelAlgoRcd>().get(sevLv);
55  severity_ = sevLv.product();
56 }
57 
58 // ------------ method called to produce the data ------------
59 void
61  const edm::EventSetup& iSetup)
62 {
63  using namespace edm;
64  using namespace std;
65 
66  // take BasicClusters
68  iEvent.getByToken(basicClustersToken_, pClusters);
69 
70  // take EcalRecHits
71  Handle<EcalRecHitCollection> recHitsHandle;
72  iEvent.getByToken(recHitsToken_,recHitsHandle);
73 
74  //Create empty output collections
75  std::vector<DetId> indexToStore;
76  indexToStore.reserve(1000);
77 
78  reco::BasicClusterCollection::const_iterator clusIt;
79 
80  std::vector<DetId> xtalsToStore;
81  xtalsToStore.reserve(50);
82  for (clusIt=pClusters->begin(); clusIt!=pClusters->end(); clusIt++) {
83  const std::vector<std::pair<DetId,float> > &clusterDetIds = clusIt->hitsAndFractions();
84  for (const auto &detidpair : clusterDetIds) {
85  indexToStore.push_back(detidpair.first);
86  }
87 
88  //below checks and additional code are only relevant for EB/EE, not for ES
89  if (clusterDetIds.front().first.subdetId()==EcalBarrel || clusterDetIds.front().first.subdetId()==EcalEndcap) {
90  std::vector<std::pair<DetId,float> >::const_iterator posCurrent;
91 
92  float eMax=0.;
93  DetId eMaxId(0);
94 
95  EcalRecHit testEcalRecHit;
96 
97  for(posCurrent = clusterDetIds.begin(); posCurrent != clusterDetIds.end(); posCurrent++)
98  {
99  EcalRecHitCollection::const_iterator itt = recHitsHandle->find((*posCurrent).first);
100  if ((!((*posCurrent).first.null())) && (itt != recHitsHandle->end()) && ((*itt).energy() > eMax) )
101  {
102  eMax = (*itt).energy();
103  eMaxId = (*itt).id();
104  }
105  }
106 
107  if (eMaxId.null())
108  continue;
109 
111 
112  xtalsToStore=topology->getWindow(eMaxId,minimalEtaSize_,minimalPhiSize_);
113 
114  for (const auto &detid : xtalsToStore) {
115  indexToStore.push_back(detid);
116  }
117 
118  }
119  }
120 
121 
123  for (EcalRecHitCollection::const_iterator it = recHitsHandle->begin(); it != recHitsHandle->end(); ++it) {
124  // also add recHits of dead TT if the corresponding TP is saturated
125  if ( it->checkFlag(EcalRecHit::kTPSaturated) ) {
126  indexToStore.push_back(it->id());
127  }
128  // add hits for severities above a threshold
129  if ( severityLevel_>=0 &&
131 
132  indexToStore.push_back(it->id());
133  }
134  if (keepNextToDead_) {
135  // also keep channels next to dead ones
136  if (EcalTools::isNextToDead(it->id(), iSetup)) {
137  indexToStore.push_back(it->id());
138  }
139  }
140 
141  if (keepNextToBoundary_){
142  // keep channels around EB/EE boundary
143  if (it->id().subdetId() == EcalBarrel){
144  EBDetId ebid(it->id());
145  if (abs(ebid.ieta())== 85)
146  indexToStore.push_back(it->id());
147  } else {
148 
149  if (EEDetId::isNextToRingBoundary(it->id()))
150  indexToStore.push_back(it->id());
151  }
152 
153  }
154 
155  }
156  }
157 
158  //unify the vector
159  std::sort(indexToStore.begin(),indexToStore.end());
160  std::unique(indexToStore.begin(),indexToStore.end());
161 
162 
163  iEvent.put(std::make_unique<DetIdCollection>(indexToStore), interestingDetIdCollection_);
164 
165 }
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
static bool isNextToDead(const DetId &id, const edm::EventSetup &es)
true if the channel is near a dead one (in the 3x3)
Definition: EcalTools.cc:52
EcalSeverityLevel::SeverityLevel severityLevel(const DetId &id) const
Evaluate status from id use channelStatus from DB.
CaloTopology const * topology(0)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
constexpr bool null() const
is this a null id ?
Definition: DetId.h:52
std::vector< EcalRecHit >::const_iterator const_iterator
int iEvent
Definition: GenABIO.cc:230
def unique(seq, keepstr=True)
Definition: tier0.py:25
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static bool isNextToRingBoundary(EEDetId id)
Definition: EEDetId.cc:375
edm::EDGetTokenT< EcalRecHitCollection > recHitsToken_
InterestingDetIdCollectionProducer(const edm::ParameterSet &)
ctor
const_iterator end() const
void produce(edm::Event &, const edm::EventSetup &) override
producer
Definition: DetId.h:18
virtual std::vector< DetId > getWindow(const DetId &id, const int &northSouthSize, const int &eastWestSize) const
const CaloSubdetectorTopology * getSubdetectorTopology(const DetId &id) const
access the subdetector Topology for the given subdetector directly
Definition: CaloTopology.cc:25
iterator find(key_type k)
HLT enums.
T get() const
Definition: EventSetup.h:68
void beginRun(edm::Run const &, const edm::EventSetup &) final
T const * product() const
Definition: ESHandle.h:84
const_iterator begin() const
Definition: Run.h:44
edm::EDGetTokenT< reco::BasicClusterCollection > basicClustersToken_
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:39