ATLAS Offline Software
Loading...
Searching...
No Matches
InDet::ClusterMakerTool Class Reference

#include <ClusterMakerTool.h>

Inheritance diagram for InDet::ClusterMakerTool:
Collaboration diagram for InDet::ClusterMakerTool:

Public Member Functions

 ClusterMakerTool (const std::string &type, const std::string &name, const IInterface *parent)
 ~ClusterMakerTool ()=default
StatusCode initialize ()
PixelCluster pixelCluster (const Identifier &clusterID, const Amg::Vector2D &localPos, std::vector< Identifier > &&rdoList, const int lvl1a, std::vector< int > &&totList, const SiWidth &width, const InDetDD::SiDetectorElement *element, bool ganged, int errorStrategy, const PixelID &pixelID, bool split, double splitProb1, double splitProb2, const PixelChargeCalibCondData *calibData, const PixelCalib::PixelOfflineCalibData *offlineCalibData, const EventContext &ctx) const
xAOD::PixelClusterxAODpixelCluster (xAOD::PixelCluster &cluster, const Amg::Vector2D &localPos, const std::vector< Identifier > &rdoList, const int lvl1a, const std::vector< int > &totList, const SiWidth &width, const InDetDD::SiDetectorElement *element, bool ganged, int errorStrategy, const PixelID &pixelID, bool split, double splitProb1, double splitProb2, const PixelChargeCalibCondData *calibData, const PixelCalib::PixelOfflineCalibData *offlineCalibData, const EventContext &ctx) const
SCT_Cluster sctCluster (const Identifier &clusterID, const Amg::Vector2D &localPos, std::vector< Identifier > &&rdoList, const SiWidth &width, const InDetDD::SiDetectorElement *element, int errorStrategy) const
ServiceHandle< StoreGateSvc > & evtStore ()
 The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
const ServiceHandle< StoreGateSvc > & detStore () const
 The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
virtual StatusCode sysInitialize () override
 Perform system initialization for an algorithm.
virtual StatusCode sysStart () override
 Handle START transition.
virtual std::vector< Gaudi::DataHandle * > inputHandles () const override
 Return this algorithm's input handles.
virtual std::vector< Gaudi::DataHandle * > outputHandles () const override
 Return this algorithm's output handles.
Gaudi::Details::PropertyBase & declareProperty (Gaudi::Property< T, V, H > &t)
void updateVHKA (Gaudi::Details::PropertyBase &)
MsgStream & msg () const
bool msgLvl (const MSG::Level lvl) const

Static Public Member Functions

static const InterfaceID & interfaceID ()

Protected Member Functions

void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce (T &h)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

template<typename ClusterType, typename IdentifierVec, typename ToTList>
ClusterType makePixelCluster (const Identifier &clusterID, const Amg::Vector2D &localPos, IdentifierVec &&rdoList, const int lvl1a, ToTList &&totList, const SiWidth &width, const InDetDD::SiDetectorElement *element, bool ganged, int errorStrategy, const PixelID &pixelID, bool split, double splitProb1, double splitProb2, const PixelChargeCalibCondData *calibData, const PixelCalib::PixelOfflineCalibData *offlineCalibData, const EventContext &ctx, xAOD::PixelCluster *cluster=nullptr) const
double getPixelCTBPhiError (int layer, int phi, int PhiClusterSize) const
Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

ToolHandle< ISiLorentzAngleToolm_pixelLorentzAngleTool {this, "PixelLorentzAngleTool", "SiLorentzAngleTool/PixelLorentzAngleTool", "Tool to retreive Lorentz angle of Pixel"}
ToolHandle< ISiLorentzAngleToolm_sctLorentzAngleTool {this, "SCTLorentzAngleTool", "SiLorentzAngleTool/SCTLorentzAngleTool", "Tool to retreive Lorentz angle of SCT"}
bool m_forceErrorStrategy1B {false}
StoreGateSvc_t m_evtStore
 Pointer to StoreGate (event store by default)
StoreGateSvc_t m_detStore
 Pointer to StoreGate (detector store by default)
std::vector< SG::VarHandleKeyArray * > m_vhka
bool m_varHandleArraysDeclared

