ATLAS Offline Software
Functions
particleGamma.h File Reference
#include <iostream>
#include <sstream>
#include <vector>
#include "CLHEP/Units/PhysicalConstants.h"
#include "GeneratorObjects/HepMcParticleLink.h"
#include "MdtCalibData/MdtFullCalibData.h"
#include "MdtCalibData/MdtTubeCalibContainer.h"
#include "MuonDigitContainer/MdtDigitContainer.h"
#include "MuonReadoutGeometry/MdtReadoutElement.h"
#include "MuonReadoutGeometry/MuonDetectorManager.h"
#include "MuonSimData/MuonSimData.h"
#include "MuonSimData/MuonSimDataCollection.h"
#include "MuonSimEvent/MDTSimHitCollection.h"
#include "PathResolver/PathResolver.h"
#include "PileUpTools/PileUpMergeSvc.h"
#include "StoreGate/StoreGateSvc.h"
#include "TrkDetDescrUtils/GeometryStatics.h"
#include "AtlasHepMC/GenParticle.h"
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

double particleGamma (const EventContext &ctx, const MDTSimHit &hit, unsigned short eventId=0)
 

Function Documentation

◆ particleGamma()

double particleGamma ( const EventContext &  ctx,
const MDTSimHit hit,
unsigned short  eventId = 0 
)

Definition at line 36 of file particleGamma.h.

36  {
37  double QGamma = -9999.;
38  const HepMcParticleLink trkParticle = HepMcParticleLink::getRedirectedLink(hit.particleLink(), eventId, ctx); // This link should now correctly resolve to the TruthEvent McEventCollection in the main StoreGateSvc.
39  HepMC::ConstGenParticlePtr genParticle = trkParticle.cptr();
40  if (genParticle) {
41  const int particleEncoding = genParticle->pdg_id();
42  if ((((int)(std::abs(particleEncoding) / 10000000) == 1) && ((int)(std::abs(particleEncoding) / 100000) == 100)) ||
43  (((int)(std::abs(particleEncoding) / 10000000) == 2) && ((int)(std::abs(particleEncoding) / 100000) == 200))) { // TODO Use functions from TruthUtils/AtlasPID.h instead?
44  const double QPx = genParticle->momentum().px();
45  const double QPy = genParticle->momentum().py();
46  const double QPz = genParticle->momentum().pz();
47  const double QE = genParticle->momentum().e();
48  const double QM2 = pow(QE, 2) - pow(QPx, 2) - pow(QPy, 2) - pow(QPz, 2);
49  if (QM2 >= 0.) {
50  const double QM = std::sqrt(QM2);
51  if (QM > 0.) {
52  QGamma = QE / QM;
53  }
54  }
55  }
56  }
57  return QGamma;
58 }
conifer::pow
constexpr int pow(int x)
Definition: conifer.h:20
HepMC::ConstGenParticlePtr
const GenParticle * ConstGenParticlePtr
Definition: GenParticle.h:38
MDTSimHit::particleLink
const HepMcParticleLink & particleLink() const
Definition: MDTSimHit.h:100