ATLAS Offline Software
Loading...
Searching...
No Matches
CaloLCOutOfClusterTool.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: CaloLCOutOfClusterTool.cxx,v 1.3 2009-01-27 09:09:16 gunal Exp $
8//
9// Description: see CaloLCOutOfClusterTool.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//-----------------------
30
32
34= default;
35
37{
38 ATH_MSG_INFO( "Initializing " << name() );
39
40 //--- check sampling names to use exclude in correction
41 std::vector<std::string>::iterator samplingIter = m_invalidSamplingNames.begin();
42 std::vector<std::string>::iterator samplingIterEnd = m_invalidSamplingNames.end();
43 for(; samplingIter!=samplingIterEnd; ++samplingIter) {
44 int theSampling(CaloSampling::Unknown);
45 for (unsigned int jsamp = 0;jsamp< CaloSampling::Unknown; jsamp++) {
47 theSampling = jsamp;
48 break;
49 }
50 }
51 if ( theSampling == CaloSampling::Unknown ) {
52 msg(MSG::ERROR) << "Calorimeter sampling "
53 << *samplingIter
54 << " is not a valid Calorimeter sampling name and will be ignored! "
55 << "Valid names are: ";
56 for (unsigned int jsamp = 0;jsamp< CaloSampling::Unknown; jsamp++) {
58 if ( jsamp < CaloSampling::Unknown-1)
59 msg() << ", ";
60 else
61 msg() << ".";
62 }
63 msg() << endmsg;
64 }
65 else {
66 m_invalidSamplings.insert(theSampling);
67 }
68 }
69
70 msg(MSG::INFO) << "Samplings to exclude from the out-of-cluster weighting:";
71 samplingIter = m_invalidSamplingNames.begin();
72 for(; samplingIter!=samplingIterEnd; ++samplingIter)
73 msg() << " " << *samplingIter;
74 msg() << endmsg;
75
76 if(m_interpolate) {
77 msg(MSG::INFO) << "Interpolation is ON, dimensions: ";
78 for(std::vector<std::string>::iterator it=m_interpolateDimensionNames.begin(); it!=m_interpolateDimensionNames.end(); ++it){
79 msg() << " " << (*it);
80 }
81 msg() << endmsg;
82 for(std::vector<std::string>::iterator it=m_interpolateDimensionNames.begin(); it!=m_interpolateDimensionNames.end(); ++it){
85 m_interpolateDimensions.push_back(int(id));
86 }else{
87 ATH_MSG_WARNING( "Dimension '" << (*it) << "' is invalid and will be excluded." );
88 }
89 }
90 }
91
92 // callback for conditions data
93 ATH_CHECK( m_key.initialize() );
94
95 return StatusCode::SUCCESS;
96}
97
98StatusCode CaloLCOutOfClusterTool::weight(CaloCluster *theCluster, const EventContext& ctx) const
99{
100 double eWeightedOrig = theCluster->e();
101 double eWeighted = theCluster->e();
102 // get ratio of weighted and em cluster
103 double ratio=1.;
104 std::vector<float> vars(5);
105
107
108 double pi0Prob = 0;
110 if (!theCluster->retrieveMoment(CaloCluster::EM_PROBABILITY,pi0Prob)) {
111 ATH_MSG_WARNING ("Cannot find cluster moment EM_PROBABILITY");
112 }
113 }
114 else if (theCluster->recoStatus().checkStatus(CaloRecoStatus::TAGGEDEM)) {
115 pi0Prob = 1.;
116 }
117 if ( pi0Prob < 0 )
118 pi0Prob = 0;
119 if ( pi0Prob > 1 )
120 pi0Prob = 1;
121
122 double classifyWeight = 1;
123 if ( m_useEmProbability ) {
124 classifyWeight = pi0Prob;
125 }
126 else if ( m_useHadProbability ) {
127 classifyWeight = 1-pi0Prob;
128 }
129 double eEM = theCluster->rawE();
130 // basic energy can be 0 if the locking tool option FixBasicEnergy is False
131 // in that case simply use the weighted energy ...
132 if ( eEM == 0. ) {
133 ATH_MSG_DEBUG("Basic Energy is 0. Use weighted energy instead");
134 eEM = eWeighted;
135 }
136
137
138
139 // subtract the samplings to ignore from eEM and eWeighted, assuming that
140 // they would not have received any cell weights in prior calibration
141 // steps
142 std::set<int>::const_iterator ivSamplingIter = m_invalidSamplings.begin();
143 std::set<int>::const_iterator ivSamplingIterEnd = m_invalidSamplings.end();
144 for(; ivSamplingIter!=ivSamplingIterEnd; ++ivSamplingIter) {
145 eEM -= theCluster->eSample((CaloCluster::CaloSample)(*ivSamplingIter));
146 eWeighted -= theCluster->eSample((CaloCluster::CaloSample)(*ivSamplingIter));
147 }
148
149 if ( eEM > 0 && eWeighted > 0 ) {
150 // ratio is em energy over weighted energy without the samplings
151 // to ignore. This means also the OOC weight tables have to be made
152 // relative to the em energy wihtout the samplings to ignore ...
153 ratio = eEM / eWeighted;
154
155 if ( ratio < 0.3 || ratio > 3 ) {
156 ATH_MSG_DEBUG("The ratio eEM/eWeighted = "
157 << ratio << " is out of normal range [0.3,3]"
158 << " - this mean we have mainly noise ... using 1 instead");
159 ratio = 1.0;
160 }
161
162 double log10cluse = log10(eEM);
163
165 const CaloLocalHadCoeff* data = *h;
166 if (!data) {
167 ATH_MSG_ERROR("Unable to access conditions object");
168 return StatusCode(CaloRecoStatus::TAGGEDUNKNOWN);
169 }
170
171 const CaloLocalHadCoeff::LocalHadDimension *logeDim = data->getArea(0)->getDimension(2);
172 double lemax = logeDim->getXmax()-0.5*logeDim->getDx();
173 double lemin = logeDim->getXmin()+0.5*logeDim->getDx();
174 if ( log10cluse > lemax )
175 log10cluse = lemax;
176 if ( log10cluse < lemin )
177 log10cluse = lemin;
178
179 double center_lambda,isolation;
180 if ( theCluster->retrieveMoment(CaloCluster::ISOLATION,isolation)
181 && theCluster->retrieveMoment(CaloCluster::CENTER_LAMBDA, center_lambda)) {
182 if ( isolation > 0 ) {
183 if ( center_lambda > 0 ) {
184 const double abseta = fabs(theCluster->eta());
185 const double log10lambda = log10(center_lambda);
186
187 vars[CaloLocalHadDefs::DIMO_SIDE] = static_cast<float> ((theCluster->eta()<0?-1.0:1.0));
188 vars[CaloLocalHadDefs::DIMO_PHI] = static_cast<float> (theCluster->phi());
189 vars[CaloLocalHadDefs::DIMO_ENER] = static_cast<float> (log10cluse);
190 vars[CaloLocalHadDefs::DIMO_ETA] = static_cast<float> (abseta);
191 vars[CaloLocalHadDefs::DIMO_LAMBDA] = static_cast<float> (log10lambda);
192
193 bool isDataOK = false;
194 double oocData(0);
195
196 // accessing coefficients (non-interpolated)
197 int iBin = data->getBin(0,vars);
198 if ( iBin >= 0 ) {
199 const CaloLocalHadCoeff::LocalHadCoeff * pData = data->getCoeff(iBin);
200 if ( pData && (*pData)[CaloLocalHadDefs::BIN_ENTRIES] > 0 ) {
201 isDataOK = true;
202 oocData = (*pData)[CaloLocalHadDefs::BIN_WEIGHT];
203 }
204 if(m_interpolate) {
205 // accesing interpolated coefficients
208 if(isa && parint[CaloLocalHadDefs::BIN_ENTRIES] > 0) {
209 isDataOK = true;
210 oocData = parint[CaloLocalHadDefs::BIN_WEIGHT];
211 }
212 }
213 }
214
215 ATH_MSG_DEBUG("interpolation_ON=" << m_interpolate
216 << "isDataOK="<<isDataOK
217 << "side = " << vars[0]
218 << ", phi = " << vars[1]
219 << ", log10cluse = " << log10cluse
220 << ", eta = " << abseta
221 << ", log10lambda = " << log10lambda
222 << ", ratio = " << ratio
223 << ", isol = " << isolation
224 << ", oocData = " << oocData);
225
226 if(isDataOK){
227 const double oocWeight = 1.+classifyWeight*isolation*ratio*oocData;
228 // loop over all cells
229 CaloCluster::cell_iterator itrCell = theCluster->cell_begin();
230 CaloCluster::cell_iterator itrCellEnd = theCluster->cell_end();
231 for (;itrCell!=itrCellEnd; ++itrCell) {
232 CaloPrefetch::nextDDE(itrCell, itrCellEnd);
233 const CaloDetDescrElement* cDDE = itrCell->caloDDE();
234 if ( cDDE && !m_invalidSamplings.contains(cDDE->getSampling())) { //Fixme ... could use a bit-pattern
235 double weight = itrCell.weight();//theCluster->getCellWeight(itrCell); // fastest!
236 weight *= oocWeight;
237 // reweight cell in cluster
238 theCluster->reweightCell(itrCell,weight);
239 } else if (cDDE) {
240 ATH_MSG_DEBUG("Exclude cell with sampling = " << cDDE->getSampling());
241 }
242 }//cell-loop
244 }
245 } // log10lambda
246 } // isolation
247 } // moments
248 }
249
250 // assume that the weighting could be called more than once. In that
251 // case eWeighted is the result of the previous step and the current
252 // e/eWeighted ratio should be multiplied with the existing
253 // OOC_WEIGHT moment
254 if ( eWeightedOrig > 0 || eWeightedOrig < 0 ) {
255 double old_weight(1);
256 if (!theCluster->retrieveMoment(CaloCluster::OOC_WEIGHT,old_weight)) {
257 ATH_MSG_ERROR("Cannot retrieve OOC_WEIGHT cluster moment." );
258 return StatusCode::FAILURE;
259 }
260 const double new_weight = old_weight*theCluster->e()/eWeightedOrig;
261 theCluster->insertMoment(CaloCluster::OOC_WEIGHT,new_weight);
262 }
263 return StatusCode::SUCCESS;
264}
#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
Header file for AthHistogramAlgorithm.
const CaloDetDescrElement * caloDDE() const
get pointer to CaloDetDescrElement (data member)
Definition CaloCell.h:321
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.
Principal data class for CaloCell clusters.
double eSample(sampling_type sampling) const
Retrieve energy in a given sampling.
bool retrieveMoment(const moment_type &momType, moment_value &momValue, bool useLink=true) const
Retrieve individual moment.
virtual double e() const
Retrieve energy independent of signal state.
void insertMoment(const moment_type &momType, const moment_value &momValue, bool useLink=true)
Set individual moment.
virtual double eta() const
Retrieve eta independent of signal state.
virtual double phi() const
Retrieve phi independent of signal state.
cell_iterator cell_end() const
Retrieve a STL-type end() iterator for the cell store.
void reweightCell(const CaloCell *pCell, double weight)
Reweight a cell with kinematic update.
cell_iterator cell_begin() const
Retrieve a STL-type begin() iterator for the cell store.
This class groups all DetDescr information related to a CaloCell.
CaloCell_ID::CaloSample getSampling() const
cell sampling
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< std::vector< std::string > > m_invalidSamplingNames
vector of names of the calorimeter samplings not to use when applying the out-of-cluster weights.
Gaudi::Property< bool > m_interpolate
interpolate correction coefficients
Gaudi::Property< bool > m_updateSamplingVars
update also sampling variables
std::vector< int > m_interpolateDimensions
actual set of dimension id's to interpolate
virtual StatusCode weight(xAOD::CaloCluster *theCluster, const EventContext &ctx) const override
Gaudi::Property< std::vector< std::string > > m_interpolateDimensionNames
vector of names of dimensions in look-up tables to interpolate
virtual StatusCode initialize() override
Gaudi::Property< bool > m_useHadProbability
look for em-probability moment and apply relative weight only
SG::ReadCondHandleKey< CaloLocalHadCoeff > m_key
name of the key for out-of-cluster weights
virtual ~CaloLCOutOfClusterTool()
Gaudi::Property< bool > m_useEmProbability
look for em-probability moment and apply relative weight only
std::set< int > m_invalidSamplings
actual set of samplings to be ignored for out-of-cluster weights
Class defines binning for user dimension.
float getDx() const
return size of bin
float getXmax() const
return maximum value for the last bin
float getXmin() const
return minimum value for the first bin
Hold binned correction data for local hadronic calibration procedure.
std::vector< float > LocalHadCoeff
Correction parameters for one general bin.
static const std::string & getSamplingName(const CaloSampling::CaloSample theSample)
Returns a string (name) for each CaloSampling.
CaloClusterCellLink::iterator cell_iterator
Iterator of the underlying CaloClusterCellLink (non-const version)
@ OOC_WEIGHT
Out-of-cluster weight (E_ooc/E_w)
@ EM_PROBABILITY
Classification probability to be em-like.
@ CENTER_LAMBDA
Shower depth at Cluster Centroid.
@ ISOLATION
Energy weighted fraction of non-clustered perimeter cells.
CaloSampling::CaloSample CaloSample
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