MuNuclearSplittingProcessXSecBias.cxx
Go to the documentation of this file.
1 //
2 // This(these) is (are) a variance reduction technique(s).
3 //
4 //
5 // Secondary biasing code adapted from EventBiasing from Jane Tinslay, SLAC. 2008.
6 // (Who's now at Cisco, I think.) XSBiasing code from http://geant4.cern.ch/collaboration/workshops/workshop2002/slides/docs/flei/biasing.pdf. (And NOT from http://www.lcsim.org/software/geant4/doxygen/html/classG4HadronicProcess.html, where I *believe* the point is to bias primaries and secondaries from sub-classes of Cross-Sections in an already chosen
7 // broader class of XSs. This doesn't help us. We need more occurrences (albeit with small
8 // weights) of a specific process.
9 //
10 //
11 // The point is to create Nsplit secondaries with each muNucl process, so
12 // that we may, e.g., create many, otherwise-rare, K0s which
13 // may charge exchange into pernicious K+s which then mimic proton dk. And,
14 // additionally now to make the XSection for the muNucl process itself higher,
15 // keeping track of appropriate weights everywhere.
16 //
17 // echurch@fnal.gov, May, 2011.
18 //
19 
20 
22 
23 #include "Geant4/G4DynamicParticle.hh"
24 #include "Geant4/G4KaonZeroShort.hh"
25 #include "Geant4/G4LorentzVector.hh"
26 #include "Geant4/G4ParticleChange.hh"
27 #include "Geant4/G4String.hh"
28 #include "Geant4/G4ThreeVector.hh"
29 #include "Geant4/G4VParticleChange.hh"
30 #include "Geant4/G4Track.hh"
31 #include "Geant4/G4Step.hh"
32 #include "Geant4/G4ParticleDefinition.hh"
33 #include "Geant4/Randomize.hh"
34 #include "Geant4/G4KaonZeroLong.hh"
35 #include "Geant4/G4ios.hh"
36 
37 #include "cetlib_except/exception.h"
38 
39 #include <TMath.h>
40 
41 namespace larg4 {
42 
43 
44 
45  G4VParticleChange* MuNuclearSplittingProcessXSecBias::PostStepDoIt(const G4Track& track, const G4Step& step)
46  {
47 
48  G4VParticleChange* pChange = pRegProcess->PostStepDoIt( track, step );
49  pChange->SetVerboseLevel(0);
50  pChange->SetSecondaryWeightByProcess(true);
51 
52  G4double secWeight = 1.0;
53  if (fNSplit>=1) secWeight = 1.0/fNSplit;
54  std::vector<G4Track*> secondaries; // Secondary store
55 
56  G4double nw=1.0;
57  // Let's go bias the underlying xsection.
58  if(xBiasMode==1)
59  {
60  G4double s = step.GetStepLength();
61  G4double x = pRegProcess->GetCurrentInteractionLength();
62  nw = pChange->GetParentWeight();
63  // need to undo the change made to the parent track in
64  // most recent step (1/wc).
65  // It seems to be possible to get here w.o. first having hit the
66  // AlongStep Method that initialied wc. So protect against division y 0.
67  if (wc < 1.0e-10) wc = 1.0;
68  nw *= 1/wc * (1. - exp(-s/x))/(1 - exp(-eFactor*s/x)) ;
69  G4cout << " MNSPXSB.PostStepDoit(): original wt = " << pChange->GetParentWeight() << " eF = " << eFactor << G4endl;
70  G4cout << " MNSPXSB.PostStepDoit(): New weight = " << nw << " wc = " << wc << " s= " << s <<" x = " << x << G4endl;
71  ((G4ParticleChange*)pChange)->ProposeParentWeight(nw) ;
72  }
73  else if (xBiasMode==2)
74  {
75  nw = pChange->GetParentWeight()*XBiasSecondaryWeight();
76  ((G4ParticleChange*)pChange)->ProposeParentWeight(nw) ;
77  // for (G4int i = 0; i < pChange->GetNumberOfSecondaries(); i++)
78  //pChange->GetSecondary(i)->SetWeight(nw);
79  if(G4UniformRand()<XBiasSurvivalProbability())
80  { // need to add the primary back in to balance the weight
81  // the position of the extra track should be moved to the next volume, but I don't know how yet.
82  // need to change number of secondary first before adding an additional track
83  pChange->SetNumberOfSecondaries(pChange->GetNumberOfSecondaries()+1);
84  G4ThreeVector mDirection = track.GetDynamicParticle()->GetMomentumDirection();
85  G4DynamicParticle* aD = new G4DynamicParticle((track.GetDynamicParticle())->GetDefinition(),
86  G4LorentzVector(mDirection,
87  (track.GetDynamicParticle())->GetKineticEnergy()+step.GetTotalEnergyDeposit()));
88  G4Track* secTrk = new G4Track(aD,step.GetDeltaTime(),step.GetDeltaPosition());
89  pChange->AddSecondary(secTrk);
90  pChange->GetSecondary(pChange->GetNumberOfSecondaries()-1)->
91  SetWeight(XBiasSurvivalProbability()*track.GetWeight());
92  }
93  else
94  {
95  nw = track.GetWeight();
96  pChange->ProposeParentWeight(nw);
97  }
98  }
99 
100  secWeight = secWeight*nw;
101 
102  // Loop over PostStepDoIt method to generate multiple secondaries.
103  // The point is for each value of ii up to Nsplit we want to re-toss
104  // all secondaries for the muNucl process. Each toss gives back numSec
105  // secondaries, which will be a different sample of species/momenta each time.
106  G4VParticleChange* particleChange = new G4VParticleChange();
107 
108  for (unsigned int ii=0; ii<(unsigned int)fNSplit; ii++) {
109  particleChange = pRegProcess->PostStepDoIt(track, step);
110  if (!particleChange)
111  throw std::runtime_error("MuNuclearSplittingProcessXSecBias::PostStepDoIt(): no particle change");
112  G4int j(0);
113  G4int numSec(particleChange->GetNumberOfSecondaries());
114  // Don't change weight of last secondary. It's the just-added primary in this mode.
115  for (j=0; j<numSec; j++)
116  {
117  G4Track* newSec = new G4Track(*(particleChange->GetSecondary(j)));
118  G4String pdgstr = newSec->GetParticleDefinition()->GetParticleName();
119  //G4double ke = newSec->GetKineticEnergy()/GeV;
120  G4int pdg = newSec->GetParticleDefinition()->GetPDGEncoding();
121  if (abs(pdg)==310 ||abs(pdg)==311 || abs(pdg)==3122 || abs(pdg)==2112)
122  {
123  // std::cout << "MuNuclSplProc: Creating " << pdgstr << " of Kinetic Energy " << ke << " GeV. numSec is " << j << std::endl;
124 
125  if (pdg==G4KaonZeroShort::KaonZeroShort()->GetPDGEncoding()&&G4UniformRand()<0.50)
126  {
127  pdg = G4KaonZeroLong::KaonZeroLong()->GetPDGEncoding();
128  pdgstr = G4KaonZeroLong::KaonZeroLong()->GetParticleName();
129  G4LorentzVector pK0L(newSec->GetMomentum(),TMath::Sqrt(TMath::Power(G4KaonZeroLong::KaonZeroLong()->GetPDGMass(),2)+TMath::Power(newSec->GetMomentum().mag(),2)));
130  G4DynamicParticle *newK0L = new G4DynamicParticle(G4KaonZeroLong::KaonZeroLong(),pK0L);
131 
132  G4Track* newSecK0L = new G4Track(newK0L,track.GetGlobalTime(),track.GetPosition());
133  secondaries.push_back(newSecK0L);
134  }
135  else
136  {
137  secondaries.push_back(newSec);
138  }
139  }
140  }
141  }
142 
143  pChange->SetNumberOfSecondaries(secondaries.size());
144  pChange->SetSecondaryWeightByProcess(true);
145  //pChange->ProposeTrackStatus(fStopAndKill); // End of the story for this muon.
146 
147  //int numSecAdd = secondaries.size();
148  std::vector<G4Track*>::iterator iter = secondaries.begin(); // Add all secondaries
149  while (iter != secondaries.end()) {
150  G4Track* myTrack = *iter;
151  G4String pdgstr = myTrack->GetParticleDefinition()->GetParticleName();
152  G4double ke = myTrack->GetKineticEnergy()/CLHEP::GeV;
153  G4int pdg = myTrack->GetParticleDefinition()->GetPDGEncoding();
154  if (abs(pdg)==130 ||abs(pdg)==310 ||abs(pdg)==311 || abs(pdg)==3122 || abs(pdg)==2112)
155  {
156  if (pdg!=2112)
157  std::cout << "MuNuclSplXSB(): Adding " << pdgstr << " of Kinetic Energy and weight " << ke << " GeV and " << secWeight << "." << std::endl;
158 
159  myTrack->SetWeight(secWeight);
160  pChange->AddSecondary(myTrack);
161  }
162  iter++;
163  }
164  return pChange;
165 
166  }
167 
169  const G4Track& track,
170  const G4Step& step
171  )
172  {
173 
174  wc = 1.0;
175 
176  if (xBiasMode==2)
177  {
178  /* I don't understand this. If I new it and Initialize(track) it works, but
179  it goes off and grinds through the event.
180  If I try to use pRegProcess->AlongStepDoIt() it crashes
181  at the GetParentWeight() just below.
182  */
183  G4VParticleChange* pC = new G4VParticleChange();
184  pC->Initialize(track);
185  //pC = pRegProcess->AlongStepDoIt( track, step);
186 
187  G4double nw = pC->GetParentWeight()*XBiasSecondaryWeight();
188  pC->ProposeParentWeight(nw) ;
189  for (G4int i = 0; i < pC->GetNumberOfSecondaries(); i++)
190  pC->GetSecondary(i)->SetWeight(nw);
191  //
192  if(G4UniformRand()<XBiasSurvivalProbability())
193  {
194  pC->SetNumberOfSecondaries(pC->GetNumberOfSecondaries()+1);
195  G4ThreeVector mDirection = track.GetDynamicParticle()->GetMomentumDirection();
196  G4DynamicParticle* aD = new G4DynamicParticle((track.GetDynamicParticle())->GetDefinition(),
197  G4LorentzVector(mDirection,
198  (track.GetDynamicParticle())->GetKineticEnergy()+step.GetTotalEnergyDeposit()));
199  G4Track* secTrk = new G4Track(aD,track.GetGlobalTime(),track.GetPosition()); // +step.GetDeltaTime(), ,+step.GetDeltaPosition()
200  secTrk->SetGoodForTrackingFlag(true);
201  pC->AddSecondary(secTrk); // Add the primary right back.
202  pC->GetSecondary(pC->GetNumberOfSecondaries()-1)->SetWeight(XBiasSurvivalProbability()*track.GetWeight());
203  //pC->ProposeTrackStatus(fStopAndKill); // 8-June-2011-21:36 MDT, from Finland, I (EDC) added this. Seems right. Makes it go, anyway. Almost too fast.
204  }
205 
206  return pC;
207  }
208  else if (xBiasMode==1)
209  {
210  fParticleChange.Initialize(track) ;
211  G4double s = step.GetStepLength();
212  G4double x = pRegProcess->GetCurrentInteractionLength();
213  //fParticleChange = &(pRegProcess->AlongStepDoIt( track, step));
214  wc = exp(-(1. - eFactor)*s/x);
215  G4double w = wc * fParticleChange.GetParentWeight();
216  fParticleChange.SetParentWeightByProcess(false);
217  fParticleChange.ProposeParentWeight(w) ;
218  if (w>100)
219  {
220  // G4cout << " MNSPXSB.AlongStepDoit(): GPW is " << w/wc << " wc = " << wc << ", CIL = "<< x << G4endl;
221  }
222 
223  return &fParticleChange;
224  }
225 
226  if (xBiasMode!=1 && xBiasMode!=2)
227  {
228  throw cet::exception("Incorrectly set Bias Mode. ")
229  << "Set XBiasMode to 0 or 1. " << xBiasMode << " not allowed!\n";
230  }
231  fParticleChange.Initialize(track) ;
232  return &fParticleChange;
233  }
234 
235 
237  {
238  G4double result = 0;
239  G4double nLTraversed = GetTotalNumberOfInteractionLengthTraversed();
240  G4double biasedProbability = 1.-std::exp(-nLTraversed);
241  G4double realProbability = 1-std::exp(-nLTraversed/eFactor);
242  result = (biasedProbability-realProbability)/biasedProbability;
243  return result;
244  }
245 
247  {
248  G4double result = 0;
249  G4double nLTraversed = GetTotalNumberOfInteractionLengthTraversed();
250  result = 1./eFactor*std::exp(-nLTraversed/eFactor*(1-1./eFactor));
251  return result;
252  }
253 
254 
256  {
258  -G4VProcess::theNumberOfInteractionLengthLeft;
259  }
260 
261 
263  G4double previousStepSize,
264  G4ForceCondition* condition )
265  {
266  // I believe this is the money line for the whole thing. By shrinking the
267  // interaction length in MuNuclear's process, it's more likely that
268  // this process will be chosen in the competition among processes.
269  if(xBiasMode)
270  {
271  G4double psgPIL = (1./eFactor) * pRegProcess->PostStepGetPhysicalInteractionLength(
272  track,
273  previousStepSize,
274  condition );
275  return psgPIL;
276  }
277 
278  return pRegProcess->PostStepGetPhysicalInteractionLength(
279  track,
280  previousStepSize,
281  condition );
282 
283 }
284 
285 
287  G4double previousStepSize,
288  G4double currentMinimumStep,
289  G4double& proposedSafety,
290  G4GPILSelection* selection )
291 {
292  if (xBiasMode==2)
293  {
294  // I am pretty sure this mode does not work because MuNuclear doesn't interact
295  // AlongStep. Below often (always?) returns a negative number, which I have to
296  // swizzle.
297  G4double intLen;
298  intLen = (1./eFactor) *
299  pRegProcess->AlongStepGetPhysicalInteractionLength( track,previousStepSize,
300  currentMinimumStep,proposedSafety,selection);
301  if (intLen<0) return currentMinimumStep;
302  else return intLen;
303  }
304  else
305  return DBL_MAX ;
306 }
307 
308 
309 } // end namespace
intermediate_table::iterator iterator
static QCString result
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
Geant4 interface.
Check of Geant4 to run the LArSoft detector simulation.
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
T abs(T value)
const double e
static constexpr double GeV
Definition: Units.h:28
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step)
list x
Definition: train.py:276
static QCString * s
Definition: config.cpp:1042
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step)
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
QTextStream & endl(QTextStream &s)