#include <TreatSecondary.h>
Public Member Functions | |
void | initTrack (const G4Track *trk) |
std::vector < math::XYZTLorentzVector > | tracks (const G4Step *step, std::string &procName, int &procID, bool &intr, double &deltaE, std::vector< int > &charges) |
TreatSecondary (const edm::ParameterSet &p) | |
virtual | ~TreatSecondary () |
Private Member Functions | |
const TreatSecondary & | operator= (const TreatSecondary &) |
TreatSecondary (const TreatSecondary &) | |
Private Attributes | |
double | eTrack |
int | killAfter |
double | minDeltaE |
int | minSec |
int | nHad |
int | nsecL |
int | step |
G4ProcessTypeEnumerator * | typeEnumerator |
int | verbosity |
Definition at line 16 of file TreatSecondary.h.
TreatSecondary::TreatSecondary | ( | const edm::ParameterSet & | p | ) |
Definition at line 16 of file TreatSecondary.cc.
References edm::ParameterSet::getUntrackedParameter(), killAfter, minDeltaE, typeEnumerator, and verbosity.
: typeEnumerator(0) { verbosity = p.getUntrackedParameter<int>("Verbosity",0); killAfter = p.getUntrackedParameter<int>("KillAfter", -1); minDeltaE = p.getUntrackedParameter<double>("MinimumDeltaE", 10.0)*MeV; edm::LogInfo("CheckSecondary") << "Instantiate CheckSecondary with Flag" << " for Killing track after "<< killAfter << " hadronic interactions\nDefine inelastic" << " if > 1 seondary or change in KE > " << minDeltaE << " MeV\n"; typeEnumerator = new G4ProcessTypeEnumerator(); }
TreatSecondary::~TreatSecondary | ( | ) | [virtual] |
Definition at line 31 of file TreatSecondary.cc.
References typeEnumerator.
{ if (typeEnumerator) delete typeEnumerator; }
TreatSecondary::TreatSecondary | ( | const TreatSecondary & | ) | [private] |
void TreatSecondary::initTrack | ( | const G4Track * | trk | ) |
Definition at line 35 of file TreatSecondary.cc.
References eTrack, LogDebug, nHad, nsecL, and step.
Referenced by StoreSecondary::update(), and CheckSecondary::update().
{ step = 0; nsecL = 0; nHad = 0; eTrack = thTk->GetKineticEnergy()/MeV; LogDebug("CheckSecondary") << "TreatSecondary::initTrack:Track: " << thTk->GetTrackID() << " Type: " << thTk->GetDefinition()->GetParticleName() << " KE " << thTk->GetKineticEnergy()/GeV << " GeV p " << thTk->GetMomentum().mag()/GeV << " GeV daughter of particle " << thTk->GetParentID(); }
const TreatSecondary& TreatSecondary::operator= | ( | const TreatSecondary & | ) | [private] |
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 50 of file TreatSecondary.cc.
References DeDxDiscriminatorTools::charge(), eTrack, i, killAfter, LogDebug, minDeltaE, nHad, nsecL, NULL, createTree::pp, proc, G4ProcessTypeEnumerator::processG4Name(), G4ProcessTypeEnumerator::processIdLong(), step, typeEnumerator, and verbosity.
Referenced by StoreSecondary::update(), and CheckSecondary::update().
{ step++; procid = -1; name = "Unknown"; hadrInt = false; deltaE = 0; std::vector<math::XYZTLorentzVector> secondaries; charges.clear(); if (aStep != NULL) { G4TrackVector* tkV = const_cast<G4TrackVector*>(aStep->GetSecondary()); G4Track* thTk = aStep->GetTrack(); const G4StepPoint* preStepPoint = aStep->GetPreStepPoint(); const G4StepPoint* postStepPoint = aStep->GetPostStepPoint(); double eTrackNew = thTk->GetKineticEnergy()/MeV; deltaE = eTrack-eTrackNew; eTrack = eTrackNew; if (tkV != 0) { int nsc = (*tkV).size(); const G4VProcess* proc = 0; if (postStepPoint) proc = postStepPoint->GetProcessDefinedStep(); procid = typeEnumerator->processIdLong(proc); G4ProcessType type = fNotDefined; if (proc) { type = proc->GetProcessType(); name = proc->GetProcessName(); } int sec = nsc - nsecL; LogDebug("CheckSecondary") << sec << " secondaries in step " << thTk->GetCurrentStepNumber() << " of track " << thTk->GetTrackID() << " from " << name << " of type " << type << " ID " << procid << " (" << typeEnumerator->processG4Name(procid) << ")"; G4TrackStatus state = thTk->GetTrackStatus(); if (state == fAlive || state == fStopButAlive) sec++; if (type == fHadronic || type == fPhotolepton_hadron || type == fDecay) { if (deltaE > minDeltaE || sec > 1) { hadrInt = true; nHad++; } LogDebug("CheckSecondary") << "Hadronic Interaction " << nHad << " of Type " << type << " with " << sec << " secondaries from process " << proc->GetProcessName() << " Delta E " << deltaE << " Flag " << hadrInt; } if (hadrInt) { math::XYZTLorentzVector secondary; if (state == fAlive || state == fStopButAlive) { G4ThreeVector pp = postStepPoint->GetMomentum(); double ee = postStepPoint->GetTotalEnergy(); secondary = math::XYZTLorentzVector(pp.x(),pp.y(),pp.z(),ee); secondaries.push_back(secondary); int charge = (int)(postStepPoint->GetCharge()); charges.push_back(charge); } for (int i=nsecL; i<nsc; i++) { G4Track* tk = (*tkV)[i]; G4ThreeVector pp = tk->GetMomentum(); double ee = tk->GetTotalEnergy(); secondary = math::XYZTLorentzVector(pp.x(),pp.y(),pp.z(),ee); secondaries.push_back(secondary); int charge = (int)(tk->GetDefinition()->GetPDGCharge()); charges.push_back(charge); } } if (killAfter >= 0 && nHad >= killAfter) thTk->SetTrackStatus(fStopAndKill); if (verbosity > 0) { sec = 0; if (state == fAlive || state == fStopButAlive) { sec++; LogDebug("CheckSecondary") << "Secondary: " << sec << " ID " << thTk->GetTrackID() << " Status " << thTk->GetTrackStatus() << " Particle " <<thTk->GetDefinition()->GetParticleName() << " Position " << postStepPoint->GetPosition() << " KE " << postStepPoint->GetKineticEnergy() << " Time " << postStepPoint->GetGlobalTime(); } for (int i=nsecL; i<nsc; i++) { sec++; G4Track* tk = (*tkV)[i]; LogDebug("CheckSecondary") << "Secondary: " << sec << " ID " << tk->GetTrackID() << " Status " << tk->GetTrackStatus() << " Particle " << tk->GetDefinition()->GetParticleName() << " Position " << tk->GetPosition() << " KE " << tk->GetKineticEnergy() << " Time " << tk->GetGlobalTime(); } } nsecL = nsc; } if (verbosity > 1) LogDebug("CheckSecondary") << "Track: " << thTk->GetTrackID() << " Status " << thTk->GetTrackStatus() << " Particle " << thTk->GetDefinition()->GetParticleName() << " at " << preStepPoint->GetPosition() << " Step: " << step << " KE " << thTk->GetKineticEnergy()/GeV << " GeV; p " << thTk->GetMomentum().mag()/GeV << " GeV/c; Step Length " << aStep->GetStepLength()<< " Energy Deposit " << aStep->GetTotalEnergyDeposit()/MeV << " MeV; Interaction " << hadrInt; } return secondaries; }
double TreatSecondary::eTrack [private] |
Definition at line 35 of file TreatSecondary.h.
Referenced by initTrack(), and tracks().
int TreatSecondary::killAfter [private] |
Definition at line 34 of file TreatSecondary.h.
Referenced by tracks(), and TreatSecondary().
double TreatSecondary::minDeltaE [private] |
Definition at line 35 of file TreatSecondary.h.
Referenced by tracks(), and TreatSecondary().
int TreatSecondary::minSec [private] |
Definition at line 34 of file TreatSecondary.h.
int TreatSecondary::nHad [private] |
Definition at line 37 of file TreatSecondary.h.
Referenced by initTrack(), and tracks().
int TreatSecondary::nsecL [private] |
Definition at line 37 of file TreatSecondary.h.
Referenced by initTrack(), and tracks().
int TreatSecondary::step [private] |
Definition at line 37 of file TreatSecondary.h.
Referenced by initTrack(), and tracks().
Definition at line 36 of file TreatSecondary.h.
Referenced by tracks(), TreatSecondary(), and ~TreatSecondary().
int TreatSecondary::verbosity [private] |
Definition at line 34 of file TreatSecondary.h.
Referenced by tracks(), and TreatSecondary().