Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronCaptureAtRest.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// G4NeutronCaptureAtRest physics process
27// Larry Felawka (TRIUMF), April 1998
28//---------------------------------------------------------------------
29
30#include <string.h>
31#include <cmath>
32#include <stdio.h>
33
35#include "G4SystemOfUnits.hh"
36#include "G4DynamicParticle.hh"
37#include "G4ParticleTypes.hh"
38#include "Randomize.hh"
41
42#define MAX_SECONDARIES 100
43
44// constructor
45
46G4NeutronCaptureAtRest::G4NeutronCaptureAtRest(const G4String& processName,
47 G4ProcessType aType ) :
48 G4VRestProcess (processName, aType), // initialization
49 massProton(G4Proton::Proton()->GetPDGMass()/GeV),
50 massNeutron(G4Neutron::Neutron()->GetPDGMass()/GeV),
51 massElectron(G4Electron::Electron()->GetPDGMass()/GeV),
52 massDeuteron(G4Deuteron::Deuteron()->GetPDGMass()/GeV),
53 massAlpha(G4Alpha::Alpha()->GetPDGMass()/GeV),
54 pdefGamma(G4Gamma::Gamma()),
55 pdefNeutron(G4Neutron::Neutron())
56{
57 G4HadronicDeprecate("G4NeutronCaptureAtRest");
58 if (verboseLevel>0) {
59 G4cout << GetProcessName() << " is created "<< G4endl;
60 }
65
67}
68
69// destructor
70
72{
74 delete [] pv;
75 delete [] eve;
76 delete [] gkin;
77}
78
80{
82}
83
85{
87}
88
89// methods.............................................................................
90
92 const G4ParticleDefinition& particle
93 )
94{
95 return ( &particle == pdefNeutron );
96
97}
98
99// Warning - this method may be optimized away if made "inline"
101{
102 return ( ngkine );
103
104}
105
106// Warning - this method may be optimized away if made "inline"
108{
109 return ( &gkin[0] );
110
111}
112
114 const G4Track& track,
116 )
117{
118 // beggining of tracking
120
121 // condition is set to "Not Forced"
123
124 // get mean life time
126
127 if ((currentInteractionLength <0.0) || (verboseLevel>2)){
128 G4cout << "G4NeutronCaptureAtRestProcess::AtRestGetPhysicalInteractionLength ";
129 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
130 track.GetDynamicParticle()->DumpInfo();
131 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
132 G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]" <<G4endl;
133 }
134
136
137}
138
140 const G4Track& track,
141 const G4Step&
142 )
143//
144// Handles Neutrons at rest; a Neutron can either create secondaries or
145// do nothing (in which case it should be sent back to decay-handling
146// section
147//
148{
149
150// Initialize ParticleChange
151// all members of G4VParticleChange are set to equal to
152// corresponding member in G4Track
153
155
156// Store some global quantities that depend on current material and particle
157
158 globalTime = track.GetGlobalTime()/s;
159 G4Material * aMaterial = track.GetMaterial();
160 const G4int numberOfElements = aMaterial->GetNumberOfElements();
161 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
162
163 const G4double* theAtomicNumberDensity = aMaterial->GetAtomicNumDensityVector();
164 G4double normalization = 0;
165 for ( G4int i1=0; i1 < numberOfElements; i1++ )
166 {
167 normalization += theAtomicNumberDensity[i1] ; // change when nucleon specific
168 // probabilities are included.
169 }
170 G4double runningSum= 0.;
171 G4double random = G4UniformRand()*normalization;
172 for ( G4int i2=0; i2 < numberOfElements; i2++ )
173 {
174 runningSum += theAtomicNumberDensity[i2]; // change when nucleon specific
175 // probabilities are included.
176 if (random<=runningSum)
177 {
178 targetCharge = G4double((*theElementVector)[i2]->GetZ());
179 targetAtomicMass = (*theElementVector)[i2]->GetN();
180 }
181 }
182 if (random>runningSum)
183 {
184 targetCharge = G4double((*theElementVector)[numberOfElements-1]->GetZ());
185 targetAtomicMass = (*theElementVector)[numberOfElements-1]->GetN();
186
187 }
188
189 if (verboseLevel>1) {
190 G4cout << "G4NeutronCaptureAtRest::AtRestDoIt is invoked " <<G4endl;
191 }
192
193 G4ParticleMomentum momentum;
194 G4float localtime;
195
197
198 GenerateSecondaries(); // Generate secondaries
199
201
202 for ( G4int isec = 0; isec < ngkine; isec++ ) {
203 G4DynamicParticle* aNewParticle = new G4DynamicParticle;
204 aNewParticle->SetDefinition( gkin[isec].GetParticleDef() );
205 aNewParticle->SetMomentum( gkin[isec].GetMomentum() * GeV );
206
207 localtime = globalTime + gkin[isec].GetTOF();
208
209 G4Track* aNewTrack = new G4Track( aNewParticle, localtime*s, position );
210 aNewTrack->SetTouchableHandle(track.GetTouchableHandle());
211 aParticleChange.AddSecondary( aNewTrack );
212
213 }
214
216
217 aParticleChange.ProposeTrackStatus(fStopAndKill); // Kill the incident Neutron
218
219// clear InteractionLengthLeft
220
222
223 return &aParticleChange;
224
225}
226
227
228void G4NeutronCaptureAtRest::GenerateSecondaries()
229{
230 static G4int index;
231 static G4int l;
232 static G4int nopt;
233 static G4int i;
234 // DHW 15 May 2011: unused: static G4ParticleDefinition* jnd;
235
236 for (i = 1; i <= MAX_SECONDARIES; ++i) {
237 pv[i].SetZero();
238 }
239
240 ngkine = 0; // number of generated secondary particles
241 ntot = 0;
242 result.SetZero();
243 result.SetMass( massNeutron );
244 result.SetKineticEnergyAndUpdate( 0. );
245 result.SetTOF( 0. );
246 result.SetParticleDef( pdefNeutron );
247
248 NeutronCapture(&nopt);
249
250 // *** CHECK WHETHER THERE ARE NEW PARTICLES GENERATED ***
251 if (ntot != 0 || result.GetParticleDef() != pdefNeutron) {
252 // *** CURRENT PARTICLE IS NOT THE SAME AS IN THE BEGINNING OR/AND ***
253 // *** ONE OR MORE SECONDARIES HAVE BEEN GENERATED ***
254
255 // --- INITIAL PARTICLE TYPE HAS BEEN CHANGED ==> PUT NEW TYPE ON ---
256 // --- THE GEANT TEMPORARY STACK ---
257
258 // --- PUT PARTICLE ON THE STACK ---
259 gkin[0] = result;
260 gkin[0].SetTOF( result.GetTOF() * 5e-11 );
261 ngkine = 1;
262
263 // --- ALL QUANTITIES ARE TAKEN FROM THE GHEISHA STACK WHERE THE ---
264 // --- CONVENTION IS THE FOLLOWING ---
265
266 // --- ONE OR MORE SECONDARIES HAVE BEEN GENERATED ---
267 for (l = 1; l <= ntot; ++l) {
268 index = l - 1;
269 // DHW 15 May 2011: unused: jnd = eve[index].GetParticleDef();
270
271 // --- ADD PARTICLE TO THE STACK IF STACK NOT YET FULL ---
272 if (ngkine < MAX_SECONDARIES) {
273 gkin[ngkine] = eve[index];
274 gkin[ngkine].SetTOF( eve[index].GetTOF() * 5e-11 );
275 ++ngkine;
276 }
277 }
278 }
279 else {
280 // --- NO SECONDARIES GENERATED AND PARTICLE IS STILL THE SAME ---
281 // --- ==> COPY EVERYTHING BACK IN THE CURRENT GEANT STACK ---
282 ngkine = 0;
283 ntot = 0;
284 globalTime += result.GetTOF() * G4float(5e-11);
285 }
286
287 // --- LIMIT THE VALUE OF NGKINE IN CASE OF OVERFLOW ---
288 ngkine = G4int(std::min(ngkine,G4int(MAX_SECONDARIES)));
289
290} // GenerateSecondaries
291
292
293void G4NeutronCaptureAtRest::Normal(G4float *ran)
294{
295 static G4int i;
296
297 // *** NVE 14-APR-1988 CERN GENEVA ***
298 // ORIGIN : H.FESEFELDT (27-OCT-1983)
299
300 *ran = G4float(-6.);
301 for (i = 1; i <= 12; ++i) {
302 *ran += G4UniformRand();
303 }
304
305} // Normal
306
307
308void G4NeutronCaptureAtRest::NeutronCapture(G4int *nopt)
309{
310 static G4int nt;
311 static G4float xp, pcm;
312 static G4float ran;
313
314 // *** ROUTINE FOR CAPTURE OF NEUTRAL BARYONS ***
315 // *** NVE 04-MAR-1988 CERN GENEVA ***
316 // ORIGIN : H.FESEFELDT (02-DEC-1986)
317
318 *nopt = 1;
319 pv[1] = result;
320 pv[2].SetZero();
321 pv[2].SetMass( AtomAs(targetAtomicMass, targetCharge) );
322 pv[2].SetMomentumAndUpdate( 0., 0., 0. );
323 pv[2].SetTOF( result.GetTOF() );
324 pv[2].SetParticleDef( NULL );
325 pv[MAX_SECONDARIES].Add( pv[1], pv[2] );
326 pv[MAX_SECONDARIES].SetMomentum( -pv[MAX_SECONDARIES].GetMomentum().x(), -pv[MAX_SECONDARIES].GetMomentum().y(), -pv[MAX_SECONDARIES].GetMomentum().z() );
328 Normal(&ran);
329 pcm = ran * G4float(.001) + G4float(.0065);
330 ran = G4UniformRand();
331 result.SetTOF( result.GetTOF() - std::log(ran) * G4float(480.) );
332 pv[3].SetZero();
333 pv[3].SetMass( 0. );
334 pv[3].SetKineticEnergyAndUpdate( pcm );
335 pv[3].SetTOF( result.GetTOF() );
336 pv[3].SetParticleDef( pdefGamma );
337 pv[3].Lor( pv[3], pv[MAX_SECONDARIES] );
338 nt = 3;
339 xp = G4float(.008) - pcm;
340 if (xp >= G4float(0.)) {
341 nt = 4;
342 pv[4].SetZero();
343 pv[4].SetMass( 0. );
344 pv[4].SetKineticEnergyAndUpdate( xp );
345 pv[4].SetTOF( result.GetTOF() );
346 pv[4].SetParticleDef( pdefGamma );
347 pv[4].Lor( pv[4], pv[MAX_SECONDARIES] );
348 }
349 result = pv[3];
350 if (nt == 4) {
351 if (ntot < MAX_SECONDARIES-1) {
352 eve[ntot++] = pv[4];
353 }
354 }
355
356} // NeutronCapture
357
358
359G4double G4NeutronCaptureAtRest::AtomAs(G4float a, G4float z)
360{
361 G4float ret_val;
362 G4double d__1, d__2;
363
364 static G4double aa;
365 static G4int ia, iz;
366 static G4double zz;
367 static G4float rma, rmd;
368 static G4int ipp;
369 static G4float rmn, rmp;
370 static G4int izz;
371 static G4float rmel;
372 static G4double mass;
373
374 // *** DETERMINATION OF THE ATOMIC MASS ***
375 // *** NVE 19-MAY-1988 CERN GENEVA ***
376 // ORIGIN : H.FESEFELDT (02-DEC-1986)
377
378 // --- GET ATOMIC (= ELECTRONS INCL.) MASSES (IN MEV) FROM RMASS ARRAY ---
379 // --- ELECTRON ---
380 rmel = massElectron * G4float(1e3);
381 // --- PROTON ---
382 rmp = massProton * G4float(1e3);
383 // --- NEUTRON ---
384 rmn = massNeutron * G4float(1e3);
385 // --- DEUTERON ---
386 rmd = massDeuteron * G4float(1e3) + rmel;
387 // --- ALPHA ---
388 rma = massAlpha * G4float(1e3) + rmel * G4float(2.);
389
390 ret_val = G4float(0.);
391 aa = a * 1.;
392 zz = z * 1.;
393 ia = G4int(a + G4float(.5));
394 if (ia < 1) {
395 return ret_val;
396 }
397 iz = G4int(z + G4float(.5));
398 if (iz < 0 || iz > ia) {
399 return ret_val;
400 }
401 mass = 0.;
402 if (ia == 1) {
403 if (iz == 0) {
404 mass = rmn;
405 }
406 else if (iz == 1) {
407 mass = rmp + rmel;
408 }
409 }
410 else if (ia == 2 && iz == 1) {
411 mass = rmd;
412 }
413 else if (ia == 4 && iz == 2) {
414 mass = rma;
415 }
416 else if ( (ia == 2 && iz != 1) || ia == 3 || (ia == 4 && iz != 2) || ia > 4) {
417 d__1 = aa / G4float(2.) - zz;
418 d__2 = zz;
419 mass = (aa - zz) * rmn + zz * rmp + zz * rmel - aa * G4float(15.67) +
420 std::pow(aa, .6666667) * G4float(17.23) + d__1 * d__1 * G4float(93.15) / aa +
421 d__2 * d__2 * G4float(.6984523) / std::pow(aa, .3333333);
422 ipp = (ia - iz) % 2;
423 izz = iz % 2;
424 if (ipp == izz) {
425 mass += (ipp + izz - 1) * G4float(12.) * std::pow(aa, -.5);
426 }
427 }
428 ret_val = mass * G4float(.001);
429 return ret_val;
430
431} // AtomAs
#define MAX_SECONDARIES
std::vector< G4Element * > G4ElementVector
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
#define G4HadronicDeprecate(name)
@ fHadronAtRest
G4ProcessType
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
float G4float
Definition: G4Types.hh:65
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
void DumpInfo(G4int mode=0) const
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetMomentum(const G4ThreeVector &momentum)
G4ParticleDefinition * GetParticleDef()
void SetParticleDef(G4ParticleDefinition *c)
void SetMomentum(G4ParticleMomentum mom)
void SetMomentumAndUpdate(G4ParticleMomentum mom)
void Lor(const G4GHEKinematicsVector &p1, const G4GHEKinematicsVector &p2)
void SetKineticEnergyAndUpdate(G4double ekin)
void Add(const G4GHEKinematicsVector &p1, const G4GHEKinematicsVector &p2)
void DeRegisterExtraProcess(G4VProcess *)
void RegisterExtraProcess(G4VProcess *)
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)
static G4HadronicProcessStore * Instance()
void PrintInfo(const G4ParticleDefinition *)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:215
const G4String & GetName() const
Definition: G4Material.hh:177
void PreparePhysicsTable(const G4ParticleDefinition &)
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4bool IsApplicable(const G4ParticleDefinition &)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
G4double GetMeanLifeTime(const G4Track &, G4ForceCondition *)
G4GHEKinematicsVector * GetSecondaryKinematics()
void AddSecondary(G4Track *aSecondary)
virtual void Initialize(const G4Track &)
Definition: G4Step.hh:78
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4double currentInteractionLength
Definition: G4VProcess.hh:297
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:92
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
#define ns
Definition: xmlparse.cc:597