Geant4  10.01
G4HadronStoppingProcess.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4HadronStoppingProcess.cc 85517 2014-10-30 15:53:26Z gcosmo $
27 //
28 //---------------------------------------------------------------------
29 //
30 // GEANT4 Class
31 //
32 // File name: G4HadronStoppingProcess
33 //
34 // Author V.Ivanchenko 21 April 2012
35 //
36 //
37 // Class Description:
38 //
39 // Base process class for nuclear capture of negatively charged particles
40 //
41 // Modifications:
42 //
43 // 20120522 M. Kelsey -- Set enableAtRestDoIt flag for G4ProcessManager
44 // 20120914 M. Kelsey -- Pass subType in base ctor, remove enable flags
45 // 20121004 K. Genser -- use G4HadronicProcessType in the constructor
46 // 20121016 K. Genser -- Reverting to use one argument c'tor
47 // 20140818 K. Genser -- Labeled tracks using G4PhysicsModelCatalog
48 //
49 //------------------------------------------------------------------------
50 
53 #include "G4HadronicProcessType.hh"
54 #include "G4EmCaptureCascade.hh"
55 #include "G4Nucleus.hh"
56 #include "G4HadFinalState.hh"
57 #include "G4HadProjectile.hh"
58 #include "G4HadSecondary.hh"
59 #include "G4Material.hh"
60 #include "G4Element.hh"
61 
62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63 
66  fElementSelector(new G4ElementSelector()),
67  fEmCascade(new G4EmCaptureCascade()), // Owned by InteractionRegistry
68  fBoundDecay(0x0),
69  emcID(G4PhysicsModelCatalog::Register(G4String((name + "_EMCascade")))),
70  ncID(G4PhysicsModelCatalog::Register(G4String((name + "_NuclearCapture")))),
71  dioID(G4PhysicsModelCatalog::Register(G4String((name + "_DIO"))))
72 {
73  // Modify G4VProcess flags to emulate G4VRest instead of G4VDiscrete
74  enableAtRestDoIt = true;
75  enablePostStepDoIt = false;
76 
78 }
79 
80 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
81 
83 {
85  delete fElementSelector;
86  // NOTE: fEmCascade and fEmBoundDecay owned by registry, not locally
87 }
88 
89 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
90 
92 {
93  return (p.GetPDGCharge() < 0.0);
94 }
95 
96 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
97 
98 void
100 {
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105 
107 {
109 }
110 
111 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
112 
115 {
116  *condition = NotForced;
117  return 0.0;
118 }
119 
120 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
121 
124 {
125  *condition = NotForced;
126  return DBL_MAX;
127 }
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
130 
132  const G4Step&)
133 {
134  // if primary is not Alive then do nothing
135  theTotalResult->Initialize(track);
136 
137  G4Nucleus* nucleus = GetTargetNucleusPointer();
138  G4Element* elm = fElementSelector->SelectZandA(track, nucleus);
139 
140  G4HadFinalState* result = 0;
141  thePro.Initialise(track);
142  G4double time0 = track.GetGlobalTime();
143  G4bool nuclearCapture = true;
144 
145  // Do the electromagnetic cascade in the nuclear field.
146  // EM cascade should keep G4HadFinalState object,
147  // because it will not be deleted at the end of this method
148  //
149  result = fEmCascade->ApplyYourself(thePro, *nucleus);
150  G4double ebound = result->GetLocalEnergyDeposit();
151  G4double edep = 0.0;
152  G4int nSecondaries = result->GetNumberOfSecondaries();
153  G4int nEmCascaceSec = nSecondaries;
154 
155  // Try decay from bound level
156  // For mu- the time of projectile should be changed.
157  // Decay should keep G4HadFinalState object,
158  // because it will not be deleted at the end of this method
159  //
160  thePro.SetBoundEnergy(ebound);
161  if(fBoundDecay) {
162  G4HadFinalState* resultDecay =
163  fBoundDecay->ApplyYourself(thePro, *nucleus);
164  G4int n = resultDecay->GetNumberOfSecondaries();
165  if(0 < n) {
166  nSecondaries += n;
167  result->AddSecondaries(resultDecay);
168  }
169  if(resultDecay->GetStatusChange() == stopAndKill) {
170  nuclearCapture = false;
171  }
172  resultDecay->Clear();
173  }
174 
175  if(nuclearCapture) {
176 
177  // select model
178  G4HadronicInteraction* model = 0;
179  try {
180  model = ChooseHadronicInteraction(thePro, *nucleus, track.GetMaterial(), elm);
181  }
182  catch(G4HadronicException & aE) {
184  ed << "Target element "<<elm->GetName()<<" Z= "
185  << nucleus->GetZ_asInt() << " A= "
186  << nucleus->GetA_asInt() << G4endl;
187  DumpState(track,"ChooseHadronicInteraction",ed);
188  ed << " No HadronicInteraction found out" << G4endl;
189  G4Exception("G4HadronStoppingProcess::AtRestDoIt", "had005",
190  FatalException, ed);
191  }
192 
193  G4HadFinalState* resultNuc = 0;
194  G4int reentryCount = 0;
195  do {
196  // sample final state
197  // nuclear interaction should keep G4HadFinalState object
198  // model should define time of each secondary particle
199  try {
200  resultNuc = model->ApplyYourself(thePro, *nucleus);
201  ++reentryCount;
202  }
203  catch(G4HadronicException aR) {
205  ed << "Call for " << model->GetModelName() << G4endl;
206  ed << "Target element "<<elm->GetName()<<" Z= "
207  << nucleus->GetZ_asInt()
208  << " A= " << nucleus->GetA_asInt() << G4endl;
209  DumpState(track,"ApplyYourself",ed);
210  ed << " ApplyYourself failed" << G4endl;
211  G4Exception("G4HadronStoppingProcess::AtRestDoIt", "had006",
212  FatalException, ed);
213  }
214 
215  // Check the result for catastrophic energy non-conservation
216  resultNuc = CheckResult(thePro, *nucleus, resultNuc);
217 
218  if(reentryCount>100) {
220  ed << "Call for " << model->GetModelName() << G4endl;
221  ed << "Target element "<<elm->GetName()<<" Z= "
222  << nucleus->GetZ_asInt()
223  << " A= " << nucleus->GetA_asInt() << G4endl;
224  DumpState(track,"ApplyYourself",ed);
225  ed << " ApplyYourself does not completed after 100 attempts" << G4endl;
226  G4Exception("G4HadronStoppingProcess::AtRestDoIt", "had006",
227  FatalException, ed);
228  }
229  }
230  while(!resultNuc);
231 
232  edep = resultNuc->GetLocalEnergyDeposit();
233  nSecondaries += resultNuc->GetNumberOfSecondaries();
234  result->AddSecondaries(resultNuc);
235  resultNuc->Clear();
236  }
237 
238  // Fill results
239  //
242  theTotalResult->SetNumberOfSecondaries(nSecondaries);
243  G4double w = track.GetWeight();
245  for(G4int i=0; i<nSecondaries; ++i) {
246  G4HadSecondary* sec = result->GetSecondary(i);
247 
248  // add track global time to the reaction time
249  G4double time = sec->GetTime();
250  if(time < 0.0) { time = 0.0; }
251  time += time0;
252 
253  // create secondary track
254  G4Track* t = new G4Track(sec->GetParticle(),
255  time,
256  track.GetPosition());
257  t->SetWeight(w*sec->GetWeight());
258 
259  // use SetCreatorModelIndex to "label" the track
260  if (i<nEmCascaceSec) {
262  } else if (nuclearCapture) {
264  } else {
266  }
267 
270  }
271  result->Clear();
272 
273  if (epReportLevel != 0) {
274  CheckEnergyMomentumConservation(track, *nucleus);
275  }
276  return theTotalResult;
277 }
278 
279 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
280 
281 void G4HadronStoppingProcess::ProcessDescription(std::ostream& outFile) const
282 {
283  outFile << "Base process for negatively charged particle capture at rest.\n";
284 }
285 
286 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double condition(const G4ErrorSymMatrix &m)
void DeRegisterExtraProcess(G4VProcess *)
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
virtual G4Element * SelectZandA(const G4Track &track, G4Nucleus *)
void AddSecondaries(const std::vector< G4HadSecondary > &addSecs)
virtual G4bool IsApplicable(const G4ParticleDefinition &)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4HadSecondary * GetSecondary(size_t i)
static G4HadronicProcessStore * Instance()
G4String name
Definition: TRTMaterials.hh:40
const G4ThreeVector & GetPosition() const
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4HadProjectile thePro
G4HadronicInteraction * fEmCascade
void CheckEnergyMomentumConservation(const G4Track &, const G4Nucleus &)
const G4String & GetModelName() const
G4double GetTime() const
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
int G4int
Definition: G4Types.hh:78
void SetWeight(G4double aValue)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetCreatorModelIndex(G4int idx)
G4HadronStoppingProcess(const G4String &name="hadronCaptureAtRest")
void ProposeWeight(G4double finalWeight)
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
void Register(T *inst)
Definition: G4AutoDelete.hh:65
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)
G4ParticleChange * theTotalResult
bool G4bool
Definition: G4Types.hh:79
Definition: G4Step.hh:76
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:352
const G4int n
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
G4double GetGlobalTime() const
const G4TouchableHandle & GetTouchableHandle() const
G4HadronicInteraction * fBoundDecay
G4Material * GetMaterial() const
void Initialise(const G4Track &aT)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4ElementSelector * fElementSelector
void RegisterExtraProcess(G4VProcess *)
virtual void Initialize(const G4Track &)
G4Nucleus * GetTargetNucleusPointer()
void SetNumberOfSecondaries(G4int totSecondaries)
G4DynamicParticle * GetParticle()
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
void SetBoundEnergy(G4double e)
virtual void ProcessDescription(std::ostream &outFile) const
void AddSecondary(G4Track *aSecondary)
G4double GetWeight() const
G4HadronicInteraction * ChooseHadronicInteraction(const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, G4Material *aMaterial, G4Element *anElement)
#define G4endl
Definition: G4ios.hh:61
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
const G4String & GetName() const
Definition: G4Element.hh:127
G4ForceCondition
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
G4int GetNumberOfSecondaries() const
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
G4double GetLocalEnergyDeposit() const
G4HadFinalState * CheckResult(const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result)
G4HadFinalStateStatus GetStatusChange() const
void PrintInfo(const G4ParticleDefinition *)
G4double GetWeight() const