Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4KaonMinusAbsorption.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// G4KaonMinusAbsorption physics process
27// Larry Felawka (TRIUMF), April 1998
28//---------------------------------------------------------------------
29
30#include <string.h>
31#include <cmath>
32#include <stdio.h>
33
35#include "G4DynamicParticle.hh"
36#include "G4ParticleTypes.hh"
37#include "Randomize.hh"
38#include "G4SystemOfUnits.hh"
41
42#define MAX_SECONDARIES 100
43
44// constructor
45
46G4KaonMinusAbsorption::G4KaonMinusAbsorption(const G4String& processName,
47 G4ProcessType aType ) :
48 G4VRestProcess (processName, aType), // initialization
49 massKaonMinus(G4KaonMinus::KaonMinus()->GetPDGMass()/GeV),
50 massGamma(G4Gamma::Gamma()->GetPDGMass()/GeV),
51 massPionZero(G4PionZero::PionZero()->GetPDGMass()/GeV),
52 massProton(G4Proton::Proton()->GetPDGMass()/GeV),
53 massLambda(G4Lambda::Lambda()->GetPDGMass()/GeV),
54 pdefKaonMinus(G4KaonMinus::KaonMinus()),
55 pdefGamma(G4Gamma::Gamma()),
56 pdefPionZero(G4PionZero::PionZero()),
57 pdefProton(G4Proton::Proton()),
58 pdefNeutron(G4Neutron::Neutron()),
59 pdefLambda(G4Lambda::Lambda()),
60 pdefDeuteron(G4Deuteron::Deuteron()),
61 pdefTriton(G4Triton::Triton()),
62 pdefAlpha(G4Alpha::Alpha())
63{
64 G4HadronicDeprecate("G4KaonMinusAbsorption");
65 if (verboseLevel>0) {
66 G4cout << GetProcessName() << " is created "<< G4endl;
67 }
72
74}
75
76// destructor
77
79{
81 delete [] pv;
82 delete [] eve;
83 delete [] gkin;
84}
85
87{
89}
90
92{
94}
95
96// methods.............................................................................
97
99 const G4ParticleDefinition& particle
100 )
101{
102 return ( &particle == pdefKaonMinus );
103
104}
105
106// Warning - this method may be optimized away if made "inline"
108{
109 return ( ngkine );
110
111}
112
113// Warning - this method may be optimized away if made "inline"
115{
116 return ( &gkin[0] );
117
118}
119
121 const G4Track& track,
123 )
124{
125 // beggining of tracking
127
128 // condition is set to "Not Forced"
130
131 // get mean life time
133
134 if ((currentInteractionLength <0.0) || (verboseLevel>2)){
135 G4cout << "G4KaonMinusAbsorptionProcess::AtRestGetPhysicalInteractionLength ";
136 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
137 track.GetDynamicParticle()->DumpInfo();
138 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
139 G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]" <<G4endl;
140 }
141
143
144}
145
147 const G4Track& track,
148 const G4Step&
149 )
150//
151// Handles KaonMinus at rest; a KaonMinus can either create secondaries or
152// do nothing (in which case it should be sent back to decay-handling
153// section
154//
155{
156
157// Initialize ParticleChange
158// all members of G4VParticleChange are set to equal to
159// corresponding member in G4Track
160
162
163// Store some global quantities that depend on current material and particle
164
165 globalTime = track.GetGlobalTime()/s;
166 G4Material * aMaterial = track.GetMaterial();
167 const G4int numberOfElements = aMaterial->GetNumberOfElements();
168 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
169
170 const G4double* theAtomicNumberDensity = aMaterial->GetAtomicNumDensityVector();
171 G4double normalization = 0;
172 for ( G4int i1=0; i1 < numberOfElements; i1++ )
173 {
174 normalization += theAtomicNumberDensity[i1] ; // change when nucleon specific
175 // probabilities are included.
176 }
177 G4double runningSum= 0.;
178 G4double random = G4UniformRand()*normalization;
179 for ( G4int i2=0; i2 < numberOfElements; i2++ )
180 {
181 runningSum += theAtomicNumberDensity[i2]; // change when nucleon specific
182 // probabilities are included.
183 if (random<=runningSum)
184 {
185 targetCharge = G4double( ((*theElementVector)[i2])->GetZ());
186 targetAtomicMass = (*theElementVector)[i2]->GetN();
187 }
188 }
189 if (random>runningSum)
190 {
191 targetCharge = G4double((*theElementVector)[numberOfElements-1]->GetZ());
192 targetAtomicMass = (*theElementVector)[numberOfElements-1]->GetN();
193
194 }
195
196 if (verboseLevel>1) {
197 G4cout << "G4KaonMinusAbsorption::AtRestDoIt is invoked " <<G4endl;
198 }
199
200 G4ParticleMomentum momentum;
201 G4float localtime;
202
204
205 GenerateSecondaries(); // Generate secondaries
206
208
209 for ( G4int isec = 0; isec < ngkine; isec++ ) {
210 G4DynamicParticle* aNewParticle = new G4DynamicParticle;
211 aNewParticle->SetDefinition( gkin[isec].GetParticleDef() );
212 aNewParticle->SetMomentum( gkin[isec].GetMomentum() * GeV );
213
214 localtime = globalTime + gkin[isec].GetTOF();
215
216 G4Track* aNewTrack = new G4Track( aNewParticle, localtime*s, position );
217 aNewTrack->SetTouchableHandle(track.GetTouchableHandle());
218 aParticleChange.AddSecondary( aNewTrack );
219
220 }
221
223
224 aParticleChange.ProposeTrackStatus(fStopAndKill); // Kill the incident KaonMinus
225
226// clear InteractionLengthLeft
227
229
230 return &aParticleChange;
231
232}
233
234
235void G4KaonMinusAbsorption::GenerateSecondaries()
236{
237 static G4int index;
238 static G4int l;
239 static G4int nopt;
240 static G4int i;
241 // DHW 15 May 2011: unused: static G4ParticleDefinition* jnd;
242
243 for (i = 1; i <= MAX_SECONDARIES; ++i) {
244 pv[i].SetZero();
245 }
246
247 ngkine = 0; // number of generated secondary particles
248 ntot = 0;
249 result.SetZero();
250 result.SetMass( massKaonMinus );
251 result.SetKineticEnergyAndUpdate( 0. );
252 result.SetTOF( 0. );
253 result.SetParticleDef( pdefKaonMinus );
254
255 KaonMinusAbsorption(&nopt);
256
257 // *** CHECK WHETHER THERE ARE NEW PARTICLES GENERATED ***
258 if (ntot != 0 || result.GetParticleDef() != pdefKaonMinus) {
259 // *** CURRENT PARTICLE IS NOT THE SAME AS IN THE BEGINNING OR/AND ***
260 // *** ONE OR MORE SECONDARIES HAVE BEEN GENERATED ***
261
262 // --- INITIAL PARTICLE TYPE HAS BEEN CHANGED ==> PUT NEW TYPE ON ---
263 // --- THE GEANT TEMPORARY STACK ---
264
265 // --- PUT PARTICLE ON THE STACK ---
266 gkin[0] = result;
267 gkin[0].SetTOF( result.GetTOF() * 5e-11 );
268 ngkine = 1;
269
270 // --- ALL QUANTITIES ARE TAKEN FROM THE GHEISHA STACK WHERE THE ---
271 // --- CONVENTION IS THE FOLLOWING ---
272
273 // --- ONE OR MORE SECONDARIES HAVE BEEN GENERATED ---
274 for (l = 1; l <= ntot; ++l) {
275 index = l - 1;
276 // DHW 15 May 2011: unused: jnd = eve[index].GetParticleDef();
277
278 // --- ADD PARTICLE TO THE STACK IF STACK NOT YET FULL ---
279 if (ngkine < MAX_SECONDARIES) {
280 gkin[ngkine] = eve[index];
281 gkin[ngkine].SetTOF( eve[index].GetTOF() * 5e-11 );
282 ++ngkine;
283 }
284 }
285 }
286 else {
287 // --- NO SECONDARIES GENERATED AND PARTICLE IS STILL THE SAME ---
288 // --- ==> COPY EVERYTHING BACK IN THE CURRENT GEANT STACK ---
289 ngkine = 0;
290 ntot = 0;
291 globalTime += result.GetTOF() * G4float(5e-11);
292 }
293
294 // --- LIMIT THE VALUE OF NGKINE IN CASE OF OVERFLOW ---
295 ngkine = G4int(std::min(ngkine,G4int(MAX_SECONDARIES)));
296
297} // GenerateSecondaries
298
299
300void G4KaonMinusAbsorption::Poisso(G4float xav, G4int *iran)
301{
302 static G4int i;
303 static G4float r, p1, p2, p3;
304 static G4int fivex;
305 static G4float rr, ran, rrr, ran1;
306
307 // *** GENERATION OF POISSON DISTRIBUTION ***
308 // *** NVE 16-MAR-1988 CERN GENEVA ***
309 // ORIGIN : H.FESEFELDT (27-OCT-1983)
310
311 // --- USE NORMAL DISTRIBUTION FOR <X> > 9.9 ---
312 if (xav > G4float(9.9)) {
313 // ** NORMAL DISTRIBUTION WITH SIGMA**2 = <X>
314 Normal(&ran1);
315 ran1 = xav + ran1 * std::sqrt(xav);
316 *iran = G4int(ran1);
317 if (*iran < 0) {
318 *iran = 0;
319 }
320 }
321 else {
322 fivex = G4int(xav * G4float(5.));
323 *iran = 0;
324 if (fivex > 0) {
325 r = std::exp(-G4double(xav));
326 ran1 = G4UniformRand();
327 if (ran1 > r) {
328 rr = r;
329 for (i = 1; i <= fivex; ++i) {
330 ++(*iran);
331 if (i <= 5) {
332 rrr = std::pow(xav, G4float(i)) / NFac(i);
333 }
334 // ** STIRLING' S FORMULA FOR LARGE NUMBERS
335 if (i > 5) {
336 rrr = std::exp(i * std::log(xav) -
337 (i + G4float(.5)) * std::log(i * G4float(1.)) +
338 i - G4float(.9189385));
339 }
340 rr += r * rrr;
341 if (ran1 <= rr) {
342 break;
343 }
344 }
345 }
346 }
347 else {
348 // ** FOR VERY SMALL XAV TRY IRAN=1,2,3
349 p1 = xav * std::exp(-G4double(xav));
350 p2 = xav * p1 / G4float(2.);
351 p3 = xav * p2 / G4float(3.);
352 ran = G4UniformRand();
353 if (ran >= p3) {
354 if (ran >= p2) {
355 if (ran >= p1) {
356 *iran = 0;
357 }
358 else {
359 *iran = 1;
360 }
361 }
362 else {
363 *iran = 2;
364 }
365 }
366 else {
367 *iran = 3;
368 }
369 }
370 }
371
372} // Poisso
373
374
375G4int G4KaonMinusAbsorption::NFac(G4int n)
376{
377 G4int ret_val;
378
379 static G4int i, j;
380
381 // *** NVE 16-MAR-1988 CERN GENEVA ***
382 // ORIGIN : H.FESEFELDT (27-OCT-1983)
383
384 ret_val = 1;
385 j = n;
386 if (j > 1) {
387 if (j > 10) {
388 j = 10;
389 }
390 for (i = 2; i <= j; ++i) {
391 ret_val *= i;
392 }
393 }
394 return ret_val;
395
396} // NFac
397
398
399void G4KaonMinusAbsorption::Normal(G4float *ran)
400{
401 static G4int i;
402
403 // *** NVE 14-APR-1988 CERN GENEVA ***
404 // ORIGIN : H.FESEFELDT (27-OCT-1983)
405
406 *ran = G4float(-6.);
407 for (i = 1; i <= 12; ++i) {
408 *ran += G4UniformRand();
409 }
410
411} // Normal
412
413
414void G4KaonMinusAbsorption::KaonMinusAbsorption(G4int *nopt)
415{
416 static G4int i;
417 static G4int nt, nbl;
418 static G4float ran, pcm;
419 static G4int isw;
420 static G4float tex;
421 static G4float ran2, tof1, ekin, ekin1, ekin2, black;
422 static G4float pnrat;
423 static G4ParticleDefinition* ipa1;
424 static G4ParticleDefinition* inve;
425
426 // *** CHARGED KAON ABSORPTION BY A NUCLEUS ***
427 // *** NVE 04-MAR-1988 CERN GENEVA ***
428 // ORIGIN : H.FESEFELDT (09-JULY-1987)
429
430 // PRODUCTION OF A HYPERFRAGMENT WITH SUBSEQUENT DECAY
431 // PANOFSKY RATIO (K- P --> LAMBDA PI0/K- P --> LAMBDA GAMMA) = 3/2
432
433 pv[1].SetZero();
434 pv[1].SetMass( massKaonMinus );
435 pv[1].SetKineticEnergyAndUpdate( 0. );
436 pv[1].SetTOF( result.GetTOF() );
437 pv[1].SetParticleDef( result.GetParticleDef() );
438 if (targetAtomicMass <= G4float(1.5)) {
439 ran = G4UniformRand();
440 tof1 = std::log(ran) * G4float(-12.5);
441 tof1 *= G4float(20.);
442 ran = G4UniformRand();
443 isw = 1;
444 if (ran < G4float(.33)) {
445 isw = 2;
446 }
447 *nopt = isw;
448 pv[3].SetZero();
449 pv[3].SetMass( massLambda );
450 pv[3].SetKineticEnergyAndUpdate( 0. );
451 pv[3].SetTOF( result.GetTOF() + tof1 );
452 pv[3].SetParticleDef( pdefLambda );
453 pcm = massKaonMinus + massProton - massLambda;
454 if (isw != 1) {
455 pv[2].SetZero();
456 pv[2].SetMass( massGamma );
457 pv[2].SetKineticEnergyAndUpdate( pcm );
458 pv[2].SetTOF( result.GetTOF() + tof1 );
459 pv[2].SetParticleDef( pdefGamma );
460 }
461 else {
462 pcm = pcm * pcm - massPionZero * massPionZero;
463 if (pcm <= G4float(0.)) {
464 pcm = G4float(0.);
465 }
466 pv[2].SetZero();
467 pv[2].SetEnergy( std::sqrt(pcm + massPionZero * massPionZero) );
468 pv[2].SetMassAndUpdate( massPionZero );
469 pv[2].SetTOF( result.GetTOF() + tof1 );
470 pv[2].SetParticleDef( pdefPionZero );
471 }
472 result = pv[2];
473 if (ntot < MAX_SECONDARIES-1) {
474 eve[ntot++] = pv[3];
475 }
476 }
477 else {
478 // **
479 // ** STAR PRODUCTION FOR PION ABSORPTION IN HEAVY ELEMENTS
480 // **
481 evapEnergy1 = G4float(.3);
482 evapEnergy3 = G4float(.15);
483 nt = 1;
484 tex = evapEnergy1;
485 black = std::log(targetAtomicMass) * G4float(.5);
486 Poisso(black, &nbl);
487 if (nt + nbl > (MAX_SECONDARIES - 2)) {
488 nbl = (MAX_SECONDARIES - 2) - nt;
489 }
490 if (nbl <= 0) {
491 nbl = 1;
492 }
493 ekin = tex / nbl;
494 ekin2 = G4float(0.);
495 for (i = 1; i <= nbl; ++i) {
496 if (nt == (MAX_SECONDARIES - 2)) {
497 continue;
498 }
499 ran2 = G4UniformRand();
500 ekin1 = -G4double(ekin) * std::log(ran2);
501 ekin2 += ekin1;
502 ipa1 = pdefNeutron;
503 pnrat = G4float(1.) - targetCharge / targetAtomicMass;
504 if (G4UniformRand() > pnrat) {
505 ipa1 = pdefProton;
506 }
507 ++nt;
508 pv[nt].SetZero();
509 pv[nt].SetMass( ipa1->GetPDGMass()/GeV );
510 pv[nt].SetKineticEnergyAndUpdate( ekin1 );
511 pv[nt].SetTOF( 2. );
512 pv[nt].SetParticleDef( ipa1 );
513 if (ekin2 > tex) {
514 break;
515 }
516 }
517 tex = evapEnergy3;
518 black = std::log(targetAtomicMass) * G4float(.5);
519 Poisso(black, &nbl);
520 if (nt + nbl > (MAX_SECONDARIES - 2)) {
521 nbl = (MAX_SECONDARIES - 2) - nt;
522 }
523 if (nbl <= 0) {
524 nbl = 1;
525 }
526 ekin = tex / nbl;
527 ekin2 = G4float(0.);
528 for (i = 1; i <= nbl; ++i) {
529 if (nt == (MAX_SECONDARIES - 2)) {
530 continue;
531 }
532 ran2 = G4UniformRand();
533 ekin1 = -G4double(ekin) * std::log(ran2);
534 ekin2 += ekin1;
535 ++nt;
536 ran = G4UniformRand();
537 inve = pdefDeuteron;
538 if (ran > G4float(.6)) {
539 inve = pdefTriton;
540 }
541 if (ran > G4float(.9)) {
542 inve = pdefAlpha;
543 }
544// PV(5,NT)=(ABS(IPA(NT))-28)*RMASS(14) <-- Wrong! (LF)
545 pv[nt].SetZero();
546 pv[nt].SetMass( inve->GetPDGMass()/GeV );
547 pv[nt].SetKineticEnergyAndUpdate( ekin1 );
548 pv[nt].SetTOF( 2. );
549 pv[nt].SetParticleDef( inve );
550 if (ekin2 > tex) {
551 break;
552 }
553 }
554 // **
555 // ** STORE ON EVENT COMMON
556 // **
557 ran = G4UniformRand();
558 tof1 = std::log(ran) * G4float(-12.5);
559 tof1 *= G4float(20.);
560 for (i = 2; i <= nt; ++i) {
561 pv[i].SetTOF( result.GetTOF() + tof1 );
562 }
563 result = pv[2];
564 for (i = 3; i <= nt; ++i) {
565 if (ntot >= MAX_SECONDARIES) {
566 break;
567 }
568 eve[ntot++] = pv[i];
569 }
570 }
571
572} // KaonMinusAbsorption
#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 SetMassAndUpdate(G4double mas)
void SetParticleDef(G4ParticleDefinition *c)
void SetKineticEnergyAndUpdate(G4double ekin)
void DeRegisterExtraProcess(G4VProcess *)
void RegisterExtraProcess(G4VProcess *)
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)
static G4HadronicProcessStore * Instance()
void PrintInfo(const G4ParticleDefinition *)
G4GHEKinematicsVector * GetSecondaryKinematics()
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetMeanLifeTime(const G4Track &, G4ForceCondition *)
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
void PreparePhysicsTable(const G4ParticleDefinition &)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
G4bool IsApplicable(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 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