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 133 of file ClusterMakerTool.cxx.

135 :
136 AthAlgTool(t,n,p)
137{
138 declareInterface<ClusterMakerTool>(this);
139}
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 581 of file ClusterMakerTool.cxx.

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

◆ initialize()

StatusCode InDet::ClusterMakerTool::initialize ( )

Definition at line 143 of file ClusterMakerTool.cxx.

143 {
144 // Code entered here will be executed once at program start.
145
146 ATH_MSG_DEBUG ( name() << " initialize()" );
147
148 if (not m_pixelLorentzAngleTool.empty()) {
150 } else {
151 m_pixelLorentzAngleTool.disable();
152 }
153 if (not m_sctLorentzAngleTool.empty()) {
155 } else {
156 m_sctLorentzAngleTool.disable();
157 }
158
159
160 return StatusCode::SUCCESS;
161}
#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 185 of file ClusterMakerTool.cxx.

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

428{
430 clusterID,
431 localPos,
432 std::move(rdoList),
433 lvl1a,
434 std::move(totList),
435 width,
436 element,
437 ganged,
438 errorStrategy,
439 pixelID,
440 split,
441 splitProb1,
442 splitProb2,
443 calibData,
444 offlineCalibData,
445 ctx);
446}
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 502 of file ClusterMakerTool.cxx.

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

465{
467 Identifier(),
468 localPos,
469 rdoList,
470 lvl1a,
471 totList,
472 width,
473 element,
474 ganged,
475 errorStrategy,
476 pixelID,
477 split,
478 splitProb1,
479 splitProb2,
480 calibData,
481 offlineCalibData,
482 ctx,
483 &cluster);
484}

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: