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  CaloSubdetectorTopology * topology_p=nullptr;
99  if (geometry_p)
100  topology_p = new EcalPreshowerTopology(geoHandle);
101 
102  // fetch the product (pSuperClusters)
103  evt.getByToken(endcapSClusterToken_, pSuperClusters);
104  const reco::SuperClusterCollection* SClusts = pSuperClusters.product();
105 
106  // fetch the product (RecHits)
107  evt.getByToken( preshHitToken_, pRecHits);
108  // pointer to the object in the product
109  const EcalRecHitCollection* rechits = pRecHits.product(); // EcalRecHitCollection hit_collection = *rhcHandle;
110 
111  LogTrace("EcalClusters") << "PreshowerPhiClusterProducerInfo: ### Total # of preshower RecHits: "<< rechits->size();
112 
113  // make the map of rechits:
114  std::map<DetId, EcalRecHit> rechits_map;
116  for (it = rechits->begin(); it != rechits->end(); it++) {
117  // remove bad ES rechits
118  if (it->recoFlag()==1 || it->recoFlag()==14 || (it->recoFlag()<=10 && it->recoFlag()>=5)) continue;
119  //Make the map of DetID, EcalRecHit pairs
120  rechits_map.insert(std::make_pair(it->id(), *it));
121  }
122  // The set of used DetID's for a given event:
123  std::set<DetId> used_strips;
124  used_strips.clear();
125  LogTrace("EcalClusters") << "PreshowerPhiClusterProducerInfo: ### rechits_map of size " << rechits_map.size() <<" was created!";
126 
127  reco::PreshowerClusterCollection clusters1, clusters2; // output collection of corrected PCs
128  reco::SuperClusterCollection new_SC; // output collection of corrected SCs
129 
130  //make cycle over super clusters
131  reco::SuperClusterCollection::const_iterator it_super;
132  int isc = 0;
133  for (it_super=SClusts->begin(); it_super!=SClusts->end(); ++it_super) {
134 
135  float e1 = 0;
136  float e2 = 0;
137  float deltaE = 0;
138 
140  ++isc;
141  LogTrace("EcalClusters")<< " superE = " << it_super->energy() << " superETA = " << it_super->eta() << " superPHI = " << it_super->phi() ;
142 
143  //cout<<"=== new SC ==="<<endl;
144  //cout<<"superE = "<<it_super->energy()<<" superETA = "<<it_super->eta()<<" superPHI = "<<it_super->phi()<<endl;
145 
146  int nBC = 0;
147  int condP1 = 1; // 0: dead channel; 1: active channel
148  int condP2 = 1;
149  float maxDeltaPhi = 0;
150  float minDeltaPhi = 0;
151  float refPhi = 0;
152 
153  reco::CaloCluster_iterator bc_iter = it_super->clustersBegin();
154  for ( ; bc_iter !=it_super->clustersEnd(); ++bc_iter) {
155  if (nBC == 0) {
156  refPhi = (*bc_iter)->phi();
157  } else {
158  if (reco::deltaPhi((*bc_iter)->phi(), refPhi) > 0 && reco::deltaPhi((*bc_iter)->phi(), refPhi) > maxDeltaPhi)
159  maxDeltaPhi = reco::deltaPhi((*bc_iter)->phi(), refPhi);
160  if (reco::deltaPhi((*bc_iter)->phi(), refPhi) < 0 && reco::deltaPhi((*bc_iter)->phi(), refPhi) < minDeltaPhi)
161  minDeltaPhi = reco::deltaPhi((*bc_iter)->phi(), refPhi);
162  //cout<<"delta phi : "<<reco::deltaPhi((*bc_iter)->phi(), refPhi)<<endl;
163  }
164  //cout<<"BC : "<<nBC<<" "<<(*bc_iter)->energy()<<" "<<(*bc_iter)->eta()<<" "<<(*bc_iter)->phi()<<endl;
165  nBC++;
166  }
167  maxDeltaPhi += esPhiClusterDeltaPhi_;
168  minDeltaPhi -= esPhiClusterDeltaPhi_;
169 
170  nBC = 0;
171  for (bc_iter = it_super->clustersBegin() ; bc_iter !=it_super->clustersEnd(); ++bc_iter) {
172  if (geometry_p) {
173 
174  // Get strip position at intersection point of the line EE - Vertex:
175  double X = (*bc_iter)->x();
176  double Y = (*bc_iter)->y();
177  double Z = (*bc_iter)->z();
178  const GlobalPoint point(X,Y,Z);
179 
180  DetId tmp1 = (dynamic_cast<const EcalPreshowerGeometry*>(geometry_p))->getClosestCellInPlane(point, 1);
181  DetId tmp2 = (dynamic_cast<const EcalPreshowerGeometry*>(geometry_p))->getClosestCellInPlane(point, 2);
182  ESDetId strip1 = (tmp1 == DetId(0)) ? ESDetId(0) : ESDetId(tmp1);
183  ESDetId strip2 = (tmp2 == DetId(0)) ? ESDetId(0) : ESDetId(tmp2);
184 
185  if (nBC == 0) {
186  if (strip1 != ESDetId(0) && strip2 != ESDetId(0)) {
187  ESChannelStatusMap::const_iterator status_p1 = channelStatus->getMap().find(strip1);
188  ESChannelStatusMap::const_iterator status_p2 = channelStatus->getMap().find(strip2);
189  if (status_p1->getStatusCode() == 1) condP1 = 0;
190  if (status_p2->getStatusCode() == 1) condP2 = 0;
191  } else if (strip1 == ESDetId(0)) {
192  condP1 = 0;
193  } else if (strip2 == ESDetId(0)) {
194  condP2 = 0;
195  }
196 
197  //cout<<"starting cluster : "<<maxDeltaPhi<<" "<<minDeltaPhi<<" "<<esPhiClusterDeltaEta_<<endl;
198  //cout<<"do plane 1 === "<<strip1.zside()<<" "<<strip1.plane()<<" "<<strip1.six()<<" "<<strip1.siy()<<" "<<strip1.strip()<<endl;
199  // Get a vector of ES clusters (found by the PreshSeeded algorithm) associated with a given EE basic cluster.
200  reco::PreshowerCluster cl1 = presh_algo->makeOneCluster(strip1,&used_strips,&rechits_map,geometry_p,topology_p,esPhiClusterDeltaEta_,minDeltaPhi,maxDeltaPhi);
201  cl1.setBCRef(*bc_iter);
202  clusters1.push_back(cl1);
203  e1 += cl1.energy();
204 
205  //cout<<"do plane 2 === "<<strip2.zside()<<" "<<strip2.plane()<<" "<<strip2.six()<<" "<<strip2.siy()<<" "<<strip2.strip()<<endl;
206  reco::PreshowerCluster cl2 = presh_algo->makeOneCluster(strip2,&used_strips,&rechits_map,geometry_p,topology_p,esPhiClusterDeltaEta_,minDeltaPhi,maxDeltaPhi);
207  cl2.setBCRef(*bc_iter);
208  clusters2.push_back(cl2);
209  e2 += cl2.energy();
210  }
211  }
212 
213  new_BC.push_back(*bc_iter);
214  nBC++;
215  } // end of cycle over BCs
216 
217  LogTrace("EcalClusters") << " For SC #" << isc-1 << ", containing "
218  << it_super->clustersSize()
219  << " basic clusters, PreshowerPhiClusterAlgo made "
220  << clusters1.size() << " in X plane and "
221  << clusters2.size()
222  << " in Y plane " << " preshower clusters " ;
223 
224  // update energy of the SuperCluster
225  if(e1+e2 > 1.0e-10) {
226 
227  e1 = e1 / mip_; // GeV to #MIPs
228  e2 = e2 / mip_;
229 
230  if (condP1 == 1 && condP2 == 1) {
231  deltaE = gamma0_*(e1 + alpha0_*e2);
232  } else if (condP1 == 1 && condP2 == 0) {
233  deltaE = gamma1_*(e1 + alpha1_*e2);
234  } else if (condP1 == 0 && condP2 == 1) {
235  deltaE = gamma2_*(e1 + alpha2_*e2);
236  } else if (condP1 == 0 && condP2 == 0) {
237  deltaE = gamma3_*(e1 + alpha3_*e2);
238  }
239  }
240 
241  //corrected Energy
242  float E = it_super->energy() + deltaE;
243 
244  LogTrace("EcalClusters") << " Creating corrected SC ";
245 
246  reco::SuperCluster sc(E, it_super->position(), it_super->seed(), new_BC, deltaE);
247  sc.setPreshowerEnergyPlane1(e1*mip_);
248  sc.setPreshowerEnergyPlane2(e2*mip_);
249  if (condP1 == 1 && condP2 == 1) sc.setPreshowerPlanesStatus(0);
250  else if (condP1 == 1 && condP2 == 0) sc.setPreshowerPlanesStatus(1);
251  else if (condP1 == 0 && condP2 == 1) sc.setPreshowerPlanesStatus(2);
252  else if (condP1 == 0 && condP2 == 0) sc.setPreshowerPlanesStatus(3);
253 
254  if (etThresh_>0){ // calling postion().theta can be expensive
255  if (sc.energy()*sin(sc.position().theta())>etThresh_ )
256  new_SC.push_back(sc);
257  } else {
258  new_SC.push_back(sc);
259  }
260 
261  } // end of cycle over SCs
262 
263  // copy the preshower clusters into collections and put in the Event:
264  clusters_p1->assign(clusters1.begin(), clusters1.end());
265  clusters_p2->assign(clusters2.begin(), clusters2.end());
266  // put collection of preshower clusters to the event
267  evt.put(std::move(clusters_p1), preshClusterCollectionX_ );
268  evt.put(std::move(clusters_p2), preshClusterCollectionY_ );
269  LogTrace("EcalClusters") << "Preshower clusters added to the event" ;
270 
271  // put collection of corrected super clusters to the event
272  superclusters_p->assign(new_SC.begin(), new_SC.end());
273  evt.put(std::move(superclusters_p), assocSClusterCollection_);
274  LogTrace("EcalClusters") << "Corrected SClusters added to the event" ;
275 
276  if (topology_p) delete topology_p;
277 }
278 
280 
281  es.get<ESGainRcd>().get(esgain_);
282  const ESGain *gain = esgain_.product();
283 
284  double ESGain = gain->getESGain();
285 
286  es.get<ESMIPToGeVConstantRcd>().get(esMIPToGeV_);
287  const ESMIPToGeVConstant *mipToGeV = esMIPToGeV_.product();
288 
289  mip_ = (ESGain == 1) ? mipToGeV->getESValueLow() : mipToGeV->getESValueHigh();
290 
291  es.get<ESChannelStatusRcd>().get(esChannelStatus_);
292 
293  es.get<ESEEIntercalibConstantsRcd>().get(esEEInterCalib_);
294  const ESEEIntercalibConstants *esEEInterCalib = esEEInterCalib_.product();
295 
296  // both planes work
297  gamma0_ = (ESGain == 1) ? 0.02 : esEEInterCalib->getGammaHigh0();
298  alpha0_ = (ESGain == 1) ? esEEInterCalib->getAlphaLow0() : esEEInterCalib->getAlphaHigh0();
299 
300  // only first plane works
301  gamma1_ = (ESGain == 1) ? (0.02 * esEEInterCalib->getGammaLow1()) : esEEInterCalib->getGammaHigh1();
302  alpha1_ = (ESGain == 1) ? esEEInterCalib->getAlphaLow1() : esEEInterCalib->getAlphaHigh1();
303 
304  // only second plane works
305  gamma2_ = (ESGain == 1) ? (0.02 * esEEInterCalib->getGammaLow2()) : esEEInterCalib->getGammaHigh2();
306  alpha2_ = (ESGain == 1) ? esEEInterCalib->getAlphaLow2() : esEEInterCalib->getAlphaHigh2();
307 
308  // both planes do not work
309  gamma3_ = (ESGain == 1) ? 0.02 : esEEInterCalib->getGammaHigh3();
310  alpha3_ = (ESGain == 1) ? esEEInterCalib->getAlphaLow3() : esEEInterCalib->getAlphaHigh3();
311 
312  es.get<ESMissingEnergyCalibrationRcd>().get(esMissingECalib_);
313  const ESMissingEnergyCalibration * esMissingECalib = esMissingECalib_.product();
314 
315  // |eta| < 1.9
316  aEta_[0] = esMissingECalib->getConstAEta0();
317  bEta_[0] = esMissingECalib->getConstBEta0();
318 
319  // 1.9 < |eta| < 2.1
320  aEta_[1] = esMissingECalib->getConstAEta1();
321  bEta_[1] = esMissingECalib->getConstBEta1();
322 
323  // 2.1 < |eta| < 2.3
324  aEta_[2] = esMissingECalib->getConstAEta2();
325  bEta_[2] = esMissingECalib->getConstBEta2();
326 
327  // 2.3 < |eta| < 2.5
328  aEta_[3] = esMissingECalib->getConstAEta3();
329  bEta_[3] = esMissingECalib->getConstBEta3();
330 
331 }
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:137
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
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:111
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:166
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:81
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:68
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