CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
EcalBasicClusterLocalContCorrection Class Reference

#include <EcalBasicClusterLocalContCorrection.h>

Public Member Functions

void checkInit () const
 
 EcalBasicClusterLocalContCorrection (edm::ConsumesCollector &&cc)
 
void init (const edm::EventSetup &es)
 
float operator() (const reco::BasicCluster &, const EcalRecHitCollection &) const
 

Private Member Functions

int getEcalModule (DetId id) const
 

Private Attributes

CaloGeometry const * caloGeometry_
 
const edm::ESGetToken< CaloGeometry, CaloGeometryRecordcaloGeometryToken_
 
EcalClusterLocalContCorrParameters const * params_
 
const edm::ESGetToken< EcalClusterLocalContCorrParameters, EcalClusterLocalContCorrParametersRcdparamsToken_
 

Detailed Description

Function to correct em object energy for energy not contained in a 5x5 crystal area in the calorimeter

$Id: EcalBasicClusterLocalContCorrection.h $Date: $Revision:

Author
Federico Ferri, CEA Saclay, November 2008

Definition at line 22 of file EcalBasicClusterLocalContCorrection.h.

Constructor & Destructor Documentation

◆ EcalBasicClusterLocalContCorrection()

EcalBasicClusterLocalContCorrection::EcalBasicClusterLocalContCorrection ( edm::ConsumesCollector &&  cc)

Definition at line 8 of file EcalBasicClusterLocalContCorrection.cc.

References gpuPixelDoublets::cc.

9  : paramsToken_{cc.esConsumes()}, caloGeometryToken_{cc.esConsumes()} {}
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
const edm::ESGetToken< EcalClusterLocalContCorrParameters, EcalClusterLocalContCorrParametersRcd > paramsToken_
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeometryToken_

Member Function Documentation

◆ checkInit()

void EcalBasicClusterLocalContCorrection::checkInit ( ) const

Definition at line 16 of file EcalBasicClusterLocalContCorrection.cc.

References Exception, and params_.

16  {
17  if (!params_) {
18  // non-initialized function parameters: throw exception
19  throw cms::Exception("EcalBasicClusterLocalContCorrection::checkInit()")
20  << "Trying to access an uninitialized crack correction function.\n"
21  "Please call `init( edm::EventSetup &)' before any use of the function.\n";
22  }
23 }
EcalClusterLocalContCorrParameters const * params_

◆ getEcalModule()

int EcalBasicClusterLocalContCorrection::getEcalModule ( DetId  id) const
private

Definition at line 157 of file EcalBasicClusterLocalContCorrection.cc.

References hcalRecHitTable_cff::ieta, and mod().

157  {
158  int mod = 0;
159  int ieta = (EBDetId(id)).ieta();
160 
161  if (fabs(ieta) <= 25)
162  mod = 0;
163  if (fabs(ieta) > 25 && fabs(ieta) <= 45)
164  mod = 1;
165  if (fabs(ieta) > 45 && fabs(ieta) <= 65)
166  mod = 2;
167  if (fabs(ieta) > 65 && fabs(ieta) <= 85)
168  mod = 3;
169 
170  return (mod);
171 }
T mod(const T &a, const T &b)
Definition: ecalDccMap.h:4

◆ init()

void EcalBasicClusterLocalContCorrection::init ( const edm::EventSetup es)

Definition at line 11 of file EcalBasicClusterLocalContCorrection.cc.

References caloGeometry_, caloGeometryToken_, edm::EventSetup::getData(), params_, and paramsToken_.

11  {
14 }
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const edm::ESGetToken< EcalClusterLocalContCorrParameters, EcalClusterLocalContCorrParametersRcd > paramsToken_
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeometryToken_
EcalClusterLocalContCorrParameters const * params_

◆ operator()()

float EcalBasicClusterLocalContCorrection::operator() ( const reco::BasicCluster basicCluster,
const EcalRecHitCollection recHit 
) const

