CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TreatSecondary.cc
Go to the documentation of this file.
2 
4 
5 #include "G4Step.hh"
6 #include "G4Track.hh"
7 #include "G4VProcess.hh"
8 #include "G4HCofThisEvent.hh"
9 #include "CLHEP/Units/GlobalSystemOfUnits.h"
10 #include "CLHEP/Units/GlobalPhysicalConstants.h"
11 
12 #include <cmath>
13 #include <iostream>
14 #include <iomanip>
15 
17 
18  verbosity = p.getUntrackedParameter<int>("Verbosity",0);
19  killAfter = p.getUntrackedParameter<int>("KillAfter", -1);
20  minDeltaE = p.getUntrackedParameter<double>("MinimumDeltaE", 10.0)*MeV;
21 
22  edm::LogInfo("CheckSecondary") << "Instantiate CheckSecondary with Flag"
23  << " for Killing track after "<< killAfter
24  << " hadronic interactions\nDefine inelastic"
25  << " if > 1 seondary or change in KE > "
26  << minDeltaE << " MeV\n";
27 
29 }
30 
32  if (typeEnumerator) delete typeEnumerator;
33 }
34 
35 void TreatSecondary::initTrack(const G4Track * thTk) {
36 
37  step = 0;
38  nsecL = 0;
39  nHad = 0;
40  eTrack = thTk->GetKineticEnergy()/MeV;
41  LogDebug("CheckSecondary") << "TreatSecondary::initTrack:Track: "
42  << thTk->GetTrackID() << " Type: "
43  << thTk->GetDefinition()->GetParticleName()
44  << " KE " << thTk->GetKineticEnergy()/GeV
45  << " GeV p " << thTk->GetMomentum().mag()/GeV
46  << " GeV daughter of particle "
47  << thTk->GetParentID();
48 }
49 
50 std::vector<math::XYZTLorentzVector> TreatSecondary::tracks(const G4Step*aStep,
51  std::string & name,
52  int & procid,
53  bool & hadrInt,
54  double & deltaE,
55  std::vector<int> & charges) {
56 
57  step++;
58  procid = -1;
59  name = "Unknown";
60  hadrInt = false;
61  deltaE = 0;
62  std::vector<math::XYZTLorentzVector> secondaries;
63  charges.clear();
64 
65  if (aStep != NULL) {
66  G4TrackVector* tkV = const_cast<G4TrackVector*>(aStep->GetSecondary());
67  G4Track* thTk = aStep->GetTrack();
68  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
69  const G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
70  double eTrackNew = thTk->GetKineticEnergy()/MeV;
71  deltaE = eTrack-eTrackNew;
72  eTrack = eTrackNew;
73  if (tkV != 0) {
74  int nsc = (*tkV).size();
75  const G4VProcess* proc = 0;
76  if (postStepPoint) proc = postStepPoint->GetProcessDefinedStep();
77  procid = typeEnumerator->processIdLong(proc);
78  G4ProcessType type = fNotDefined;
79  if (proc) {
80  type = proc->GetProcessType();
81  name = proc->GetProcessName();
82  }
83  int sec = nsc - nsecL;
84  LogDebug("CheckSecondary") << sec << " secondaries in step "
85  << thTk->GetCurrentStepNumber()
86  << " of track " << thTk->GetTrackID()
87  << " from " << name << " of type "
88  << type << " ID " << procid << " ("
89  << typeEnumerator->processG4Name(procid)
90  << ")";
91 
92  G4TrackStatus state = thTk->GetTrackStatus();
93  if (state == fAlive || state == fStopButAlive) sec++;
94 
95  if (type == fHadronic || type == fPhotolepton_hadron || type == fDecay) {
96  if (deltaE > minDeltaE || sec > 1) {
97  hadrInt = true;
98  nHad++;
99  }
100  LogDebug("CheckSecondary") << "Hadronic Interaction " << nHad
101  << " of Type " << type << " with "
102  << sec << " secondaries from process "
103  << proc->GetProcessName() << " Delta E "
104  << deltaE << " Flag " << hadrInt;
105  }
106  if (hadrInt) {
107  math::XYZTLorentzVector secondary;
108  if (state == fAlive || state == fStopButAlive) {
109  G4ThreeVector pp = postStepPoint->GetMomentum();
110  double ee = postStepPoint->GetTotalEnergy();
111  secondary = math::XYZTLorentzVector(pp.x(),pp.y(),pp.z(),ee);
112  secondaries.push_back(secondary);
113  int charge = (int)(postStepPoint->GetCharge());
114  charges.push_back(charge);
115  }
116  for (int i=nsecL; i<nsc; i++) {
117  G4Track* tk = (*tkV)[i];
118  G4ThreeVector pp = tk->GetMomentum();
119  double ee = tk->GetTotalEnergy();
120  secondary = math::XYZTLorentzVector(pp.x(),pp.y(),pp.z(),ee);
121  secondaries.push_back(secondary);
122  int charge = (int)(tk->GetDefinition()->GetPDGCharge());
123  charges.push_back(charge);
124  }
125  }
126 
127  if (killAfter >= 0 && nHad >= killAfter)
128  thTk->SetTrackStatus(fStopAndKill);
129 
130  if (verbosity > 0) {
131  sec = 0;
132  if (state == fAlive || state == fStopButAlive) {
133  sec++;
134  LogDebug("CheckSecondary") << "Secondary: " << sec << " ID "
135  << thTk->GetTrackID() << " Status "
136  << thTk->GetTrackStatus() << " Particle "
137  <<thTk->GetDefinition()->GetParticleName()
138  << " Position "
139  << postStepPoint->GetPosition() << " KE "
140  << postStepPoint->GetKineticEnergy()
141  << " Time "
142  << postStepPoint->GetGlobalTime();
143  }
144  for (int i=nsecL; i<nsc; i++) {
145  sec++;
146  G4Track* tk = (*tkV)[i];
147  LogDebug("CheckSecondary") << "Secondary: " << sec << " ID "
148  << tk->GetTrackID() << " Status "
149  << tk->GetTrackStatus() << " Particle "
150  << tk->GetDefinition()->GetParticleName()
151  << " Position " << tk->GetPosition()
152  << " KE " << tk->GetKineticEnergy()
153  << " Time " << tk->GetGlobalTime();
154  }
155  }
156  nsecL = nsc;
157  }
158 
159  if (verbosity > 1)
160  LogDebug("CheckSecondary") << "Track: " << thTk->GetTrackID()
161  << " Status " << thTk->GetTrackStatus()
162  << " Particle "
163  << thTk->GetDefinition()->GetParticleName()
164  << " at " << preStepPoint->GetPosition()
165  << " Step: " << step << " KE "
166  << thTk->GetKineticEnergy()/GeV << " GeV; p "
167  << thTk->GetMomentum().mag()/GeV
168  << " GeV/c; Step Length "
169  << aStep->GetStepLength()<< " Energy Deposit "
170  << aStep->GetTotalEnergyDeposit()/MeV
171  << " MeV; Interaction " << hadrInt;
172  }
173  return secondaries;
174 }
#define LogDebug(id)
type
Definition: HCALResponse.h:22
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
tuple pp
Definition: createTree.py:15
TrainProcessor *const proc
Definition: MVATrainer.cc:101
virtual ~TreatSecondary()
#define NULL
Definition: scimark2.h:8
std::vector< math::XYZTLorentzVector > tracks(const G4Step *step, std::string &procName, int &procID, bool &intr, double &deltaE, std::vector< int > &charges)
G4ProcessTypeEnumerator * typeEnumerator
double charge(const std::vector< uint8_t > &Ampls)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
TreatSecondary(const edm::ParameterSet &p)
int processIdLong(const G4VProcess *)
void initTrack(const G4Track *trk)
char state
Definition: procUtils.cc:75
std::vector< std::string > processG4Name(std::string)