Public Member Functions | Private Member Functions | Private Attributes | List of all members
gar::garg4::MuNuclearSplittingProcess Class Reference

#include <MuNuclearSplittingProcess.h>

Inheritance diagram for gar::garg4::MuNuclearSplittingProcess:

Public Member Functions

 MuNuclearSplittingProcess ()
 
 ~MuNuclearSplittingProcess ()
 
void SetNSplit (G4int nTrx)
 
void SetIsActive (G4bool doIt)
 

Private Member Functions

G4VParticleChange * PostStepDoIt (const G4Track &track, const G4Step &step)
 

Private Attributes

G4int fNSplit
 
G4bool fActive
 

Detailed Description

Definition at line 31 of file MuNuclearSplittingProcess.h.

Constructor & Destructor Documentation

gar::garg4::MuNuclearSplittingProcess::MuNuclearSplittingProcess ( )
inline

Definition at line 34 of file MuNuclearSplittingProcess.h.

34 {};
gar::garg4::MuNuclearSplittingProcess::~MuNuclearSplittingProcess ( )
inline

Definition at line 35 of file MuNuclearSplittingProcess.h.

35 {};

Member Function Documentation

G4VParticleChange * gar::garg4::MuNuclearSplittingProcess::PostStepDoIt ( const G4Track &  track,
const G4Step &  step 
)
private

Definition at line 40 of file MuNuclearSplittingProcess.cxx.

41  {
42 
43  G4VParticleChange* particleChange = new G4VParticleChange();
44  G4double weight = track.GetWeight()/fNSplit;
45  std::vector<G4Track*> secondaries; // Secondary store
46 
47  // Loop over PostStepDoIt method to generate multiple secondaries.
48  // The point is for each value of ii up to Nsplit we want to re-toss
49  // all secondaries for the muNucl process. Each toss gives back numSec
50  // secondaries, which will be a different sample of species/momenta each time.
51  for (unsigned int ii=0; ii<(unsigned int)fNSplit; ii++) {
52  particleChange = pRegProcess->PostStepDoIt(track, step);
53  if (!particleChange)
54  throw std::runtime_error("MuNuclearSplittingProcess::PostStepDoIt(): no particle change");
55  G4int j(0);
56  G4int numSec(particleChange->GetNumberOfSecondaries());
57  for (j=0; j<numSec; j++)
58  {
59  G4Track* newSec = new G4Track(*(particleChange->GetSecondary(j)));
60  G4String pdgstr = newSec->GetParticleDefinition()->GetParticleName();
61  G4double ke = newSec->GetKineticEnergy()/CLHEP::GeV;
62  G4int pdg = newSec->GetParticleDefinition()->GetPDGEncoding();
63  if (abs(pdg)==310 ||abs(pdg)==311 || abs(pdg)==3122 || abs(pdg)==2112)
64  {
65  // std::cout << "MuNuclSplProc: Creating " << pdgstr << " of Kinetic Energy " << ke << " GeV. numSec is " << j << std::endl;
66 
67  if (pdg==G4KaonZeroShort::KaonZeroShort()->GetPDGEncoding()&&G4UniformRand()<0.50)
68  {
69  pdg = G4KaonZeroLong::KaonZeroLong()->GetPDGEncoding();
70  pdgstr = G4KaonZeroLong::KaonZeroLong()->GetParticleName();
71  G4LorentzVector pK0L(newSec->GetMomentum(),TMath::Sqrt(TMath::Power(G4KaonZeroLong::KaonZeroLong()->GetPDGMass(),2)+TMath::Power(newSec->GetMomentum().mag(),2)));
72  G4DynamicParticle *newK0L = new G4DynamicParticle(G4KaonZeroLong::KaonZeroLong(),pK0L);
73 
74  G4Track* newSecK0L = new G4Track(newK0L,track.GetGlobalTime(),track.GetPosition());
75  secondaries.push_back(newSecK0L);
76  }
77  else
78  {
79  secondaries.push_back(newSec);
80  }
81 
82  if (abs(pdg)==130 ||abs(pdg)==310 ||abs(pdg)==311 || abs(pdg)==3122)
83  {
84  // WrappedmuNuclear always produces K0s if it produces a K0 at all.
85  // Let's make half of these K0Ls.
86  std::cout << "MuNuclSplProc: Creating " << pdgstr << " of Kinetic Energy " << ke << " GeV. numSec is " << j << std::endl;
87  }
88 
89  }
90  }
91  }
92 
93  particleChange->SetNumberOfSecondaries(secondaries.size());
94  particleChange->SetSecondaryWeightByProcess(true);
95  //int numSecAdd = secondaries.size();
96  std::vector<G4Track*>::iterator iter = secondaries.begin(); // Add all secondaries
97  while (iter != secondaries.end()) {
98  G4Track* myTrack = *iter;
99  G4String pdgstr = myTrack->GetParticleDefinition()->GetParticleName();
100  //G4double ke = myTrack->GetKineticEnergy()/GeV;
101  G4int pdg = myTrack->GetParticleDefinition()->GetPDGEncoding();
102  if (abs(pdg)==130 ||abs(pdg)==310 ||abs(pdg)==311 || abs(pdg)==3122 || abs(pdg)==2112)
103  // std::cout << "MuNuclSplProc: Adding " << pdgstr << " of Kinetic Energy " << ke << " GeV." << std::endl;
104  // Will ask for PdgCode and only Add K0s.
105  myTrack->SetWeight(weight);
106  particleChange->AddSecondary(myTrack);
107  iter++;
108  }
109  return particleChange;
110  }
intermediate_table::iterator iterator
weight
Definition: test.py:257
T abs(T value)
static constexpr double GeV
Definition: Units.h:28
QTextStream & endl(QTextStream &s)
void gar::garg4::MuNuclearSplittingProcess::SetIsActive ( G4bool  doIt)
inline
void gar::garg4::MuNuclearSplittingProcess::SetNSplit ( G4int  nTrx)
inline

Member Data Documentation

G4bool gar::garg4::MuNuclearSplittingProcess::fActive
private

Definition at line 43 of file MuNuclearSplittingProcess.h.

G4int gar::garg4::MuNuclearSplittingProcess::fNSplit
private

Definition at line 38 of file MuNuclearSplittingProcess.h.


The documentation for this class was generated from the following files: