LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
larg4::MuNuclearSplittingProcess Class Reference

#include "MuNuclearSplittingProcess.h"

Inheritance diagram for larg4::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 29 of file MuNuclearSplittingProcess.h.

Constructor & Destructor Documentation

larg4::MuNuclearSplittingProcess::MuNuclearSplittingProcess ( )
inline

Definition at line 32 of file MuNuclearSplittingProcess.h.

32 {};
larg4::MuNuclearSplittingProcess::~MuNuclearSplittingProcess ( )
inline

Definition at line 33 of file MuNuclearSplittingProcess.h.

33 {};

Member Function Documentation

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

Definition at line 39 of file MuNuclearSplittingProcess.cxx.

References fNSplit, and weight.

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

Definition at line 36 of file MuNuclearSplittingProcess.h.

References fActive, and fNSplit.

void larg4::MuNuclearSplittingProcess::SetNSplit ( G4int  nTrx)
inline

Definition at line 35 of file MuNuclearSplittingProcess.h.

References fNSplit.

Member Data Documentation

G4bool larg4::MuNuclearSplittingProcess::fActive
private

Definition at line 41 of file MuNuclearSplittingProcess.h.

Referenced by SetIsActive().

G4int larg4::MuNuclearSplittingProcess::fNSplit
private

Definition at line 36 of file MuNuclearSplittingProcess.h.

Referenced by PostStepDoIt(), SetIsActive(), and SetNSplit().


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