00001 #include <iostream>
00002 #include <sstream>
00003 #include <istream>
00004 #include <fstream>
00005 #include <iomanip>
00006 #include <string>
00007 #include <cmath>
00008 #include <functional>
00009 #include <stdlib.h>
00010 #include <string.h>
00011
00012 #include "HLTrigger/HLTanalyzers/interface/HLTAlCa.h"
00013
00014 HLTAlCa::HLTAlCa() {
00015 evtCounter=0;
00016
00017
00018 _Monte=false;
00019 _Debug=false;
00020
00021 TheMapping = new EcalElectronicsMapping();
00022 first_ = true;
00023 }
00024
00025
00026 void HLTAlCa::setup(const edm::ParameterSet& pSet, TTree* HltTree) {
00027
00028 edm::ParameterSet myEmParams = pSet.getParameter<edm::ParameterSet>("RunParameters") ;
00029 std::vector<std::string> parameterNames = myEmParams.getParameterNames() ;
00030
00031 for ( std::vector<std::string>::iterator iParam = parameterNames.begin();
00032 iParam != parameterNames.end(); iParam++ ){
00033 if ( (*iParam) == "Monte" ) _Monte = myEmParams.getParameter<bool>( *iParam );
00034 else if ( (*iParam) == "Debug" ) _Debug = myEmParams.getParameter<bool>( *iParam );
00035 }
00036
00037
00038 clusSeedThr_ = pSet.getParameter<double> ("clusSeedThr");
00039 clusSeedThrEndCap_ = pSet.getParameter<double> ("clusSeedThrEndCap");
00040 clusEtaSize_ = pSet.getParameter<int> ("clusEtaSize");
00041 clusPhiSize_ = pSet.getParameter<int> ("clusPhiSize");
00042 seleXtalMinEnergy_ = pSet.getParameter<double> ("seleXtalMinEnergy");
00043 RegionalMatch_ = pSet.getUntrackedParameter<bool>("RegionalMatch",true);
00044 ptMinEMObj_ = pSet.getParameter<double>("ptMinEMObj");
00045 EMregionEtaMargin_ = pSet.getParameter<double>("EMregionEtaMargin");
00046 EMregionPhiMargin_ = pSet.getParameter<double>("EMregionPhiMargin");
00047 Jets_ = pSet.getUntrackedParameter<bool>("Jets",false);
00048 JETSdoCentral_ = pSet.getUntrackedParameter<bool>("JETSdoCentral",true);
00049 JETSdoForward_ = pSet.getUntrackedParameter<bool>("JETSdoForward",true);
00050 JETSdoTau_ = pSet.getUntrackedParameter<bool>("JETSdoTau",true);
00051 JETSregionEtaMargin_ = pSet.getUntrackedParameter<double>("JETS_regionEtaMargin",1.0);
00052 JETSregionPhiMargin_ = pSet.getUntrackedParameter<double>("JETS_regionPhiMargin",1.0);
00053 Ptmin_jets_ = pSet.getUntrackedParameter<double>("Ptmin_jets",0.);
00054
00055
00056
00057 edm::ParameterSet posCalcParameters =
00058 pSet.getParameter<edm::ParameterSet>("posCalcParameters");
00059 posCalculator_ = PositionCalc(posCalcParameters);
00060
00061
00062 const int kMaxClus = 2000;
00063 ptClusAll = new float[kMaxClus];
00064 etaClusAll = new float[kMaxClus];
00065 phiClusAll = new float[kMaxClus];
00066 s4s9ClusAll = new float[kMaxClus];
00067
00068
00069 HltTree->Branch("ohHighestEnergyEERecHit",&ohHighestEnergyEERecHit,"ohHighestEnergyEERecHit/F");
00070 HltTree->Branch("ohHighestEnergyEBRecHit",&ohHighestEnergyEBRecHit,"ohHighestEnergyEBRecHit/F");
00071 HltTree->Branch("ohHighestEnergyHBHERecHit",&ohHighestEnergyHBHERecHit,"ohHighestEnergyHBHERecHit/F");
00072 HltTree->Branch("ohHighestEnergyHORecHit",&ohHighestEnergyHORecHit,"ohHighestEnergyHORecHit/F");
00073 HltTree->Branch("ohHighestEnergyHFRecHit",&ohHighestEnergyHFRecHit,"ohHighestEnergyHFRecHit/F");
00074 HltTree->Branch("Nalcapi0clusters",&Nalcapi0clusters,"Nalcapi0clusters/I");
00075 HltTree->Branch("ohAlcapi0ptClusAll",ptClusAll,"ohAlcapi0ptClusAll[Nalcapi0clusters]/F");
00076 HltTree->Branch("ohAlcapi0etaClusAll",etaClusAll,"ohAlcapi0etaClusAll[Nalcapi0clusters]/F");
00077 HltTree->Branch("ohAlcapi0phiClusAll",phiClusAll,"ohAlcapi0phiClusAll[Nalcapi0clusters]/F");
00078 HltTree->Branch("ohAlcapi0s4s9ClusAll",s4s9ClusAll,"ohAlcapi0s4s9ClusAll[Nalcapi0clusters]/F");
00079 }
00080
00081
00082 void HLTAlCa::analyze(const edm::Handle<EBRecHitCollection> & ebrechits,
00083 const edm::Handle<EERecHitCollection> & eerechits,
00084 const edm::Handle<HBHERecHitCollection> & hbherechits,
00085 const edm::Handle<HORecHitCollection> & horechits,
00086 const edm::Handle<HFRecHitCollection> & hfrechits,
00087 const edm::Handle<EBRecHitCollection> & pi0ebrechits,
00088 const edm::Handle<EERecHitCollection> & pi0eerechits,
00089 const edm::Handle<l1extra::L1EmParticleCollection> & l1EGIso,
00090 const edm::Handle<l1extra::L1EmParticleCollection> & l1EGNonIso,
00091 const edm::Handle<l1extra::L1JetParticleCollection> & l1extjetc,
00092 const edm::Handle<l1extra::L1JetParticleCollection> & l1extjetf,
00093 const edm::Handle<l1extra::L1JetParticleCollection> & l1exttaujet,
00094 const edm::ESHandle< EcalElectronicsMapping > & ecalmapping,
00095 const edm::ESHandle<CaloGeometry> & geoHandle,
00096 const edm::ESHandle<CaloTopology> & pTopology,
00097 const edm::ESHandle<L1CaloGeometry> & l1CaloGeom,
00098 TTree* HltTree) {
00099
00100
00101
00102 ohHighestEnergyEERecHit = -1.0;
00103 ohHighestEnergyEBRecHit = -1.0;
00104 ohHighestEnergyHBHERecHit = -1.0;
00105 ohHighestEnergyHORecHit = -1.0;
00106 ohHighestEnergyHFRecHit = -1.0;
00107
00108 int neerechits = 0;
00109 int nebrechits = 0;
00110 int nhbherechits = 0;
00111 int nhorechits = 0;
00112 int nhfrechits = 0;
00113
00114 if (ebrechits.isValid()) {
00115 EBRecHitCollection myebrechits;
00116 myebrechits = * ebrechits;
00117
00118 nebrechits = myebrechits.size();
00119 float ebrechitenergy = -1.0;
00120
00121 typedef EBRecHitCollection::const_iterator ebrechititer;
00122
00123 for (ebrechititer i=myebrechits.begin(); i!=myebrechits.end(); i++) {
00124 ebrechitenergy = i->energy();
00125 if(ebrechitenergy > ohHighestEnergyEBRecHit)
00126 ohHighestEnergyEBRecHit = ebrechitenergy;
00127 }
00128 }
00129 else {nebrechits = 0;}
00130
00131 if (eerechits.isValid()) {
00132 EERecHitCollection myeerechits;
00133 myeerechits = * eerechits;
00134
00135 neerechits = myeerechits.size();
00136 float eerechitenergy = -1.0;
00137
00138 typedef EERecHitCollection::const_iterator eerechititer;
00139
00140 for (eerechititer i=myeerechits.begin(); i!=myeerechits.end(); i++) {
00141 eerechitenergy = i->energy();
00142 if(eerechitenergy > ohHighestEnergyEERecHit)
00143 ohHighestEnergyEERecHit = eerechitenergy;
00144 }
00145 }
00146 else {neerechits = 0;}
00147
00148 if (hbherechits.isValid()) {
00149 HBHERecHitCollection myhbherechits;
00150 myhbherechits = * hbherechits;
00151
00152 nhbherechits = myhbherechits.size();
00153 float hbherechitenergy = -1.0;
00154
00155 typedef HBHERecHitCollection::const_iterator hbherechititer;
00156
00157 for (hbherechititer i=myhbherechits.begin(); i!=myhbherechits.end(); i++) {
00158 hbherechitenergy = i->energy();
00159 if(hbherechitenergy > ohHighestEnergyHBHERecHit)
00160 ohHighestEnergyHBHERecHit = hbherechitenergy;
00161 }
00162 }
00163 else {nhbherechits = 0;}
00164
00165 if (horechits.isValid()) {
00166 HORecHitCollection myhorechits;
00167 myhorechits = * horechits;
00168
00169 nhorechits = myhorechits.size();
00170 float horechitenergy = -1.0;
00171
00172 typedef HORecHitCollection::const_iterator horechititer;
00173
00174 for (horechititer i=myhorechits.begin(); i!=myhorechits.end(); i++) {
00175 horechitenergy = i->energy();
00176 if(horechitenergy > ohHighestEnergyHORecHit)
00177 ohHighestEnergyHORecHit = horechitenergy;
00178 }
00179 }
00180 else {nhorechits = 0;}
00181
00182 if (hfrechits.isValid()) {
00183 HFRecHitCollection myhfrechits;
00184 myhfrechits = * hfrechits;
00185
00186 nhfrechits = myhfrechits.size();
00187 float hfrechitenergy = -1.0;
00188
00189 typedef HFRecHitCollection::const_iterator hfrechititer;
00190
00191 for (hfrechititer i=myhfrechits.begin(); i!=myhfrechits.end(); i++) {
00192 hfrechitenergy = i->energy();
00193 if(hfrechitenergy > ohHighestEnergyHFRecHit)
00194 ohHighestEnergyHFRecHit = hfrechitenergy;
00195 }
00196 }
00197 else {nhfrechits = 0;}
00198
00200
00202
00203 if (first_) {
00204
00205 const EcalElectronicsMapping* TheMapping_ = ecalmapping.product();
00206 *TheMapping = *TheMapping_;
00207 first_ = false;
00208
00209 geometry_eb = geoHandle->getSubdetectorGeometry(DetId::Ecal,EcalBarrel);
00210 geometry_ee = geoHandle->getSubdetectorGeometry(DetId::Ecal,EcalEndcap);
00211 geometry_es = geoHandle->getSubdetectorGeometry(DetId::Ecal, EcalPreshower);
00212
00213 topology_eb = pTopology->getSubdetectorTopology(DetId::Ecal,EcalBarrel);
00214 topology_ee = pTopology->getSubdetectorTopology(DetId::Ecal,EcalEndcap);
00215 }
00216
00217
00218
00219
00220 FEDListUsed.clear();
00221 std::vector<int>::iterator it;
00222 if( RegionalMatch_){
00223
00224 for( l1extra::L1EmParticleCollection::const_iterator emItr = l1EGIso->begin();
00225 emItr != l1EGIso->end() ;++emItr ){
00226
00227 float pt = emItr -> pt();
00228
00229 if( pt< ptMinEMObj_ ) continue;
00230
00231 int etaIndex = emItr->gctEmCand()->etaIndex() ;
00232 int phiIndex = emItr->gctEmCand()->phiIndex() ;
00233 double etaLow = l1CaloGeom->etaBinLowEdge( etaIndex ) ;
00234 double etaHigh = l1CaloGeom->etaBinHighEdge( etaIndex ) ;
00235 double phiLow = l1CaloGeom->emJetPhiBinLowEdge( phiIndex ) ;
00236 double phiHigh = l1CaloGeom->emJetPhiBinHighEdge( phiIndex ) ;
00237
00238 std::vector<int> feds = ListOfFEDS(etaLow, etaHigh, phiLow, phiHigh, EMregionEtaMargin_, EMregionPhiMargin_);
00239 for (int n=0; n < (int)feds.size(); n++) {
00240 int fed = feds[n];
00241 it = find(FEDListUsed.begin(),FEDListUsed.end(),fed);
00242 if( it == FEDListUsed.end()){
00243 FEDListUsed.push_back(fed);
00244 }
00245 }
00246 }
00247
00248 for( l1extra::L1EmParticleCollection::const_iterator emItr = l1EGNonIso->begin();
00249 emItr != l1EGNonIso->end() ;++emItr ){
00250
00251 float pt = emItr -> pt();
00252
00253 if( pt< ptMinEMObj_ ) continue;
00254
00255 int etaIndex = emItr->gctEmCand()->etaIndex() ;
00256 int phiIndex = emItr->gctEmCand()->phiIndex() ;
00257 double etaLow = l1CaloGeom->etaBinLowEdge( etaIndex ) ;
00258 double etaHigh = l1CaloGeom->etaBinHighEdge( etaIndex ) ;
00259 double phiLow = l1CaloGeom->emJetPhiBinLowEdge( phiIndex ) ;
00260 double phiHigh = l1CaloGeom->emJetPhiBinHighEdge( phiIndex ) ;
00261
00262 std::vector<int> feds = ListOfFEDS(etaLow, etaHigh, phiLow, phiHigh, EMregionEtaMargin_, EMregionPhiMargin_);
00263 for (int n=0; n < (int)feds.size(); n++) {
00264 int fed = feds[n];
00265 it = find(FEDListUsed.begin(),FEDListUsed.end(),fed);
00266 if( it == FEDListUsed.end()){
00267 FEDListUsed.push_back(fed);
00268 }
00269 }
00270 }
00271
00272 if( Jets_ ){
00273
00274 double epsilon = 0.01;
00275
00276 if (JETSdoCentral_) {
00277
00278 for (l1extra::L1JetParticleCollection::const_iterator jetItr=l1extjetc->begin(); jetItr != l1extjetc->end(); jetItr++) {
00279
00280 double pt = jetItr-> pt();
00281 double eta = jetItr-> eta();
00282 double phi = jetItr-> phi();
00283
00284 if (pt < Ptmin_jets_ ) continue;
00285
00286 std::vector<int> feds = ListOfFEDS(eta, eta, phi-epsilon, phi+epsilon, JETSregionEtaMargin_, JETSregionPhiMargin_);
00287 for (int n=0; n < (int)feds.size(); n++) {
00288 int fed = feds[n];
00289 it = find(FEDListUsed.begin(),FEDListUsed.end(),fed);
00290 if( it == FEDListUsed.end()){
00291 FEDListUsed.push_back(fed);
00292 }
00293 }
00294 }
00295
00296 }
00297
00298 if (JETSdoForward_) {
00299
00300 for (l1extra::L1JetParticleCollection::const_iterator jetItr=l1extjetf->begin(); jetItr != l1extjetf->end(); jetItr++) {
00301
00302 double pt = jetItr -> pt();
00303 double eta = jetItr -> eta();
00304 double phi = jetItr -> phi();
00305
00306 if (pt < Ptmin_jets_ ) continue;
00307
00308 std::vector<int> feds = ListOfFEDS(eta, eta, phi-epsilon, phi+epsilon, JETSregionEtaMargin_, JETSregionPhiMargin_);
00309
00310 for (int n=0; n < (int)feds.size(); n++) {
00311 int fed = feds[n];
00312 it = find(FEDListUsed.begin(),FEDListUsed.end(),fed);
00313 if( it == FEDListUsed.end()){
00314 FEDListUsed.push_back(fed);
00315 }
00316 }
00317
00318 }
00319 }
00320
00321 if (JETSdoTau_) {
00322
00323 for (l1extra::L1JetParticleCollection::const_iterator jetItr=l1exttaujet->begin(); jetItr != l1exttaujet->end(); jetItr++) {
00324
00325 double pt = jetItr -> pt();
00326 double eta = jetItr -> eta();
00327 double phi = jetItr -> phi();
00328
00329 if (pt < Ptmin_jets_ ) continue;
00330
00331 std::vector<int> feds = ListOfFEDS(eta, eta, phi-epsilon, phi+epsilon, JETSregionEtaMargin_, JETSregionPhiMargin_);
00332 for (int n=0; n < (int)feds.size(); n++) {
00333 int fed = feds[n];
00334 it = find(FEDListUsed.begin(),FEDListUsed.end(),fed);
00335 if( it == FEDListUsed.end()){
00336 FEDListUsed.push_back(fed);
00337 }
00338 }
00339
00340 }
00341 }
00342 }
00343 }
00344
00346
00347 FEDListUsedBarrel.clear();
00348 FEDListUsedEndcap.clear();
00349 for( int j=0; j< int(FEDListUsed.size());j++){
00350 int fed = FEDListUsed[j];
00351
00352 if( fed >= 10 && fed <= 45){
00353 FEDListUsedBarrel.push_back(fed);
00354 }else FEDListUsedEndcap.push_back(fed);
00355 }
00356
00357
00358
00359 if (pi0ebrechits.isValid()) {
00360
00361 const EcalRecHitCollection *hitCollection_p = pi0ebrechits.product();
00362
00363 std::vector<EcalRecHit> seeds;
00364 seeds.clear();
00365
00366 std::vector<EBDetId> usedXtals;
00367 usedXtals.clear();
00368
00369 detIdEBRecHits.clear();
00370 EBRecHits.clear();
00371
00372 detIdEBRecHits.clear();
00373 EBRecHits.clear();
00374
00376 EBRecHitCollection::const_iterator itb;
00377
00378 for (itb=pi0ebrechits->begin(); itb!=pi0ebrechits->end(); itb++) {
00379 double energy = itb->energy();
00380 if( energy < seleXtalMinEnergy_) continue;
00381
00382 EBDetId det = itb->id();
00383
00384 if (RegionalMatch_){
00385 int fed = TheMapping->DCCid(det);
00386 it = find(FEDListUsedBarrel.begin(),FEDListUsedBarrel.end(),fed);
00387 if(it == FEDListUsedBarrel.end()) continue;
00388 }
00389
00390 detIdEBRecHits.push_back(det);
00391 EBRecHits.push_back(*itb);
00392
00393 if (energy > clusSeedThr_) seeds.push_back(*itb);
00394 }
00395
00396 int nClus;
00397 std::vector<float> eClus;
00398 std::vector<float> etClus;
00399 std::vector<float> etaClus;
00400 std::vector<float> phiClus;
00401 std::vector<EBDetId> max_hit;
00402 std::vector< std::vector<EcalRecHit> > RecHitsCluster;
00403 std::vector<float> s4s9Clus;
00404
00405 nClus=0;
00406
00407
00408 sort(seeds.begin(), seeds.end(), eecalRecHitLess());
00409
00410 Nalcapi0clusters = 0;
00411 nClusAll = 0;
00412
00413 for (std::vector<EcalRecHit>::iterator itseed=seeds.begin(); itseed!=seeds.end(); itseed++) {
00414 EBDetId seed_id = itseed->id();
00415 std::vector<EBDetId>::const_iterator usedIds;
00416
00417 bool seedAlreadyUsed=false;
00418 for(usedIds=usedXtals.begin(); usedIds!=usedXtals.end(); usedIds++){
00419 if(*usedIds==seed_id){
00420 seedAlreadyUsed=true;
00421 break;
00422 }
00423 }
00424
00425 if(seedAlreadyUsed)continue;
00426
00427 std::vector<DetId> clus_v = topology_eb->getWindow(seed_id,clusEtaSize_,clusPhiSize_);
00428 std::vector< std::pair<DetId, float> > clus_used;
00429
00430 std::vector<EcalRecHit> RecHitsInWindow;
00431
00432 double simple_energy = 0;
00433
00434 for (std::vector<DetId>::iterator det=clus_v.begin(); det!=clus_v.end(); det++) {
00435 EBDetId EBdet = *det;
00436
00437 bool HitAlreadyUsed=false;
00438 for(usedIds=usedXtals.begin(); usedIds!=usedXtals.end(); usedIds++){
00439 if(*usedIds==*det){
00440 HitAlreadyUsed=true;
00441 break;
00442 }
00443 }
00444
00445 if(HitAlreadyUsed)continue;
00446
00448 if (RegionalMatch_){
00449 int fed = TheMapping->DCCid(EBdet);
00450 it = find(FEDListUsedBarrel.begin(),FEDListUsedBarrel.end(),fed);
00451 if(it == FEDListUsedBarrel.end()) continue;
00452 }
00453
00454
00455 std::vector<EBDetId>::iterator itdet = find( detIdEBRecHits.begin(),detIdEBRecHits.end(),EBdet);
00456 if(itdet == detIdEBRecHits.end()) continue;
00457
00458 int nn = int(itdet - detIdEBRecHits.begin());
00459 usedXtals.push_back(*det);
00460 RecHitsInWindow.push_back(EBRecHits[nn]);
00461 clus_used.push_back( std::pair<DetId, float>(*det, 1) );
00462 simple_energy = simple_energy + EBRecHits[nn].energy();
00463 }
00464
00465 math::XYZPoint clus_pos = posCalculator_.Calculate_Location(clus_used,hitCollection_p,geometry_eb,geometry_es);
00466
00467 float theta_s = 2. * atan(exp(-clus_pos.eta()));
00468 float p0x_s = simple_energy * sin(theta_s) * cos(clus_pos.phi());
00469 float p0y_s = simple_energy * sin(theta_s) * sin(clus_pos.phi());
00470 float et_s = sqrt( p0x_s*p0x_s + p0y_s*p0y_s);
00471
00472 eClus.push_back(simple_energy);
00473 etClus.push_back(et_s);
00474 etaClus.push_back(clus_pos.eta());
00475 phiClus.push_back(clus_pos.phi());
00476 max_hit.push_back(seed_id);
00477 RecHitsCluster.push_back(RecHitsInWindow);
00478
00479
00480
00481
00482
00483 float s4s9_tmp[4];
00484 for(int i=0;i<4;i++)s4s9_tmp[i]= 0;
00485
00486 int seed_ieta = seed_id.ieta();
00487 int seed_iphi = seed_id.iphi();
00488
00489 convxtalid( seed_iphi,seed_ieta);
00490
00491 for(unsigned int j=0; j<RecHitsInWindow.size();j++){
00492 EBDetId det = (EBDetId)RecHitsInWindow[j].id();
00493
00494 int ieta = det.ieta();
00495 int iphi = det.iphi();
00496
00497 convxtalid(iphi,ieta);
00498
00499 float en = RecHitsInWindow[j].energy();
00500
00501 int dx = diff_neta_s(seed_ieta,ieta);
00502 int dy = diff_nphi_s(seed_iphi,iphi);
00503
00504 if(dx <= 0 && dy <=0) s4s9_tmp[0] += en;
00505 if(dx >= 0 && dy <=0) s4s9_tmp[1] += en;
00506 if(dx <= 0 && dy >=0) s4s9_tmp[2] += en;
00507 if(dx >= 0 && dy >=0) s4s9_tmp[3] += en;
00508 }
00509
00510 float s4s9_max = *std::max_element( s4s9_tmp,s4s9_tmp+4)/simple_energy;
00511 s4s9Clus.push_back(s4s9_max);
00512
00513 if(nClusAll<MAXCLUS){
00514 ptClusAll[nClusAll] = et_s;
00515 etaClusAll[nClusAll] = clus_pos.eta();
00516 phiClusAll[nClusAll] = clus_pos.phi();
00517 s4s9ClusAll[nClusAll] = s4s9_max;
00518 nClusAll++;
00519 Nalcapi0clusters++;
00520 }
00521
00522 nClus++;
00523 }
00524 }
00525
00526
00527
00528 if (pi0eerechits.isValid()) {
00529
00530 const EcalRecHitCollection *hitCollection_e = pi0eerechits.product();
00531
00532 detIdEERecHits.clear();
00533 EERecHits.clear();
00534
00535 std::vector<EcalRecHit> seedsEndCap;
00536 seedsEndCap.clear();
00537
00538 std::vector<EEDetId> usedXtalsEndCap;
00539 usedXtalsEndCap.clear();
00540
00542 EERecHitCollection::const_iterator ite;
00543 for (ite=pi0eerechits->begin(); ite!=pi0eerechits->end(); ite++) {
00544 double energy = ite->energy();
00545 if( energy < seleXtalMinEnergy_) continue;
00546
00547 EEDetId det = ite->id();
00548 if (RegionalMatch_){
00549 EcalElectronicsId elid = TheMapping->getElectronicsId(det);
00550 int fed = elid.dccId();
00551 it = find(FEDListUsedEndcap.begin(),FEDListUsedEndcap.end(),fed);
00552 if(it == FEDListUsedEndcap.end()) continue;
00553 }
00554
00555 detIdEERecHits.push_back(det);
00556 EERecHits.push_back(*ite);
00557
00558
00559 if (energy > clusSeedThrEndCap_) seedsEndCap.push_back(*ite);
00560
00561 }
00562
00563
00564 std::auto_ptr< EERecHitCollection > pi0EERecHitCollection( new EERecHitCollection );
00565
00566 int nClusEndCap;
00567 std::vector<float> eClusEndCap;
00568 std::vector<float> etClusEndCap;
00569 std::vector<float> etaClusEndCap;
00570 std::vector<float> phiClusEndCap;
00571 std::vector< std::vector<EcalRecHit> > RecHitsClusterEndCap;
00572 std::vector<float> s4s9ClusEndCap;
00573 nClusEndCap=0;
00574
00575
00576 sort(seedsEndCap.begin(), seedsEndCap.end(), eecalRecHitLess());
00577
00578 for (std::vector<EcalRecHit>::iterator itseed=seedsEndCap.begin(); itseed!=seedsEndCap.end(); itseed++) {
00579 EEDetId seed_id = itseed->id();
00580 std::vector<EEDetId>::const_iterator usedIds;
00581
00582 bool seedAlreadyUsed=false;
00583 for(usedIds=usedXtalsEndCap.begin(); usedIds!=usedXtalsEndCap.end(); usedIds++){
00584 if(*usedIds==seed_id){
00585 seedAlreadyUsed=true;
00586 break;
00587 }
00588 }
00589
00590 if(seedAlreadyUsed)continue;
00591
00592 std::vector<DetId> clus_v = topology_ee->getWindow(seed_id,clusEtaSize_,clusPhiSize_);
00593 std::vector< std::pair<DetId, float> > clus_used;
00594
00595 std::vector<EcalRecHit> RecHitsInWindow;
00596
00597 double simple_energy = 0;
00598
00599 for (std::vector<DetId>::iterator det=clus_v.begin(); det!=clus_v.end(); det++) {
00600 EEDetId EEdet = *det;
00601
00602 bool HitAlreadyUsed=false;
00603 for(usedIds=usedXtalsEndCap.begin(); usedIds!=usedXtalsEndCap.end(); usedIds++){
00604 if(*usedIds==*det){
00605 HitAlreadyUsed=true;
00606 break;
00607 }
00608 }
00609
00610 if(HitAlreadyUsed)continue;
00611
00612
00613 if (RegionalMatch_){
00614 EcalElectronicsId elid = TheMapping->getElectronicsId(EEdet);
00615 int fed = elid.dccId();
00616 it = find(FEDListUsedEndcap.begin(),FEDListUsedEndcap.end(),fed);
00617 if(it == FEDListUsedEndcap.end()) continue;
00618 }
00619
00620 std::vector<EEDetId>::iterator itdet = find( detIdEERecHits.begin(),detIdEERecHits.end(),EEdet);
00621 if(itdet == detIdEERecHits.end()) continue;
00622
00623 int nn = int(itdet - detIdEERecHits.begin());
00624 usedXtalsEndCap.push_back(*det);
00625 RecHitsInWindow.push_back(EERecHits[nn]);
00626 clus_used.push_back( std::pair<DetId, float>(*det, 1) );
00627 simple_energy = simple_energy + EERecHits[nn].energy();
00628 }
00629
00630 math::XYZPoint clus_pos = posCalculator_.Calculate_Location(clus_used,hitCollection_e,geometry_ee,geometry_es);
00631
00632 float theta_s = 2. * atan(exp(-clus_pos.eta()));
00633 float p0x_s = simple_energy * sin(theta_s) * cos(clus_pos.phi());
00634 float p0y_s = simple_energy * sin(theta_s) * sin(clus_pos.phi());
00635 float et_s = sqrt( p0x_s*p0x_s + p0y_s*p0y_s);
00636
00637 eClusEndCap.push_back(simple_energy);
00638 etClusEndCap.push_back(et_s);
00639 etaClusEndCap.push_back(clus_pos.eta());
00640 phiClusEndCap.push_back(clus_pos.phi());
00641 RecHitsClusterEndCap.push_back(RecHitsInWindow);
00642
00643
00644
00645 float s4s9_tmp[4];
00646 for(int i=0;i<4;i++) s4s9_tmp[i]= 0;
00647
00648 int ixSeed = seed_id.ix();
00649 int iySeed = seed_id.iy();
00650 for(unsigned int j=0; j<RecHitsInWindow.size();j++){
00651 EEDetId det_this = (EEDetId)RecHitsInWindow[j].id();
00652 int dx = ixSeed - det_this.ix();
00653 int dy = iySeed - det_this.iy();
00654
00655 float en = RecHitsInWindow[j].energy();
00656
00657 if(dx <= 0 && dy <=0) s4s9_tmp[0] += en;
00658 if(dx >= 0 && dy <=0) s4s9_tmp[1] += en;
00659 if(dx <= 0 && dy >=0) s4s9_tmp[2] += en;
00660 if(dx >= 0 && dy >=0) s4s9_tmp[3] += en;
00661 }
00662
00663 float s4s9_max = *std::max_element( s4s9_tmp,s4s9_tmp+4)/simple_energy;
00664 s4s9ClusEndCap.push_back(s4s9_max);
00665
00666 if(nClusAll<MAXCLUS){
00667 ptClusAll[nClusAll] = et_s;
00668 etaClusAll[nClusAll] = clus_pos.eta();
00669 phiClusAll[nClusAll] = clus_pos.phi();
00670 s4s9ClusAll[nClusAll] = s4s9_max;
00671 nClusAll++;
00672 Nalcapi0clusters++;
00673 }
00674
00675 nClusEndCap++;
00676 }
00677 }
00678
00679
00680 if(!(pi0eerechits.isValid()) && !(pi0ebrechits.isValid()))
00681 Nalcapi0clusters = 0;
00682 }
00683
00684
00686
00688
00689
00691 std::vector<int> HLTAlCa::ListOfFEDS(double etaLow, double etaHigh, double phiLow,
00692 double phiHigh, double etamargin, double phimargin)
00693 {
00694
00695 std::vector<int> FEDs;
00696
00697 if (phimargin > Geom::pi()) phimargin = Geom::pi() ;
00698
00699 etaLow -= etamargin;
00700 etaHigh += etamargin;
00701 double phiMinus = phiLow - phimargin;
00702 double phiPlus = phiHigh + phimargin;
00703
00704 bool all = false;
00705 double dd = fabs(phiPlus-phiMinus);
00706
00707 if (dd > 2.*Geom::pi() ) all = true;
00708
00709 while (phiPlus > Geom::pi()) { phiPlus -= 2.*Geom::pi() ; }
00710 while (phiMinus < 0) { phiMinus += 2.*Geom::pi() ; }
00711 if ( phiMinus > Geom::pi()) phiMinus -= 2.*Geom::pi() ;
00712
00713 double dphi = phiPlus - phiMinus;
00714 if (dphi < 0) dphi += 2.*Geom::pi() ;
00715
00716 if (dphi > Geom::pi()) {
00717 int fed_low1 = TheMapping->GetFED(etaLow,phiMinus*180./Geom::pi());
00718 int fed_low2 = TheMapping->GetFED(etaLow,phiPlus*180./Geom::pi());
00719
00720 if (fed_low1 == fed_low2) all = true;
00721 int fed_hi1 = TheMapping->GetFED(etaHigh,phiMinus*180./Geom::pi());
00722 int fed_hi2 = TheMapping->GetFED(etaHigh,phiPlus*180./Geom::pi());
00723
00724 if (fed_hi1 == fed_hi2) all = true;
00725 }
00726
00727 if (all) {
00728
00729 phiMinus = -20 * Geom::pi() / 180.;
00730 phiPlus = -40 * Geom::pi() / 180.;
00731 }
00732
00733 const EcalEtaPhiRegion ecalregion(etaLow,etaHigh,phiMinus,phiPlus);
00734
00735 FEDs = TheMapping->GetListofFEDs(ecalregion);
00736
00737 return FEDs;
00738
00739 }
00740
00741
00743 int HLTAlCa::convertSmToFedNumbBarrel(int ieta, int smId){
00744
00745 if( ieta<=-1) return smId - 9;
00746 else return smId + 27;
00747
00748
00749 }
00750
00751
00752 void HLTAlCa::convxtalid(Int_t &nphi,Int_t &neta)
00753 {
00754
00755
00756
00757
00758 if(neta > 0) neta -= 1;
00759 if(nphi > 359) nphi=nphi-360;
00760
00761
00762 if(nphi >359 || nphi <0 || neta< -85 || neta > 84)
00763 {
00764 std::cout <<" unexpected fatal error in HLTPi0::convxtalid "<< nphi << " " << neta << " " <<std::endl;
00765
00766 }
00767 }
00768
00769
00770
00771
00772 int HLTAlCa::diff_neta_s(Int_t neta1, Int_t neta2){
00773 Int_t mdiff;
00774 mdiff=(neta1-neta2);
00775 return mdiff;
00776 }
00777
00778
00779
00780 int HLTAlCa::diff_nphi_s(Int_t nphi1,Int_t nphi2) {
00781 Int_t mdiff;
00782 if(std::abs(nphi1-nphi2) < (360-std::abs(nphi1-nphi2))) {
00783 mdiff=nphi1-nphi2;
00784 }
00785 else {
00786 mdiff=360-std::abs(nphi1-nphi2);
00787 if(nphi1>nphi2) mdiff=-mdiff;
00788 }
00789 return mdiff;
00790 }
00791