CMS 3D CMS Logo

EcalDigiSelector.cc
Go to the documentation of this file.
1 #include <vector>
17 #include <TMath.h>
18 #include <iostream>
19 #include <set>
20 
21 using namespace reco;
23  selectedEcalEBDigiCollection_ = ps.getParameter<std::string>("selectedEcalEBDigiCollection");
24  selectedEcalEEDigiCollection_ = ps.getParameter<std::string>("selectedEcalEEDigiCollection");
25 
26  barrelSuperClusterProducer_ =
27  consumes<reco::SuperClusterCollection>(ps.getParameter<edm::InputTag>("barrelSuperClusterProducer"));
28  endcapSuperClusterProducer_ =
29  consumes<reco::SuperClusterCollection>(ps.getParameter<edm::InputTag>("endcapSuperClusterProducer"));
30 
31  EcalEBRecHitToken_ = consumes<EcalRecHitCollection>(ps.getParameter<edm::InputTag>("EcalEBRecHitTag"));
32  EcalEERecHitToken_ = consumes<EcalRecHitCollection>(ps.getParameter<edm::InputTag>("EcalEERecHitTag"));
33 
34  EcalEBDigiToken_ = consumes<EBDigiCollection>(ps.getParameter<edm::InputTag>("EcalEBDigiTag"));
35  EcalEEDigiToken_ = consumes<EEDigiCollection>(ps.getParameter<edm::InputTag>("EcalEEDigiTag"));
36 
37  cluster_pt_thresh_ = ps.getParameter<double>("cluster_pt_thresh");
38  single_cluster_thresh_ = ps.getParameter<double>("single_cluster_thresh");
39 
40  nclus_sel_ = ps.getParameter<int>("nclus_sel");
41  produces<EBDigiCollection>(selectedEcalEBDigiCollection_);
42  produces<EEDigiCollection>(selectedEcalEEDigiCollection_);
43 }
44 
46  //Get BarrelSuperClusters to start.
47  edm::Handle<reco::SuperClusterCollection> pBarrelSuperClusters;
48 
49  evt.getByToken(barrelSuperClusterProducer_, pBarrelSuperClusters);
50 
51  const reco::SuperClusterCollection& BarrelSuperClusters = *pBarrelSuperClusters;
52  //Got BarrelSuperClusters
53 
54  //Get BarrelSuperClusters to start.
55  edm::Handle<reco::SuperClusterCollection> pEndcapSuperClusters;
56 
57  evt.getByToken(endcapSuperClusterProducer_, pEndcapSuperClusters);
58 
59  const reco::SuperClusterCollection& EndcapSuperClusters = *pEndcapSuperClusters;
60  //Got EndcapSuperClusters
61 
62  reco::SuperClusterCollection saveBarrelSuperClusters;
63  reco::SuperClusterCollection saveEndcapSuperClusters;
64  bool meet_single_thresh = false;
65  //Loop over barrel superclusters, and apply threshold
66  for (int loop = 0; loop < int(BarrelSuperClusters.size()); loop++) {
67  SuperCluster clus1 = BarrelSuperClusters[loop];
68  float eta1 = clus1.eta();
69  float energy1 = clus1.energy();
70  float theta1 = 2 * atan(exp(-1. * eta1));
71  float cluspt1 = energy1 * sin(theta1);
72  if (cluspt1 > cluster_pt_thresh_) {
73  saveBarrelSuperClusters.push_back(clus1);
74  if (cluspt1 > single_cluster_thresh_)
75  meet_single_thresh = true;
76  }
77  }
78 
79  //Loop over endcap superclusters, and apply threshold
80  for (int loop = 0; loop < int(EndcapSuperClusters.size()); loop++) {
81  SuperCluster clus1 = EndcapSuperClusters[loop];
82  float eta1 = clus1.eta();
83  float energy1 = clus1.energy();
84  float theta1 = 2 * atan(exp(-1. * eta1));
85  float cluspt1 = energy1 * sin(theta1);
86  if (cluspt1 > cluster_pt_thresh_) {
87  saveEndcapSuperClusters.push_back(clus1);
88  if (cluspt1 > single_cluster_thresh_)
89  meet_single_thresh = true;
90  }
91  }
92 
93  auto SEBDigiCol = std::make_unique<EBDigiCollection>();
94  auto SEEDigiCol = std::make_unique<EEDigiCollection>();
95  int TotClus = saveBarrelSuperClusters.size() + saveEndcapSuperClusters.size();
96 
97  if (TotClus >= nclus_sel_ || meet_single_thresh) {
98  if (!saveBarrelSuperClusters.empty()) {
100  es.get<CaloTopologyRecord>().get(pTopology);
101  const CaloTopology* topology = pTopology.product();
102 
103  //get barrel digi collection
105  const EBDigiCollection* digis = nullptr;
106  evt.getByToken(EcalEBDigiToken_, pdigis);
107  digis = pdigis.product(); // get a ptr to the product
108 
110  const EcalRecHitCollection* rechits = nullptr;
111  evt.getByToken(EcalEBRecHitToken_, prechits);
112  rechits = prechits.product(); // get a ptr to the product
113 
114  if (digis) {
115  std::vector<DetId> saveTheseDetIds;
116  //pick out the detids for the 3x3 in each of the selected superclusters
117  for (int loop = 0; loop < int(saveBarrelSuperClusters.size()); loop++) {
118  SuperCluster clus1 = saveBarrelSuperClusters[loop];
119  const CaloClusterPtr& bcref = clus1.seed();
120  const BasicCluster* bc = bcref.get();
121  //Get the maximum detid
122  DetId maxDetId = EcalClusterTools::getMaximum(*bc, rechits).first;
123  // Loop over the 3x3 array centered on maximum detid
124  for (DetId detId : CaloRectangleRange(1, maxDetId, *topology))
125  saveTheseDetIds.push_back(detId);
126  }
127  for (int detloop = 0; detloop < int(saveTheseDetIds.size()); ++detloop) {
128  EBDetId detL = EBDetId(saveTheseDetIds[detloop]);
129 
130  for (EBDigiCollection::const_iterator blah = digis->begin(); blah != digis->end(); blah++) {
131  if (detL == blah->id()) {
132  EBDataFrame myDigi = (*blah);
133  SEBDigiCol->push_back(detL);
134 
135  EBDataFrame df(SEBDigiCol->back());
136  for (int iq = 0; iq < myDigi.size(); ++iq) {
137  df.setSample(iq, myDigi.sample(iq).raw());
138  }
139  //ebcounter++;
140  }
141  }
142  //if (ebcounter >= int(saveTheseDetIds.size())) break;
143  } //loop over dets
144  }
145 
146  } //If barrel superclusters need saving.
147 
148  if (!saveEndcapSuperClusters.empty()) {
149  edm::ESHandle<CaloTopology> pTopology;
150  es.get<CaloTopologyRecord>().get(pTopology);
151  const CaloTopology* topology = pTopology.product();
152 
153  //Get endcap rec hit collection
154  //get endcap digi collection
156  const EEDigiCollection* digis = nullptr;
157  evt.getByToken(EcalEEDigiToken_, pdigis);
158  digis = pdigis.product(); // get a ptr to the product
159 
161  const EcalRecHitCollection* rechits = nullptr;
162  evt.getByToken(EcalEERecHitToken_, prechits);
163  rechits = prechits.product(); // get a ptr to the product
164 
165  if (digis) {
166  //std::vector<DetId> saveTheseDetIds;
167  std::set<DetId> saveTheseDetIds;
168  //pick out the digis for the 3x3 in each of the selected superclusters
169  for (int loop = 0; loop < int(saveEndcapSuperClusters.size()); loop++) {
170  SuperCluster clus1 = saveEndcapSuperClusters[loop];
171  const CaloClusterPtr& bcref = clus1.seed();
172  const BasicCluster* bc = bcref.get();
173  //Get the maximum detid
174  DetId maxDetId = EcalClusterTools::getMaximum(*bc, rechits).first;
175  // Loop over the 3x3 array centered on maximum detid
176  for (DetId detId : CaloRectangleRange(1, maxDetId, *topology))
177  saveTheseDetIds.insert(detId);
178  }
179  int eecounter = 0;
180  for (EEDigiCollection::const_iterator blah = digis->begin(); blah != digis->end(); blah++) {
181  std::set<DetId>::const_iterator finder = saveTheseDetIds.find(blah->id());
182  if (finder != saveTheseDetIds.end()) {
183  EEDetId detL = EEDetId(*finder);
184 
185  if (detL == blah->id()) {
186  EEDataFrame myDigi = (*blah);
187  SEEDigiCol->push_back(detL);
188  EEDataFrame df(SEEDigiCol->back());
189  for (int iq = 0; iq < myDigi.size(); ++iq) {
190  df.setSample(iq, myDigi.sample(iq).raw());
191  }
192  eecounter++;
193  }
194  }
195  if (eecounter >= int(saveTheseDetIds.size()))
196  break;
197  } //loop over digis
198  }
199  } //If endcap superclusters need saving.
200 
201  } //If we're actually saving stuff
202 
203  //Okay, either my collections have been filled with the requisite Digis, or they haven't.
204 
205  //Empty collection, or full, still put in event.
206  SEBDigiCol->sort();
207  SEEDigiCol->sort();
208  evt.put(std::move(SEBDigiCol), selectedEcalEBDigiCollection_);
209  evt.put(std::move(SEEDigiCol), selectedEcalEEDigiCollection_);
210 }
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
CaloTopology const * topology(0)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:139
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
EcalMGPASample sample(int i) const
Definition: EcalDataFrame.h:29
const_iterator begin() const
uint16_t raw() const
get the raw word
int size() const
Definition: EcalDataFrame.h:26
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:180
EcalDigiSelector(const edm::ParameterSet &ps)
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
double energy() const
cluster energy
Definition: CaloCluster.h:148
Definition: DetId.h:17
T const * product() const
Definition: Handle.h:69
boost::transform_iterator< IterHelp, boost::counting_iterator< int > > const_iterator
const_iterator end() const
fixed size matrix
void setSample(int i, EcalMGPASample sam)
Definition: EcalDataFrame.h:43
T get() const
Definition: EventSetup.h:73
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:77
void produce(edm::Event &, const edm::EventSetup &) override
T const * product() const
Definition: ESHandle.h:86
def move(src, dest)
Definition: eostools.py:511