Detailed Description

Definition at line 55 of file ClusterMakerTool.h.

Member Typedef Documentation

◆ StoreGateSvc_t

typedef ServiceHandle<StoreGateSvc> AthCommonDataStore< AthCommonMsg< AlgTool > >::StoreGateSvc_t
privateinherited

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ ClusterMakerTool()

InDet::ClusterMakerTool::ClusterMakerTool ( const std::string & type,
const std::string & name,
const IInterface * parent )

Definition at line 139 of file ClusterMakerTool.cxx.

141 :
142 AthAlgTool(t,n,p)
143{
144 declareInterface<ClusterMakerTool>(this);
145}
AthAlgTool()
Default constructor:

◆ ~ClusterMakerTool()

InDet::ClusterMakerTool::~ClusterMakerTool ( )
default

Member Function Documentation

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareGaudiProperty ( Gaudi::Property< T, V, H > & hndl,
const SG::VarHandleKeyType &  )
inlineprivateinherited

specialization for handling Gaudi::Property<SG::VarHandleKey>

Definition at line 156 of file AthCommonDataStore.h.

158 {
160 hndl.value(),
161 hndl.documentation());
162
163 }
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)

◆ declareProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareProperty ( Gaudi::Property< T, V, H > & t)
inlineinherited

Definition at line 145 of file AthCommonDataStore.h.

