00001 #include "Calibration/EcalCalibAlgos/interface/Pi0FixedMassWindowCalibration.h"
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include "CondFormats/DataRecord/interface/EcalIntercalibConstantsRcd.h"
00011 #include "Calibration/Tools/interface/Pi0CalibXMLwriter.h"
00012
00013
00014 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00015 #include "DataFormats/EgammaReco/interface/BasicClusterFwd.h"
00016
00017 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00018 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00019 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00020 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00021 #include "Geometry/CaloTopology/interface/EcalEndcapTopology.h"
00022 #include "Geometry/CaloTopology/interface/EcalBarrelTopology.h"
00023
00024
00025 #include "DataFormats/EgammaReco/interface/ClusterShape.h"
00026 #include "DataFormats/EgammaReco/interface/ClusterShapeFwd.h"
00027
00028
00029
00030 using namespace std;
00031
00032
00033
00034 Pi0FixedMassWindowCalibration::Pi0FixedMassWindowCalibration(const edm::ParameterSet& iConfig) :
00035 theMaxLoops( iConfig.getUntrackedParameter<unsigned int>("maxLoops",0) ),
00036 ecalHitsProducer_( iConfig.getParameter< std::string > ("ecalRecHitsProducer") ),
00037 barrelHits_( iConfig.getParameter< std::string > ("barrelHitCollection") )
00038 {
00039
00040 std::cout << "[Pi0FixedMassWindowCalibration] Constructor "<<std::endl;
00041
00042 std::string verbosityString = iConfig.getParameter<std::string>("VerbosityLevel");
00043 if (verbosityString == "DEBUG") verbosity = IslandClusterAlgo::pDEBUG;
00044 else if (verbosityString == "WARNING") verbosity = IslandClusterAlgo::pWARNING;
00045 else if (verbosityString == "INFO") verbosity = IslandClusterAlgo::pINFO;
00046 else verbosity = IslandClusterAlgo::pERROR;
00047
00048
00049 barrelClusterCollection_ = iConfig.getParameter<std::string>("barrelClusterCollection");
00050
00051
00052 double barrelSeedThreshold = iConfig.getParameter<double>("IslandBarrelSeedThr");
00053 double endcapSeedThreshold = iConfig.getParameter<double>("IslandEndcapSeedThr");
00054
00055
00056
00057 selePi0PtGammaOneMin_ = iConfig.getParameter<double>("selePi0PtGammaOneMin");
00058 selePi0PtGammaTwoMin_ = iConfig.getParameter<double>("selePi0PtGammaTwoMin");
00059
00060 selePi0DRBelt_ = iConfig.getParameter<double>("selePi0DRBelt");
00061 selePi0DetaBelt_ = iConfig.getParameter<double>("selePi0DetaBelt");
00062
00063 selePi0PtPi0Min_ = iConfig.getParameter<double>("selePi0PtPi0Min");
00064
00065 selePi0S4S9GammaOneMin_ = iConfig.getParameter<double>("selePi0S4S9GammaOneMin");
00066 selePi0S4S9GammaTwoMin_ = iConfig.getParameter<double>("selePi0S4S9GammaTwoMin");
00067 selePi0S9S25GammaOneMin_ = iConfig.getParameter<double>("selePi0S9S25GammaOneMin");
00068 selePi0S9S25GammaTwoMin_ = iConfig.getParameter<double>("selePi0S9S25GammaTwoMin");
00069
00070 selePi0EtBeltIsoRatioMax_ = iConfig.getParameter<double>("selePi0EtBeltIsoRatioMax");
00071
00072 selePi0MinvMeanFixed_ = iConfig.getParameter<double>("selePi0MinvMeanFixed");
00073 selePi0MinvSigmaFixed_ = iConfig.getParameter<double>("selePi0MinvSigmaFixed");
00074
00075
00076
00077
00078
00079 edm::ParameterSet posCalcParameters =
00080 iConfig.getParameter<edm::ParameterSet>("posCalcParameters");
00081 posCalculator_ = PositionCalc(posCalcParameters);
00082 shapeAlgo_ = ClusterShapeAlgo(posCalcParameters);
00083 clustershapecollectionEB_ = iConfig.getParameter<std::string>("clustershapecollectionEB");
00084
00085
00086 barrelClusterShapeAssociation_ = iConfig.getParameter<std::string>("barrelShapeAssociation");
00087
00088 island_p = new IslandClusterAlgo(barrelSeedThreshold, endcapSeedThreshold, posCalculator_,verbosity);
00089
00090 theParameterSet=iConfig;
00091
00092
00093 }
00094
00095
00096
00097
00098
00099 Pi0FixedMassWindowCalibration::~Pi0FixedMassWindowCalibration()
00100 {
00101
00102 }
00103
00104
00105
00106
00107 void Pi0FixedMassWindowCalibration::beginOfJob( )
00108 {
00109
00110
00111
00112 isfirstcall_=true;
00113
00114
00115
00116 }
00117
00118
00119 void Pi0FixedMassWindowCalibration::endOfJob()
00120 {
00121
00122 std::cout << "[Pi0FixedMassWindowCalibration] endOfJob"<<endl;
00123
00124
00125
00126 Pi0CalibXMLwriter barrelWriter(EcalBarrel,99);
00127
00128 std::vector<DetId>::const_iterator barrelIt=barrelCells.begin();
00129 for (; barrelIt!=barrelCells.end(); barrelIt++) {
00130 EBDetId eb(*barrelIt);
00131 int ieta = eb.ieta();
00132 int iphi = eb.iphi();
00133 int sign = eb.zside()>0 ? 1 : 0;
00134 barrelWriter.writeLine(eb,newCalibs_barl[abs(ieta)-1][iphi-1][sign]);
00135 }
00136
00137 }
00138
00139
00140
00141 void Pi0FixedMassWindowCalibration::startingNewLoop(unsigned int iLoop )
00142 {
00143
00144 for (int sign=0; sign<2; sign++) {
00145 for (int ieta=0; ieta<85; ieta++) {
00146 for (int iphi=0; iphi<360; iphi++) {
00147 wxtals[ieta][iphi][sign]=0.;
00148 mwxtals[ieta][iphi][sign]=0.;
00149 }
00150 }
00151 }
00152 std::cout << "[Pi0FixedMassWindowCalibration] Starting loop number " << iLoop<<std::endl;
00153
00154 }
00155
00156
00157
00158
00159
00160 edm::EDLooper::Status
00161 Pi0FixedMassWindowCalibration::endOfLoop(const edm::EventSetup& iSetup, unsigned int iLoop)
00162 {
00163 std::cout << "[Pi0FixedMassWindowCalibration] Ending loop " << iLoop<<std::endl;
00164
00165 for (int sign=0; sign<2; sign++) {
00166 for (int ieta=0; ieta<85; ieta++) {
00167 for (int iphi=0; iphi<360; iphi++) {
00168
00169 if (wxtals[ieta][iphi][sign] == 0 )
00170 {
00171 newCalibs_barl[ieta][iphi][sign]=oldCalibs_barl[ieta][iphi][sign];
00172 } else {
00173 newCalibs_barl[ieta][iphi][sign]=oldCalibs_barl[ieta][iphi][sign]*(mwxtals[ieta][iphi][sign]/wxtals[ieta][iphi][sign]);
00174 }
00175 cout<< " New calibration constant: ieta iphi sign - old,mwxtals,wxtals,new: "<<ieta<<" "<<iphi<<" "<<sign<<" - "<<oldCalibs_barl[ieta][iphi][sign]<<" "<<mwxtals[ieta][iphi][sign]<<" "<<wxtals[ieta][iphi][sign]<<" "<<newCalibs_barl[ieta][iphi][sign]<<endl;
00176
00177 }
00178 }
00179 }
00180
00181 Pi0CalibXMLwriter barrelWriter(EcalBarrel,iLoop+1);
00182
00183 std::vector<DetId>::const_iterator barrelIt=barrelCells.begin();
00184 for (; barrelIt!=barrelCells.end(); barrelIt++) {
00185 EBDetId eb(*barrelIt);
00186 int ieta = eb.ieta();
00187 int iphi = eb.iphi();
00188 int sign = eb.zside()>0 ? 1 : 0;
00189 barrelWriter.writeLine(eb,newCalibs_barl[abs(ieta)-1][iphi-1][sign]);
00190 if (iphi==1) {
00191 std::cout << "Calib constant for barrel crystal "
00192 << " (" << ieta << "," << iphi << ") changed from "
00193 << oldCalibs_barl[abs(ieta)-1][iphi-1][sign] << " to "
00194 << newCalibs_barl[abs(ieta)-1][iphi-1][sign] << std::endl;
00195 }
00196 }
00197
00198
00199
00200 for (int sign=0; sign<2; sign++) {
00201 for (int ieta=0; ieta<85; ieta++) {
00202 for (int iphi=0; iphi<360; iphi++) {
00203 oldCalibs_barl[ieta][iphi][sign]=newCalibs_barl[ieta][iphi][sign];
00204 newCalibs_barl[ieta][iphi][sign]=0;
00205 }
00206 }
00207 }
00208
00209
00210 if ( iLoop == theMaxLoops-1 || iLoop >= theMaxLoops ) return kStop;
00211 else return kContinue;
00212 }
00213
00214
00215
00216
00217 edm::EDLooper::Status
00218 Pi0FixedMassWindowCalibration::duringLoop(const edm::Event& event,
00219 const edm::EventSetup& setup )
00220 {
00221 using namespace edm;
00222 using namespace std;
00223
00224
00225
00226
00227 if (isfirstcall_){
00228
00229
00230 for (int sign=0; sign<2; sign++) {
00231 for (int ieta=0; ieta<85; ieta++) {
00232 for (int iphi=0; iphi<360; iphi++) {
00233 oldCalibs_barl[ieta][iphi][sign]=0.;
00234 newCalibs_barl[ieta][iphi][sign]=0.;
00235 wxtals[ieta][iphi][sign]=0.;
00236 mwxtals[ieta][iphi][sign]=0.;
00237 }
00238 }
00239 }
00240
00241
00242
00243 edm::ESHandle<EcalIntercalibConstants> pIcal;
00244 EcalIntercalibConstantMap imap;
00245
00246 try {
00247 setup.get<EcalIntercalibConstantsRcd>().get(pIcal);
00248 std::cout << "Taken EcalIntercalibConstants" << std::endl;
00249 imap = pIcal.product()->getMap();
00250 std::cout << "imap.size() = " << imap.size() << std::endl;
00251 } catch ( std::exception& ex ) {
00252 std::cerr << "Error! can't get EcalIntercalibConstants " << std::endl;
00253 }
00254
00255
00256 edm::ESHandle<CaloGeometry> geoHandle;
00257 setup.get<CaloGeometryRecord>().get(geoHandle);
00258 const CaloGeometry& geometry = *geoHandle;
00259
00260
00261
00262 barrelCells = geometry.getValidDetIds(DetId::Ecal, EcalBarrel);
00263 std::vector<DetId>::const_iterator barrelIt;
00264 for (barrelIt=barrelCells.begin(); barrelIt!=barrelCells.end(); barrelIt++) {
00265 EBDetId eb(*barrelIt);
00266
00267
00268 EcalIntercalibConstantMap::const_iterator itcalib = imap.find(eb.rawId());
00269 if ( itcalib == imap.end() ) {
00270
00271 }
00272 EcalIntercalibConstant calib = (*itcalib);
00273 int sign = eb.zside()>0 ? 1 : 0;
00274 oldCalibs_barl[abs(eb.ieta())-1][eb.iphi()-1][sign] = calib;
00275 if (eb.iphi()==1) std::cout << "Read old constant for crystal "
00276 << " (" << eb.ieta() << "," << eb.iphi()
00277 << ") : " << calib << std::endl;
00278
00279
00280
00281 }
00282 isfirstcall_=false;
00283 }
00284
00285
00286
00287 nevent++;
00288
00289 if ((nevent<100 && nevent%10==0)
00290 ||(nevent<1000 && nevent%100==0)
00291 ||(nevent<10000 && nevent%100==0)
00292 ||(nevent<100000 && nevent%1000==0)
00293 ||(nevent<10000000 && nevent%1000==0))
00294 std::cout << "[Pi0FixedMassWindowCalibration] Events processed: "<<nevent<<std::endl;
00295
00296 recHitsEB_map = new std::map<DetId, EcalRecHit>();
00297
00298 EcalRecHitCollection* recalibEcalRecHitCollection( new EcalRecHitCollection );
00299
00300 int nRecHitsEB=0;
00301 Handle<EcalRecHitCollection> pEcalRecHitBarrelCollection;
00302 event.getByLabel(ecalHitsProducer_, barrelHits_, pEcalRecHitBarrelCollection);
00303 const EcalRecHitCollection* ecalRecHitBarrelCollection = pEcalRecHitBarrelCollection.product();
00304 cout << " ECAL Barrel RecHits # "<< ecalRecHitBarrelCollection->size() <<endl;
00305 for(EcalRecHitCollection::const_iterator aRecHitEB = ecalRecHitBarrelCollection->begin(); aRecHitEB != ecalRecHitBarrelCollection->end(); aRecHitEB++) {
00306
00307
00308
00309 EBDetId ebrhdetid = aRecHitEB->detid();
00310
00311
00312
00313
00314 int sign = ebrhdetid.zside()>0 ? 1 : 0;
00315 EcalRecHit aHit(aRecHitEB->id(),aRecHitEB->energy()*oldCalibs_barl[abs(ebrhdetid.ieta())-1][ebrhdetid.iphi()-1][sign],aRecHitEB->time());
00316 recalibEcalRecHitCollection->push_back(aHit);
00317
00318 nRecHitsEB++;
00319
00320 }
00321
00322
00323 int irecalib=0;
00324 for(EcalRecHitCollection::const_iterator aRecHitEB = recalibEcalRecHitCollection->begin(); aRecHitEB != recalibEcalRecHitCollection->end(); aRecHitEB++) {
00325
00326
00327
00328
00329
00330
00331
00332 std::pair<DetId, EcalRecHit> map_entry(aRecHitEB->id(), *aRecHitEB);
00333 recHitsEB_map->insert(map_entry);
00334
00335
00336 irecalib++;
00337
00338 }
00339
00340
00341
00342 edm::ESHandle<CaloGeometry> geoHandle;
00343 setup.get<CaloGeometryRecord>().get(geoHandle);
00344
00345 const CaloSubdetectorGeometry *geometry_p;
00346 CaloSubdetectorTopology *topology_p;
00347
00348 std::string clustershapetag;
00349 geometry_p = geoHandle->getSubdetectorGeometry(DetId::Ecal, EcalBarrel);
00350 topology_p = new EcalBarrelTopology(geoHandle);
00351
00352 const CaloSubdetectorGeometry *geometryES_p;
00353 geometryES_p = geoHandle->getSubdetectorGeometry(DetId::Ecal, EcalPreshower);
00354
00355
00356
00357
00358
00359
00360
00361
00362
00363
00364
00365
00366
00367
00368
00369
00370
00371
00372
00373
00374
00375
00376
00377
00378
00379
00380
00381 reco::BasicClusterCollection clusters_recalib;
00382 clusters_recalib = island_p->makeClusters(recalibEcalRecHitCollection,geometry_p,topology_p,geometryES_p,IslandClusterAlgo::barrel);
00383
00384
00385 std::vector <reco::ClusterShape> ClusVec_recalib;
00386 for (int erg=0;erg<int(clusters_recalib.size());++erg){
00387 reco::ClusterShape TestShape_recalib = shapeAlgo_.Calculate(clusters_recalib[erg],recalibEcalRecHitCollection,geometry_p,topology_p);
00388 ClusVec_recalib.push_back(TestShape_recalib);
00389 }
00390
00391
00392 std::auto_ptr<reco::ClusterShapeCollection> clustersshapes_p_recalib(new reco::ClusterShapeCollection);
00393 clustersshapes_p_recalib->assign(ClusVec_recalib.begin(), ClusVec_recalib.end());
00394
00395 cout<<"[Pi0FixedMassWindowCalibration][recalibration] Basic Cluster collection size: "<<clusters_recalib.size()<<endl;
00396 cout<<"[Pi0FixedMassWindowCalibration][recalibration] Basic Cluster Shape Collection size: "<<clustersshapes_p_recalib->size()<<endl;
00397
00398
00399
00400
00401
00402
00403 static const int MAXBCEB = 200;
00404 static const int MAXBCEBRH = 200;
00405 int nIslandBCEB;
00406 float eIslandBCEB[MAXBCEB];
00407 float etIslandBCEB[MAXBCEB];
00408 float etaIslandBCEB[MAXBCEB];
00409 float phiIslandBCEB[MAXBCEB];
00410 float e2x2IslandBCEB[MAXBCEB];
00411 float e3x3IslandBCEB[MAXBCEB];
00412 float e5x5IslandBCEB[MAXBCEB];
00413
00414
00415 int nIslandBCEBRecHits[MAXBCEB];
00416
00417 int ietaIslandBCEBRecHits[MAXBCEB][MAXBCEBRH];
00418 int iphiIslandBCEBRecHits[MAXBCEB][MAXBCEBRH];
00419 int zsideIslandBCEBRecHits[MAXBCEB][MAXBCEBRH];
00420 float eIslandBCEBRecHits[MAXBCEB][MAXBCEBRH];
00421
00422 nIslandBCEB=0;
00423 for(int i=0; i<MAXBCEB; i++){
00424 eIslandBCEB[i] = 0;
00425 etIslandBCEB[i] = 0;
00426 etaIslandBCEB[i] = 0;
00427 phiIslandBCEB[i] = 0;
00428 e2x2IslandBCEB[i] = 0;
00429 e3x3IslandBCEB[i] = 0;
00430 e5x5IslandBCEB[i] = 0;
00431 nIslandBCEBRecHits[i] = 0;
00432 for(int j=0;j<MAXBCEBRH;j++){
00433
00434 ietaIslandBCEBRecHits[i][j] = 0;
00435 iphiIslandBCEBRecHits[i][j] = 0;
00436 zsideIslandBCEBRecHits[i][j] = 0;
00437 eIslandBCEBRecHits[i][j] = 0;
00438 }
00439 }
00440
00441
00442 int iClus_recalib=0;
00443 for(reco::BasicClusterCollection::const_iterator aClus = clusters_recalib.begin(); aClus != clusters_recalib.end(); aClus++) {
00444 cout<<" CLUSTER [recalibration] : #,NHits,e,et,eta,phi,e2x2,e3x3,e5x5: "<<iClus_recalib<<" "<<aClus->size()<<" "<<aClus->energy()<<" "<<aClus->energy()*sin(aClus->position().theta())<<" "<<aClus->position().eta()<<" "<<aClus->position().phi()<<" "<<(*clustersshapes_p_recalib)[iClus_recalib].e2x2()<<" "<<(*clustersshapes_p_recalib)[iClus_recalib].e3x3()<<" "<<(*clustersshapes_p_recalib)[iClus_recalib].e5x5()<<endl;
00445
00446 eIslandBCEB[nIslandBCEB] = aClus->energy();
00447 etIslandBCEB[nIslandBCEB] = aClus->energy()*sin(aClus->position().theta());
00448 etaIslandBCEB[nIslandBCEB] = aClus->position().eta();
00449 phiIslandBCEB[nIslandBCEB] = aClus->position().phi();
00450
00451 e2x2IslandBCEB[nIslandBCEB] = (*clustersshapes_p_recalib)[nIslandBCEB].e2x2();
00452 e3x3IslandBCEB[nIslandBCEB] = (*clustersshapes_p_recalib)[nIslandBCEB].e3x3();
00453 e5x5IslandBCEB[nIslandBCEB] = (*clustersshapes_p_recalib)[nIslandBCEB].e5x5();
00454
00455 nIslandBCEBRecHits[nIslandBCEB] = aClus->size();
00456
00457 std::vector<std::pair< DetId,float> > hits = aClus->hitsAndFractions();
00458 std::vector<std::pair< DetId,float> >::iterator hit;
00459 std::map<DetId, EcalRecHit>::iterator aHit;
00460
00461 int irhcount=0;
00462 for(hit = hits.begin(); hit != hits.end(); hit++)
00463 {
00464
00465 aHit = recHitsEB_map->find((*hit).first);
00466
00467
00468 EBDetId sel_rh = aHit->second.detid();
00469
00470
00471
00472
00473 ietaIslandBCEBRecHits[nIslandBCEB][irhcount] = sel_rh.ieta();
00474 iphiIslandBCEBRecHits[nIslandBCEB][irhcount] = sel_rh.iphi();
00475 zsideIslandBCEBRecHits[nIslandBCEB][irhcount] = sel_rh.zside();
00476 eIslandBCEBRecHits[nIslandBCEB][irhcount] = aHit->second.energy();
00477
00478 irhcount++;
00479 }
00480 nIslandBCEB++;
00481 iClus_recalib++;
00482
00483
00484 }
00485
00486
00487
00488 if (nIslandBCEB > 1)
00489 {
00490 for(int i=0 ; i<nIslandBCEB ; i++)
00491 {
00492 for(int j=i+1 ; j<nIslandBCEB ; j++)
00493 {
00494
00495 if( etIslandBCEB[i]>selePi0PtGammaOneMin_ && etIslandBCEB[j]>selePi0PtGammaOneMin_)
00496 {
00497
00498 float theta_0 = 2. * atan(exp(-etaIslandBCEB[i]));
00499 float theta_1 = 2. * atan(exp(-etaIslandBCEB[j]));
00500
00501 float p0x = eIslandBCEB[i] * sin(theta_0) * cos(phiIslandBCEB[i]);
00502 float p1x = eIslandBCEB[j] * sin(theta_1) * cos(phiIslandBCEB[j]);
00503
00504 float p0y = eIslandBCEB[i] * sin(theta_0) * sin(phiIslandBCEB[i]);
00505 float p1y = eIslandBCEB[j] * sin(theta_1) * sin(phiIslandBCEB[j]);
00506
00507 float p0z = eIslandBCEB[i] * cos(theta_0);
00508 float p1z = eIslandBCEB[j] * cos(theta_1);
00509
00510 float pi0_px = p0x + p1x;
00511 float pi0_py = p0y + p1y;
00512 float pi0_pz = p0z + p1z;
00513
00514 float pi0_ptot = sqrt (pi0_px*pi0_px + pi0_py*pi0_py + pi0_pz*pi0_pz);
00515
00516 float pi0_theta = acos(pi0_pz/pi0_ptot);
00517 float pi0_eta = -log(tan(pi0_theta/2));
00518 float pi0_phi = atan(pi0_py/pi0_px);
00519
00520
00521
00522
00523 float et_belt=0;
00524 for(Int_t k=0 ; k<nIslandBCEB ; k++)
00525 {
00526 if( (k != i) && (k != j) )
00527 {
00528 float dr_pi0_k = sqrt( (etaIslandBCEB[k]-pi0_eta)*(etaIslandBCEB[k]-pi0_eta) + (phiIslandBCEB[k]-pi0_phi)*(phiIslandBCEB[k]-pi0_phi) );
00529 float deta_pi0_k = fabs(etaIslandBCEB[k]-pi0_eta);
00530 if ( (dr_pi0_k<selePi0DRBelt_) && (deta_pi0_k<selePi0DetaBelt_) ) et_belt = et_belt + etIslandBCEB[k];
00531 }
00532 }
00533
00534
00535 float pt_pi0 = sqrt( (p0x+p1x)*(p0x+p1x) + (p0y+p1y)*(p0y+p1y));
00536
00537
00538
00539 if (pt_pi0 > selePi0PtPi0Min_)
00540 {
00541 float m_inv = sqrt ( (eIslandBCEB[i] + eIslandBCEB[j])*(eIslandBCEB[i] + eIslandBCEB[j]) - (p0x+p1x)*(p0x+p1x) - (p0y+p1y)*(p0y+p1y) - (p0z+p1z)*(p0z+p1z) );
00542 cout <<" pi0 (pt>2.5 GeV) m_inv = "<<m_inv<<endl;
00543
00544 float s4s9_1 = e2x2IslandBCEB[i]/e3x3IslandBCEB[i];
00545 float s4s9_2 = e2x2IslandBCEB[j]/e3x3IslandBCEB[j];
00546
00547 float s9s25_1 = e3x3IslandBCEB[i]/e5x5IslandBCEB[i];
00548 float s9s25_2 = e3x3IslandBCEB[j]/e5x5IslandBCEB[j];
00549
00550
00551
00552
00553 if (s4s9_1>selePi0S4S9GammaOneMin_ && s4s9_2>selePi0S4S9GammaTwoMin_ && s9s25_1>selePi0S9S25GammaOneMin_ && s9s25_2>selePi0S9S25GammaTwoMin_ && (et_belt/pt_pi0 < selePi0EtBeltIsoRatioMax_) )
00554 {
00555
00556 if ( m_inv>(selePi0MinvMeanFixed_ - 2*selePi0MinvSigmaFixed_) && m_inv<(selePi0MinvMeanFixed_ + 2*selePi0MinvSigmaFixed_) )
00557 {
00558
00559 cout<<" Pi0 Good candidate : minv = "<<m_inv<<endl;
00560 for(int kk=0 ; kk<nIslandBCEBRecHits[i] ; kk++)
00561 {
00562 int ieta_xtal = ietaIslandBCEBRecHits[i][kk];
00563 int iphi_xtal = iphiIslandBCEBRecHits[i][kk];
00564 int sign = zsideIslandBCEBRecHits[i][kk]>0 ? 1 : 0;
00565 wxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign] = wxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign] + eIslandBCEBRecHits[i][kk]/e3x3IslandBCEB[i];
00566 mwxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign] = mwxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign] + (selePi0MinvMeanFixed_/m_inv)*(selePi0MinvMeanFixed_/m_inv)*(eIslandBCEBRecHits[i][kk]/e3x3IslandBCEB[i]);
00567 cout<< "[Pi0FixedMassWindowCalibration] eta, phi, sign, e, e3x3, wxtals and mwxtals: "<<ieta_xtal<<" "<<iphi_xtal<<" "<<sign<<" "<<eIslandBCEBRecHits[i][kk]<<" "<<e3x3IslandBCEB[i]<<" "<<wxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign]<<" "<<mwxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign]<<endl;
00568 }
00569
00570 for(int kk=0 ; kk<nIslandBCEBRecHits[j] ; kk++)
00571 {
00572 int ieta_xtal = ietaIslandBCEBRecHits[j][kk];
00573 int iphi_xtal = iphiIslandBCEBRecHits[j][kk];
00574 int sign = zsideIslandBCEBRecHits[j][kk]>0 ? 1 : 0;
00575 wxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign] = wxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign] + eIslandBCEBRecHits[j][kk]/e3x3IslandBCEB[j];
00576 mwxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign] = mwxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign] + (selePi0MinvMeanFixed_/m_inv)*(selePi0MinvMeanFixed_/m_inv)*(eIslandBCEBRecHits[j][kk]/e3x3IslandBCEB[j]);
00577 cout<< "[Pi0FixedMassWindowCalibration] eta, phi, sign, e, e3x3, wxtals and mwxtals: "<<ieta_xtal<<" "<<iphi_xtal<<" "<<sign<<" "<<eIslandBCEBRecHits[j][kk]<<" "<<e3x3IslandBCEB[j]<<" "<<wxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign]<<" "<<mwxtals[abs(ieta_xtal)-1][iphi_xtal-1][sign]<<endl;
00578 }
00579
00580 }
00581 }
00582 }
00583 }
00584
00585
00586 }
00587 }
00588
00589
00590 } else {
00591 cout<< " Not enough ECAL Barrel Basic Clusters: "<<nIslandBCEB<<endl;
00592 }
00593
00594
00595
00596
00597
00598
00599 return kContinue;
00600 }
00601
00602