CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
TreatSecondary Class Reference

#include <TreatSecondary.h>

Public Member Functions

void initTrack (const G4Track *trk)
 
const TreatSecondaryoperator= (const TreatSecondary &)=delete
 
std::vector< math::XYZTLorentzVectortracks (const G4Step *step, std::string &procName, int &procID, bool &intr, double &deltaE, std::vector< int > &charges)
 
 TreatSecondary (const edm::ParameterSet &p)
 
 TreatSecondary (const TreatSecondary &)=delete
 
virtual ~TreatSecondary ()
 

Private Attributes

double eTrack
 
int killAfter
 
double minDeltaE
 
int minSec
 
int nHad
 
int nsecL
 
int step
 
G4ProcessTypeEnumeratortypeEnumerator
 
int verbosity
 

Detailed Description

Definition at line 16 of file TreatSecondary.h.

Constructor & Destructor Documentation

◆ TreatSecondary() [1/2]

TreatSecondary::TreatSecondary ( const edm::ParameterSet p)

Definition at line 20 of file TreatSecondary.cc.

References killAfter, minDeltaE, AlCaHLTBitMon_ParallelJobs::p, typeEnumerator, and verbosity.

20  : typeEnumerator(nullptr) {
21  verbosity = p.getUntrackedParameter<int>("Verbosity", 0);
22  killAfter = p.getUntrackedParameter<int>("KillAfter", -1);
23  minDeltaE = p.getUntrackedParameter<double>("MinimumDeltaE", 10.0) * MeV;
24 
25  edm::LogInfo("CheckSecondary") << "Instantiate CheckSecondary with Flag"
26  << " for Killing track after " << killAfter
27  << " hadronic interactions\nDefine inelastic"
28  << " if > 1 seondary or change in KE > " << minDeltaE << " MeV\n";
29 
31 }
G4ProcessTypeEnumerator * typeEnumerator
Log< level::Info, false > LogInfo

◆ TreatSecondary() [2/2]

TreatSecondary::TreatSecondary ( const TreatSecondary )
delete

◆ ~TreatSecondary()

TreatSecondary::~TreatSecondary ( )
virtual

Definition at line 33 of file TreatSecondary.cc.

References typeEnumerator.

33  {
34  if (typeEnumerator)
35  delete typeEnumerator;
36 }
G4ProcessTypeEnumerator * typeEnumerator

Member Function Documentation

◆ initTrack()

void TreatSecondary::initTrack ( const G4Track *  trk)

Definition at line 38 of file TreatSecondary.cc.

References eTrack, LogDebug, nHad, and nsecL.

Referenced by StoreSecondary::update(), and CheckSecondary::update().

38  {
39  step = 0;
40  nsecL = 0;
41  nHad = 0;
42  eTrack = thTk->GetKineticEnergy() / MeV;
43  LogDebug("CheckSecondary") << "TreatSecondary::initTrack:Track: " << thTk->GetTrackID()
44  << " Type: " << thTk->GetDefinition()->GetParticleName() << " KE "
45  << thTk->GetKineticEnergy() / GeV << " GeV p " << thTk->GetMomentum().mag() / GeV
46  << " GeV daughter of particle " << thTk->GetParentID();
47 }
step
Definition: StallMonitor.cc:83
#define LogDebug(id)

◆ operator=()

const TreatSecondary& TreatSecondary::operator= ( const TreatSecondary )
delete

◆ tracks()

std::vector< math::XYZTLorentzVector > TreatSecondary::tracks ( const G4Step *  step,
std::string &  procName,
int &  procID,
bool &  intr,
double &  deltaE,
std::vector< int > &  charges 
)

Definition at line 49 of file TreatSecondary.cc.

References ALCARECOTkAlJpsiMuMu_cff::charge, CosmicGenFilterHelix_cfi::charges, EgHLTOffHistBins_cfi::deltaE, eTrack, mps_fire::i, createfilelist::int, LogDebug, Skims_PA_cff::name, nHad, nsecL, createTree::pp, ValidateTausOnZEEFastSim_cff::proc, G4ProcessTypeEnumerator::processG4Name(), G4ProcessTypeEnumerator::processIdLong(), fileinputsource_cfi::sec, typeEnumerator, and verbosity.

Referenced by StoreSecondary::update(), and CheckSecondary::update().

