CMS 3D CMS Logo

PiZeroAnalyzer.cc
Go to the documentation of this file.
1 #include <iostream>
2 
4 
5 //#define TWOPI 6.283185308
6 //
7 
18 using namespace std;
19 
21  fName_ = pset.getUntrackedParameter<std::string>("Name");
22  prescaleFactor_ = pset.getUntrackedParameter<int>("prescaleFactor", 1);
23 
24  barrelEcalHits_token_ = consumes<edm::SortedCollection<EcalRecHit, edm::StrictWeakOrdering<EcalRecHit> > >(
25  pset.getParameter<edm::InputTag>("barrelEcalHits"));
26  endcapEcalHits_token_ = consumes<edm::SortedCollection<EcalRecHit, edm::StrictWeakOrdering<EcalRecHit> > >(
27  pset.getParameter<edm::InputTag>("endcapEcalHits"));
28 
29  nEvt_ = 0;
30 
31  // parameters for Pizero finding
32  seleXtalMinEnergy_ = pset.getParameter<double>("seleXtalMinEnergy");
33  clusSeedThr_ = pset.getParameter<double>("clusSeedThr");
34  clusEtaSize_ = pset.getParameter<int>("clusEtaSize");
35  clusPhiSize_ = pset.getParameter<int>("clusPhiSize");
36  ParameterLogWeighted_ = pset.getParameter<bool>("ParameterLogWeighted");
37  ParameterX0_ = pset.getParameter<double>("ParameterX0");
38  ParameterT0_barl_ = pset.getParameter<double>("ParameterT0_barl");
39  ParameterW0_ = pset.getParameter<double>("ParameterW0");
40 
41  selePtGammaOne_ = pset.getParameter<double>("selePtGammaOne");
42  selePtGammaTwo_ = pset.getParameter<double>("selePtGammaTwo");
43  seleS4S9GammaOne_ = pset.getParameter<double>("seleS4S9GammaOne");
44  seleS4S9GammaTwo_ = pset.getParameter<double>("seleS4S9GammaTwo");
45  selePtPi0_ = pset.getParameter<double>("selePtPi0");
46  selePi0Iso_ = pset.getParameter<double>("selePi0Iso");
47  selePi0BeltDR_ = pset.getParameter<double>("selePi0BeltDR");
48  selePi0BeltDeta_ = pset.getParameter<double>("selePi0BeltDeta");
49  seleMinvMaxPi0_ = pset.getParameter<double>("seleMinvMaxPi0");
50  seleMinvMinPi0_ = pset.getParameter<double>("seleMinvMinPi0");
51 
52  posCalcParameters_ = pset.getParameter<edm::ParameterSet>("posCalcParameters");
53 }
54 
56 
58  edm::Run const& /* iRun */,
59  edm::EventSetup const& /* iSetup */) {
60  currentFolder_.str("");
61  currentFolder_ << "Egamma/PiZeroAnalyzer/";
62  ibooker.setCurrentFolder(currentFolder_.str());
63 
64  hMinvPi0EB_ = ibooker.book1D("Pi0InvmassEB", "Pi0 Invariant Mass in EB", 100, 0., 0.5);
65  hMinvPi0EB_->setAxisTitle("Inv Mass [GeV] ", 1);
66 
67  hPt1Pi0EB_ = ibooker.book1D("Pt1Pi0EB", "Pt 1st most energetic Pi0 photon in EB", 100, 0., 20.);
68  hPt1Pi0EB_->setAxisTitle("1st photon Pt [GeV] ", 1);
69 
70  hPt2Pi0EB_ = ibooker.book1D("Pt2Pi0EB", "Pt 2nd most energetic Pi0 photon in EB", 100, 0., 20.);
71  hPt2Pi0EB_->setAxisTitle("2nd photon Pt [GeV] ", 1);
72 
73  hPtPi0EB_ = ibooker.book1D("PtPi0EB", "Pi0 Pt in EB", 100, 0., 20.);
74  hPtPi0EB_->setAxisTitle("Pi0 Pt [GeV] ", 1);
75 
76  hIsoPi0EB_ = ibooker.book1D("IsoPi0EB", "Pi0 Iso in EB", 50, 0., 1.);
77  hIsoPi0EB_->setAxisTitle("Pi0 Iso", 1);
78 }
79 
81  using namespace edm;
82 
83  if (nEvt_ % prescaleFactor_)
84  return;
85  nEvt_++;
86  LogInfo("PiZeroAnalyzer") << "PiZeroAnalyzer Analyzing event number: " << e.id() << " Global Counter " << nEvt_
87  << "\n";
88 
89  // Get EcalRecHits
90  bool validEcalRecHits = true;
91  Handle<EcalRecHitCollection> barrelHitHandle;
93  e.getByToken(barrelEcalHits_token_, barrelHitHandle);
94  if (!barrelHitHandle.isValid()) {
95  edm::LogError("PhotonProducer") << "Error! Can't get the product: barrelEcalHits_token_";
96  validEcalRecHits = false;
97  }
98 
99  Handle<EcalRecHitCollection> endcapHitHandle;
100  e.getByToken(endcapEcalHits_token_, endcapHitHandle);
102  if (!endcapHitHandle.isValid()) {
103  edm::LogError("PhotonProducer") << "Error! Can't get the product: endcapEcalHits_token_";
104  validEcalRecHits = false;
105  }
106 
107  if (validEcalRecHits)
108  makePizero(esup, barrelHitHandle, endcapHitHandle);
109 }
110 
114  const EcalRecHitCollection* hitCollection_p = rhEB.product();
115 
116  edm::ESHandle<CaloGeometry> geoHandle;
117  es.get<CaloGeometryRecord>().get(geoHandle);
118 
119  edm::ESHandle<CaloTopology> theCaloTopology;
120  es.get<CaloTopologyRecord>().get(theCaloTopology);
121 
122  const CaloSubdetectorTopology* topology_p;
123  const CaloSubdetectorGeometry* geometry_p = geoHandle->getSubdetectorGeometry(DetId::Ecal, EcalBarrel);
124  const CaloSubdetectorGeometry* geometryES_p = geoHandle->getSubdetectorGeometry(DetId::Ecal, EcalPreshower);
125 
126  // Parameters for the position calculation:
127  PositionCalc posCalculator_ = PositionCalc(posCalcParameters_);
128  //
129  std::map<DetId, EcalRecHit> recHitsEB_map;
130  //
131  std::vector<EcalRecHit> seeds;
132 
133  seeds.clear();
134  //
135  vector<EBDetId> usedXtals;
136  usedXtals.clear();
137  //
139  //
140  static const int MAXCLUS = 2000;
141  int nClus = 0;
142  vector<float> eClus;
143  vector<float> etClus;
144  vector<float> etaClus;
145  vector<float> phiClus;
146  vector<EBDetId> max_hit;
147  vector<vector<EcalRecHit> > RecHitsCluster;
148  vector<float> s4s9Clus;
149 
150  // find cluster seeds in EB
151  for (itb = rhEB->begin(); itb != rhEB->end(); ++itb) {
152  EBDetId id(itb->id());
153  double energy = itb->energy();
154  if (energy > seleXtalMinEnergy_) {
155  std::pair<DetId, EcalRecHit> map_entry(itb->id(), *itb);
156  recHitsEB_map.insert(map_entry);
157  }
158  if (energy > clusSeedThr_)
159  seeds.push_back(*itb);
160  } // Eb rechits
161 
162  sort(seeds.begin(), seeds.end(), [](auto& x, auto& y) { return (x.energy() > y.energy()); });
163  for (std::vector<EcalRecHit>::iterator itseed = seeds.begin(); itseed != seeds.end(); itseed++) {
164  EBDetId seed_id = itseed->id();
165  std::vector<EBDetId>::const_iterator usedIds;
166 
167  bool seedAlreadyUsed = false;
168  for (usedIds = usedXtals.begin(); usedIds != usedXtals.end(); usedIds++) {
169  if (*usedIds == seed_id) {
170  seedAlreadyUsed = true;
171  //cout<< " Seed with energy "<<itseed->energy()<<" was used !"<<endl;
172  break;
173  }
174  }
175  if (seedAlreadyUsed)
176  continue;
177  topology_p = theCaloTopology->getSubdetectorTopology(DetId::Ecal, EcalBarrel);
178  std::vector<DetId> clus_v = topology_p->getWindow(seed_id, clusEtaSize_, clusPhiSize_);
179  //std::vector<DetId> clus_used;
180  std::vector<std::pair<DetId, float> > clus_used;
181 
182  vector<EcalRecHit> RecHitsInWindow;
183 
184  double simple_energy = 0;
185 
186  for (std::vector<DetId>::iterator det = clus_v.begin(); det != clus_v.end(); det++) {
187  // EBDetId EBdet = *det;
188  // cout<<" det "<< EBdet<<" ieta "<<EBdet.ieta()<<" iphi "<<EBdet.iphi()<<endl;
189  bool HitAlreadyUsed = false;
190  for (usedIds = usedXtals.begin(); usedIds != usedXtals.end(); usedIds++) {
191  if (*usedIds == *det) {
192  HitAlreadyUsed = true;
193  break;
194  }
195  }
196  if (HitAlreadyUsed)
197  continue;
198  if (recHitsEB_map.find(*det) != recHitsEB_map.end()) {
199  // cout<<" Used det "<< EBdet<<endl;
200  std::map<DetId, EcalRecHit>::iterator aHit;
201  aHit = recHitsEB_map.find(*det);
202  usedXtals.push_back(*det);
203  RecHitsInWindow.push_back(aHit->second);
204  clus_used.push_back(std::pair<DetId, float>(*det, 1.));
205  simple_energy = simple_energy + aHit->second.energy();
206  }
207  }
208 
209  math::XYZPoint clus_pos = posCalculator_.Calculate_Location(clus_used, hitCollection_p, geometry_p, geometryES_p);
210  float theta_s = 2. * atan(exp(-clus_pos.eta()));
211  float p0x_s = simple_energy * sin(theta_s) * cos(clus_pos.phi());
212  float p0y_s = simple_energy * sin(theta_s) * sin(clus_pos.phi());
213  // float p0z_s = simple_energy * cos(theta_s);
214  float et_s = sqrt(p0x_s * p0x_s + p0y_s * p0y_s);
215 
216  //cout << " Simple Clustering: E,Et,px,py,pz: "<<simple_energy<<" "<<et_s<<" "<<p0x_s<<" "<<p0y_s<<" "<<endl;
217 
218  eClus.push_back(simple_energy);
219  etClus.push_back(et_s);
220  etaClus.push_back(clus_pos.eta());
221  phiClus.push_back(clus_pos.phi());
222  max_hit.push_back(seed_id);
223  RecHitsCluster.push_back(RecHitsInWindow);
224  //Compute S4/S9 variable
225  //We are not sure to have 9 RecHits so need to check eta and phi:
226  float s4s9_[4];
227  for (int i = 0; i < 4; i++)
228  s4s9_[i] = itseed->energy();
229  for (unsigned int j = 0; j < RecHitsInWindow.size(); j++) {
230  //cout << " Simple cluster rh, ieta, iphi : "<<((EBDetId)RecHitsInWindow[j].id()).ieta()<<" "<<((EBDetId)RecHitsInWindow[j].id()).iphi()<<endl;
231  if ((((EBDetId)RecHitsInWindow[j].id()).ieta() == seed_id.ieta() - 1 && seed_id.ieta() != 1) ||
232  (seed_id.ieta() == 1 && (((EBDetId)RecHitsInWindow[j].id()).ieta() == seed_id.ieta() - 2))) {
233  if (((EBDetId)RecHitsInWindow[j].id()).iphi() == seed_id.iphi() - 1 ||
234  ((EBDetId)RecHitsInWindow[j].id()).iphi() - 360 == seed_id.iphi() - 1) {
235  s4s9_[0] += RecHitsInWindow[j].energy();
236  } else {
237  if (((EBDetId)RecHitsInWindow[j].id()).iphi() == seed_id.iphi()) {
238  s4s9_[0] += RecHitsInWindow[j].energy();
239  s4s9_[1] += RecHitsInWindow[j].energy();
240  } else {
241  if (((EBDetId)RecHitsInWindow[j].id()).iphi() == seed_id.iphi() + 1 ||
242  ((EBDetId)RecHitsInWindow[j].id()).iphi() - 360 == seed_id.iphi() + 1) {
243  s4s9_[1] += RecHitsInWindow[j].energy();
244  }
245  }
246  }
247  } else {
248  if (((EBDetId)RecHitsInWindow[j].id()).ieta() == seed_id.ieta()) {
249  if (((EBDetId)RecHitsInWindow[j].id()).iphi() == seed_id.iphi() - 1 ||
250  ((EBDetId)RecHitsInWindow[j].id()).iphi() - 360 == seed_id.iphi() - 1) {
251  s4s9_[0] += RecHitsInWindow[j].energy();
252  s4s9_[3] += RecHitsInWindow[j].energy();
253  } else {
254  if (((EBDetId)RecHitsInWindow[j].id()).iphi() == seed_id.iphi() + 1 ||
255  ((EBDetId)RecHitsInWindow[j].id()).iphi() - 360 == seed_id.iphi() + 1) {
256  s4s9_[1] += RecHitsInWindow[j].energy();
257  s4s9_[2] += RecHitsInWindow[j].energy();
258  }
259  }
260  } else {
261  if ((((EBDetId)RecHitsInWindow[j].id()).ieta() == seed_id.ieta() + 1 && seed_id.ieta() != -1) ||
262  (seed_id.ieta() == -1 && (((EBDetId)RecHitsInWindow[j].id()).ieta() == seed_id.ieta() + 2))) {
263  if (((EBDetId)RecHitsInWindow[j].id()).iphi() == seed_id.iphi() - 1 ||
264  ((EBDetId)RecHitsInWindow[j].id()).iphi() - 360 == seed_id.iphi() - 1) {
265  s4s9_[3] += RecHitsInWindow[j].energy();
266  } else {
267  if (((EBDetId)RecHitsInWindow[j].id()).iphi() == seed_id.iphi()) {
268  s4s9_[2] += RecHitsInWindow[j].energy();
269  s4s9_[3] += RecHitsInWindow[j].energy();
270  } else {
271  if (((EBDetId)RecHitsInWindow[j].id()).iphi() == seed_id.iphi() + 1 ||
272  ((EBDetId)RecHitsInWindow[j].id()).iphi() - 360 == seed_id.iphi() + 1) {
273  s4s9_[2] += RecHitsInWindow[j].energy();
274  }
275  }
276  }
277  } else {
278  cout << " (EBDetId)RecHitsInWindow[j].id()).ieta() " << ((EBDetId)RecHitsInWindow[j].id()).ieta()
279  << " seed_id.ieta() " << seed_id.ieta() << endl;
280  cout << " Problem with S4 calculation " << endl;
281  return;
282  }
283  }
284  }
285  }
286  s4s9Clus.push_back(*max_element(s4s9_, s4s9_ + 4) / simple_energy);
287  // cout<<" s4s9Clus[0] "<<s4s9_[0]/simple_energy<<" s4s9Clus[1] "<<s4s9_[1]/simple_energy<<" s4s9Clus[2] "<<s4s9_[2]/simple_energy<<" s4s9Clus[3] "<<s4s9_[3]/simple_energy<<endl;
288  // cout<<" Max "<<*max_element( s4s9_,s4s9_+4)/simple_energy<<endl;
289  nClus++;
290  if (nClus == MAXCLUS)
291  return;
292  } // End loop over seed clusters
293 
294  // cout<< " Pi0 clusters: "<<nClus<<endl;
295 
296  // Selection, based on Simple clustering
297  //pi0 candidates
298  static const int MAXPI0S = 200;
299  int npi0_s = 0;
300 
301  vector<EBDetId> scXtals;
302  scXtals.clear();
303 
304  if (nClus <= 1)
305  return;
306  for (Int_t i = 0; i < nClus; i++) {
307  for (Int_t j = i + 1; j < nClus; j++) {
308  // cout<<" i "<<i<<" etClus[i] "<<etClus[i]<<" j "<<j<<" etClus[j] "<<etClus[j]<<endl;
309  if (etClus[i] > selePtGammaOne_ && etClus[j] > selePtGammaTwo_ && s4s9Clus[i] > seleS4S9GammaOne_ &&
310  s4s9Clus[j] > seleS4S9GammaTwo_) {
311  float theta_0 = 2. * atan(exp(-etaClus[i]));
312  float theta_1 = 2. * atan(exp(-etaClus[j]));
313 
314  float p0x = eClus[i] * sin(theta_0) * cos(phiClus[i]);
315  float p1x = eClus[j] * sin(theta_1) * cos(phiClus[j]);
316  float p0y = eClus[i] * sin(theta_0) * sin(phiClus[i]);
317  float p1y = eClus[j] * sin(theta_1) * sin(phiClus[j]);
318  float p0z = eClus[i] * cos(theta_0);
319  float p1z = eClus[j] * cos(theta_1);
320 
321  float pt_pi0 = sqrt((p0x + p1x) * (p0x + p1x) + (p0y + p1y) * (p0y + p1y));
322  // cout<<" pt_pi0 "<<pt_pi0<<endl;
323  if (pt_pi0 < selePtPi0_)
324  continue;
325  float m_inv = sqrt((eClus[i] + eClus[j]) * (eClus[i] + eClus[j]) - (p0x + p1x) * (p0x + p1x) -
326  (p0y + p1y) * (p0y + p1y) - (p0z + p1z) * (p0z + p1z));
327  if ((m_inv < seleMinvMaxPi0_) && (m_inv > seleMinvMinPi0_)) {
328  //New Loop on cluster to measure isolation:
329  vector<int> IsoClus;
330  IsoClus.clear();
331  float Iso = 0;
332  TVector3 pi0vect = TVector3((p0x + p1x), (p0y + p1y), (p0z + p1z));
333  for (Int_t k = 0; k < nClus; k++) {
334  if (k == i || k == j)
335  continue;
336  TVector3 Clusvect = TVector3(eClus[k] * sin(2. * atan(exp(-etaClus[k]))) * cos(phiClus[k]),
337  eClus[k] * sin(2. * atan(exp(-etaClus[k]))) * sin(phiClus[k]),
338  eClus[k] * cos(2. * atan(exp(-etaClus[k]))));
339  float dretaclpi0 = fabs(etaClus[k] - pi0vect.Eta());
340  float drclpi0 = Clusvect.DeltaR(pi0vect);
341 
342  if ((drclpi0 < selePi0BeltDR_) && (dretaclpi0 < selePi0BeltDeta_)) {
343  Iso = Iso + etClus[k];
344  IsoClus.push_back(k);
345  }
346  }
347 
348  if (Iso / pt_pi0 < selePi0Iso_) {
349  hMinvPi0EB_->Fill(m_inv);
350  hPt1Pi0EB_->Fill(etClus[i]);
351  hPt2Pi0EB_->Fill(etClus[j]);
352  hPtPi0EB_->Fill(pt_pi0);
353  hIsoPi0EB_->Fill(Iso / pt_pi0);
354 
355  npi0_s++;
356  }
357 
358  if (npi0_s == MAXPI0S)
359  return;
360  }
361  }
362  }
363  }
364 }
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX)
Definition: DQMStore.cc:239
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:34
T getUntrackedParameter(std::string const &, T const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:418
PiZeroAnalyzer(const edm::ParameterSet &)
void makePizero(const edm::EventSetup &es, const edm::Handle< EcalRecHitCollection > eb, const edm::Handle< EcalRecHitCollection > ee)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< EcalRecHit >::const_iterator const_iterator
int iphi() const
get the crystal iphi
Definition: EBDetId.h:51
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
int ieta() const
get the crystal ieta
Definition: EBDetId.h:49
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
bool isValid() const
Definition: HandleBase.h:70
void analyze(const edm::Event &, const edm::EventSetup &) override
const_iterator end() const
virtual std::vector< DetId > getWindow(const DetId &id, const int &northSouthSize, const int &eastWestSize) const
T const * product() const
Definition: Handle.h:69
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const CaloSubdetectorTopology * getSubdetectorTopology(const DetId &id) const
access the subdetector Topology for the given subdetector directly
Definition: CaloTopology.cc:17
edm::EventID id() const
Definition: EventBase.h:59
HLT enums.
math::XYZPoint Calculate_Location(const HitsAndFractions &iDetIds, const edm::SortedCollection< HitType > *iRecHits, const CaloSubdetectorGeometry *iSubGeom, const CaloSubdetectorGeometry *iESGeom=0)
Definition: PositionCalc.h:65
T get() const
Definition: EventSetup.h:73
~PiZeroAnalyzer() override
const_iterator begin() const
Definition: Run.h:45
virtual void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)