ATLAS Offline Software
Loading...
Searching...
No Matches
SingleTrackValidation Class Reference

#include <SingleTrackValidation.h>

Inheritance diagram for SingleTrackValidation:
Collaboration diagram for SingleTrackValidation:

Classes

class  Clockwork

Public Member Functions

 SingleTrackValidation (const std::string &name, ISvcLocator *pSvcLocator)
 ~SingleTrackValidation ()
StatusCode initialize () override
StatusCode execute () override
StatusCode finalize () override
virtual StatusCode sysInitialize () override
 Override sysInitialize.
virtual const DataObjIDColl & extraOutputDeps () const override
 Return the list of extra output dependencies.
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 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

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

 SingleTrackValidation (const SingleTrackValidation &)
SingleTrackValidationoperator= (const SingleTrackValidation &)
Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

SG::ReadHandleKey< McEventCollectionm_truthKey { this, "TruthKey", "TruthEvent" }
SG::ReadCondHandleKey< CaloDetDescrManagerm_caloMgrKey
SG::ReadCondHandleKey< AtlasFieldCacheCondObjm_fieldCacheCondObjInputKey
ServiceHandle< IPartPropSvc > m_ppSvc {this, "PartPropSvc", "PartPropSvc"}
ServiceHandle< ITHistSvc > m_histSvc { this, "THistSvc", "THistSvc", "Histogramming svc" }
Clockworkm_c
TH1F * m_histos [162] {}
DataObjIDColl m_extendedExtraObjects
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 20 of file SingleTrackValidation.h.

Member Typedef Documentation

◆ StoreGateSvc_t

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

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ SingleTrackValidation() [1/2]

SingleTrackValidation::SingleTrackValidation ( const std::string & name,
ISvcLocator * pSvcLocator )

Definition at line 108 of file SingleTrackValidation.cxx.

108 :
109 AthAlgorithm(name,pSvcLocator),m_c(new Clockwork())
110{
111 for (unsigned int i=0;i<162;++i) m_histos[i] = nullptr;
112}
AthAlgorithm()
Default constructor:

◆ ~SingleTrackValidation()

SingleTrackValidation::~SingleTrackValidation ( )

Definition at line 114 of file SingleTrackValidation.cxx.

114 {
115 if (m_c!=nullptr){ delete m_c; m_c=nullptr; }
116}

◆ SingleTrackValidation() [2/2]

SingleTrackValidation::SingleTrackValidation ( const SingleTrackValidation & )
private

Member Function Documentation

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Algorithm > >::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< Algorithm > >::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< Algorithm > >::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< Algorithm > >::evtStore ( )
inlineinherited

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

Definition at line 85 of file AthCommonDataStore.h.

◆ execute()

StatusCode SingleTrackValidation::execute ( )
override

Definition at line 257 of file SingleTrackValidation.cxx.