145 {
146 typedef typename SG::HandleClassifier<T>::type htype;
148 }
Gaudi::Details::PropertyBase & declareGaudiProperty(Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
specialization for handling Gaudi::Property<SG::VarHandleKey>

◆ detStore()

const ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::detStore ( ) const
inlineinherited

The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 95 of file AthCommonDataStore.h.

◆ evtStore()

ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::evtStore ( )
inlineinherited

The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 85 of file AthCommonDataStore.h.

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::extraDeps_update_handler ( Gaudi::Details::PropertyBase & ExtraDeps)
protectedinherited

Add StoreName to extra input/output deps as needed.

use the logic of the VarHandleKey to parse the DataObjID keys supplied via the ExtraInputs and ExtraOuputs Properties to add the StoreName if it's not explicitly given

◆ getPixelCTBPhiError()

double InDet::ClusterMakerTool::getPixelCTBPhiError ( int layer,
int phi,
int PhiClusterSize ) const
private

Definition at line 587 of file ClusterMakerTool.cxx.

588 {
589
590 double sigmaL0Phi1[3] = { 8.2*micrometer, 9.7*micrometer, 14.6*micrometer};
591 double sigmaL1Phi1[3] = {14.6*micrometer, 9.3*micrometer, 14.6*micrometer};
592 double sigmaL2Phi1[3] = {14.6*micrometer, 8.6*micrometer, 14.6*micrometer};
593 double sigmaL0Phi0[3] = {14.6*micrometer, 13.4*micrometer, 13.0*micrometer};
594 double sigmaL1Phi0[3] = {14.6*micrometer, 8.5*micrometer, 11.0*micrometer};
595 double sigmaL2Phi0[3] = {14.6*micrometer, 11.6*micrometer, 9.3*micrometer};
596
597 if(phiClusterSize > 3) return 14.6*micrometer;
598
599 if(layer == 0 && phi == 0) return sigmaL0Phi0[phiClusterSize-1];
600 if(layer == 1 && phi == 0) return sigmaL1Phi0[phiClusterSize-1];
601 if(layer == 2 && phi == 0) return sigmaL2Phi0[phiClusterSize-1];
602 if(layer == 0 && phi == 1) return sigmaL0Phi1[phiClusterSize-1];
603 if(layer == 1 && phi == 1) return sigmaL1Phi1[phiClusterSize-1];
604 if(layer == 2 && phi == 1) return sigmaL2Phi1[phiClusterSize-1];
605
606 // shouldn't really happen...
607 ATH_MSG_WARNING("Unexpected layer and phi numbers: layer = "
608 << layer << " and phi = " << phi);
609 return 14.6*micrometer;
610
611}
Scalar phi() const
phi method
#define ATH_MSG_WARNING(x)

◆ initialize()

StatusCode InDet::ClusterMakerTool::initialize ( )

Definition at line 149 of file ClusterMakerTool.cxx.

149 {
150 // Code entered here will be executed once at program start.
151
152 ATH_MSG_DEBUG ( name() << " initialize()" );
153
154 if (not m_pixelLorentzAngleTool.empty()) {
156 } else {
157 m_pixelLorentzAngleTool.disable();
158 }
159 if (not m_sctLorentzAngleTool.empty()) {
161 } else {
162 m_sctLorentzAngleTool.disable();
163 }
164
165
166 return StatusCode::SUCCESS;
167}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_DEBUG(x)
ToolHandle< ISiLorentzAngleTool > m_sctLorentzAngleTool
ToolHandle< ISiLorentzAngleTool > m_pixelLorentzAngleTool

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::inputHandles ( ) const
overridevirtualinherited

Return this algorithm's input handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ interfaceID()

const InterfaceID & InDet::ClusterMakerTool::interfaceID ( )
inlinestatic

Definition at line 67 of file ClusterMakerTool.h.

67{ return IID_ClusterMakerTool; };
static const InterfaceID IID_ClusterMakerTool("InDet::ClusterMakerTool", 1, 0)

◆ makePixelCluster()

template<typename ClusterType, typename IdentifierVec, typename ToTList>
ClusterType InDet::ClusterMakerTool::makePixelCluster ( const Identifier & clusterID,
const Amg::Vector2D & localPos,
IdentifierVec && rdoList,
const int lvl1a,
ToTList && totList,
const SiWidth & width,
const InDetDD::SiDetectorElement * element,
bool ganged,
int errorStrategy,
const PixelID & pixelID,
bool split,
double splitProb1,
double splitProb2,
const PixelChargeCalibCondData * calibData,
const PixelCalib::PixelOfflineCalibData * offlineCalibData,
const EventContext & ctx,
xAOD::PixelCluster * cluster = nullptr ) const
private

Definition at line 191 of file ClusterMakerTool.cxx.

208 {
209
210 ATH_MSG_VERBOSE("ClusterMakerTool called, number ");
211 if ( errorStrategy==2 && m_forceErrorStrategy1B ) errorStrategy=1;
212
213 // Fill vector of charges and compute charge balance
214 const InDetDD::PixelModuleDesign* design = (dynamic_cast<const InDetDD::PixelModuleDesign*>(&element->design()));
215 if (not design){
216 throw std::runtime_error( "Dynamic cast failed for design in ClusterMakerTool.cxx");
217 }
218 int rowMin = design->rows();
219 int rowMax = 0;
220 int colMin = design->columns();
221 int colMax = 0;
222 float qRowMin = 0; float qRowMax = 0;
223 float qColMin = 0; float qColMax = 0;
224 std::vector<float> chargeList;
225 int nRDO=rdoList.size();
226 if (calibData) {
227 chargeList.reserve(nRDO);
228 IdentifierHash moduleHash = element->identifyHash(); // wafer hash
229 for (int i=0; i<nRDO; i++) {
230 Identifier pixid=rdoList[i];
231 int ToT=totList[i];
232
233 float charge = ToT;
234 assert( element->identifyHash() == pixelID.wafer_hash(pixelID.wafer_id(pixid)));
235 std::array<InDetDD::PixelDiodeTree::CellIndexType,2> diode_idx
237 pixelID.eta_index(pixid));
238 InDetDD::PixelDiodeTree::DiodeProxy si_param ( design->diodeProxyFromIdx(diode_idx));
239 std::uint32_t feValue = design->getFE(si_param);
240 auto diode_type = design->getDiodeType(si_param);
241 if ( design->getReadoutTechnology() == InDetDD::PixelReadoutTechnology::FEI3
242 && design->numberOfConnectedCells( design->readoutIdOfCell(InDetDD::SiCellId(diode_idx[0],diode_idx[1])))>1) {
244 }
245
246 charge = calibData->getCharge(diode_type, moduleHash, feValue, ToT);
247 if (design->getReadoutTechnology() != InDetDD::PixelReadoutTechnology::RD53 && (moduleHash<12 || moduleHash>2035)) {
248 charge = ToT/8.0*(8000.0-1200.0)+1200.0;
249 }
250 chargeList.push_back(charge);
251 }
252 }
253
254 for (int i=0; i<nRDO; i++) {
255 Identifier pixid=rdoList[i];
256 int ToT=totList[i];
257
258 float charge = ToT;
259 if (calibData) { charge=chargeList[i]; }
260
261 // std::cout << "tot, charge = " << ToT << " " << charge << std::endl;
262 int row = pixelID.phi_index(pixid);
263 int col = pixelID.eta_index(pixid);
264 if (row == rowMin) qRowMin += charge;
265 if (row < rowMin){
266 rowMin = row;
267 qRowMin = charge;
268 }
269
270 if (row == rowMax) qRowMax += charge;
271 if (row > rowMax){
272 rowMax = row;
273 qRowMax = charge;
274 }
275 if (col == colMin) qColMin += charge;
276 if (col < colMin){
277 colMin = col;
278 qColMin = charge;
279 }
280
281 if (col == colMax) qColMax += charge;
282 if (col > colMax){
283 colMax = col;
284 qColMax = charge;
285 }
286 }
287
288 Identifier newClusterID = pixelID.pixel_id(pixelID.wafer_id(clusterID),rowMin,colMin);
289 // Compute omega for charge interpolation correction (if required)
290 // Two pixels may have charge=0 (very rarely, hopefully)
291 float omegax = -1;
292 float omegay = -1;
293 if(qRowMin+qRowMax > 0) omegax = qRowMax/float(qRowMin+qRowMax);
294 if(qColMin+qColMax > 0) omegay = qColMax/float(qColMin+qColMax);
295
296 ATH_MSG_VERBOSE("omega = " << omegax << " " << omegay);
297
298// ask for Lorentz correction, get global position
299 double shift = m_pixelLorentzAngleTool->getLorentzShift(element->identifyHash(), ctx);
300 Amg::Vector2D locpos(localPos[Trk::locX]+shift, localPos[Trk::locY]);
301// find global position of element
302 const Amg::Transform3D& T = element->surface().transform();
303 double Ax[3] = {T(0,0),T(1,0),T(2,0)};
304 double Ay[3] = {T(0,1),T(1,1),T(2,1)};
305 double R [3] = {T(0,3),T(1,3),T(2,3)};
306
307 const Amg::Vector2D& M = locpos;
308 Amg::Vector3D globalPos(M[0]*Ax[0]+M[1]*Ay[0]+R[0],M[0]*Ax[1]+M[1]*Ay[1]+R[1],M[0]*Ax[2]+M[1]*Ay[2]+R[2]);
309
310 // error matrix
311 const Amg::Vector2D& colRow = width.colRow();// made ref to avoid
312 // unnecessary copy EJWM
313 auto errorMatrix = Amg::MatrixX(2,2);
314 errorMatrix.setIdentity();
315
316 // switches are more readable **OPT**
317 // actually they're slower as well (so I'm told) so perhaps
318 // this should be re-written at some point EJWM
319 double eta = std::abs(globalPos.eta());
320 double zPitch = width.z()/colRow.y();
321
322 const AtlasDetectorID* aid = element->getIdHelper();
323
325 throw std::runtime_error( "Wrong helper type in ClusterMakerTool.cxx.");
326 }
327 const PixelID* pid = static_cast<const PixelID*>(aid);
328 int layer = pid->layer_disk(clusterID);
329 int phimod = pid->phi_module(clusterID);
330 switch (errorStrategy){
331 case 0:
332 errorMatrix.fillSymmetric(0,0,square(width.phiR())*ONE_TWELFTH);
333 errorMatrix.fillSymmetric(1,1,square(width.z())*ONE_TWELFTH);
334 break;
335 case 1:
336 errorMatrix.fillSymmetric(0,0,square(width.phiR()/colRow.x())*ONE_TWELFTH);
337 errorMatrix.fillSymmetric(1,1,square(width.z()/colRow.y())*ONE_TWELFTH);
338 break;
339 case 2:
340 // use parameterization only if the cluster does not
341 // contain long pixels or ganged pixels
342 // Also require calibration service is available....
343 if (!ganged && zPitch>399*micrometer && zPitch<401*micrometer) {
344 if (offlineCalibData) {
345 if (element->isBarrel()) {
346 int ibin = offlineCalibData->getPixelClusterErrorData()->getBarrelBin(eta,int(colRow.y()),int(colRow.x()));
347 double phiError = offlineCalibData->getPixelClusterErrorData()->getPixelBarrelPhiError(ibin);
348 double etaError = offlineCalibData->getPixelClusterErrorData()->getPixelBarrelEtaError(ibin);
349 errorMatrix.fillSymmetric(0,0,pow(phiError,2));
350 errorMatrix.fillSymmetric(1,1,pow(etaError,2));
351 }
352 else {
353 int ibin = offlineCalibData->getPixelClusterErrorData()->getEndcapBin(int(colRow.y()),int(colRow.x()));
354 double phiError = offlineCalibData->getPixelClusterErrorData()->getPixelEndcapPhiError(ibin);
355 double etaError = offlineCalibData->getPixelClusterErrorData()->getPixelEndcapRError(ibin);
356 errorMatrix.fillSymmetric(0,0,square(phiError));
357 errorMatrix.fillSymmetric(1,1,square(etaError));
358 }
359 }
360 }else{// cluster with ganged and/or long pixels
361 errorMatrix.fillSymmetric(0,0,square(width.phiR()/colRow.x())*ONE_TWELFTH);
362 errorMatrix.fillSymmetric(1,1,square(zPitch)*ONE_TWELFTH);
363 }
364 break;
365
366 case 10:
367 errorMatrix.fillSymmetric(0,0,square( getPixelCTBPhiError(layer,phimod,int(colRow.x()))));
368 errorMatrix.fillSymmetric(1,1,square(width.z()/colRow.y())*ONE_TWELFTH);
369 break;
370
371 default:
372 errorMatrix.fillSymmetric(0,0,square(width.phiR()/colRow.x())*ONE_TWELFTH);
373 errorMatrix.fillSymmetric(1,1,square(width.z()/colRow.y())*ONE_TWELFTH);
374 break;
375 }
376
377 //1) We want to move always for the Trk::PixelCluster
378 //2) We forward the rdoList and chargeList that could have be passed
379 // different ways.
380 //
381 static_assert(std::is_same_v<ClusterType, PixelCluster> ||
382 std::is_same_v<std::remove_pointer_t<ClusterType>,
384 "Not an InDet::PixelCluster or xAOD::PixelCluster");
385 if constexpr (std::is_same<ClusterType, InDet::PixelCluster>::value) {
386 return newInDetpixelCluster(newClusterID, locpos,
387 globalPos,
388 std::forward<IdentifierVec>(rdoList),
389 lvl1a,
390 std::forward<ToTList>(totList),
391 std::move(chargeList),
392 width,
393 element,
394 std::move(errorMatrix),
395 omegax,
396 omegay,
397 split,
398 splitProb1,
399 splitProb2);
400 } else{
401 return AddNewxAODpixelCluster(*cluster)(newClusterID,
402 locpos,
403 globalPos,
404 std::forward<IdentifierVec>(rdoList),
405 lvl1a,
406 std::forward<ToTList>(totList),
407 chargeList,
408 width,
409 element,
410 errorMatrix,
411 split,
412 splitProb1,
413 splitProb2);
414 }
415}
Scalar eta() const
pseudorapidity method
#define ATH_MSG_VERBOSE(x)
double charge(const T &p)
Definition AtlasPID.h:997
const double width
constexpr int pow(int base, int exp) noexcept
virtual HelperType helper() const
Type of helper, defaulted to 'Unimplemented'.
static constexpr std::array< PixelDiodeTree::CellIndexType, 2 > makeCellIndex(T local_x_idx, T local_y_idx)
Create a 2D cell index from the indices in local-x (phi, row) and local-y (eta, column) direction.
virtual const SiDetectorDesign & design() const override final
access to the local description (inline):
virtual IdentifierHash identifyHash() const override final
identifier hash (inline)
const AtlasDetectorID * getIdHelper() const
Returns the id helper (inline)
Trk::Surface & surface()
Element Surface.
double getPixelCTBPhiError(int layer, int phi, int PhiClusterSize) const
int getBarrelBin(double eta, int etaClusterSize, int phiClusterSize) const
int getEndcapBin(int etaClusterSize, int phiClusterSize) const
PixelClusterErrorData * getPixelClusterErrorData()
float getCharge(InDetDD::PixelDiodeType type, unsigned int moduleHash, unsigned int FE, float ToT) const
int eta_index(const Identifier &id) const
Definition PixelID.h:645
Identifier wafer_id(int barrel_ec, int layer_disk, int phi_module, int eta_module) const
For a single crystal.
Definition PixelID.h:360
Identifier pixel_id(int barrel_ec, int layer_disk, int phi_module, int eta_module, int phi_index, int eta_index) const
For an individual pixel.
Definition PixelID.h:428
IdentifierHash wafer_hash(Identifier wafer_id) const
wafer hash from id
Definition PixelID.h:383
int phi_index(const Identifier &id) const
Definition PixelID.h:639
const Amg::Transform3D & transform() const
Returns HepGeom::Transform3D by reference.
std::vector< std::string > split(const std::string &s, const std::string &t=":")
Definition hcg.cxx:177
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 2, 1 > Vector2D
Eigen::Matrix< double, 3, 1 > Vector3D
double R(const INavigable4Momentum *p1, const double v_eta, const double v_phi)
phimod(flags, cells_name, *args, **kw)
float etaError(const U &p)
row
Appending html table to final .html summary file.
@ layer
Definition HitInfo.h:79
unsigned long long T
@ locY
local cartesian
Definition ParamDefs.h:38
@ locX
Definition ParamDefs.h:37
PixelCluster_v1 PixelCluster
Define the version of the pixel cluster class.

