CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
InterestingDetIdCollectionProducer.cc
Go to the documentation of this file.
2 
6 
8 
11 
15 
19 
23 
24 // $Id: InterestingDetIdCollectionProducer.cc,v 1.10 2011/05/19 14:29:48 argiro Exp $
25 
27 {
28 
29  recHitsLabel_ = iConfig.getParameter< edm::InputTag > ("recHitsLabel");
30  basicClusters_ = iConfig.getParameter< edm::InputTag > ("basicClustersLabel");
31 
32  interestingDetIdCollection_ = iConfig.getParameter<std::string>("interestingDetIdCollection");
33 
34  minimalEtaSize_ = iConfig.getParameter<int> ("etaSize");
35  minimalPhiSize_ = iConfig.getParameter<int> ("phiSize");
36  if ( minimalPhiSize_ % 2 == 0 || minimalEtaSize_ % 2 == 0)
37  edm::LogError("InterestingDetIdCollectionProducerError") << "Size of eta/phi should be odd numbers";
38 
39  //register your products
40  produces< DetIdCollection > (interestingDetIdCollection_) ;
41 
42  severityLevel_ = iConfig.getParameter<int>("severityLevel");
43  keepNextToDead_ = iConfig.getParameter<bool>("keepNextToDead");
44  keepNextToBoundary_ = iConfig.getParameter<bool>("keepNextToBoundary");
45 }
46 
47 
49 {}
50 
52 {
53  edm::ESHandle<CaloTopology> theCaloTopology;
54  iSetup.get<CaloTopologyRecord>().get(theCaloTopology);
55  caloTopology_ = &(*theCaloTopology);
56 
58  iSetup.get<EcalSeverityLevelAlgoRcd>().get(sevLv);
59  severity_ = sevLv.product();
60 }
61 
62 // ------------ method called to produce the data ------------
63 void
65  const edm::EventSetup& iSetup)
66 {
67  using namespace edm;
68  using namespace std;
69 
70  // take BasicClusters
72  iEvent.getByLabel(basicClusters_, pClusters);
73 
74  // take EcalRecHits
75  Handle<EcalRecHitCollection> recHitsHandle;
76  iEvent.getByLabel(recHitsLabel_,recHitsHandle);
77 
78  //Create empty output collections
79  std::vector<DetId> indexToStore;
80  indexToStore.reserve(1000);
81 
82  reco::BasicClusterCollection::const_iterator clusIt;
83 
84  std::vector<DetId> xtalsToStore;
85  xtalsToStore.reserve(50);
86  for (clusIt=pClusters->begin(); clusIt!=pClusters->end(); clusIt++) {
87  //PG barrel
88 
89  float eMax=0.;
90  DetId eMaxId(0);
91 
92  std::vector<std::pair<DetId,float> > clusterDetIds = (*clusIt).hitsAndFractions();
93  std::vector<std::pair<DetId,float> >::iterator posCurrent;
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 
110  const CaloSubdetectorTopology* topology = caloTopology_->getSubdetectorTopology(eMaxId.det(),eMaxId.subdetId());
111 
112  xtalsToStore=topology->getWindow(eMaxId,minimalEtaSize_,minimalPhiSize_);
113  std::vector<std::pair<DetId,float > > xtalsInClus=(*clusIt).hitsAndFractions();
114 
115  for (unsigned int ii=0;ii<xtalsInClus.size();ii++)
116  {
117  xtalsToStore.push_back(xtalsInClus[ii].first);
118  }
119 
120  indexToStore.insert(indexToStore.end(),xtalsToStore.begin(),xtalsToStore.end());
121  }
122 
123 
124 
125  for (EcalRecHitCollection::const_iterator it = recHitsHandle->begin(); it != recHitsHandle->end(); ++it) {
126  // also add recHits of dead TT if the corresponding TP is saturated
127  if ( it->checkFlag(EcalRecHit::kTPSaturated) ) {
128  indexToStore.push_back(it->id());
129  }
130  // add hits for severities above a threshold
131  if ( severityLevel_>=0 &&
133 
134  indexToStore.push_back(it->id());
135  }
136  if (keepNextToDead_) {
137  // also keep channels next to dead ones
138  if (EcalTools::isNextToDead(it->id(), iSetup)) {
139  indexToStore.push_back(it->id());
140  }
141  }
142 
143  if (keepNextToBoundary_){
144  // keep channels around EB/EE boundary
145  if (it->id().subdetId() == EcalBarrel){
146  EBDetId ebid(it->id());
147  if (abs(ebid.ieta())== 85)
148  indexToStore.push_back(it->id());
149  } else {
150 
151  if (EEDetId::isNextToRingBoundary(it->id()))
152  indexToStore.push_back(it->id());
153  }
154 
155  }
156 
157  }
158 
159  //unify the vector
160  std::sort(indexToStore.begin(),indexToStore.end());
161  std::unique(indexToStore.begin(),indexToStore.end());
162 
163  std::auto_ptr< DetIdCollection > detIdCollection (new DetIdCollection(indexToStore) ) ;
164 
165 
166  iEvent.put( detIdCollection, interestingDetIdCollection_ );
167 
168 }
T getParameter(std::string const &) const
virtual void produce(edm::Event &, const edm::EventSetup &)
producer
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:53
std::vector< EcalRecHit >::const_iterator const_iterator
#define abs(x)
Definition: mlp_lapack.h:159
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
float energy() const
Definition: CaloRecHit.h:19
static bool isNextToRingBoundary(EEDetId id)
Definition: EEDetId.cc:391
InterestingDetIdCollectionProducer(const edm::ParameterSet &)
ctor
bool first
Definition: L1TdeRCT.cc:94
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
Definition: DetId.h:20
virtual std::vector< DetId > getWindow(const DetId &id, const int &northSouthSize, const int &eastWestSize) const
void beginRun(edm::Run &, const edm::EventSetup &)
const T & get() const
Definition: EventSetup.h:55
bool null() const
is this a null id ?
Definition: DetId.h:47
T const * product() const
Definition: ESHandle.h:62
const CaloSubdetectorTopology * getSubdetectorTopology(const DetId &id) const
access the subdetector Topology for the given subdetector directly
Definition: CaloTopology.cc:26
edm::EDCollection< DetId > DetIdCollection
Detector det() const
get the detector field from this detid
Definition: DetId.h:37
EcalSeverityLevel::SeverityLevel severityLevel(const DetId &id, const EcalRecHitCollection &rhs) const
Evaluate status from id.
Definition: Run.h:33