CMS 3D CMS Logo

TreatSecondary.cc
Go to the documentation of this file.
3 
5 
6 #include <CLHEP/Units/GlobalPhysicalConstants.h>
7 #include <CLHEP/Units/SystemOfUnits.h>
8 #include "G4HCofThisEvent.hh"
9 #include "G4Step.hh"
10 #include "G4Track.hh"
11 #include "G4VProcess.hh"
12 
13 #include <cmath>
14 #include <iomanip>
15 #include <iostream>
16 
17 using CLHEP::GeV;
18 using CLHEP::MeV;
19 
20 TreatSecondary::TreatSecondary(const edm::ParameterSet &p) : 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 }
32 
34  if (typeEnumerator)
35  delete typeEnumerator;
36 }
37 
38 void TreatSecondary::initTrack(const G4Track *thTk) {
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 }
48 
49 std::vector<math::XYZTLorentzVector> TreatSecondary::tracks(
50  const G4Step *aStep, std::string &name, int &procid, bool &hadrInt, double &deltaE, std::vector<int> &charges) {
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 }
virtual ~TreatSecondary()
std::vector< math::XYZTLorentzVector > tracks(const G4Step *step, std::string &procName, int &procID, bool &intr, double &deltaE, std::vector< int > &charges)
G4ProcessTypeEnumerator * typeEnumerator
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
TreatSecondary(const edm::ParameterSet &p)
int processIdLong(const G4VProcess *p) const
Log< level::Info, false > LogInfo
void initTrack(const G4Track *trk)
std::string processG4Name(int) const
step
Definition: StallMonitor.cc:83
charges
only generated particles of these IDs are considered
#define LogDebug(id)