Go to the documentation of this file.00001 #include "RecoEcal/EgammaCoreTools/plugins/EcalClusterCrackCorrection.h"
00002 #include "TVector2.h"
00003 #include "TMath.h"
00004 #include "DataFormats/EgammaReco/interface/BasicClusterFwd.h"
00005 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00006 #include "RecoEcal/EgammaCoreTools/interface/PositionCalc.h"
00007 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00008 #include "Geometry/CaloTopology/interface/CaloSubdetectorTopology.h"
00009 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00010 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00011 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00012 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00013 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00014 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
00015 #include "Geometry/CaloGeometry/interface/TruncatedPyramid.h"
00016
00018 #include "FWCore/Framework/interface/EDAnalyzer.h"
00019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00020 #include "FWCore/Utilities/interface/InputTag.h"
00021
00022 #include "FWCore/Utilities/interface/EDMException.h"
00023 #include <iostream>
00024
00025
00026 float EcalClusterCrackCorrection::getValue( const reco::BasicCluster & basicCluster, const EcalRecHitCollection & recHit) const
00027 {
00028
00029 checkInit();
00030
00031
00032
00033 EcalFunctionParameters::const_iterator it;
00034 std::cout << "[[EcalClusterCrackCorrectionBaseClass::getValue]] "
00035 << params_->params().size() << " parameters:";
00036 for ( it = params_->params().begin(); it != params_->params().end(); ++it ) {
00037 std::cout << " " << *it;
00038 }
00039 std::cout << "\n";
00040 return 1;
00041 }
00042
00043
00044 float EcalClusterCrackCorrection::getValue( const reco::CaloCluster & seedbclus) const
00045 {
00046 checkInit();
00047
00048
00049 double correction_factor=1.;
00050 double fetacor=1.;
00051 double fphicor=1.;
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062 if (TMath::Abs(seedbclus.eta()) >1.4442 ) return 1.;
00063
00064 edm::ESHandle<CaloGeometry> pG;
00065 es_->get<CaloGeometryRecord>().get(pG);
00066
00067 const CaloSubdetectorGeometry* geom=pG->getSubdetectorGeometry(DetId::Ecal,EcalBarrel);
00068
00069 const math::XYZPoint position_ = seedbclus.position();
00070 double Theta = -position_.theta()+0.5*TMath::Pi();
00071 double Eta = position_.eta();
00072 double Phi = TVector2::Phi_mpi_pi(position_.phi());
00073
00074
00075
00076 const float X0 = 0.89; const float T0 = 7.4;
00077 double depth = X0 * (T0 + log(seedbclus.energy()));
00078
00079
00080
00081
00082 std::vector< std::pair<DetId, float> > crystals_vector = seedbclus.hitsAndFractions();
00083 float dphimin=999.;
00084 float detamin=999.;
00085 int ietaclosest = 0;
00086 int iphiclosest = 0;
00087 for (unsigned int icry=0; icry!=crystals_vector.size(); ++icry) {
00088 EBDetId crystal(crystals_vector[icry].first);
00089 const CaloCellGeometry* cell=geom->getGeometry(crystal);
00090 GlobalPoint center_pos = (dynamic_cast<const TruncatedPyramid*>(cell))->getPosition(depth);
00091 double EtaCentr = center_pos.eta();
00092 double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());
00093 if (TMath::Abs(EtaCentr-Eta) < detamin) {
00094 detamin = TMath::Abs(EtaCentr-Eta);
00095 ietaclosest = crystal.ieta();
00096 }
00097 if (TMath::Abs(TVector2::Phi_mpi_pi(PhiCentr-Phi)) < dphimin) {
00098 dphimin = TMath::Abs(TVector2::Phi_mpi_pi(PhiCentr-Phi));
00099 iphiclosest = crystal.iphi();
00100 }
00101 }
00102 EBDetId crystalseed(ietaclosest, iphiclosest);
00103
00104
00105 const CaloCellGeometry* cell=geom->getGeometry(crystalseed);
00106 GlobalPoint center_pos = (dynamic_cast<const TruncatedPyramid*>(cell))->getPosition(depth);
00107
00108
00109 if (ietaclosest<0) iphiclosest = 361 - iphiclosest;
00110 int iphimod20 = iphiclosest%20;
00111 if ( iphimod20 >1 ) fphicor=1.;
00112
00113 else{
00114 double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());
00115 double PhiWidth = (TMath::Pi()/180.);
00116 double PhiCry = (TVector2::Phi_mpi_pi(Phi-PhiCentr))/PhiWidth;
00117 if (PhiCry>0.5) PhiCry=0.5;
00118 if (PhiCry<-0.5) PhiCry=-0.5;
00119
00120 if (ietaclosest<0) PhiCry *= -1.;
00121
00122
00123 double g[5];
00124 int offset = iphimod20==0 ?
00125 10
00126 : 15;
00127 for (int k=0; k!=5; ++k) g[k] = (params_->params())[k+offset];
00128
00129 fphicor=0.;
00130 for (int k=0; k!=5; ++k) fphicor += g[k]*std::pow(PhiCry,k);
00131 }
00132
00133
00134 int ietamod20 = ietaclosest%20;
00135 if (TMath::Abs(ietaclosest) <25 || (TMath::Abs(ietamod20)!=5 && TMath::Abs(ietamod20)!=6) ) fetacor = 1.;
00136
00137 else
00138 {
00139 double ThetaCentr = -center_pos.theta()+0.5*TMath::Pi();
00140 double ThetaWidth = (TMath::Pi()/180.)*TMath::Cos(ThetaCentr);
00141 double EtaCry = (Theta-ThetaCentr)/ThetaWidth;
00142 if (EtaCry>0.5) EtaCry=0.5;
00143 if (EtaCry<-0.5) EtaCry=-0.5;
00144
00145 if (ietaclosest<0) EtaCry *= -1.;
00146
00147
00148 double f[5];
00149 int offset = TMath::Abs(ietamod20)==5 ?
00150 0
00151 : 5;
00152 for (int k=0; k!=5; ++k) f[k] = (params_->params())[k+offset];
00153
00154 fetacor=0.;
00155 for (int k=0; k!=5; ++k) fetacor += f[k]*std::pow(EtaCry,k);
00156 }
00157
00158 correction_factor = 1./(fetacor*fphicor);
00159
00160
00161
00162 return correction_factor;
00163 }
00164
00165
00166
00167
00168
00169
00170
00171
00172
00173
00174 float EcalClusterCrackCorrection::getValue( const reco::SuperCluster & superCluster, const int mode ) const
00175 {
00176 checkInit();
00177
00178
00179
00180
00181
00182
00183 return getValue(*(superCluster.seed()));
00184
00185
00186 }
00187
00188
00189
00190 #include "RecoEcal/EgammaCoreTools/interface/EcalClusterFunctionFactory.h"
00191 DEFINE_EDM_PLUGIN( EcalClusterFunctionFactory, EcalClusterCrackCorrection, "EcalClusterCrackCorrection");