47 seleXtalMinEnergy_ = pset.
getParameter<
double> (
"seleXtalMinEnergy");
48 clusSeedThr_ = pset.
getParameter<
double> (
"clusSeedThr");
51 ParameterLogWeighted_ = pset.
getParameter<
bool> (
"ParameterLogWeighted");
52 ParameterX0_ = pset.
getParameter<
double> (
"ParameterX0");
53 ParameterT0_barl_ = pset.
getParameter<
double> (
"ParameterT0_barl");
54 ParameterW0_ = pset.
getParameter<
double> (
"ParameterW0");
56 selePtGammaOne_ = pset.
getParameter<
double> (
"selePtGammaOne");
57 selePtGammaTwo_ = pset.
getParameter<
double> (
"selePtGammaTwo");
58 seleS4S9GammaOne_ = pset.
getParameter<
double> (
"seleS4S9GammaOne");
59 seleS4S9GammaTwo_ = pset.
getParameter<
double> (
"seleS4S9GammaTwo");
62 selePi0BeltDR_ = pset.
getParameter<
double> (
"selePi0BeltDR");
63 selePi0BeltDeta_ = pset.
getParameter<
double> (
"selePi0BeltDeta");
64 seleMinvMaxPi0_ = pset.
getParameter<
double> (
"seleMinvMaxPi0");
65 seleMinvMinPi0_ = pset.
getParameter<
double> (
"seleMinvMinPi0");
95 if (verbosity_ > 0 ) {
112 currentFolder_.str(
"");
113 currentFolder_ <<
"Egamma/PiZeroAnalyzer/";
119 hMinvPi0EB_ =
dbe_->
book1D(
"Pi0InvmassEB",
"Pi0 Invariant Mass in EB",100,0.,0.5);
122 hPt1Pi0EB_ =
dbe_->
book1D(
"Pt1Pi0EB",
"Pt 1st most energetic Pi0 photon in EB",100,0.,20.);
125 hPt2Pi0EB_ =
dbe_->
book1D(
"Pt2Pi0EB",
"Pt 2nd most energetic Pi0 photon in EB",100,0.,20.);
128 hPtPi0EB_ =
dbe_->
book1D(
"PtPi0EB",
"Pi0 Pt in EB",100,0.,20.);
131 hIsoPi0EB_ =
dbe_->
book1D(
"IsoPi0EB",
"Pi0 Iso in EB",50,0.,1.);
149 if (nEvt_% prescaleFactor_ )
return;
151 LogInfo(
"PiZeroAnalyzer") <<
"PiZeroAnalyzer Analyzing event number: " << e.
id() <<
" Global Counter " << nEvt_ <<
"\n";
155 bool validEcalRecHits=
true;
158 e.
getByLabel(barrelEcalHits_, barrelHitHandle);
159 if (!barrelHitHandle.isValid()) {
160 edm::LogError(
"PhotonProducer") <<
"Error! Can't get the product "<<barrelEcalHits_.label();
161 validEcalRecHits=
false;
165 e.
getByLabel(endcapEcalHits_, endcapHitHandle);
167 if (!endcapHitHandle.isValid()) {
168 edm::LogError(
"PhotonProducer") <<
"Error! Can't get the product "<<endcapEcalHits_.label();
169 validEcalRecHits=
false;
172 if (validEcalRecHits) makePizero(esup, barrelHitHandle, endcapHitHandle);
200 std::map<DetId, EcalRecHit> recHitsEB_map;
202 std::vector<EcalRecHit> seeds;
206 vector<EBDetId> usedXtals;
211 static const int MAXCLUS = 2000;
214 vector<float> etClus;
215 vector<float> etaClus;
216 vector<float> phiClus;
217 vector<EBDetId> max_hit;
218 vector< vector<EcalRecHit> > RecHitsCluster;
219 vector<float> s4s9Clus;
222 for(itb=rhEB->begin(); itb!=rhEB->end(); ++itb){
224 double energy = itb->energy();
225 if (energy > seleXtalMinEnergy_) {
226 std::pair<DetId, EcalRecHit> map_entry(itb->id(), *itb);
227 recHitsEB_map.insert(map_entry);
229 if (energy > clusSeedThr_) seeds.push_back(*itb);
233 for (std::vector<EcalRecHit>::iterator itseed=seeds.begin(); itseed!=seeds.end(); itseed++) {
234 EBDetId seed_id = itseed->id();
235 std::vector<EBDetId>::const_iterator usedIds;
237 bool seedAlreadyUsed=
false;
238 for(usedIds=usedXtals.begin(); usedIds!=usedXtals.end(); usedIds++){
239 if(*usedIds==seed_id){
240 seedAlreadyUsed=
true;
245 if(seedAlreadyUsed)
continue;
247 std::vector<DetId> clus_v = topology_p->
getWindow(seed_id,clusEtaSize_,clusPhiSize_);
249 std::vector<std::pair<DetId, float> > clus_used;
251 vector<EcalRecHit> RecHitsInWindow;
253 double simple_energy = 0;
255 for (std::vector<DetId>::iterator det=clus_v.begin(); det!=clus_v.end(); det++) {
258 bool HitAlreadyUsed=
false;
259 for(usedIds=usedXtals.begin(); usedIds!=usedXtals.end(); usedIds++){
265 if(HitAlreadyUsed)
continue;
266 if (recHitsEB_map.find(*det) != recHitsEB_map.end()){
268 std::map<DetId, EcalRecHit>::iterator aHit;
269 aHit = recHitsEB_map.find(*det);
270 usedXtals.push_back(*det);
271 RecHitsInWindow.push_back(aHit->second);
272 clus_used.push_back( std::pair<DetId, float>(*det, 1.) );
273 simple_energy = simple_energy + aHit->second.energy();
277 math::XYZPoint clus_pos = posCalculator_.Calculate_Location(clus_used,hitCollection_p,geometry_p,geometryES_p);
278 float theta_s = 2. * atan(
exp(-clus_pos.eta()));
279 float p0x_s = simple_energy *
sin(theta_s) *
cos(clus_pos.phi());
280 float p0y_s = simple_energy *
sin(theta_s) *
sin(clus_pos.phi());
282 float et_s =
sqrt( p0x_s*p0x_s + p0y_s*p0y_s);
286 eClus.push_back(simple_energy);
287 etClus.push_back(et_s);
288 etaClus.push_back(clus_pos.eta());
289 phiClus.push_back(clus_pos.phi());
290 max_hit.push_back(seed_id);
291 RecHitsCluster.push_back(RecHitsInWindow);
295 for(
int i=0;
i<4;
i++)s4s9_[
i]= itseed->energy();
296 for(
unsigned int j=0;
j<RecHitsInWindow.size();
j++){
298 if((((
EBDetId)RecHitsInWindow[
j].
id()).ieta() == seed_id.
ieta()-1 && seed_id.
ieta()!=1 ) || ( seed_id.
ieta()==1 && (((
EBDetId)RecHitsInWindow[
j].
id()).ieta() == seed_id.
ieta()-2))){
299 if(((
EBDetId)RecHitsInWindow[
j].id()).iphi() == seed_id.
iphi()-1 ||((
EBDetId)RecHitsInWindow[
j].
id()).iphi()-360 == seed_id.
iphi()-1 ){
300 s4s9_[0]+=RecHitsInWindow[
j].energy();
302 if(((
EBDetId)RecHitsInWindow[
j].
id()).iphi() == seed_id.
iphi()){
303 s4s9_[0]+=RecHitsInWindow[
j].energy();
304 s4s9_[1]+=RecHitsInWindow[
j].energy();
306 if(((
EBDetId)RecHitsInWindow[
j].
id()).iphi() == seed_id.
iphi()+1 ||((
EBDetId)RecHitsInWindow[
j].
id()).iphi()-360 == seed_id.
iphi()+1 ){
307 s4s9_[1]+=RecHitsInWindow[
j].energy();
312 if(((
EBDetId)RecHitsInWindow[
j].
id()).ieta() == seed_id.
ieta()){
313 if(((
EBDetId)RecHitsInWindow[
j].
id()).iphi() == seed_id.
iphi()-1 ||((
EBDetId)RecHitsInWindow[
j].
id()).iphi()-360 == seed_id.
iphi()-1 ){
314 s4s9_[0]+=RecHitsInWindow[
j].energy();
315 s4s9_[3]+=RecHitsInWindow[
j].energy();
317 if(((
EBDetId)RecHitsInWindow[
j].
id()).iphi() == seed_id.
iphi()+1 ||((
EBDetId)RecHitsInWindow[
j].
id()).iphi()-360 == seed_id.
iphi()+1 ){
318 s4s9_[1]+=RecHitsInWindow[
j].energy();
319 s4s9_[2]+=RecHitsInWindow[
j].energy();
323 if((((
EBDetId)RecHitsInWindow[
j].
id()).ieta() == seed_id.
ieta()+1 && seed_id.
ieta()!=-1 ) || ( seed_id.
ieta()==-1 && (((
EBDetId)RecHitsInWindow[
j].
id()).ieta() == seed_id.
ieta()+2))){
324 if(((
EBDetId)RecHitsInWindow[
j].id()).iphi() == seed_id.
iphi()-1 ||((
EBDetId)RecHitsInWindow[
j].
id()).iphi()-360 == seed_id.
iphi()-1 ){
325 s4s9_[3]+=RecHitsInWindow[
j].energy();
327 if(((
EBDetId)RecHitsInWindow[
j].
id()).iphi() == seed_id.
iphi()){
328 s4s9_[2]+=RecHitsInWindow[
j].energy();
329 s4s9_[3]+=RecHitsInWindow[
j].energy();
331 if(((
EBDetId)RecHitsInWindow[
j].
id()).iphi() == seed_id.
iphi()+1 ||((
EBDetId)RecHitsInWindow[
j].
id()).iphi()-360 == seed_id.
iphi()+1 ){
332 s4s9_[2]+=RecHitsInWindow[
j].energy();
337 cout<<
" (EBDetId)RecHitsInWindow[j].id()).ieta() "<<((
EBDetId)RecHitsInWindow[
j].
id()).ieta()<<
" seed_id.ieta() "<<seed_id.
ieta()<<endl;
338 cout<<
" Problem with S4 calculation "<<endl;
return;
343 s4s9Clus.push_back(*max_element( s4s9_,s4s9_+4)/simple_energy);
347 if (nClus == MAXCLUS)
return;
354 static const int MAXPI0S = 200;
357 vector<EBDetId> scXtals;
360 if (nClus <= 1)
return;
361 for(Int_t
i=0 ;
i<nClus ;
i++){
362 for(Int_t
j=
i+1 ;
j<nClus ;
j++){
364 if( etClus[
i]>selePtGammaOne_ && etClus[
j]>selePtGammaTwo_ && s4s9Clus[
i]>seleS4S9GammaOne_ && s4s9Clus[
j]>seleS4S9GammaTwo_){
365 float theta_0 = 2. * atan(
exp(-etaClus[
i]));
366 float theta_1 = 2. * atan(
exp(-etaClus[
j]));
368 float p0x = eClus[
i] *
sin(theta_0) *
cos(phiClus[i]);
369 float p1x = eClus[
j] *
sin(theta_1) *
cos(phiClus[j]);
370 float p0y = eClus[
i] *
sin(theta_0) *
sin(phiClus[i]);
371 float p1y = eClus[
j] *
sin(theta_1) *
sin(phiClus[j]);
372 float p0z = eClus[
i] *
cos(theta_0);
373 float p1z = eClus[
j] *
cos(theta_1);
375 float pt_pi0 =
sqrt( (p0x+p1x)*(p0x+p1x) + (p0y+p1y)*(p0y+p1y));
377 if (pt_pi0 < selePtPi0_)
continue;
378 float m_inv =
sqrt ( (eClus[i] + eClus[j])*(eClus[i] + eClus[j]) - (p0x+p1x)*(p0x+p1x) - (p0y+p1y)*(p0y+p1y) - (p0z+p1z)*(p0z+p1z) );
379 if ( (m_inv<seleMinvMaxPi0_) && (m_inv>seleMinvMinPi0_) ){
385 TVector3 pi0vect = TVector3((p0x+p1x), (p0y+p1y), (p0z+p1z));
386 for(Int_t
k=0 ;
k<nClus ;
k++){
387 if(
k==i ||
k==j)
continue;
388 TVector3 Clusvect = TVector3(eClus[
k] *
sin(2. * atan(
exp(-etaClus[
k]))) *
cos(phiClus[k]), eClus[k] *
sin(2. * atan(
exp(-etaClus[k]))) *
sin(phiClus[k]) , eClus[k] *
cos(2. * atan(
exp(-etaClus[k]))));
389 float dretaclpi0 = fabs(etaClus[k] - pi0vect.Eta());
390 float drclpi0 = Clusvect.DeltaR(pi0vect);
392 if((drclpi0<selePi0BeltDR_) && (dretaclpi0<selePi0BeltDeta_) ){
394 Iso = Iso + etClus[
k];
395 IsoClus.push_back(k);
400 if(Iso/pt_pi0<selePi0Iso_){
402 hMinvPi0EB_->Fill(m_inv);
403 hPt1Pi0EB_->Fill(etClus[i]);
404 hPt2Pi0EB_->Fill(etClus[j]);
405 hPtPi0EB_->Fill(pt_pi0);
406 hIsoPi0EB_->Fill(Iso/pt_pi0);
412 if(npi0_s == MAXPI0S)
return;
427 bool outputMEsInRootFile = parameters_.getParameter<
bool>(
"OutputMEsInRootFile");
428 std::string
outputFileName = parameters_.getParameter<std::string>(
"OutputFileName");
429 if(outputMEsInRootFile){
433 edm::LogInfo(
"PiZeroAnalyzer") <<
"Analyzed " << nEvt_ <<
"\n";
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
PiZeroAnalyzer(const edm::ParameterSet &)
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE")
void makePizero(const edm::EventSetup &es, const edm::Handle< EcalRecHitCollection > eb, const edm::Handle< EcalRecHitCollection > ee)
Sin< T >::type sin(const T &t)
std::vector< T >::const_iterator const_iterator
Exp< T >::type exp(const T &t)
int iphi() const
get the crystal iphi
Cos< T >::type cos(const T &t)
void setVerbose(unsigned level)
int ieta() const
get the crystal ieta
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
virtual std::vector< DetId > getWindow(const DetId &id, const int &northSouthSize, const int &eastWestSize) const
XYZPointD XYZPoint
point in space with cartesian internal representation
virtual void analyze(const edm::Event &, const edm::EventSetup &)
T const * product() const
virtual ~PiZeroAnalyzer()
void showDirStructure(void) const
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
void setCurrentFolder(const std::string &fullpath)