◆ msg()

MsgStream & AthCommonMsg< AlgTool >::msg ( ) const
inlineinherited

Definition at line 24 of file AthCommonMsg.h.

24 {
25 return this->msgStream();
26 }

◆ msgLvl()

bool AthCommonMsg< AlgTool >::msgLvl ( const MSG::Level lvl) const
inlineinherited

Definition at line 30 of file AthCommonMsg.h.

30 {
31 return this->msgLevel(lvl);
32 }

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::outputHandles ( ) const
overridevirtualinherited

Return this algorithm's output handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ pixelCluster()

PixelCluster InDet::ClusterMakerTool::pixelCluster ( const Identifier & clusterID,
const Amg::Vector2D & localPos,
std::vector< Identifier > && rdoList,
const int lvl1a,
std::vector< int > && totList,
const SiWidth & width,
const InDetDD::SiDetectorElement * element,
bool ganged,
int errorStrategy,
const PixelID & pixelID,
bool split,
double splitProb1,
double splitProb2,
const PixelChargeCalibCondData * calibData,
const PixelCalib::PixelOfflineCalibData * offlineCalibData,
const EventContext & ctx ) const

Definition at line 417 of file ClusterMakerTool.cxx.

434{
436 clusterID,
437 localPos,
438 std::move(rdoList),
439 lvl1a,
440 std::move(totList),
441 width,
442 element,
443 ganged,
444 errorStrategy,
445 pixelID,
446 split,
447 splitProb1,
448 splitProb2,
449 calibData,
450 offlineCalibData,
451 ctx);
452}
ClusterType makePixelCluster(const Identifier &clusterID, const Amg::Vector2D &localPos, IdentifierVec &&rdoList, const int lvl1a, ToTList &&totList, const SiWidth &width, const InDetDD::SiDetectorElement *element, bool ganged, int errorStrategy, const PixelID &pixelID, bool split, double splitProb1, double splitProb2, const PixelChargeCalibCondData *calibData, const PixelCalib::PixelOfflineCalibData *offlineCalibData, const EventContext &ctx, xAOD::PixelCluster *cluster=nullptr) const