Definition at line 28 of file EcalBasicClusterLocalContCorrection.cc.

References funct::abs(), ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), funct::cos(), hcalRecHitTable_cff::depth, DetId::Ecal, EcalBarrel, reco::CaloCluster::energy(), PV3DBase< T, PVType, FrameType >::eta(), dqmdumpme::first, relativeConstraints::geom, reco::CaloCluster::hitsAndFractions(), EBDetId::ieta(), EBDetId::iphi(), dqm-mbProfile::log, M_PI, PV3DBase< T, PVType, FrameType >::phi(), VtxSmearedParameters_cfi::Phi, reco::CaloCluster::position(), reco::CaloCluster::seed(), PV3DBase< T, PVType, FrameType >::theta(), and ecalPiZeroTask_cfi::X0.

29  {
30  checkInit();
31 
32  // number of parameters needed by this parametrization
33  constexpr size_t nparams = 24;
34 
35  //correction factor to be returned, and to be calculated in this present function:
36  double correction_factor = 1.;
37  double fetacor = 1.; //eta dependent part of the correction factor
38  double fphicor = 1.; //phi dependent part of the correction factor
39 
40  //--------------if barrel calculate local position wrt xtal center -------------------
43 
44  const math::XYZPoint &position_ = basicCluster.position();
45  double Theta = -position_.theta() + 0.5 * M_PI;
46  double Eta = position_.eta();
47  double Phi = TVector2::Phi_mpi_pi(position_.phi());
48 
49  //Calculate expected depth of the maximum shower from energy (like in PositionCalc::Calculate_Location()):
50  // The parameters X0 and T0 are hardcoded here because these values were used to calculate the corrections:
51  constexpr float X0 = 0.89;
52  constexpr float T0 = 7.4;
53  double depth = X0 * (T0 + log(basicCluster.energy()));
54 
55  //search which crystal is closest to the cluster position and call it crystalseed:
56  //std::vector<DetId> crystals_vector = *scRef.getHitsByDetId(); //deprecated
57  std::vector<std::pair<DetId, float> > crystals_vector = basicCluster.hitsAndFractions();
58  float dphimin = 999.;
59  float detamin = 999.;
60  int ietaclosest = 0;
61  int iphiclosest = 0;
62 
63  for (unsigned int icry = 0; icry != crystals_vector.size(); ++icry) {
64  EBDetId crystal(crystals_vector[icry].first);
65  auto cell = geom->getGeometry(crystal); // problema qui
66  GlobalPoint center_pos = cell->getPosition(depth);
67  double EtaCentr = center_pos.eta();
68  double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());
69  if (std::abs(EtaCentr - Eta) < detamin) {
70  detamin = std::abs(EtaCentr - Eta);
71  ietaclosest = crystal.ieta();
72  }
73  if (std::abs(TVector2::Phi_mpi_pi(PhiCentr - Phi)) < dphimin) {
74  dphimin = std::abs(TVector2::Phi_mpi_pi(PhiCentr - Phi));
75  iphiclosest = crystal.iphi();
76  }
77  }
78 
79  EBDetId crystalseed(ietaclosest, iphiclosest);
80 
81  // Get center cell position from shower depth
82  auto cell = geom->getGeometry(crystalseed);
83  GlobalPoint center_pos = cell->getPosition(depth);
84 
85  //PHI
86  double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());
87  double PhiWidth = (M_PI / 180.);
88  double PhiCry = (TVector2::Phi_mpi_pi(Phi - PhiCentr)) / PhiWidth;
89  if (PhiCry > 0.5)
90  PhiCry = 0.5;
91  if (PhiCry < -0.5)
92  PhiCry = -0.5;
93  //flip to take into account ECAL barrel symmetries:
94  if (ietaclosest < 0)
95  PhiCry *= -1.;
96 
97  //ETA
98  double ThetaCentr = -center_pos.theta() + 0.5 * M_PI;
99  double ThetaWidth = (M_PI / 180.) * std::cos(ThetaCentr);
100  double EtaCry = (Theta - ThetaCentr) / ThetaWidth;
101  if (EtaCry > 0.5)
102  EtaCry = 0.5;
103  if (EtaCry < -0.5)
104  EtaCry = -0.5;
105  //flip to take into account ECAL barrel symmetries:
106  if (ietaclosest < 0)
107  EtaCry *= -1.;
108 
109  //-------------- end calculate local position -------------
110 
111  size_t payloadsize = params_->params().size();
112 
113  if (payloadsize < nparams)
114  edm::LogError("Invalid Payload") << "Parametrization requires " << nparams << " parameters but only " << payloadsize
115  << " are found in DB. Perhaps incompatible Global Tag" << std::endl;
116 
117  if (payloadsize > nparams)
118  edm::LogWarning("Size mismatch ") << "Parametrization requires " << nparams << " parameters but " << payloadsize
119  << " are found in DB. Perhaps incompatible Global Tag" << std::endl;
120 
121  std::pair<double, double> localPosition(EtaCry, PhiCry);
122 
123  //--- local cluster coordinates
124  float localEta = localPosition.first;
125  float localPhi = localPosition.second;
126 
127  //--- ecal module
128  int imod = getEcalModule(basicCluster.seed());
129 
130  //-- corrections parameters
131  float pe[3], pp[3];
132  pe[0] = (params_->params())[0 + imod * 3];
133  pe[1] = (params_->params())[1 + imod * 3];
134  pe[2] = (params_->params())[2 + imod * 3];
135  pp[0] = (params_->params())[12 + imod * 3];
136  pp[1] = (params_->params())[13 + imod * 3];
137  pp[2] = (params_->params())[14 + imod * 3];
138 
139  //--- correction vs local eta
140  fetacor = pe[0] + pe[1] * localEta + pe[2] * localEta * localEta;
141 
142  //--- correction vs local phi
143  fphicor = pp[0] + pp[1] * localPhi + pp[2] * localPhi * localPhi;
144 
145  //if the seed crystal is neighbourgh of a supermodule border, don't apply the phi dependent containment corrections, but use the larger crack corrections instead.
146  int iphimod20 = std::abs(iphiclosest % 20);
147  if (iphimod20 <= 1)
148  fphicor = 1.;
149 
150  correction_factor = (1. / fetacor) * (1. / fphicor);
151 
152  //return the correction factor. Use it to multiply the cluster energy.
153  return correction_factor;
154 }
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:153
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:209
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:218
T eta() const
Definition: PV3DBase.h:73
Log< level::Error, false > LogError
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
EcalFunctionParameters & params()
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define M_PI
double energy() const
cluster energy
Definition: CaloCluster.h:148
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
EcalClusterLocalContCorrParameters const * params_
Log< level::Warning, false > LogWarning
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72

Member Data Documentation

◆ caloGeometry_

CaloGeometry const* EcalBasicClusterLocalContCorrection::caloGeometry_
private

Definition at line 42 of file EcalBasicClusterLocalContCorrection.h.

Referenced by init().

◆ caloGeometryToken_

const edm::ESGetToken<CaloGeometry, CaloGeometryRecord> EcalBasicClusterLocalContCorrection::caloGeometryToken_
private

Definition at line 39 of file EcalBasicClusterLocalContCorrection.h.

Referenced by init().

◆ params_

EcalClusterLocalContCorrParameters const* EcalBasicClusterLocalContCorrection::params_
private

Definition at line 41 of file EcalBasicClusterLocalContCorrection.h.

Referenced by checkInit(), and init().

◆ paramsToken_

const edm::ESGetToken<EcalClusterLocalContCorrParameters, EcalClusterLocalContCorrParametersRcd> EcalBasicClusterLocalContCorrection::paramsToken_
private

Definition at line 38 of file EcalBasicClusterLocalContCorrection.h.

Referenced by init().