ATLAS Offline Software
Loading...
Searching...
No Matches
CaloLCWeightTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
5//-----------------------------------------------------------------------
6// File and Version Information:
7// $Id: CaloLCWeightTool.cxx,v 1.9 2009-01-27 09:09:15 gunal Exp $
8//
9// Description: see CaloLCWeightTool.h
10//
11// Environment:
12// Software developed for the ATLAS Detector at CERN LHC
13//
14// Author List:
15// Sven Menke
16//
17//-----------------------------------------------------------------------
18
19//-----------------------
20// This Class's Header --
21//-----------------------
22#include "CaloLCWeightTool.h"
28
31
33
35{
36 if(m_interpolate) {
37 msg(MSG::INFO) << "Interpolation is ON, dimensions: ";
38 for(std::vector<std::string>::iterator it=m_interpolateDimensionNames.begin(); it!=m_interpolateDimensionNames.end(); ++it){
39 msg(MSG::INFO) << " " << (*it);
40 }
41 msg() << endmsg;
42 for(std::vector<std::string>::iterator it=m_interpolateDimensionNames.begin(); it!=m_interpolateDimensionNames.end(); ++it){
45 m_interpolateDimensions.push_back(int(id));
46 }else{
47 ATH_MSG_WARNING( "Dimension '" << (*it) << "' is invalid and will be excluded." );
48 }
49 }
50 }
51
52 ATH_MSG_INFO( "Initializing " << name() );
53
54 ATH_CHECK( m_key.initialize() );
55
56 ATH_CHECK( detStore()->retrieve( m_calo_id, "CaloCell_ID") );
57
58 m_sampnames.reserve(CaloSampling::Unknown);
59 for (int iSamp=0;iSamp<CaloSampling::Unknown;iSamp++) {
61 }
62
63 ATH_CHECK(m_noiseCDOKey.initialize());
64
65 return StatusCode::SUCCESS;
66}
67
68StatusCode CaloLCWeightTool::weight(xAOD::CaloCluster *theCluster, const EventContext& ctx) const
69{
70 const CaloLocalHadCoeff* data(nullptr);
72
74 const CaloNoise* noiseCDO=*noiseHdl;
75
76 data = *rch;
77 if(data==nullptr) {
78 ATH_MSG_ERROR("Unable to access conditions object");
79 return StatusCode::FAILURE;
80 }
81 // this is not super effective, but try to put is here first, and optimize later
82 std::vector<int> isAmpMap(CaloSampling::Unknown,-1);
83 for (int iArea=0;iArea<data->getSizeAreaSet();iArea++) {
84 for (int iSamp=0;iSamp<CaloSampling::Unknown;iSamp++) {
85 if ( m_sampnames[iSamp] == data->getArea(iArea)->getTitle() ) {
87 isAmpMap[iSamp] = iArea;
88 break;
89 }
90 }
91 }
92
93 double eEM = theCluster->e();
94
95 std::vector<float> vars(5);
96
99
100 double pi0Prob = 0;
101 if ( m_useHadProbability) {
102 if (!theCluster->retrieveMoment(xAOD::CaloCluster::EM_PROBABILITY,pi0Prob)) {
103 ATH_MSG_WARNING ("Cannot find cluster moment EM_PROBABILITY");
104 }
105 }
106 else if (theCluster->recoStatus().checkStatus(CaloRecoStatus::TAGGEDEM)) {
107 pi0Prob = 1.;
108 }
109 if ( pi0Prob < 0 )
110 pi0Prob = 0;
111 if ( pi0Prob > 1 )
112 pi0Prob = 1;
113
114 if ( eEM > 0 ) {
115 // loop over all cells
116
117 xAOD::CaloCluster::cell_iterator itrCell = theCluster->cell_begin();
118 xAOD::CaloCluster::cell_iterator itrCellEnd = theCluster->cell_end();
119 for (;itrCell!=itrCellEnd; ++itrCell) {
120 CaloPrefetch::nextDDE(itrCell, itrCellEnd);
121 // check calo and sampling index for current cell
122 Identifier myId = itrCell->ID();
123 CaloCell_ID::CaloSample theSample = CaloCell_ID::CaloSample(m_calo_id->calo_sample(myId));
124 if ( isAmpMap[theSample] >= 0 ) {
125 double sigma = noiseCDO->getNoise(itrCell->ID(),itrCell->gain());
126 double energy = fabs(itrCell->e());
127 double ratio = 0;
128 if ( std::isfinite(sigma) && sigma > 0 )
129 ratio = energy/sigma;
130 if ( ratio > m_signalOverNoiseCut ) {
131 double volume = 0;
132 double density = 0;
133 const CaloDetDescrElement* myCDDE = itrCell->caloDDE();
134 if ( myCDDE ) {
135 volume = myCDDE->volume();
136 }
137 if ( volume > 0 )
138 density = energy/volume;
139 if ( density > 0 ) {
140 double abseta = fabs(itrCell->eta());
141 double log10edens = log10(density);
142 double log10cluse = log10(eEM);
143 const CaloLocalHadCoeff::LocalHadDimension *logeDim = data->getArea(isAmpMap[theSample])->getDimension(3);
144 double lemax = logeDim->getXmax()-0.5*logeDim->getDx();
145 if ( log10cluse > lemax ) log10cluse = lemax;
146
147 vars[CaloLocalHadDefs::DIMW_SIDE] = static_cast<float> ((itrCell->eta()<0?-1.0:1.0));
148 vars[CaloLocalHadDefs::DIMW_ETA] = static_cast<float> (abseta);
149 vars[CaloLocalHadDefs::DIMW_PHI] = static_cast<float> (itrCell->phi());
150 vars[CaloLocalHadDefs::DIMW_ENER] = static_cast<float> (log10cluse);
151 vars[CaloLocalHadDefs::DIMW_EDENS] = static_cast<float> (log10edens);
152
153 bool isDataOK = false;
154 double wData(0);
155
156 // accessing coefficients (non-interpolated)
157 int iBin = data->getBin(isAmpMap[theSample],vars);
158 if ( iBin >= 0 ) {
159 const CaloLocalHadCoeff::LocalHadCoeff * pData = data->getCoeff(iBin);
160 if ( pData && (*pData)[CaloLocalHadDefs::BIN_ENTRIES] > 10 ) {
161 isDataOK = true;
162 wData = (*pData)[CaloLocalHadDefs::BIN_WEIGHT];
163 }
164 if(m_interpolate) {
165 // accesing interpolated coefficients
166 bool isa = CaloLCCoeffHelper::Interpolate(data, isAmpMap[theSample], vars, parint, m_interpolateDimensions);
167 if(isa && parint[CaloLocalHadDefs::BIN_ENTRIES] > 10) {
168 isDataOK = true;
169 wData = parint[CaloLocalHadDefs::BIN_WEIGHT];
170 }
171 }
172 }
173
174 if(isDataOK) {
175 ATH_MSG_DEBUG(" weight("
176 << theSample << ", "
177 << vars[0] << ", "
178 << vars[1] << ", "
179 << vars[2] << ", "
180 << vars[3] << ", "
181 << vars[4] << ") = "
182 << wData);
183 double weight = itrCell.weight();//theCluster->getCellWeight(itrCell); // fastest!
184 weight *= (pi0Prob + (1-pi0Prob)*wData);
185 // reweight cell in cluster
186 theCluster->reweightCell(itrCell,weight);
187 }
188 } // density
189 } // noise cut
190 } // sampling
191 } // itrCell
193 } // eEM
194
195 // assume that the weighting could be called more than once. In that case
196 // eEM is the result of the previous step and the current e/eEM ratio
197 // should be multiplied with the existing HAD_WEIGHT moment
198 double new_weight (1);
199 if (!theCluster->retrieveMoment(xAOD::CaloCluster::HAD_WEIGHT,new_weight)) {
200 ATH_MSG_ERROR("Cannot retrieve HAD_WEIGHT cluster moment." );
201 return StatusCode::FAILURE;
202 }
203
204 if ( eEM > 0 || eEM < 0 ) {
205 new_weight *= theCluster->e()/eEM;
206 }
207 theCluster->insertMoment(xAOD::CaloCluster::HAD_WEIGHT,new_weight);
208
209 return StatusCode::SUCCESS;
210}
#define endmsg
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
char data[hepevt_bytes_allocation_ATLAS]
Definition HepEvt.cxx:11
CaloSampling::CaloSample CaloSample
Definition CaloCell_ID.h:53
virtual double e() const override final
get energy (data member) (synonym to method energy()
Definition CaloCell.h:333
virtual double phi() const override final
get phi (through CaloDetDescrElement)
Definition CaloCell.h:375
const CaloDetDescrElement * caloDDE() const
get pointer to CaloDetDescrElement (data member)
Definition CaloCell.h:321
virtual double eta() const override final
get eta (through CaloDetDescrElement)
Definition CaloCell.h:382
CaloGain::CaloGain gain() const
get gain (data member )
Definition CaloCell.h:361
Identifier ID() const
get ID (from cached data member) non-virtual and inline for fast access
Definition CaloCell.h:295
static void calculateKine(xAOD::CaloCluster *clu, const bool useweight=true, const bool updateLayers=true, const bool useGPUCriteria=false)
Helper class to calculate cluster kinematics based on cells.
This class groups all DetDescr information related to a CaloCell.
static bool Interpolate(const CaloLocalHadCoeff *m_data, const unsigned int n_area, std::vector< float > &x, CaloLocalHadCoeff::LocalHadCoeff &pars, const std::vector< int > &dim, double xfit=0.)
static CaloLocalHadDefs::LocalHadDimensionId getDimensionId(const std::string &dimensionName)
Gaudi::Property< bool > m_useHadProbability
look for em-probability moment and apply relative weight only
Gaudi::Property< bool > m_interpolate
interpolate correction coefficients
Gaudi::Property< double > m_signalOverNoiseCut
minimal signal/elec_noise ratio for a cell to be weighted
Gaudi::Property< bool > m_updateSamplingVars
update also sampling variables
const CaloCell_ID * m_calo_id
SG::ReadCondHandleKey< CaloNoise > m_noiseCDOKey
SG::ReadCondHandleKey< CaloLocalHadCoeff > m_key
name of the key for had cell weights
std::vector< int > m_interpolateDimensions
actual set of dimension id's to interpolate
std::vector< std::string > m_sampnames
vector of names of individual samplings
virtual StatusCode initialize() override
Gaudi::Property< std::vector< std::string > > m_interpolateDimensionNames
vector of names of dimensions in look-up tables to interpolate
virtual StatusCode weight(xAOD::CaloCluster *theCluster, const EventContext &ctx) const override
virtual ~CaloLCWeightTool()
Class defines binning for user dimension.
float getDx() const
return size of bin
float getXmax() const
return maximum value for the last bin
Hold binned correction data for local hadronic calibration procedure.
std::vector< float > LocalHadCoeff
Correction parameters for one general bin.
float getNoise(const IdentifierHash h, const int gain) const
Accessor by IdentifierHash and gain.
Definition CaloNoise.h:35
virtual bool checkStatus(const StatusIndicator &statusIndicator) const
Check status.
static const std::string & getSamplingName(const CaloSampling::CaloSample theSample)
Returns a string (name) for each CaloSampling.
bool retrieveMoment(MomentType type, double &value) const
Retrieve individual moment.
void reweightCell(cell_iterator it, const double weight)
Method to reweight a cell in the cluster (Beware: Kinematics not updated!)
CaloClusterCellLink::iterator cell_iterator
Iterator of the underlying CaloClusterCellLink (non-const version)
virtual double e() const
The total energy of the particle.
void insertMoment(MomentType type, double value)
const_cell_iterator cell_end() const
@ EM_PROBABILITY
Classification probability to be em-like.
@ HAD_WEIGHT
Hadronic weight (E_w/E_em)
const_cell_iterator cell_begin() const
Iterator of the underlying CaloClusterCellLink (const version)
CaloRecoStatus & recoStatus()
Accesssor to CaloRecoStatus (non-const)
void nextDDE(Iter iter, Iter endIter)
Prefetch next CaloDDE.
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
LocalHadDimensionId
enums to identify user dimensions id number DIMC_* - classification, DIMW_*-weighting,...
MsgStream & msg
Definition testRead.cxx:32