CMS 3D CMS Logo

PreshowerPhiClusterProducer.cc
Go to the documentation of this file.
1 #include <vector>
2 #include <memory>
3 
37 #include <fstream>
38 
40 
41 using namespace edm;
42 using namespace std;
43 
45 
46  // use configuration file to setup input/output collection names
47  preshHitToken_ =
48  consumes<EcalRecHitCollection>(ps.getParameter<edm::InputTag>("preshRecHitProducer"));
49 
50  // Name of a SuperClusterCollection to make associations:
51  endcapSClusterToken_ =
52  consumes<reco::SuperClusterCollection>(ps.getParameter<edm::InputTag>("endcapSClusterProducer"));
53 
54  // Output collections:
55  preshClusterCollectionX_ = ps.getParameter<std::string>("preshClusterCollectionX");
56  preshClusterCollectionY_ = ps.getParameter<std::string>("preshClusterCollectionY");
57 
58  assocSClusterCollection_ = ps.getParameter<std::string>("assocSClusterCollection");
59 
60  produces< reco::PreshowerClusterCollection >(preshClusterCollectionX_);
61  produces< reco::PreshowerClusterCollection >(preshClusterCollectionY_);
62  produces< reco::SuperClusterCollection >(assocSClusterCollection_);
63 
64  float esStripECut = ps.getParameter<double>("esStripEnergyCut");
65  esPhiClusterDeltaEta_ = ps.getParameter<double>("esPhiClusterDeltaEta");
66  esPhiClusterDeltaPhi_ = ps.getParameter<double>("esPhiClusterDeltaPhi");
67 
68  etThresh_ = ps.getParameter<double>("etThresh");
69 
70  presh_algo = new PreshowerPhiClusterAlgo(esStripECut);
71 }
72 
74  delete presh_algo;
75 }
76 
78 
81 
82  // get the ECAL geometry:
84  es.get<CaloGeometryRecord>().get(geoHandle);
85 
86  // retrieve ES-EE intercalibration constants and channel status
87  set(es);
88  const ESChannelStatus *channelStatus = esChannelStatus_.product();
89 
91 
92  // create unique_ptr to a PreshowerClusterCollection
93  auto clusters_p1 = std::make_unique<reco::PreshowerClusterCollection>();
94  auto clusters_p2 = std::make_unique<reco::PreshowerClusterCollection>();
95  // create new collection of corrected super clusters
96  auto superclusters_p = std::make_unique<reco::SuperClusterCollection>();
97 
98  // fetch the product (pSuperClusters)
99  evt.getByToken(endcapSClusterToken_, pSuperClusters);
100  const reco::SuperClusterCollection* SClusts = pSuperClusters.product();
101 
102  // fetch the product (RecHits)
103  evt.getByToken( preshHitToken_, pRecHits);
104  // pointer to the object in the product
105  const EcalRecHitCollection* rechits = pRecHits.product(); // EcalRecHitCollection hit_collection = *rhcHandle;
106 
107  LogTrace("EcalClusters") << "PreshowerPhiClusterProducerInfo: ### Total # of preshower RecHits: "<< rechits->size();
108 
109  // make the map of rechits:
110  std::map<DetId, EcalRecHit> rechits_map;
112  for (it = rechits->begin(); it != rechits->end(); it++) {
113  // remove bad ES rechits
114  if (it->recoFlag()==1 || it->recoFlag()==14 || (it->recoFlag()<=10 && it->recoFlag()>=5)) continue;
115  //Make the map of DetID, EcalRecHit pairs
116  rechits_map.insert(std::make_pair(it->id(), *it));
117  }
118  // The set of used DetID's for a given event:
119  std::set<DetId> used_strips;
120  used_strips.clear();
121  LogTrace("EcalClusters") << "PreshowerPhiClusterProducerInfo: ### rechits_map of size " << rechits_map.size() <<" was created!";
122 
123  reco::PreshowerClusterCollection clusters1, clusters2; // output collection of corrected PCs
124  reco::SuperClusterCollection new_SC; // output collection of corrected SCs
125 
126  //make cycle over super clusters
127  reco::SuperClusterCollection::const_iterator it_super;
128  int isc = 0;
129  for (it_super=SClusts->begin(); it_super!=SClusts->end(); ++it_super) {
130 
131  float e1 = 0;
132  float e2 = 0;
133  float deltaE = 0;
134 
136  ++isc;
137  LogTrace("EcalClusters")<< " superE = " << it_super->energy() << " superETA = " << it_super->eta() << " superPHI = " << it_super->phi() ;
138 
139  //cout<<"=== new SC ==="<<endl;
140  //cout<<"superE = "<<it_super->energy()<<" superETA = "<<it_super->eta()<<" superPHI = "<<it_super->phi()<<endl;
141 
142  int nBC = 0;
143  int condP1 = 1; // 0: dead channel; 1: active channel
144  int condP2 = 1;
145  float maxDeltaPhi = 0;
146  float minDeltaPhi = 0;
147  float refPhi = 0;
148 
149  reco::CaloCluster_iterator bc_iter = it_super->clustersBegin();
150  for ( ; bc_iter !=it_super->clustersEnd(); ++bc_iter) {
151  if (nBC == 0) {
152  refPhi = (*bc_iter)->phi();
153  } else {
154  if (reco::deltaPhi((*bc_iter)->phi(), refPhi) > 0 && reco::deltaPhi((*bc_iter)->phi(), refPhi) > maxDeltaPhi)
155  maxDeltaPhi = reco::deltaPhi((*bc_iter)->phi(), refPhi);
156  if (reco::deltaPhi((*bc_iter)->phi(), refPhi) < 0 && reco::deltaPhi((*bc_iter)->phi(), refPhi) < minDeltaPhi)
157  minDeltaPhi = reco::deltaPhi((*bc_iter)->phi(), refPhi);
158  //cout<<"delta phi : "<<reco::deltaPhi((*bc_iter)->phi(), refPhi)<<endl;
159  }
160  //cout<<"BC : "<<nBC<<" "<<(*bc_iter)->energy()<<" "<<(*bc_iter)->eta()<<" "<<(*bc_iter)->phi()<<endl;
161  nBC++;
162  }
163  maxDeltaPhi += esPhiClusterDeltaPhi_;
164  minDeltaPhi -= esPhiClusterDeltaPhi_;
165 
166  nBC = 0;
167  for (bc_iter = it_super->clustersBegin() ; bc_iter !=it_super->clustersEnd(); ++bc_iter) {
168  if (geometry_p) {
169 
170  // Get strip position at intersection point of the line EE - Vertex:
171  double X = (*bc_iter)->x();
172  double Y = (*bc_iter)->y();
173  double Z = (*bc_iter)->z();
174  const GlobalPoint point(X,Y,Z);
175 
176  DetId tmp1 = (dynamic_cast<const EcalPreshowerGeometry*>(geometry_p))->getClosestCellInPlane(point, 1);
177  DetId tmp2 = (dynamic_cast<const EcalPreshowerGeometry*>(geometry_p))->getClosestCellInPlane(point, 2);
178  ESDetId strip1 = (tmp1 == DetId(0)) ? ESDetId(0) : ESDetId(tmp1);
179  ESDetId strip2 = (tmp2 == DetId(0)) ? ESDetId(0) : ESDetId(tmp2);
180 
181  if (nBC == 0) {
182  if (strip1 != ESDetId(0) && strip2 != ESDetId(0)) {
183  ESChannelStatusMap::const_iterator status_p1 = channelStatus->getMap().find(strip1);
184  ESChannelStatusMap::const_iterator status_p2 = channelStatus->getMap().find(strip2);
185  if (status_p1->getStatusCode() == 1) condP1 = 0;
186  if (status_p2->getStatusCode() == 1) condP2 = 0;
187  } else if (strip1 == ESDetId(0)) {
188  condP1 = 0;
189  } else if (strip2 == ESDetId(0)) {
190  condP2 = 0;
191  }
192 
193  //cout<<"starting cluster : "<<maxDeltaPhi<<" "<<minDeltaPhi<<" "<<esPhiClusterDeltaEta_<<endl;
194  //cout<<"do plane 1 === "<<strip1.zside()<<" "<<strip1.plane()<<" "<<strip1.six()<<" "<<strip1.siy()<<" "<<strip1.strip()<<endl;
195  // Get a vector of ES clusters (found by the PreshSeeded algorithm) associated with a given EE basic cluster.
196  reco::PreshowerCluster cl1 = presh_algo->makeOneCluster(strip1,&used_strips,&rechits_map,geometry_p,esPhiClusterDeltaEta_,minDeltaPhi,maxDeltaPhi);
197  cl1.setBCRef(*bc_iter);
198  clusters1.push_back(cl1);
199  e1 += cl1.energy();
200 
201  //cout<<"do plane 2 === "<<strip2.zside()<<" "<<strip2.plane()<<" "<<strip2.six()<<" "<<strip2.siy()<<" "<<strip2.strip()<<endl;
202  reco::PreshowerCluster cl2 = presh_algo->makeOneCluster(strip2,&used_strips,&rechits_map,geometry_p,esPhiClusterDeltaEta_,minDeltaPhi,maxDeltaPhi);
203  cl2.setBCRef(*bc_iter);
204  clusters2.push_back(cl2);
205  e2 += cl2.energy();
206  }
207  }
208 
209  new_BC.push_back(*bc_iter);
210  nBC++;
211  } // end of cycle over BCs
212 
213  LogTrace("EcalClusters") << " For SC #" << isc-1 << ", containing "
214  << it_super->clustersSize()
215  << " basic clusters, PreshowerPhiClusterAlgo made "
216  << clusters1.size() << " in X plane and "
217  << clusters2.size()
218  << " in Y plane " << " preshower clusters " ;
219 
220  // update energy of the SuperCluster
221  if(e1+e2 > 1.0e-10) {
222 
223  e1 = e1 / mip_; // GeV to #MIPs
224  e2 = e2 / mip_;
225 
226  if (condP1 == 1 && condP2 == 1) {
227  deltaE = gamma0_*(e1 + alpha0_*e2);
228  } else if (condP1 == 1 && condP2 == 0) {
229  deltaE = gamma1_*(e1 + alpha1_*e2);
230  } else if (condP1 == 0 && condP2 == 1) {
231  deltaE = gamma2_*(e1 + alpha2_*e2);
232  } else if (condP1 == 0 && condP2 == 0) {
233  deltaE = gamma3_*(e1 + alpha3_*e2);
234  }
235  }
236 
237  //corrected Energy
238  float E = it_super->energy() + deltaE;
239 
240  LogTrace("EcalClusters") << " Creating corrected SC ";
241 
242  reco::SuperCluster sc(E, it_super->position(), it_super->seed(), new_BC, deltaE);
243  sc.setPreshowerEnergyPlane1(e1*mip_);
244  sc.setPreshowerEnergyPlane2(e2*mip_);
245  if (condP1 == 1 && condP2 == 1) sc.setPreshowerPlanesStatus(0);
246  else if (condP1 == 1 && condP2 == 0) sc.setPreshowerPlanesStatus(1);
247  else if (condP1 == 0 && condP2 == 1) sc.setPreshowerPlanesStatus(2);
248  else if (condP1 == 0 && condP2 == 0) sc.setPreshowerPlanesStatus(3);
249 
250  if (etThresh_>0){ // calling postion().theta can be expensive
251  if (sc.energy()*sin(sc.position().theta())>etThresh_ )
252  new_SC.push_back(sc);
253  } else {
254  new_SC.push_back(sc);
255  }
256 
257  } // end of cycle over SCs
258 
259  // copy the preshower clusters into collections and put in the Event:
260  clusters_p1->assign(clusters1.begin(), clusters1.end());
261  clusters_p2->assign(clusters2.begin(), clusters2.end());
262  // put collection of preshower clusters to the event
263  evt.put(std::move(clusters_p1), preshClusterCollectionX_ );
264  evt.put(std::move(clusters_p2), preshClusterCollectionY_ );
265  LogTrace("EcalClusters") << "Preshower clusters added to the event" ;
266 
267  // put collection of corrected super clusters to the event
268  superclusters_p->assign(new_SC.begin(), new_SC.end());
269  evt.put(std::move(superclusters_p), assocSClusterCollection_);
270  LogTrace("EcalClusters") << "Corrected SClusters added to the event" ;
271 }
272 
274 
275  es.get<ESGainRcd>().get(esgain_);
276  const ESGain *gain = esgain_.product();
277 
278  double ESGain = gain->getESGain();
279 
280  es.get<ESMIPToGeVConstantRcd>().get(esMIPToGeV_);
281  const ESMIPToGeVConstant *mipToGeV = esMIPToGeV_.product();
282 
283  mip_ = (ESGain == 1) ? mipToGeV->getESValueLow() : mipToGeV->getESValueHigh();
284 
285  es.get<ESChannelStatusRcd>().get(esChannelStatus_);
286 
287  es.get<ESEEIntercalibConstantsRcd>().get(esEEInterCalib_);
288  const ESEEIntercalibConstants *esEEInterCalib = esEEInterCalib_.product();
289 
290  // both planes work
291  gamma0_ = (ESGain == 1) ? 0.02 : esEEInterCalib->getGammaHigh0();
292  alpha0_ = (ESGain == 1) ? esEEInterCalib->getAlphaLow0() : esEEInterCalib->getAlphaHigh0();
293 
294  // only first plane works
295  gamma1_ = (ESGain == 1) ? (0.02 * esEEInterCalib->getGammaLow1()) : esEEInterCalib->getGammaHigh1();
296  alpha1_ = (ESGain == 1) ? esEEInterCalib->getAlphaLow1() : esEEInterCalib->getAlphaHigh1();
297 
298  // only second plane works
299  gamma2_ = (ESGain == 1) ? (0.02 * esEEInterCalib->getGammaLow2()) : esEEInterCalib->getGammaHigh2();
300  alpha2_ = (ESGain == 1) ? esEEInterCalib->getAlphaLow2() : esEEInterCalib->getAlphaHigh2();
301 
302  // both planes do not work
303  gamma3_ = (ESGain == 1) ? 0.02 : esEEInterCalib->getGammaHigh3();
304  alpha3_ = (ESGain == 1) ? esEEInterCalib->getAlphaLow3() : esEEInterCalib->getAlphaHigh3();
305 
306  es.get<ESMissingEnergyCalibrationRcd>().get(esMissingECalib_);
307  const ESMissingEnergyCalibration * esMissingECalib = esMissingECalib_.product();
308 
309  // |eta| < 1.9
310  aEta_[0] = esMissingECalib->getConstAEta0();
311  bEta_[0] = esMissingECalib->getConstBEta0();
312 
313  // 1.9 < |eta| < 2.1
314  aEta_[1] = esMissingECalib->getConstAEta1();
315  bEta_[1] = esMissingECalib->getConstBEta1();
316 
317  // 2.1 < |eta| < 2.3
318  aEta_[2] = esMissingECalib->getConstAEta2();
319  bEta_[2] = esMissingECalib->getConstBEta2();
320 
321  // 2.3 < |eta| < 2.5
322  aEta_[3] = esMissingECalib->getConstAEta3();
323  bEta_[3] = esMissingECalib->getConstBEta3();
324 
325 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
T getParameter(std::string const &) const
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:49
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
Definition: ESGain.h:7
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:140
const self & getMap() const
float getESValueLow() const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< EcalRecHit >::const_iterator const_iterator
#define X(str)
Definition: MuonsGrabber.cc:48
void set(const edm::EventSetup &es)
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
void setBCRef(const CaloClusterPtr &r)
DetIds of component RecHits – now inherited from CaloCluster.
const_iterator find(uint32_t rawId) const
float energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:104
std::vector< PreshowerCluster > PreshowerClusterCollection
collection of PreshowerCluster objects
double energy() const
cluster energy
Definition: CaloCluster.h:126
#define LogTrace(id)
float theta() const
Momentum polar angle. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:144
const_iterator end() const
void produce(edm::Event &evt, const edm::EventSetup &es) override
Definition: DetId.h:18
T const * product() const
Definition: Handle.h:74
float getESGain() const
Definition: ESGain.h:13
std::vector< Item >::const_iterator const_iterator
HLT enums.
size_type size() const
T get() const
Definition: EventSetup.h:71
float getESValueHigh() const
def move(src, dest)
Definition: eostools.py:511
PreshowerPhiClusterProducer(const edm::ParameterSet &ps)
*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
const_iterator begin() const