◆ renounce()

std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > AthCommonDataStore< AthCommonMsg< AlgTool > >::renounce ( T & h)
inlineprotectedinherited

Definition at line 380 of file AthCommonDataStore.h.

381 {
382 h.renounce();
384 }
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)

◆ renounceArray()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ sctCluster()

SCT_Cluster InDet::ClusterMakerTool::sctCluster ( const Identifier & clusterID,
const Amg::Vector2D & localPos,
std::vector< Identifier > && rdoList,
const SiWidth & width,
const InDetDD::SiDetectorElement * element,
int errorStrategy ) const

Definition at line 508 of file ClusterMakerTool.cxx.

514{
515
516 double shift =
517 m_sctLorentzAngleTool->getLorentzShift(element->identifyHash(), Gaudi::Hive::currentContext());
518 Amg::Vector2D locpos(localPos[Trk::locX] + shift, localPos[Trk::locY]);
519
520 // error matrix
521 const Amg::Vector2D& colRow = width.colRow(); // made ref to avoid
522 // unnecessary copy EJWM
523
524 auto errorMatrix = Amg::MatrixX(2,2);
525 errorMatrix.setIdentity();
526
527 // switches are more readable **OPT**
528 // actually they're slower as well (so I'm told) so perhaps
529 // this should be re-written at some point EJWM
530
531 switch (errorStrategy){
532 case 0:
533 errorMatrix.fillSymmetric(0,0,square(width.phiR())*ONE_TWELFTH);
534 errorMatrix.fillSymmetric(1,1,square(width.z())*ONE_TWELFTH);
535 break;
536 case 1:
537 // mat(1,1) = pow(width.phiR()/colRow.x(),2)/12;
538 // single strip - resolution close to pitch/sqrt(12)
539 // two-strip hits: better resolution, approx. 40% lower
540 if(colRow.x() == 1){
541 errorMatrix.fillSymmetric(0,0,square(1.05*width.phiR())*ONE_TWELFTH);
542 }
543 else if(colRow.x() == 2){
544 errorMatrix.fillSymmetric(0,0,square(0.27*width.phiR())*ONE_TWELFTH);
545 }
546 else{
547 errorMatrix.fillSymmetric(0,0,square(width.phiR())*ONE_TWELFTH);
548 }
549 errorMatrix.fillSymmetric(1,1,square(width.z()/colRow.y())*ONE_TWELFTH);
550 break;
551 default:
552 // single strip - resolution close to pitch/sqrt(12)
553 // two-strip hits: better resolution, approx. 40% lower
554 if(colRow.x() == 1){
555 errorMatrix.fillSymmetric(0,0,square(width.phiR())*ONE_TWELFTH);
556 }
557 else if(colRow.x() == 2){
558 errorMatrix.fillSymmetric(0,0,square(0.27*width.phiR())*ONE_TWELFTH);
559 }
560 else{
561 errorMatrix.fillSymmetric(0,0,square(width.phiR())*ONE_TWELFTH);
562 }
563 errorMatrix.fillSymmetric(1,1,square(width.z()/colRow.y())*ONE_TWELFTH);
564 break;
565 }
566
567 auto designShape = element->design().shape();
568 // rotation for endcap SCT
569 if(designShape == InDetDD::Trapezoid || designShape == InDetDD::Annulus) {
570 double sn = element->sinStereoLocal(localPos);
571 double sn2 = sn*sn;
572 double cs2 = 1.-sn2;
573 double w = element->phiPitch(localPos)/element->phiPitch();
574 double v0 = (errorMatrix)(0,0)*w*w;
575 double v1 = (errorMatrix)(1,1);
576 errorMatrix.fillSymmetric(0,0,cs2*v0+sn2*v1);
577 errorMatrix.fillSymmetric(0,1,sn*sqrt(cs2)*(v0-v1));
578 errorMatrix.fillSymmetric(1,1,sn2*v0+cs2*v1);
579 } //else if (designShape == InDetDD::PolarAnnulus) {// Polar rotation for endcap}
580
581 return SCT_Cluster(clusterID, locpos, std::move(rdoList), width, element, std::move(errorMatrix));
582
583}
virtual DetectorShape shape() const
Shape of element.
double phiPitch() const
Pitch (inline methods)
double sinStereoLocal(const Amg::Vector2D &localPos) const
Angle of strip in local frame with respect to the etaAxis.