50  {
51  step++;
52  procid = -1;
53  name = "Unknown";
54  hadrInt = false;
55  deltaE = 0;
56  std::vector<math::XYZTLorentzVector> secondaries;
57  charges.clear();
58 
59  if (aStep != nullptr) {
60  const G4TrackVector *tkV = aStep->GetSecondary();
61  G4Track *thTk = aStep->GetTrack();
62  const G4StepPoint *preStepPoint = aStep->GetPreStepPoint();
63  const G4StepPoint *postStepPoint = aStep->GetPostStepPoint();
64  double eTrackNew = thTk->GetKineticEnergy() / MeV;
65  deltaE = eTrack - eTrackNew;
66  eTrack = eTrackNew;
67  if (tkV != nullptr && postStepPoint != nullptr) {
68  int nsc = (*tkV).size();
69  const G4VProcess *proc = postStepPoint->GetProcessDefinedStep();
70  if (proc != nullptr) {
71  G4ProcessType type = proc->GetProcessType();
73  name = proc->GetProcessName();
74  int sec = nsc - nsecL;
75  LogDebug("CheckSecondary") << sec << " secondaries in step " << thTk->GetCurrentStepNumber() << " of track "
76  << thTk->GetTrackID() << " from " << name << " of type " << type << " ID " << procid
77  << " (" << typeEnumerator->processG4Name(procid) << ")";
78 
79  // hadronic interaction
80  if (procid >= 121 && procid <= 151) {
81  LogDebug("CheckSecondary") << "Hadronic Interaction " << nHad << " of Type " << procid << " with " << sec
82  << " secondaries from process " << proc->GetProcessName() << " Delta E " << deltaE
83  << " Flag " << hadrInt;
84  math::XYZTLorentzVector secondary;
85  for (int i = nsecL; i < nsc; ++i) {
86  G4Track *tk = (*tkV)[i];
87  G4ThreeVector pp = tk->GetMomentum();
88  double ee = tk->GetTotalEnergy();
89  secondary = math::XYZTLorentzVector(pp.x(), pp.y(), pp.z(), ee);
90  secondaries.push_back(secondary);
91  int charge = (int)(tk->GetDefinition()->GetPDGCharge());
92  charges.push_back(charge);
93  }
94  if (verbosity > 0) {
95  for (int i = nsecL; i < nsc; i++) {
96  G4Track *tk = (*tkV)[i];
97  LogDebug("CheckSecondary") << "Secondary: " << sec << " ID " << tk->GetTrackID() << " Status "
98  << tk->GetTrackStatus() << " Particle "
99  << tk->GetDefinition()->GetParticleName() << " Position " << tk->GetPosition()
100  << " KE " << tk->GetKineticEnergy() << " Time " << tk->GetGlobalTime();
101  }
102  }
103  }
104  nsecL = nsc;
105  }
106  }
107  if (verbosity > 1) {
108  LogDebug("CheckSecondary") << "Track: " << thTk->GetTrackID() << " Status " << thTk->GetTrackStatus()
109  << " Particle " << thTk->GetDefinition()->GetParticleName() << " at "
110  << preStepPoint->GetPosition() << " Step: " << step << " KE "
111  << thTk->GetKineticEnergy() / GeV << " GeV; p " << thTk->GetMomentum().mag() / GeV
112  << " GeV/c; Step Length " << aStep->GetStepLength() << " Energy Deposit "
113  << aStep->GetTotalEnergyDeposit() / MeV << " MeV; Interaction " << hadrInt;
114  }
115  }
116  return secondaries;
117 }
G4ProcessTypeEnumerator * typeEnumerator
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int processIdLong(const G4VProcess *p) const
std::string processG4Name(int) const
step
Definition: StallMonitor.cc:83
charges
only generated particles of these IDs are considered
#define LogDebug(id)

Member Data Documentation

◆ eTrack

double TreatSecondary::eTrack
private

Definition at line 29 of file TreatSecondary.h.

Referenced by initTrack(), and tracks().

◆ killAfter

int TreatSecondary::killAfter
private

Definition at line 28 of file TreatSecondary.h.

Referenced by TreatSecondary().

◆ minDeltaE

double TreatSecondary::minDeltaE
private

Definition at line 29 of file TreatSecondary.h.

Referenced by TreatSecondary().

◆ minSec

int TreatSecondary::minSec
private

Definition at line 28 of file TreatSecondary.h.

◆ nHad

int TreatSecondary::nHad
private

Definition at line 31 of file TreatSecondary.h.

Referenced by initTrack(), and tracks().

◆ nsecL

int TreatSecondary::nsecL
private

Definition at line 31 of file TreatSecondary.h.

Referenced by initTrack(), and tracks().

◆ step

int TreatSecondary::step
private

Definition at line 31 of file TreatSecondary.h.

◆ typeEnumerator

G4ProcessTypeEnumerator* TreatSecondary::typeEnumerator
private

Definition at line 30 of file TreatSecondary.h.

Referenced by tracks(), TreatSecondary(), and ~TreatSecondary().

◆ verbosity

int TreatSecondary::verbosity
private

Definition at line 28 of file TreatSecondary.h.

Referenced by tracks(), and TreatSecondary().