ATLAS Offline Software
Loading...
Searching...
No Matches
LArOFCCondAlg.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3 */
4
5#include <cmath>
6
7//FIXME
8//are these threadsafe?
10#include "CaloDetDescr/CaloDetDescrElement.h"
12
13#include <Eigen/Dense>
14
15#include "LArOFCCondAlg.h"
16
18
21
22LArOFCCondAlg::LArOFCCondAlg(const std::string &name,
23 ISvcLocator *pSvcLocator)
24 : ::AthCondAlgorithm(name, pSvcLocator),
25 m_LArOnOffIdMappingObjKey("LArOnOffIdMap"),
26 m_LArShapeObjKey("LArShapeSym"),
27 m_LArNoiseObjKey("LArNoiseSym"),
28 m_LArPedestalObjKey("LArPedestal"),
29 m_LArAutoCorrTotalObjKey("LArAutoCorrTotal"),
30 m_LArOFCObjKey("LArOFC"),
31 m_Nminbias(0), m_isMC(true),
32 m_isSuperCell(false), m_firstSample(0),
33 m_useHighestGainAutoCorr(false), m_Dump(false) {
34 declareProperty("LArOnOffIdMappingObjKey", m_LArOnOffIdMappingObjKey,
35 "Key to read LArOnOffIdMapping object");
36 declareProperty("LArShapeObjKey", m_LArShapeObjKey,
37 "Key to read LArShape object");
38 declareProperty("LArNoiseObjKey", m_LArNoiseObjKey,
39 "Key to read LArNoise object");
40 declareProperty("LArPedestalObjKey", m_LArPedestalObjKey,
41 "Key to read LArPedestal object");
42 declareProperty("LArAutoCorrTotalObjKey", m_LArAutoCorrTotalObjKey,
43 "Key to read LArAutoCorrTotal object");
44 declareProperty("LArOFCObjKey", m_LArOFCObjKey,
45 "Key to write LArOFC object");
46 declareProperty("Nminbias", m_Nminbias);
47 declareProperty("isMC", m_isMC);
48 declareProperty("isSuperCell", m_isSuperCell);
50 "firstSample", m_firstSample,
51 "First sample to use for in-time event on the full pulse shape");
52 declareProperty("useHighestGainAutoCorr",m_useHighestGainAutoCorr);
53 declareProperty("DumpOFCCondAlg",m_Dump);
54 }
55
57
59 ATH_MSG_DEBUG("initialize " << name());
60
61 // ReadCondHandle initialization
62 ATH_CHECK(m_LArShapeObjKey.initialize());
63
66
67 //WriteHandle initialization
68 ATH_CHECK(m_LArOFCObjKey.initialize());
69
72
73 // Number of gains (does this have to be in initialize now b/c of AthenaMT?)
74 // Copied from LArADC2MeVCondAlg.cxx
75 if (m_isSuperCell) {
76 m_nGains = 1;
77 } else {
78 m_nGains = 3;
79 }
80
81 return StatusCode::SUCCESS;
82}
83
84StatusCode LArOFCCondAlg::execute(const EventContext& ctx) const {
85
86 // WriteHandle setup
88 // So the following should not be called usually?!
89 if (writeHandle.isValid()) {
91 "CondHandle "
92 << writeHandle.fullKey() << " is already valid.");
93 return StatusCode::SUCCESS;
94 }
95
96 // Identifier helper
97 // Copied from LArADC2MeVCondAlg.cxx
98 const LArOnlineID_Base *larOnlineID = nullptr;
99 if (m_isSuperCell) {
100 const LArOnline_SuperCellID *scidhelper;
101 ATH_CHECK(detStore()->retrieve(scidhelper, "LArOnline_SuperCellID"));
102 larOnlineID = scidhelper; // cast to base-class
103 } else { // regular cells
104 const LArOnlineID *idhelper;
105 ATH_CHECK(detStore()->retrieve(idhelper, "LArOnlineID"));
106 larOnlineID = idhelper; // cast to base-class
107 }
108
109 // Mapping helper
110 const LArOnOffIdMapping *larOnOffIdMapping = nullptr;
111 SG::ReadCondHandle<LArOnOffIdMapping> larOnOffIdMappingHdl{
113 };
114 larOnOffIdMapping = *larOnOffIdMappingHdl;
115 if (larOnOffIdMapping == nullptr) {
116 ATH_MSG_ERROR("Failed to retrieve LArOnOffIdMapping object");
117 return StatusCode::FAILURE;
118 }
119 writeHandle.addDependency(larOnOffIdMappingHdl);
120
121 // Get pointers to inputs
122 // Retrieve validity ranges and determine their intersection
124 // FIXME: should check if handle is properly created and/or check if handle is
125 // properly retrieved
126 // operator star of a ReadCondHandle returns a const pointer to type T
127 const ILArShape *larShape{ *ShapeHdl };
128 if (larShape == nullptr) {
129 ATH_MSG_ERROR("Failed to retrieve LArShape object");
130 return StatusCode::FAILURE;
131 }
132 writeHandle.addDependency(ShapeHdl);
133
135 const LArAutoCorrTotal *larAutoCorrTotal = nullptr;
136 larAutoCorrTotal= *AutoCorrTotalHdl;
137 if (larAutoCorrTotal == nullptr) {
138 ATH_MSG_ERROR("Failed to retrieve LArADC2MeV object");
139 return StatusCode::FAILURE;
140 }
141 writeHandle.addDependency(AutoCorrTotalHdl);
142
143 // Consider the determinstic objects
144 const ILArNoise *larNoise = nullptr;
145 const ILArPedestal *larPedestal = nullptr;
146
147 if (m_isMC) {
149 larNoise = *NoiseHdl;
150 if (larNoise == nullptr) {
151 ATH_MSG_ERROR("Failed to retrieve object LArNoise");
152 return StatusCode::FAILURE;
153 }
154 writeHandle.addDependency(NoiseHdl);
155 } else {
157 larPedestal = *PedestalHdl;
158 if (larPedestal == nullptr) {
159 ATH_MSG_ERROR("Failed to retrieve object LArPedestal");
160 return StatusCode::FAILURE;
161 }
162 writeHandle.addDependency(PedestalHdl);
163 }
164
165 ATH_MSG_INFO("IOV found from intersection for LArOFCCondObj object: "
166 << writeHandle.getRange());
167
168 // make output object
169 // dimensions: number of gains x number of channel IDs x elements of
170 // OFC
171 std::unique_ptr<LArOFC> larOFC =
172 std::make_unique<LArOFC>(larOnlineID, larOnOffIdMapping, m_nGains);
173
174 std::vector<float> OFCa_tmp, OFCb_tmp;
175
176 for (const HWIdentifier chid : larOnlineID->channel_range()) {
177 const IdentifierHash hid = larOnlineID->channel_Hash(chid);
178
179 //if (!(larOnOffIdMapping->isOnlineConnected(chid))) continue;
180 if (larOnOffIdMapping->isOnlineConnected(chid)) {
181 for (size_t igain = 0; igain < m_nGains; igain++) {
182
183
184 //:::::::::::::::::::::::::::::::
185 //retrieve the data
186 //:::::::::::::::::::::::::::::::
187 ILArShape::ShapeRef_t Shape = larShape->Shape(chid,igain);
188 unsigned int nsamples_shape = Shape.size();
189 ILArShape::ShapeRef_t ShapeDer = larShape->ShapeDer(chid,igain);
190 //:::::::::::::::::::::::::::::::
191
192 // get Noise autocorrelation for gain
193 int igain_autocorr = 0;
194 // to use only Autocorr fro highest gain in optimization: HEC/FCAL=> medium gain EM=>high gain
196 &&(larOnlineID->isHECchannel(chid) || larOnlineID->isFCALchannel(chid)) ) igain_autocorr=1;
197
198 const std::vector<double> AutoCorr =
199 larAutoCorrTotal->autoCorrTotal(chid,igain_autocorr,m_Nminbias);
200 //unsigned int nsamples_AC_OFC=AutoCorr.size()+1;
201 unsigned int nsamples_AC_OFC = (1+((int)(sqrt(1+8*AutoCorr.size()))))/2;
202
203 const std::vector<double>& rmsSampl =
204 larAutoCorrTotal->samplRMS(chid,igain_autocorr,m_Nminbias);
205 unsigned int nsamples2 = rmsSampl.size();
206 if (nsamples2 != nsamples_AC_OFC) {
207 ATH_MSG_WARNING( " bad size for rmsSampl " );
208 //return (m_OFCtmp); // return empty vector
209 }
210 //:::::::::::::::::::::::::::::::
211 //unsigned int iBeginOfNSamples=findTheNSamples(Shape,
212 // nsamples_AC_OFC,
213 // nsamples_shape);
214 unsigned int firstSample = m_firstSample;
215 if(larOnlineID->isHECchannel(chid) && m_firstSample == 0 && nsamples_AC_OFC==4){
216 firstSample=1;
217 }
218 unsigned int iBeginOfNSamples = firstSample;
219 if(nsamples_AC_OFC + iBeginOfNSamples > nsamples_shape)
220 iBeginOfNSamples=0;
221 //:::::::::::::::::::::::::::::::
222
223 if(m_isMC) {
224 }
225 else
226 {
227 float RMSpedestal = larPedestal->pedestalRMS(chid,igain);
228 if(RMSpedestal>= (1.0+LArElecCalib::ERRORCODE))
229 ;
230 else
231 {
232 ATH_MSG_WARNING(" PedestalRMS vector empty for "
233 <<chid<<" at gain "<<igain );
234 }
235 }
236 //:::::::::::::::::::::::::::::::
237 //protection against missing data
238 //:::::::::::::::::::::::::::::::
239 if(Shape.size()==0 || ShapeDer.size()==0 || AutoCorr.empty())
240 {
241 ATH_MSG_WARNING("Some data are missing -> OFC will be empty for "
242 <<chid<<" at gain "<<igain );
243 //return (m_OFCtmp);
244 //returns an empty vector
245 }
246 if (Shape.size()!=ShapeDer.size()) {
247 ATH_MSG_ERROR("Got invalid shape object: Size=" << Shape.size() << ", DerSize=" << ShapeDer.size());
248 return StatusCode::SUCCESS;
249 }
250 //:::::::::::::::::::::::::::::::
251 unsigned int l,c,i;
252 //:::::::::::::::::::::::::::::::
253 //calculations
254 //:::::::::::::::::::::::::::::::
255 // fill and inverrt AC matrix
256 //HepMatrix AC(nsamples_AC_OFC,nsamples_AC_OFC),
257 //ACinv(nsamples_AC_OFC,nsamples_AC_OFC);
258 Eigen::MatrixXf AC = Eigen::MatrixXf::Zero(nsamples_AC_OFC,nsamples_AC_OFC);
259 Eigen::MatrixXf ACinv = Eigen::MatrixXf::Zero(nsamples_AC_OFC,nsamples_AC_OFC);
260 for(l=0;l<nsamples_AC_OFC;++l) { // l=line c=column
261 for(c=0;c<nsamples_AC_OFC;++c) {
262 if (l==c) {
263 AC(l,c)=1.;
264 }
265 else {
266 int i1=std::min(l,c);
267 int i2=std::max(l,c);
268 int index = i1*nsamples_AC_OFC - i1*(i1+1)/2 -(i1+1) + i2;
269 AC(l,c)=AutoCorr[index];
270 }
271 AC(l,c) = AC(l,c)*rmsSampl[l]*rmsSampl[c];
272 }
273 }
274 ACinv=AC.inverse();
275 //:::::::::::::::::::::::::::::::
276
277
278 float ACinv_PS[32];//ACinv_PS
279 float ACinv_PSD[32]; //ACinv_PSD
280 //Q1 Q2 Q3 DELTA
281 float Q1=0.;
282 float Q2=0.;
283 float Q3=0.;
284
285 for(l=0;l<nsamples_AC_OFC;++l)
286 {
287 ACinv_PS[l]=0.; ACinv_PSD[l]=0.;
288 for(c=0;c<nsamples_AC_OFC;++c){
289 ACinv_PS[l]+=ACinv(l,c)*Shape[c+iBeginOfNSamples];
290 ACinv_PSD[l]+=ACinv(l,c)*ShapeDer[c+iBeginOfNSamples];
291 }
292 Q1+=Shape[l+iBeginOfNSamples]*ACinv_PS[l];
293 Q2+=ShapeDer[l+iBeginOfNSamples]*ACinv_PSD[l];
294 Q3+=ShapeDer[l+iBeginOfNSamples]*ACinv_PS[l];
295 }
296 float DELTA=Q1*Q2-Q3*Q3;
297 //:::::::::::::::::::::::::::::::
298 //OFCa
299 OFCa_tmp.resize(nsamples_AC_OFC);
300 for(i=0;i<nsamples_AC_OFC;++i)
301 OFCa_tmp[i]=(ACinv_PS[i]*Q2-ACinv_PSD[i]*Q3)/DELTA;
302 //OFCb
303 OFCb_tmp.resize(nsamples_AC_OFC);
304 for(i=0;i<nsamples_AC_OFC;++i)
305 OFCb_tmp[i]=(ACinv_PS[i]*Q3-ACinv_PSD[i]*Q1)/DELTA;
306
307 //for debugging only
308 if(m_Dump)
309 {
310 std::cout<<larOnlineID
311 ->show_to_string(larOnOffIdMapping->cnvToIdentifier(chid))
312 <<" gain="<<igain<<" Nminbias="<<m_Nminbias<<std::endl;
313 std::cout<<"Shape: ";
314 for(c=0;c<nsamples_shape;++c)
315 std::cout<<Shape[c]<<" ";
316 std::cout<<std::endl;
317 std::cout<<"ShapeDer: ";
318 for(c=0;c<nsamples_shape;++c)
319 std::cout<<ShapeDer[c]<<" ";
320 std::cout<<std::endl;
321 for(c=0;c<nsamples_AC_OFC;++c)
322 std::cout<<Shape[c+iBeginOfNSamples]<<" ";
323 std::cout<<" <- "<<iBeginOfNSamples<<std::endl;
324 for(i=0;i<nsamples_AC_OFC;++i) std::cout<<ACinv_PS[i]<<" ";
325 std::cout<<std::endl;
326 for(i=0;i<nsamples_AC_OFC;++i) std::cout<<ACinv_PSD[i]<<" ";
327 std::cout<<std::endl;
328 std::cout<<" Q1="<<Q1<<" Q2="<<Q2<<" Q3="<<Q3
329 <<" DELTA="<<DELTA<<std::endl;
330 std::cout << " OFCa: ";
331 for(i=0;i<nsamples_AC_OFC;++i)
332 std::cout<<(ACinv_PS[i]*Q2-ACinv_PSD[i]*Q3)/DELTA<<" ";
333 std::cout<<std::endl;
334 }
335 bool stat = larOFC->setOFC(hid, igain, std::make_pair(OFCa_tmp, OFCb_tmp));
336 if (!stat) {
337 msg(MSG::ERROR) << "LArOFC::setOFC fails for gain " << igain << ", hash " << hid << endmsg;
338 }
339 } // end loop over gains
340 } else { // end loop over connected channels -- now, set empty for disc. chnanels
341 for (unsigned igain=0;igain<m_nGains;++igain) {
342 std::vector<float> empty;
343 bool stat = larOFC->setOFC(hid,igain, std::make_pair(empty, empty));
344 if (!stat) {
345 msg(MSG::ERROR) << "LArOFC::setOFC fails for gain " << igain << ", hash " << hid << endmsg;
346 }
347 } // end loop over gains of disconnected channels
348 } // end loop over disconnected channels
349
350 } // end loop over all channels
351
352
353 ATH_CHECK(writeHandle.record(std::move(larOFC)));
354 ATH_MSG_INFO("Wrote LArOFC obj to CondStore");
355 return StatusCode::SUCCESS;
356}
357
359 unsigned int nsamples_AC_OFC,
360 unsigned int nsamples_shape) const
361{
362 unsigned int i_ShapeMax=0;
363 double ShapeMax=0;
364 for(unsigned int i=0;i<nsamples_shape;++i)
365 {
366 double value=Shape[i];
367 if(value>ShapeMax) { ShapeMax=value; i_ShapeMax=i; }
368 else if(value<0 && i>3) break;//after the peak
369 }
370
371 unsigned int tmp=int(nsamples_AC_OFC/2.);
372 if(tmp>i_ShapeMax) return 0;
373 else return i_ShapeMax-tmp;
374}
#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)
Definition of CaloDetDescrManager.
static const double DELTA
static const Attributes_t empty
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
const ServiceHandle< StoreGateSvc > & detStore() const
Base class for conditions algorithms.
std::string show_to_string(Identifier id, const IdContext *context=0, char sep='.') const
or provide the printout in string form
virtual float pedestalRMS(const HWIdentifier &id, int gain) const =0
access to RMS of Pedestal index by Identifier, and gain setting
LArVectorProxy ShapeRef_t
This class defines the interface for accessing Shape (Nsample variable, Dt = 25 ns fixed) @stereotype...
Definition ILArShape.h:26
virtual ShapeRef_t Shape(const HWIdentifier &id, int gain, int tbin=0, int mode=0) const =0
virtual ShapeRef_t ShapeDer(const HWIdentifier &id, int gain, int tbin=0, int mode=0) const =0
This is a "hash" representation of an Identifier.
const std::vector< double > autoCorrTotal(const IdentifierHash &hid, int gain, float Nminbias) const
const std::vector< double > samplRMS(const IdentifierHash &hid, int gain, float Nminbias) const
unsigned int m_firstSample
virtual StatusCode execute(const EventContext &ctx) const override
virtual StatusCode initialize() override
SG::WriteCondHandleKey< LArOFC > m_LArOFCObjKey
unsigned int findTheNSamples(ILArShape::ShapeRef_t Shape, unsigned int nsamples_AC_OFC, unsigned int nsamples_shape) const
SG::ReadCondHandleKey< LArOnOffIdMapping > m_LArOnOffIdMappingObjKey
bool m_useHighestGainAutoCorr
SG::ReadCondHandleKey< ILArShape > m_LArShapeObjKey
SG::ReadCondHandleKey< LArAutoCorrTotal > m_LArAutoCorrTotalObjKey
LArOFCCondAlg(const std::string &name, ISvcLocator *pSvcLocator)
SG::ReadCondHandleKey< ILArNoise > m_LArNoiseObjKey
virtual ~LArOFCCondAlg() override
SG::ReadCondHandleKey< ILArPedestal > m_LArPedestalObjKey
Identifier cnvToIdentifier(const HWIdentifier &sid) const
create an Identifier from a HWIdentifier (inline)
bool isOnlineConnected(const HWIdentifier &sid) const
Test whether a HWIdentifier is connected of not (inline)
Helper for the Liquid Argon Calorimeter cell identifiers.
id_range channel_range() const
IdentifierHash channel_Hash(HWIdentifier channelId) const
Create channel_hash from channel_Id.
bool isFCALchannel(const HWIdentifier id) const
virtual bool isHECchannel(const HWIdentifier id) const =0
void addDependency(const EventIDRange &range)
const EventIDRange & getRange() const
StatusCode record(const EventIDRange &range, T *t)
record handle, with explicit range DEPRECATED
const DataObjID & fullKey() const
Definition index.py:1