◆ sysInitialize()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysInitialize ( )
overridevirtualinherited

Perform system initialization for an algorithm.

We override this to declare all the elements of handle key arrays at the end of initialization. See comments on updateVHKA.

Reimplemented in asg::AsgMetadataTool, AthCheckedComponent< AthAlgTool >, AthCheckedComponent<::AthAlgTool >, and DerivationFramework::CfAthAlgTool.

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysStart ( )
overridevirtualinherited

Handle START transition.

We override this in order to make sure that conditions handle keys can cache a pointer to the conditions container.

◆ updateVHKA()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::updateVHKA ( Gaudi::Details::PropertyBase & )
inlineinherited

Definition at line 308 of file AthCommonDataStore.h.

308 {
309 // debug() << "updateVHKA for property " << p.name() << " " << p.toString()
310 // << " size: " << m_vhka.size() << endmsg;
311 for (auto &a : m_vhka) {
313 for (auto k : keys) {
314 k->setOwner(this);
315 }
316 }
317 }
std::vector< SG::VarHandleKeyArray * > m_vhka

◆ xAODpixelCluster()

xAOD::PixelCluster * InDet::ClusterMakerTool::xAODpixelCluster ( xAOD::PixelCluster & cluster,
const Amg::Vector2D & localPos,
const std::vector< Identifier > & rdoList,
const int lvl1a,
const std::vector< int > & totList,
const SiWidth & width,
const InDetDD::SiDetectorElement * element,
bool ganged,
int errorStrategy,
const PixelID & pixelID,
bool split,
double splitProb1,
double splitProb2,
const PixelChargeCalibCondData * calibData,
const PixelCalib::PixelOfflineCalibData * offlineCalibData,
const EventContext & ctx ) const