257 {
258
259 if (m_c->cpuTime==0) {
260 m_c->cpuTime=getCpu();
261 m_c->cpuTime+=1; // Fill the histogram with -1 for the first event
262 }
263 m_c->cpuTime= getCpu()-m_c->cpuTime;
264 m_histos[156]->Fill( m_c->cpuTime/100. , 1. );
265
266 const EventContext& context = getContext();
267 int RunNum=context.eventID().run_number();
268 int EvtNum=context.eventID().event_number();
269 double RunStr=double(RunNum);
270 double EvtStr=double(EvtNum);
271 m_c->EventNo=EvtStr;
272 m_c->RunNo=RunStr;
273 m_histos[160]->Fill(EvtStr);
274 m_histos[159]->Fill(RunNum);
275
276 MagField::AtlasFieldCache fieldCache;
277 // Get field cache object
278 SG::ReadCondHandle<AtlasFieldCacheCondObj> readHandle{m_fieldCacheCondObjInputKey};
279 const AtlasFieldCacheCondObj* fieldCondObj{*readHandle};
280 if (fieldCondObj == nullptr) {
281 ATH_MSG_ERROR("Failed to retrieve AtlasFieldCacheCondObj with key " << m_fieldCacheCondObjInputKey.key());
282 return StatusCode::FAILURE;
283 }
284 fieldCondObj->getInitializedCache (fieldCache);
285
286 // Get the MC Truth Information
287 SG::ReadHandle<McEventCollection> mcEvent{m_truthKey};
288 for (const HepMC::GenEvent* e : *mcEvent) {
289
290 // Get just the primary, call it "theParticle"
291 auto theParticle = *HepMC::begin(*e);
292
293 // Fetch whatever particle properties will be used in the following:
294 const HepPDT::ParticleDataTable * dataTable = m_c->partPropSvc->PDT();
295 const HepPDT::ParticleData * particleData = dataTable->particle(iabs(theParticle->pdg_id()));
296
297 // Get the kinematic variables:
298 HepLorentzVector momentum(theParticle->momentum().px(),
299 theParticle->momentum().py(),
300 theParticle->momentum().pz(),
301 theParticle->momentum().e());
302 Point3D<double> origin(theParticle->production_vertex()->position().x(),
303 theParticle->production_vertex()->position().y(),
304 theParticle->production_vertex()->position().z());
305 double charge = theParticle->pdg_id() > 0 ? particleData->charge() : - particleData->charge();
306 // Put Eta and Phi into the Ntuple
307 m_c->phi = theParticle->momentum().phi();
308 m_c->eta = -log(tan(theParticle->momentum().theta()/2));
309 if (!finite(m_c->eta)) m_c->eta=0;
310 m_c->pt = theParticle->momentum().perp();
311 int partID = theParticle->pdg_id();
312 double pID = double(partID);
313 m_c->PDG = pID;
314 m_c->Energy = theParticle->momentum().e();
315 m_histos[2]->Fill( theParticle->momentum().phi() );
316 double myEta = -log(tan(theParticle->momentum().theta()/2));
317 if (!finite(myEta)) m_histos[0]->Fill(0);
318 else m_histos[0]->Fill( myEta );
319 m_histos[158]->Fill( pID );
320 m_histos[1]->Fill( theParticle->momentum().perp()/Units::GeV );
321 m_histos[157]->Fill( m_c->Energy = theParticle->momentum().e()/Units::GeV );
322
323 // Make an extrapolator:
324 const Genfun::AtlasBComponent Bx(0,&fieldCache);
325 const Genfun::AtlasBComponent By(1,&fieldCache);
326 const Genfun::AtlasBComponent Bz(2,&fieldCache);
327 GeoXPEngine extrapolator(Bx, By, Bz, origin, momentum, charge);
328
329 // Extrapolate to the first layer in the barrel:
330 m_c->x = 0;
331 m_c->y = 0;
332 m_c->z = 0;
333 double x=0,y=0,z=0;
334 bool hitsBarrel=false;
335 //bool hitsEndcap=false;
336 for (double t = 0; t< 50; t += 0.1) {
337 x = extrapolator.x()(t);
338 y = extrapolator.y()(t);
339 z = extrapolator.z()(t);
340 double magicZ=3640.0*mm;
341 double magicR=1500.0*mm;
342 if (x*x+y*y > magicR*magicR) {
343 m_c->x = x;
344 m_c->y = y;
345 m_c->z = z;
346 hitsBarrel=true;
347 break;
348 }
349 else if (z*z > magicZ*magicZ) {
350 m_c->x = x;
351 m_c->y = y;
352 m_c->z = z;
353 //hitsEndcap=true;
354 break;
355 }
356 }
357
358 m_histos[3]->Fill( x );
359 m_histos[4]->Fill( y );
360 m_histos[5]->Fill( z );
361
362 // You have an x,y, and z position. Now go and get the Element corresponding to
363 // that hit position. There are four, one for each sampling layer:
364 double radImpact = std::sqrt(x*x+y*y+z*z);
365 double phiImpact = std::atan2(y,x);
366 double thetaImpact = std::acos(z/radImpact);
367 double etaImpact = -std::log(std::tan(thetaImpact/2));
368
369 SG::ReadCondHandle<CaloDetDescrManager> caloMgrHandle{m_caloMgrKey};
370 ATH_CHECK(caloMgrHandle.isValid());
371 const CaloDetDescrManager* caloMgr = *caloMgrHandle;
372 const CaloDetDescrElement *element[15]={nullptr};
373
374 for (int i=0;i<4;i++) {
375 try {
376 element[i] = caloMgr->get_element(CaloCell_ID::LAREM,i, hitsBarrel, etaImpact, phiImpact);
377 }
378 catch (const LArID_Exception & e) {
379 std::cerr << "SingleTrackValidation EXCEPTION (LAREM)" << e.message() << std::endl;
380 }
381 }
382 for (int i=0;i<4;i++) {
383 try {
384 element[i+4] = caloMgr->get_element(CaloCell_ID::LAREM,i, hitsBarrel, etaImpact, phiImpact);
385 }
386 catch (const LArID_Exception & e) {
387 std::cerr << "SingleTrackValidation EXCEPTION (LAREM)" << e.message() << std::endl;
388 }
389 }
390 for (int i=0;i<4;i++) {
391 try {
392 element[i+8] = caloMgr->get_element(CaloCell_ID::LARHEC,i, hitsBarrel, etaImpact, phiImpact);
393 }
394 catch (const LArID_Exception & e) {
395 std::cerr << "SingleTrackValidation EXCEPTION in (LARHEC)" << e.message() << std::endl;
396 }
397 }
398 for (int i=1;i<4;i++) {
399 try {
400 element[i+11] = caloMgr->get_element(CaloCell_ID::LARFCAL,i, hitsBarrel, etaImpact, phiImpact);
401 }
402 catch (const LArID_Exception & e) {
403 std::cerr << "SingleTrackValidation EXCEPTIONin LARFCAL" << e.message() << std::endl;
404 }
405 }
406
407
408 // Now go and find out how much energy is there:
409 for (int z=0;z<15;z++){
410 m_c->s_c00[z]=0;
411 m_c->s_t00[z]=0;
412 }
413
414 std::string lArKey = hitsBarrel ? "LArHitEMB" : "LArHitEMEC" ;
415
416 double eSum [15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
417 double eEta [15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
418 double eEtaEta [15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
419 double ePhi [15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
420 double ePhiPhi [15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
421 int hit_count [15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
422 double eX [15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
423 double eXX [15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
424 double eY [15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
425 double eYY [15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
426 double c00 [15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
427 double t00 [15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
428 // double width [15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
429 double e_dep = 0;
430
431 for (int i=0;i<4;i++) { // Loop over the four LAr collections
432 if (i==0) {
433 lArKey="LArHitEMB";
434 } else if (i==1) {
435 lArKey="LArHitEMEC";
436 } else if (i==2) {
437 lArKey="LArHitHEC";
438 } else if (i==3) {
439 lArKey="LArHitFCAL";
440 }
441
442 SG::ReadHandle<LArHitContainer> larHitContainer(lArKey, "StoreGateSvc");
443 if (larHitContainer.isValid()) {
444 for (const LArHit* larHit : *larHitContainer) {
445 const CaloDetDescrElement *hitElement = caloMgr->get_element(larHit->cellID());
446 int samplingLayer = m_c->cellId->sampling(larHit->cellID());
447 double energy = larHit->energy();
448
449 for (int j=0;j<15;j++) {
450 if (hitElement==element[j]) {
451 c00[j] += energy;
452 }
453 }
454
455 if (lArKey=="LArHitEMEC") samplingLayer += 4;
456 else if (lArKey=="LArHitHEC") samplingLayer += 8;
457 else if (lArKey=="LArHitFCAL") samplingLayer += 11;
458
459 // Calculate phi of hit w.r.t phiImpact
460 // to avoid problems due to the 2pi modulus of phi
461 // we are thus calculating deltaphi directly, so don't subtract phiImpact again later
462 double hitPhi = hitElement->phi() - phiImpact;
463 if (hitPhi < -pi) hitPhi += twopi;
464 if (hitPhi > pi) hitPhi -= twopi;
465
466 eSum [samplingLayer]+=energy;
467 eEta [samplingLayer]+=energy*hitElement->eta();
468 eEtaEta [samplingLayer]+=energy*hitElement->eta()*hitElement->eta();
469 ePhi [samplingLayer]+=energy*hitPhi;
470 ePhiPhi [samplingLayer]+=energy*hitPhi*hitPhi;
471 eX [samplingLayer]+=energy*hitElement->x();
472 eXX [samplingLayer]+=energy*hitElement->x()*hitElement->x();
473 eY [samplingLayer]+=energy*hitElement->y();
474 eYY [samplingLayer]+=energy*hitElement->y()*hitElement->y();
475 t00 [samplingLayer]+=energy*larHit->time();
476 hit_count[samplingLayer]+=1;
477 }
478 }
479 }
480
481 for (int i=0;i<15;i++) {
482 if (eSum[i]!=0) eEta[i] /= eSum[i];
483 if (eSum[i]!=0) eEtaEta[i] /= eSum[i];
484 if (eSum[i]!=0) ePhi[i] /= eSum[i];
485 if (eSum[i]!=0) ePhiPhi[i] /= eSum[i];
486 if (eSum[i]!=0) eY[i] /= eSum[i];
487 if (eSum[i]!=0) eYY[i] /= eSum[i];
488 if (eSum[i]!=0) eX[i] /= eSum[i];
489 if (eSum[i]!=0) eXX[i] /= eSum[i];
490 if (eSum[i]!=0) t00[i] /= eSum[i];
491 e_dep+=eSum[i];
492 }
493
494 double dThetaDEta = -1.0/cosh(etaImpact);
495
496 //Fill the E Sum, center cell E, hit count fields:
497 m_c->E_Deposit=e_dep;
498 for (int z=0;z<15;z++){
499 m_c->s_sumE[z]=eSum[z];
500 m_c->s_c00[z]=c00[z];
501 m_c->s_t00[z]=t00[z];
502 m_c->s_hits[z]=hit_count[z];
503 if (z<12){
504 m_c->s_deltaPhi[z]=radImpact*std::sin(thetaImpact)*(ePhi[z]);
505 m_c->s_sigmaPhi[z]=radImpact*std::sin(thetaImpact)*std::sqrt(ePhiPhi[z]- ePhi[z]*ePhi[z]);
506 m_c->s_deltaEta[z]=radImpact*dThetaDEta*(eEta[z]-etaImpact);
507 m_c->s_sigmaEta[z]=radImpact*std::fabs(dThetaDEta)*std::sqrt(eEtaEta[z]- eEta[z]*eEta[z]);
508 } else {
509 m_c->s_deltaPhi[z]=(eX[z]-x);
510 m_c->s_sigmaPhi[z]=std::sqrt(eXX[z]- eX[z]*eX[z]);
511 m_c->s_deltaEta[z]=(eY[z]-y);
512 m_c->s_sigmaEta[z]=std::sqrt(eYY[z]-eY[z]*eY[z]);
513 }
514 m_c->s_widthX[z]=std::sqrt(eXX[z]-eX[z]*eX[z]);
515 m_c->s_widthY[z]=std::sqrt(eYY[z]-eY[z]*eY[z]);
516 }
517
518 m_histos[161]->Fill(e_dep/Units::GeV);
519 for (int i=0;i<15;i++){
520 m_histos[6+i]->Fill( c00[i]/Units::GeV );
521 m_histos[21+i]->Fill( hit_count[i] );
522 m_histos[36+i]->Fill( eSum[i]/Units::GeV );
523 m_histos[111+i]->Fill( t00[i] );
524 m_histos[126+i]->Fill( sqrt(eXX[i]-eX[i]*eX[i]) );
525 m_histos[141+i]->Fill( sqrt(eYY[i]-eY[i]*eY[i]) );
526 if (i<8){
527 m_histos[51+i]->Fill( radImpact*std::sin(thetaImpact)*ePhi[i] );
528 m_histos[66+i]->Fill( radImpact*std::sin(thetaImpact)*std::sqrt(ePhiPhi[i]-ePhi[i]*ePhi[i]) );
529 m_histos[81+i]->Fill( radImpact*dThetaDEta*(eEta[i]-etaImpact) );
530 m_histos[96+i]->Fill( radImpact*std::fabs(dThetaDEta)*std::sqrt(eEtaEta[i]-eEta[i]*eEta[i]) );
531 } else {
532 m_histos[51+i]->Fill( eX[i]-x );
533 m_histos[66+i]->Fill( std::sqrt(eXX[i]-eX[i]*eX[i]) );
534 m_histos[81+i]->Fill( eY[i]-y );
535 m_histos[96+i]->Fill( std::sqrt(eYY[i]-eY[i]*eY[i]) );
536 }
537 }
538
539 ATH_CHECK(ntupleSvc()->writeRecord(m_c->nt));
540
541 }
542
543 m_c->cpuTime= getCpu();
544
545 return StatusCode::SUCCESS;
546
547}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
double charge(const T &p)
Definition AtlasPID.h:997
int iabs(int a)
INTupleSvc * ntupleSvc()
int getCpu()
#define pi
#define y
#define x
#define z
constexpr double twopi
void getInitializedCache(MagField::AtlasFieldCache &cache) const
get B field cache for evaluation as a function of 2-d or 3-d position.
const CaloDetDescrElement * get_element(const Identifier &cellId) const
get element by its identifier
SG::ReadCondHandleKey< CaloDetDescrManager > m_caloMgrKey
SG::ReadHandleKey< McEventCollection > m_truthKey
SG::ReadCondHandleKey< AtlasFieldCacheCondObj > m_fieldCacheCondObjInputKey
GenEvent::particle_iterator begin(HepMC::GenEvent &e)
Definition GenEvent.h:622

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::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

◆ extraOutputDeps()

const DataObjIDColl & AthAlgorithm::extraOutputDeps ( ) const
overridevirtualinherited

Return the list of extra output dependencies.

This list is extended to include symlinks implied by inheritance relations.

Definition at line 50 of file AthAlgorithm.cxx.

51{
52 // If we didn't find any symlinks to add, just return the collection
53 // from the base class. Otherwise, return the extended collection.
54 if (!m_extendedExtraObjects.empty()) {
56 }
57 return Algorithm::extraOutputDeps();
58}
DataObjIDColl m_extendedExtraObjects

◆ finalize()

StatusCode SingleTrackValidation::finalize ( )
override

Definition at line 549 of file SingleTrackValidation.cxx.

549 {
550 return StatusCode::SUCCESS;
551}

◆ initialize()

StatusCode SingleTrackValidation::initialize ( )
override

Definition at line 118 of file SingleTrackValidation.cxx.

118 {
119 std::string names[162] = { "eta", "pt", "phi", "pos_x", "pos_y", "pos_z",
120 "emb0_cell", "emb1_cell", "emb2_cell", "emb3_cell", "emec0_cell", "emec1_cell", "emec2_cell", "emec3_cell",
121 "hec0_cell", "hec1_cell", "hec2_cell", "hec3_cell", "fc1_cell", "fc2_cell", "fc3_cell",
122 "emb0_hits", "emb1_hits", "emb2_hits", "emb3_hits", "emec0_hits", "emec1_hits", "emec2_hits", "emec3_hits",
123 "hec0_hits", "hec1_hits", "hec2_hits", "hec3_hits", "fc1_hits", "fc2_hits", "fc3_hits",
124 "emb0_sumE", "emb1_sumE", "emb2_sumE", "emb3_sumE", "emec0_sumE", "emec1_sumE", "emec2_sumE", "emec3_sumE",
125 "hec0_sumE", "hec1_sumE", "hec2_sumE", "hec3_sumE", "fc1_sumE", "fc2_sumE", "fc3_sumE",
126 "emb0_dPhi", "emb1_dPhi", "emb2_dPhi", "emb3_dPhi", "emec0_dPhi", "emec1_dPhi", "emec2_dPhi", "emec3_dPhi",
127 "hec0_dPhi", "hec1_dPhi", "hec2_dPhi", "hec3_dPhi", "fc1_dX", "fc2_dX", "fc3_dX",
128 "emb0_sPhi", "emb1_sPhi", "emb2_sPhi", "emb3_sPhi", "emec0_sPhi", "emec1_sPhi", "emec2_sPhi", "emec3_sPhi",
129 "hec0_sPhi", "hec1_sPhi", "hec2_sPhi", "hec3_sPhi", "fc1_sX", "fc2_sX", "fc3_sX",
130 "emb0_dEta", "emb1_dEta", "emb2_dEta", "emb3_dEta", "emec0_dEta", "emec1_dEta", "emec2_dEta", "emec3_dEta",
131 "hec0_dEta", "hec1_dEta", "hec2_dEta", "hec3_dEta", "fc1_dY", "fc2_dY", "fc3_dY",
132 "emb0_sEta", "emb1_sEta", "emb2_sEta", "emb3_sEta", "emec0_sEta", "emec1_sEta", "emec2_sEta", "emec3_sEta",
133 "hec0_sEta", "hec1_sEta", "hec2_sEta", "hec3_sEta", "fc1_sY", "fc2_sY", "fc3_sY",
134 "emb0_time", "emb1_time", "emb2_time", "emb3_time", "emec0_time", "emec1_time", "emec2_time", "emec3_time",
135 "hec0_time", "hec1_time", "hec2_time", "hec3_time", "fc1_time", "fc2_time", "fc3_time",
136 "emb0_widthX", "emb1_widthX", "emb2_widthX", "emb3_widthX", "emec0_widthX", "emec1_widthX", "emec2_widthX", "emec3_widthX",
137 "hec0_widthX", "hec1_widthX", "hec2_widthX", "hec3_widthX", "fc1_widthX", "fc2_widthX", "fc3_widthX",
138 "emb0_widthY", "emb1_widthY", "emb2_widthY", "emb3_widthY", "emec0_widthY", "emec1_widthY", "emec2_widthY", "emec3_widthY",
139 "hec0_widthY", "hec1_widthY", "hec2_widthY", "hec3_widthY", "fc1_widthY", "fc2_widthY", "fc3_widthY",
140 "cpuTime", "Energy", "PDG_ID", "RunNo", "EventNo", "E_Dep" };
141
142 double lim[162][2] = { {0,5}, {0,100}, {-4,4}, {-1600,1600}, {-1600,1600}, {-4000,4000},
143 {0,0.25}, {0,3}, {0,6}, {0,0.1}, {0,0.25}, {0,2}, {0,3}, {0,0.1}, {0,1}, {0,10}, {0,10}, {0,10}, {0,0.1}, {0,0.1}, {0,0.1}, //cell
144 {0,100}, {0,600}, {0,600}, {0,50}, {0,50}, {0,400}, {0,500}, {0,200}, {0,100}, {0,1000}, {0,1000}, {0,100}, {0,150}, {0,50}, {0,10}, //hits
145 {0,0.4}, {0,5}, {0,10}, {0,0.2}, {0,0.1}, {0,3}, {0,5}, {0,0.2}, {0,1}, {0,10}, {0,10}, {0,0.1}, {0,1}, {0,0.1}, {0,0.1}, //sumE
146 {-500,500}, {-100,100}, {-15,15}, {-200,200}, {-3000,3000}, {-60,60}, {-25,25}, {-200,200}, {-50,50}, {-50,50}, {-50,50}, {-50,50}, {-60,60}, {-500,500}, {-200,200}, //dPhi
147 {0,1000}, {0,500}, {0,200}, {0,500}, {0,2000}, {0,250}, {0,300}, {0,500}, {0,200}, {0,200}, {0,200}, {0,200}, {0,100}, {0,100}, {0,50}, //sPhi
148 {-150,150}, {-15,15}, {-20,20}, {-200,200}, {-0,2500}, {-50,20}, {-15,15}, {-150,150}, {-50,50}, {-50,50}, {-50,50}, {-50,50}, {-60,60}, {-500,500}, {-200,200}, //dEta
149 {0,500}, {0,100}, {0,100}, {0,400}, {0,1000}, {0,150}, {0,100}, {0,400}, {0,200}, {0,200}, {0,200}, {0,200}, {0,60}, {0,100}, {0,50}, //sPhi
150 {0,750}, {0,25}, {0,20}, {0,1000}, {0,10000}, {0,40}, {0,30}, {0,1000}, {0,1000}, {0,100}, {0,100}, {0,200}, {0,10}, {0,500}, {0,100}, //time
151 {-150,150}, {-15,15}, {-20,20}, {-200,200}, {-0,2500}, {-50,20}, {-15,15}, {-150,150}, {-50,50}, {-50,50}, {-50,50}, {-50,50}, {-60,60}, {-500,500}, {-200,200}, //widthX
152 {-150,150}, {-15,15}, {-20,20}, {-200,200}, {-0,2500}, {-50,20}, {-15,15}, {-150,150}, {-50,50}, {-50,50}, {-50,50}, {-50,50}, {-60,60}, {-500,500}, {-200,200}, //widthY
153 {0,50}, {0,100}, {-25,25}, {0,10}, {0,1000}, {0,10}};
154
155 //-------------------------------------------------------------------------//
156 // //
157 // Initialize the Particle Property Service. This is necessary in order //
158 // to obtain charge & type & other properties of the primary particle and //
159 // other particles that may turn up in the debris. //
160 // //
161 ATH_CHECK(m_ppSvc.retrieve());
162 m_c->partPropSvc = m_ppSvc.get();
163 ATH_CHECK(m_histSvc.retrieve());
164 m_c->histSvc = m_histSvc.get();
165 ATH_CHECK(detStore()->retrieve(m_c->cellId, "CaloCell_ID"));
166 ATH_CHECK(m_truthKey.initialize());
167 ATH_CHECK(m_caloMgrKey.initialize());
169
170 //----------------Now initialize the ntuples ----------------------//
171 // //
172 //==~ ~ ~==//
173 std::string filename="";
174 for (int i=0;i<162;i++){
175 m_histos[i] = new TH1F( names[i].data(), names[i].data(),50,lim[i][0],lim[i][1]);
176 filename = "/file1/Electron/";
177 filename.append(names[i]);
178 if (m_c->histSvc->regHist( filename.data() , m_histos[i] ).isFailure()){
179 ATH_MSG_WARNING( "Failed to register historam " << names[i] << ". Not sure what will happen now..." );
180 }
181 }
182
183
184 NTupleFilePtr file(ntupleSvc(),"/NTUPLES/FILE");
185 if (!file) throw std::runtime_error ("Ntuple MGR not open");
186 NTuple::Directory *col=ntupleSvc()->createDirectory("/NTUPLES/FILE/COL");
187 NTuplePtr nt(ntupleSvc(),"/NTUPLES/FILE/COL/SingleTrackValidation");
188 if (!nt) nt=ntupleSvc()->book(col, 1, CLID_ColumnWiseTuple, "SingleTrackValidation");
189
190 if (nt->addItem("Eta", m_c->eta ).isFailure() ||
191 nt->addItem("Pt", m_c->pt ).isFailure() ||
192 nt->addItem("BarrelX", m_c->x ).isFailure() ||
193 nt->addItem("BarrelY", m_c->y ).isFailure() ||
194 nt->addItem("BarrelZ", m_c->z ).isFailure() ||
195 nt->addItem("Phi", m_c->phi ).isFailure() ){
196 ATH_MSG_WARNING( "Registration of some of the ntuple branches failed. No idea what will happen next..." );
197 }
198
199
200 char title[80];
201
202 // Indices from 0 to 11: 4 EMB, 4 EMEC, and 4 FCAL sampling layers
203 // Handling FCAL layers separately for different titles (could do
204 // something more complicated if we wanted)
205 for (int i=0;i<15;i++){
206 if (i<12) sprintf(title,"S%i_C00",i);
207 else sprintf(title,"FC%i_C00",i-11);
208 if (nt->addItem(title,m_c->s_c00[i]).isFailure()) ATH_MSG_INFO( "Registration of a branch failed in the ntupler..." );
209 if (i<12) sprintf(title,"S%i_SumE",i);
210 else sprintf(title,"FC%i_SumE",i-11);
211 if (nt->addItem(title,m_c->s_sumE[i]).isFailure()) ATH_MSG_INFO( "Registration of a branch failed in the ntupler..." );
212 if (i<12) sprintf(title,"S%i_Hits",i);
213 else sprintf(title,"FC%i_Hits",i-11);
214 if (nt->addItem(title,m_c->s_hits[i]).isFailure()) ATH_MSG_INFO( "Registration of a branch failed in the ntupler..." );
215 if (i<12) sprintf(title,"S%i_DeltaPhi",i);
216 else sprintf(title,"FC%i_DeltaX",i-11);
217 if (nt->addItem(title,m_c->s_deltaPhi[i]).isFailure()) ATH_MSG_INFO( "Registration of a branch failed in the ntupler..." );
218 if (i<12) sprintf(title,"S%i_SigmaPhi",i);
219 else sprintf(title,"FC%i_SigmaX",i-11);
220 if (nt->addItem(title,m_c->s_sigmaPhi[i]).isFailure()) ATH_MSG_INFO( "Registration of a branch failed in the ntupler..." );
221 if (i<12) sprintf(title,"S%i_DeltaEta",i);
222 else sprintf(title,"FC%i_DeltaY",i-11);
223 if (nt->addItem(title,m_c->s_deltaEta[i]).isFailure()) ATH_MSG_INFO( "Registration of a branch failed in the ntupler..." );
224 if (i<12) sprintf(title,"S%i_SigmaEta",i);
225 else sprintf(title,"FC%i_SigmaY",i-11);
226 if (nt->addItem(title,m_c->s_sigmaEta[i]).isFailure()) ATH_MSG_INFO( "Registration of a branch failed in the ntupler..." );
227 if (i<12) sprintf(title,"S%i_Time",i);
228 else sprintf(title,"FC%i_Time",i-11);
229 if (nt->addItem(title,m_c->s_t00[i]).isFailure()) ATH_MSG_INFO( "Registration of a branch failed in the ntupler..." );
230 if (i<12) sprintf(title,"S%i_WidthX",i);
231 else sprintf(title,"FC%i_WidthX",i-11);
232 if (nt->addItem(title,m_c->s_widthX[i]).isFailure()) ATH_MSG_INFO( "Registration of a branch failed in the ntupler..." );
233 if (i<12) sprintf(title,"S%i_WidthY",i);
234 else sprintf(title,"FC%i_WidthY",i-11);
235 if (nt->addItem(title,m_c->s_widthY[i]).isFailure()) ATH_MSG_INFO( "Registration of a branch failed in the ntupler..." );
236 }
237
238 if (nt->addItem("CPU" , m_c->cpuTime ).isFailure() ||
239 nt->addItem("TrackEnergy" , m_c->Energy ).isFailure() ||
240 nt->addItem("ParticleID" , m_c->PDG ).isFailure() ||
241 nt->addItem("Run#" , m_c->RunNo ).isFailure() ||
242 nt->addItem("Event#" , m_c->EventNo ).isFailure() ||
243 nt->addItem("DepositedEnergy", m_c->E_Deposit ).isFailure() ){
244 ATH_MSG_WARNING( "Registration of some of the ntuple branches failed. No idea what will happen next..." );
245 }
246
247 m_c->cpuTime=0.0;
248 m_c->nt = nt;
249
250 //==~ ~ ~==//
251 // //
252 //------------------------Done with initializations------------------------//
253
254 return StatusCode::SUCCESS;
255}
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
char data[hepevt_bytes_allocation_ATLAS]
Definition HepEvt.cxx:11
const ServiceHandle< StoreGateSvc > & detStore() const
ServiceHandle< IPartPropSvc > m_ppSvc
ServiceHandle< ITHistSvc > m_histSvc
retrieve(aClass, aKey=None)
Definition PyKernel.py:110
TH1F(name, title, nxbins, bins_par2, bins_par3=None, path='', **kwargs)
TFile * file

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Algorithm > >::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.

◆ msg()

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

Definition at line 24 of file AthCommonMsg.h.

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

◆ msgLvl()

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

Definition at line 30 of file AthCommonMsg.h.

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

◆ operator=()

SingleTrackValidation & SingleTrackValidation::operator= ( const SingleTrackValidation & )
private

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Algorithm > >::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.

◆ 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< Algorithm > >::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< Algorithm > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ sysInitialize()

StatusCode AthAlgorithm::sysInitialize ( )
overridevirtualinherited

Override sysInitialize.

Override sysInitialize from the base class.

Loop through all output handles, and if they're WriteCondHandles, automatically register them and this Algorithm with the CondSvc

Scan through all outputHandles, and if they're WriteCondHandles, register them with the CondSvc

Reimplemented from AthCommonDataStore< AthCommonMsg< Algorithm > >.

Reimplemented in AthAnalysisAlgorithm, AthFilterAlgorithm, AthHistogramAlgorithm, and PyAthena::Alg.

Definition at line 66 of file AthAlgorithm.cxx.

66 {
68
69 if (sc.isFailure()) {
70 return sc;
71 }
72 ServiceHandle<ICondSvc> cs("CondSvc",name());
73 for (auto h : outputHandles()) {
74 if (h->isCondition() && h->mode() == Gaudi::DataHandle::Writer) {
75 // do this inside the loop so we don't create the CondSvc until needed
76 if ( cs.retrieve().isFailure() ) {
77 ATH_MSG_WARNING("no CondSvc found: won't autoreg WriteCondHandles");
78 return StatusCode::SUCCESS;
79 }
80 if (cs->regHandle(this,*h).isFailure()) {
81 sc = StatusCode::FAILURE;
82 ATH_MSG_ERROR("unable to register WriteCondHandle " << h->fullKey()
83 << " with CondSvc");
84 }
85 }
86 }
87 return sc;
88}
static Double_t sc
virtual StatusCode sysInitialize() override
Override sysInitialize.
AthCommonDataStore(const std::string &name, T... args)
virtual std::vector< Gaudi::DataHandle * > outputHandles() const override
::StatusCode StatusCode
StatusCode definition for legacy code.

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< Algorithm > >::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< Algorithm > >::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

Member Data Documentation

◆ m_c

Clockwork* SingleTrackValidation::m_c
private

Definition at line 48 of file SingleTrackValidation.h.

◆ m_caloMgrKey

SG::ReadCondHandleKey<CaloDetDescrManager> SingleTrackValidation::m_caloMgrKey
private
Initial value:
{ this
, "CaloDetDescrManager"
, "CaloDetDescrManager"
, "SG Key for CaloDetDescrManager in the Condition Store" }

Definition at line 34 of file SingleTrackValidation.h.

34 { this
35 , "CaloDetDescrManager"
36 , "CaloDetDescrManager"
37 , "SG Key for CaloDetDescrManager in the Condition Store" };

◆ m_detStore

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

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_evtStore

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

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_extendedExtraObjects

DataObjIDColl AthAlgorithm::m_extendedExtraObjects
privateinherited

Definition at line 79 of file AthAlgorithm.h.

◆ m_fieldCacheCondObjInputKey

SG::ReadCondHandleKey<AtlasFieldCacheCondObj> SingleTrackValidation::m_fieldCacheCondObjInputKey
private
Initial value:
{this
, "AtlasFieldCacheCondObj"
, "fieldCondObj"
, "Name of the Magnetic Field conditions object key"}

Definition at line 40 of file SingleTrackValidation.h.

40 {this
41 , "AtlasFieldCacheCondObj"
42 , "fieldCondObj"
43 , "Name of the Magnetic Field conditions object key"};

◆ m_histos

TH1F* SingleTrackValidation::m_histos[162] {}
private

Definition at line 50 of file SingleTrackValidation.h.

50{};

◆ m_histSvc

ServiceHandle<ITHistSvc> SingleTrackValidation::m_histSvc { this, "THistSvc", "THistSvc", "Histogramming svc" }
private

Definition at line 46 of file SingleTrackValidation.h.

46{ this, "THistSvc", "THistSvc", "Histogramming svc" };

◆ m_ppSvc

ServiceHandle<IPartPropSvc> SingleTrackValidation::m_ppSvc {this, "PartPropSvc", "PartPropSvc"}
private

Definition at line 45 of file SingleTrackValidation.h.

45{this, "PartPropSvc", "PartPropSvc"};

◆ m_truthKey

SG::ReadHandleKey<McEventCollection> SingleTrackValidation::m_truthKey { this, "TruthKey", "TruthEvent" }
private

Definition at line 32 of file SingleTrackValidation.h.

32{ this, "TruthKey", "TruthEvent" };

◆ m_varHandleArraysDeclared

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

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

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

Definition at line 398 of file AthCommonDataStore.h.


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