All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MuNuclearSplittingProcess.cxx
Go to the documentation of this file.
1 //
2 // Code adapted from EventBiasing from Jane Tinslay, SLAC. 2008. Now at Cisco,
3 // I think.
4 //
5 // This is a variance reduction technique.
6 //
7 // The point is to create Nsplit secondaries with each muNucl process, so
8 // that we may, e.g., create many, otherwise-rare, K0s which
9 // may charge exchange into pernicious K+s which then mimic proton dk.
10 //
11 // echurch@fnal.gov, May, 2011.
12 //
13 
14 
16 
17 #include <stdexcept> // std::runtime_error
18 
19 #include "Geant4/G4VParticleChange.hh"
20 
21 #include "Geant4/G4Track.hh"
22 #include "Geant4/G4ParticleDefinition.hh"
23 #include "Geant4/Randomize.hh"
24 #include "Geant4/G4KaonZeroLong.hh"
25 #include "Geant4/G4DynamicParticle.hh"
26 #include "Geant4/G4KaonZeroShort.hh"
27 #include "Geant4/G4LorentzVector.hh"
28 #include "Geant4/G4String.hh"
29 #include "Geant4/G4ThreeVector.hh"
30 #include "Geant4/G4VProcess.hh"
31 
32 #include <TMath.h>
33 
34 namespace larg4 {
35 
36 G4VParticleChange* MuNuclearSplittingProcess::PostStepDoIt(const G4Track& track, const G4Step& step)
37 {
38 
39  G4VParticleChange* particleChange = new G4VParticleChange();
40  G4double weight = track.GetWeight()/fNSplit;
41  std::vector<G4Track*> secondaries; // Secondary store
42 
43  // Loop over PostStepDoIt method to generate multiple secondaries.
44  // The point is for each value of ii up to Nsplit we want to re-toss
45  // all secondaries for the muNucl process. Each toss gives back numSec
46  // secondaries, which will be a different sample of species/momenta each time.
47  for (unsigned int ii=0; ii<(unsigned int)fNSplit; ii++) {
48  particleChange = pRegProcess->PostStepDoIt(track, step);
49  if (!particleChange)
50  throw std::runtime_error("MuNuclearSplittingProcess::PostStepDoIt(): no particle change");
51  G4int j(0);
52  G4int numSec(particleChange->GetNumberOfSecondaries());
53  for (j=0; j<numSec; j++)
54  {
55  G4Track* newSec = new G4Track(*(particleChange->GetSecondary(j)));
56  G4String pdgstr = newSec->GetParticleDefinition()->GetParticleName();
57  G4double ke = newSec->GetKineticEnergy()/CLHEP::GeV;
58  G4int pdg = newSec->GetParticleDefinition()->GetPDGEncoding();
59  if (abs(pdg)==310 ||abs(pdg)==311 || abs(pdg)==3122 || abs(pdg)==2112)
60  {
61  // std::cout << "MuNuclSplProc: Creating " << pdgstr << " of Kinetic Energy " << ke << " GeV. numSec is " << j << std::endl;
62 
63  if (pdg==G4KaonZeroShort::KaonZeroShort()->GetPDGEncoding()&&G4UniformRand()<0.50)
64  {
65  pdg = G4KaonZeroLong::KaonZeroLong()->GetPDGEncoding();
66  pdgstr = G4KaonZeroLong::KaonZeroLong()->GetParticleName();
67  G4LorentzVector pK0L(newSec->GetMomentum(),TMath::Sqrt(TMath::Power(G4KaonZeroLong::KaonZeroLong()->GetPDGMass(),2)+TMath::Power(newSec->GetMomentum().mag(),2)));
68  G4DynamicParticle *newK0L = new G4DynamicParticle(G4KaonZeroLong::KaonZeroLong(),pK0L);
69 
70  G4Track* newSecK0L = new G4Track(newK0L,track.GetGlobalTime(),track.GetPosition());
71  secondaries.push_back(newSecK0L);
72  }
73  else
74  {
75  secondaries.push_back(newSec);
76  }
77 
78  if (abs(pdg)==130 ||abs(pdg)==310 ||abs(pdg)==311 || abs(pdg)==3122)
79  {
80  // WrappedmuNuclear always produces K0s if it produces a K0 at all.
81  // Let's make half of these K0Ls.
82  std::cout << "MuNuclSplProc: Creating " << pdgstr << " of Kinetic Energy " << ke << " GeV. numSec is " << j << std::endl;
83  }
84 
85  }
86  }
87  }
88 
89  particleChange->SetNumberOfSecondaries(secondaries.size());
90  particleChange->SetSecondaryWeightByProcess(true);
91  //int numSecAdd = secondaries.size();
92  std::vector<G4Track*>::iterator iter = secondaries.begin(); // Add all secondaries
93  while (iter != secondaries.end()) {
94  G4Track* myTrack = *iter;
95  G4String pdgstr = myTrack->GetParticleDefinition()->GetParticleName();
96  //G4double ke = myTrack->GetKineticEnergy()/GeV;
97  G4int pdg = myTrack->GetParticleDefinition()->GetPDGEncoding();
98  if (abs(pdg)==130 ||abs(pdg)==310 ||abs(pdg)==311 || abs(pdg)==3122 || abs(pdg)==2112)
99  // std::cout << "MuNuclSplProc: Adding " << pdgstr << " of Kinetic Energy " << ke << " GeV." << std::endl;
100  // Will ask for PdgCode and only Add K0s.
101  myTrack->SetWeight(weight);
102  particleChange->AddSecondary(myTrack);
103  iter++;
104  }
105  return particleChange;
106 }
107 
108 } // end namespace
util::quantities::gigaelectronvolt GeV
var pdg
Definition: selectors.fcl:14
process_name use argoneut_mc_hitfinder track
T abs(T value)
Check of Geant4 to run the LArSoft detector simulation.
Var Sqrt(const Var &v)
Use to take sqrt of a var.
Definition: Var.cxx:172
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step)
BEGIN_PROLOG could also be cout