Definition at line 454 of file ClusterMakerTool.cxx.

471{
473 Identifier(),
474 localPos,
475 rdoList,
476 lvl1a,
477 totList,
478 width,
479 element,
480 ganged,
481 errorStrategy,
482 pixelID,
483 split,
484 splitProb1,
485 splitProb2,
486 calibData,
487 offlineCalibData,
488 ctx,
489 &cluster);
490}

Member Data Documentation

◆ m_detStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_evtStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_evtStore
privateinherited

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_forceErrorStrategy1B

bool InDet::ClusterMakerTool::m_forceErrorStrategy1B {false}
private

Definition at line 172 of file ClusterMakerTool.h.

172{false};

◆ m_pixelLorentzAngleTool

ToolHandle<ISiLorentzAngleTool> InDet::ClusterMakerTool::m_pixelLorentzAngleTool {this, "PixelLorentzAngleTool", "SiLorentzAngleTool/PixelLorentzAngleTool", "Tool to retreive Lorentz angle of Pixel"}
private

Definition at line 166 of file ClusterMakerTool.h.

167{this, "PixelLorentzAngleTool", "SiLorentzAngleTool/PixelLorentzAngleTool", "Tool to retreive Lorentz angle of Pixel"};

◆ m_sctLorentzAngleTool

ToolHandle<ISiLorentzAngleTool> InDet::ClusterMakerTool::m_sctLorentzAngleTool {this, "SCTLorentzAngleTool", "SiLorentzAngleTool/SCTLorentzAngleTool", "Tool to retreive Lorentz angle of SCT"}
private

Definition at line 169 of file ClusterMakerTool.h.

170{this, "SCTLorentzAngleTool", "SiLorentzAngleTool/SCTLorentzAngleTool", "Tool to retreive Lorentz angle of SCT"};

◆ m_varHandleArraysDeclared

bool AthCommonDataStore< AthCommonMsg< AlgTool > >::m_varHandleArraysDeclared
privateinherited

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

std::vector<SG::VarHandleKeyArray*> AthCommonDataStore< AthCommonMsg< AlgTool > >::m_vhka
privateinherited

Definition at line 398 of file AthCommonDataStore.h.


The documentation for this class was generated from the following files: