Geant4-11
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes | Private Member Functions | Private Attributes
G4PhaseSpaceDecayChannel Class Reference

#include <G4PhaseSpaceDecayChannel.hh>

Inheritance diagram for G4PhaseSpaceDecayChannel:
G4VDecayChannel

Public Types

enum  { MAX_N_DAUGHTERS =5 }
 

Public Member Functions

virtual G4DecayProductsDecayIt (G4double)
 
void DumpInfo ()
 
 G4PhaseSpaceDecayChannel (const G4String &theParentName, G4double theBR, G4int theNumberOfDaughters, const G4String &theDaughterName1, const G4String &theDaughterName2="", const G4String &theDaughterName3="", const G4String &theDaughterName4="", const G4String &theDaughterName5="")
 
 G4PhaseSpaceDecayChannel (G4int Verbose=1)
 
G4int GetAngularMomentum ()
 
G4double GetBR () const
 
G4ParticleDefinitionGetDaughter (G4int anIndex)
 
G4double GetDaughterMass (G4int anIndex) const
 
const G4StringGetDaughterName (G4int anIndex) const
 
const G4StringGetKinematicsName () const
 
G4int GetNumberOfDaughters () const
 
G4ParticleDefinitionGetParent ()
 
G4double GetParentMass () const
 
const G4StringGetParentName () const
 
const G4ThreeVectorGetPolarization () const
 
G4double GetRangeMass () const
 
G4int GetVerboseLevel () const
 
G4bool IsOKWithParentMass (G4double parentMass)
 
G4bool operator!= (const G4VDecayChannel &r) const
 
G4bool operator< (const G4VDecayChannel &right) const
 
G4bool operator== (const G4VDecayChannel &r) const
 
G4bool SampleDaughterMasses ()
 
void SetBR (G4double value)
 
void SetDaughter (G4int anIndex, const G4ParticleDefinition *particle_type)
 
void SetDaughter (G4int anIndex, const G4String &particle_name)
 
G4bool SetDaughterMasses (G4double masses[])
 
void SetNumberOfDaughters (G4int value)
 
void SetParent (const G4ParticleDefinition *particle_type)
 
void SetParent (const G4String &particle_name)
 
void SetPolarization (const G4ThreeVector &)
 
void SetRangeMass (G4double val)
 
void SetVerboseLevel (G4int value)
 
virtual ~G4PhaseSpaceDecayChannel ()
 

Static Public Member Functions

static G4double Pmx (G4double e, G4double p1, G4double p2)
 

Protected Member Functions

void CheckAndFillDaughters ()
 
void CheckAndFillParent ()
 
void ClearDaughtersName ()
 
G4double DynamicalMass (G4double massPDG, G4double width, G4double maxDev=1.0) const
 

Protected Attributes

G4String ** daughters_name = nullptr
 
G4Mutex daughtersMutex
 
G4ParticleDefinition ** G4MT_daughters = nullptr
 
G4doubleG4MT_daughters_mass = nullptr
 
G4doubleG4MT_daughters_width = nullptr
 
G4ParticleDefinitionG4MT_parent = nullptr
 
G4double G4MT_parent_mass = 0.0
 
G4String kinematics_name = ""
 
G4int numberOfDaughters = 0
 
G4Stringparent_name = nullptr
 
G4ThreeVector parent_polarization
 
G4Mutex parentMutex
 
G4ParticleTableparticletable = nullptr
 
G4double rangeMass = 2.5
 
G4double rbranch = 0.0
 
G4int verboseLevel = 1
 

Static Protected Attributes

static const G4String noName = " "
 

Private Member Functions

void FillDaughters ()
 
void FillParent ()
 
const G4StringGetNoName () const
 
G4DecayProductsManyBodyDecayIt ()
 
G4DecayProductsOneBodyDecayIt ()
 
G4DecayProductsThreeBodyDecayIt ()
 
G4DecayProductsTwoBodyDecayIt ()
 

Private Attributes

G4Cache< G4doublecurrent_parent_mass
 
G4double givenDaughterMasses [MAX_N_DAUGHTERS]
 
G4bool useGivenDaughterMass = false
 

Detailed Description

Definition at line 42 of file G4PhaseSpaceDecayChannel.hh.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
MAX_N_DAUGHTERS 

Definition at line 46 of file G4PhaseSpaceDecayChannel.hh.

Constructor & Destructor Documentation

◆ G4PhaseSpaceDecayChannel() [1/2]

G4PhaseSpaceDecayChannel::G4PhaseSpaceDecayChannel ( G4int  Verbose = 1)

Definition at line 42 of file G4PhaseSpaceDecayChannel.cc.

43 : G4VDecayChannel("Phase Space", Verbose)
44{
45}

◆ G4PhaseSpaceDecayChannel() [2/2]

G4PhaseSpaceDecayChannel::G4PhaseSpaceDecayChannel ( const G4String theParentName,
G4double  theBR,
G4int  theNumberOfDaughters,
const G4String theDaughterName1,
const G4String theDaughterName2 = "",
const G4String theDaughterName3 = "",
const G4String theDaughterName4 = "",
const G4String theDaughterName5 = "" 
)

Definition at line 47 of file G4PhaseSpaceDecayChannel.cc.

56 : G4VDecayChannel("Phase Space", theParentName,theBR, theNumberOfDaughters,
57 theDaughterName1, theDaughterName2,
58 theDaughterName3, theDaughterName4, theDaughterName5)
59{
60}

◆ ~G4PhaseSpaceDecayChannel()

G4PhaseSpaceDecayChannel::~G4PhaseSpaceDecayChannel ( )
virtual

Definition at line 63 of file G4PhaseSpaceDecayChannel.cc.

64{
65}

Member Function Documentation

◆ CheckAndFillDaughters()

void G4VDecayChannel::CheckAndFillDaughters ( )
inlineprotectedinherited

◆ CheckAndFillParent()

void G4VDecayChannel::CheckAndFillParent ( )
inlineprotectedinherited

◆ ClearDaughtersName()

void G4VDecayChannel::ClearDaughtersName ( )
protectedinherited

Definition at line 183 of file G4VDecayChannel.cc.

184{
186 if ( daughters_name != nullptr)
187 {
188 if (numberOfDaughters>0)
189 {
190#ifdef G4VERBOSE
191 if (verboseLevel>1)
192 {
193 G4cout << "G4VDecayChannel::ClearDaughtersName() "
194 << " for " << *parent_name << G4endl;
195 }
196#endif
197 for (G4int index=0; index<numberOfDaughters; ++index)
198 {
199 delete daughters_name[index];
200 }
201 }
202 delete [] daughters_name;
203 daughters_name = nullptr;
204 }
205
206 delete [] G4MT_daughters;
207 delete [] G4MT_daughters_mass;
208 delete [] G4MT_daughters_width;
209 G4MT_daughters_width = nullptr;
210 G4MT_daughters = nullptr;
211 G4MT_daughters_mass = nullptr;
212
214}
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4String * parent_name
G4String ** daughters_name
G4double * G4MT_daughters_mass
G4double * G4MT_daughters_width

References G4VDecayChannel::daughters_name, G4VDecayChannel::daughtersMutex, G4cout, G4endl, G4VDecayChannel::G4MT_daughters, G4VDecayChannel::G4MT_daughters_mass, G4VDecayChannel::G4MT_daughters_width, G4VDecayChannel::numberOfDaughters, G4VDecayChannel::parent_name, and G4VDecayChannel::verboseLevel.

Referenced by G4DalitzDecayChannel::operator=(), G4KL3DecayChannel::operator=(), G4MuonDecayChannel::operator=(), G4MuonDecayChannelWithSpin::operator=(), G4MuonRadiativeDecayChannelWithSpin::operator=(), G4NeutronBetaDecayChannel::operator=(), G4PionRadiativeDecayChannel::operator=(), G4TauLeptonicDecayChannel::operator=(), G4VDecayChannel::operator=(), G4VDecayChannel::SetNumberOfDaughters(), and G4VDecayChannel::~G4VDecayChannel().

◆ DecayIt()

G4DecayProducts * G4PhaseSpaceDecayChannel::DecayIt ( G4double  parentMass)
virtual

Implements G4VDecayChannel.

Definition at line 68 of file G4PhaseSpaceDecayChannel.cc.

69{
70#ifdef G4VERBOSE
71 if (GetVerboseLevel()>1)
72 G4cout << "G4PhaseSpaceDecayChannel::DecayIt()" << G4endl;
73#endif
74
75 G4DecayProducts* products = nullptr;
76
79
80 if (parentMass >0.0) current_parent_mass.Put( parentMass );
82
83 switch (numberOfDaughters)
84 {
85 case 0:
86#ifdef G4VERBOSE
87 if (GetVerboseLevel()>0)
88 {
89 G4cout << "G4PhaseSpaceDecayChannel::DecayIt() -";
90 G4cout << " daughters not defined " << G4endl;
91 }
92#endif
93 break;
94 case 1:
95 products = OneBodyDecayIt();
96 break;
97 case 2:
98 products = TwoBodyDecayIt();
99 break;
100 case 3:
101 products = ThreeBodyDecayIt();
102 break;
103 default:
104 products = ManyBodyDecayIt();
105 break;
106 }
107#ifdef G4VERBOSE
108 if ((products == nullptr) && (GetVerboseLevel()>0))
109 {
110 G4cout << "G4PhaseSpaceDecayChannel::DecayIt() - ";
111 G4cout << *parent_name << " cannot decay " << G4endl;
112 DumpInfo();
113 }
114#endif
115 return products;
116}
void Put(const value_type &val) const
Definition: G4Cache.hh:321
G4Cache< G4double > current_parent_mass
G4int GetVerboseLevel() const
G4double G4MT_parent_mass
void CheckAndFillDaughters()

References G4VDecayChannel::CheckAndFillDaughters(), G4VDecayChannel::CheckAndFillParent(), current_parent_mass, G4VDecayChannel::DumpInfo(), G4cout, G4endl, G4VDecayChannel::G4MT_parent_mass, G4VDecayChannel::GetVerboseLevel(), ManyBodyDecayIt(), G4VDecayChannel::numberOfDaughters, OneBodyDecayIt(), G4VDecayChannel::parent_name, G4Cache< VALTYPE >::Put(), ThreeBodyDecayIt(), and TwoBodyDecayIt().

◆ DumpInfo()

void G4VDecayChannel::DumpInfo ( )
inherited

Definition at line 560 of file G4VDecayChannel.cc.

561{
562 G4cout << " BR: " << rbranch << " [" << kinematics_name << "]";
563 G4cout << " : " ;
564 for (G4int index=0; index<numberOfDaughters; ++index)
565 {
566 if(daughters_name[index] != nullptr)
567 {
568 G4cout << " " << *(daughters_name[index]);
569 }
570 else
571 {
572 G4cout << " not defined ";
573 }
574 }
575 G4cout << G4endl;
576}
G4String kinematics_name

References G4VDecayChannel::daughters_name, G4cout, G4endl, G4VDecayChannel::kinematics_name, G4VDecayChannel::numberOfDaughters, and G4VDecayChannel::rbranch.

Referenced by G4GeneralPhaseSpaceDecay::DecayIt(), DecayIt(), and G4KL3DecayChannel::G4KL3DecayChannel().

◆ DynamicalMass()

G4double G4VDecayChannel::DynamicalMass ( G4double  massPDG,
G4double  width,
G4double  maxDev = 1.0 
) const
protectedinherited

Definition at line 585 of file G4VDecayChannel.cc.

587{
588 if (width<=0.0) return massPDG;
589 if (maxDev >rangeMass) maxDev = rangeMass;
590 if (maxDev <=-1.*rangeMass) return massPDG; // cannot calculate
591
592 G4double x = G4UniformRand()*(maxDev+rangeMass) - rangeMass;
594 const std::size_t MAX_LOOP=10000;
595 for (std::size_t loop_counter=0; loop_counter <MAX_LOOP; ++loop_counter)
596 {
597 if ( y * (width*width*x*x + massPDG*massPDG*width*width)
598 <= massPDG*massPDG*width*width ) break;
599 x = G4UniformRand()*(maxDev+rangeMass) - rangeMass;
600 y = G4UniformRand();
601 }
602 G4double mass = massPDG + x*width;
603 return mass;
604}
double G4double
Definition: G4Types.hh:83
#define G4UniformRand()
Definition: Randomize.hh:52

References G4UniformRand, and G4VDecayChannel::rangeMass.

Referenced by ThreeBodyDecayIt(), and TwoBodyDecayIt().

◆ FillDaughters()

void G4VDecayChannel::FillDaughters ( )
privateinherited

Definition at line 308 of file G4VDecayChannel.cc.

309{
311
312 // Double check, check again if another thread has already filled this, in
313 // case do not need to do anything
314 if ( G4MT_daughters != nullptr ) return;
315
316 G4int index;
317
318#ifdef G4VERBOSE
319 if (verboseLevel>1) G4cout << "G4VDecayChannel::FillDaughters()" << G4endl;
320#endif
321 if (G4MT_daughters != nullptr)
322 {
323 delete [] G4MT_daughters;
324 G4MT_daughters = nullptr;
325 }
326
327 // parent mass
329 G4double parentmass = G4MT_parent->GetPDGMass();
330
331 //
332 G4double sumofdaughtermass = 0.0;
333 G4double sumofdaughterwidthsq = 0.0;
334
335 if ((numberOfDaughters <=0) || (daughters_name == nullptr) )
336 {
337#ifdef G4VERBOSE
338 if (verboseLevel>0)
339 {
340 G4cout << "G4VDecayChannel::FillDaughters() - "
341 << "[ " << G4MT_parent->GetParticleName() << " ]"
342 << "numberOfDaughters is not defined yet";
343 }
344#endif
345 G4MT_daughters = nullptr;
346 G4Exception("G4VDecayChannel::FillDaughters()",
347 "PART011", FatalException,
348 "Cannot fill daughters: numberOfDaughters is not defined yet");
349 }
350
351 // create and set the array of pointers to daughter particles
353 if (G4MT_daughters_mass != nullptr) delete [] G4MT_daughters_mass;
354 if (G4MT_daughters_width != nullptr) delete [] G4MT_daughters_width;
357 // loop over all daughters
358 for (index=0; index<numberOfDaughters; ++index)
359 {
360 if (daughters_name[index] == nullptr)
361 {
362 // daughter name is not defined
363#ifdef G4VERBOSE
364 if (verboseLevel>0)
365 {
366 G4cout << "G4VDecayChannel::FillDaughters() - "
367 << "[ " << G4MT_parent->GetParticleName() << " ]"
368 << index << "-th daughter is not defined yet" << G4endl;
369 }
370#endif
371 G4MT_daughters[index] = nullptr;
372 G4Exception("G4VDecayChannel::FillDaughters()",
373 "PART011", FatalException,
374 "Cannot fill daughters: name of daughter is not defined yet");
375 }
376 // search daughter particles in the particle table
378 if (G4MT_daughters[index] == nullptr )
379 {
380 // cannot find the daughter particle
381#ifdef G4VERBOSE
382 if (verboseLevel>0)
383 {
384 G4cout << "G4VDecayChannel::FillDaughters() - "
385 << "[ " << G4MT_parent->GetParticleName() << " ]"
386 << index << ":" << *daughters_name[index]
387 << " is not defined !!" << G4endl;
388 G4cout << " The BR of this decay mode is set to zero." << G4endl;
389 }
390#endif
391 SetBR(0.0);
392 return;
393 }
394#ifdef G4VERBOSE
395 if (verboseLevel>1)
396 {
397 G4cout << index << ":" << *daughters_name[index];
398 G4cout << ":" << G4MT_daughters[index] << G4endl;
399 }
400#endif
402 G4double d_width = G4MT_daughters[index]->GetPDGWidth();
403 G4MT_daughters_width[index] = d_width;
404 sumofdaughtermass += G4MT_daughters[index]->GetPDGMass();
405 sumofdaughterwidthsq += d_width*d_width;
406 } // end loop over all daughters
407
408 // check sum of daghter mass
409 G4double widthMass = std::sqrt(G4MT_parent->GetPDGWidth()
411 +sumofdaughterwidthsq);
412 if ( (G4MT_parent->GetParticleType() != "nucleus")
413 && (numberOfDaughters !=1)
414 && (sumofdaughtermass > parentmass + rangeMass*widthMass) )
415 {
416 // !!! illegal mass !!!
417#ifdef G4VERBOSE
418 if (GetVerboseLevel()>0)
419 {
420 G4cout << "G4VDecayChannel::FillDaughters() - "
421 << "[ " << G4MT_parent->GetParticleName() << " ]"
422 << " Energy/Momentum conserevation breaks " << G4endl;
423 if (GetVerboseLevel()>1)
424 {
425 G4cout << " parent:" << *parent_name
426 << " mass:" << parentmass/GeV << "[GeV/c/c]" << G4endl;
427 for (index=0; index<numberOfDaughters; ++index)
428 {
429 G4cout << " daughter " << index << ":" << *daughters_name[index]
430 << " mass:" << G4MT_daughters[index]->GetPDGMass()/GeV
431 << "[GeV/c/c]" << G4endl;
432 }
433 }
434 }
435#endif
436 }
437}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static constexpr double GeV
Definition: G4SIunits.hh:203
const G4String & GetParticleType() const
G4double GetPDGWidth() const
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
void SetBR(G4double value)
G4ParticleTable * particletable

References G4VDecayChannel::CheckAndFillParent(), G4VDecayChannel::daughters_name, G4VDecayChannel::daughtersMutex, FatalException, G4ParticleTable::FindParticle(), G4cout, G4endl, G4Exception(), G4VDecayChannel::G4MT_daughters, G4VDecayChannel::G4MT_daughters_mass, G4VDecayChannel::G4MT_daughters_width, G4VDecayChannel::G4MT_parent, G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetParticleType(), G4ParticleDefinition::GetPDGMass(), G4ParticleDefinition::GetPDGWidth(), G4VDecayChannel::GetVerboseLevel(), GeV, G4VDecayChannel::numberOfDaughters, G4VDecayChannel::parent_name, G4VDecayChannel::particletable, G4VDecayChannel::rangeMass, G4VDecayChannel::SetBR(), and G4VDecayChannel::verboseLevel.

Referenced by G4VDecayChannel::CheckAndFillDaughters().

◆ FillParent()

void G4VDecayChannel::FillParent ( )
privateinherited

Definition at line 440 of file G4VDecayChannel.cc.

441{
442 G4AutoLock lock(&parentMutex);
443
444 // Double check, check again if another thread has already filled this, in
445 // case do not need to do anything
446 if ( G4MT_parent != nullptr ) return;
447
448 if (parent_name == nullptr)
449 {
450 // parent name is not defined
451#ifdef G4VERBOSE
452 if (verboseLevel>0)
453 {
454 G4cout << "G4VDecayChannel::FillParent() - "
455 << "parent name is not defined !!" << G4endl;
456 }
457#endif
458 G4MT_parent = nullptr;
459 G4Exception("G4VDecayChannel::FillParent()",
460 "PART012", FatalException,
461 "Cannot fill parent: parent name is not defined yet");
462 return;
463 }
464
465 // search parent particle in the particle table
467 if (G4MT_parent == nullptr)
468 {
469 // parent particle does not exist
470#ifdef G4VERBOSE
471 if (verboseLevel>0)
472 {
473 G4cout << "G4VDecayChannel::FillParent() - "
474 << *parent_name << " does not exist !!" << G4endl;
475 }
476#endif
477 G4Exception("G4VDecayChannel::FillParent()",
478 "PART012", FatalException,
479 "Cannot fill parent: parent does not exist");
480 return;
481 }
483}

References FatalException, G4ParticleTable::FindParticle(), G4cout, G4endl, G4Exception(), G4VDecayChannel::G4MT_parent, G4VDecayChannel::G4MT_parent_mass, G4ParticleDefinition::GetPDGMass(), G4VDecayChannel::parent_name, G4VDecayChannel::parentMutex, G4VDecayChannel::particletable, and G4VDecayChannel::verboseLevel.

Referenced by G4VDecayChannel::CheckAndFillParent().

◆ GetAngularMomentum()

G4int G4VDecayChannel::GetAngularMomentum ( )
inherited

Definition at line 492 of file G4VDecayChannel.cc.

493{
494 // determine angular momentum
495
496 // fill pointers to daughter particles if not yet set
498
499 const G4int PiSpin = G4MT_parent->GetPDGiSpin();
500 const G4int PParity = G4MT_parent->GetPDGiParity();
501 if (2==numberOfDaughters) // up to now we can only handle two particle decays
502 {
503 const G4int D1iSpin = G4MT_daughters[0]->GetPDGiSpin();
504 const G4int D1Parity = G4MT_daughters[0]->GetPDGiParity();
505 const G4int D2iSpin = G4MT_daughters[1]->GetPDGiSpin();
506 const G4int D2Parity = G4MT_daughters[1]->GetPDGiParity();
507 const G4int MiniSpin = std::abs (D1iSpin - D2iSpin);
508 const G4int MaxiSpin = D1iSpin + D2iSpin;
509 const G4int lMax = (PiSpin+D1iSpin+D2iSpin)/2; // l is always int
510 G4int lMin;
511#ifdef G4VERBOSE
512 if (verboseLevel>1)
513 {
514 G4cout << "iSpin: " << PiSpin << " -> " << D1iSpin
515 << " + " << D2iSpin << G4endl;
516 G4cout << "2*jmin, 2*jmax, lmax " << MiniSpin << " " << MaxiSpin
517 << " " << lMax << G4endl;
518 }
519#endif
520 for (G4int j=MiniSpin; j<=MaxiSpin; j+=2) // loop over all possible
521 { // spin couplings
522 lMin = std::abs(PiSpin-j)/2;
523#ifdef G4VERBOSE
524 if (verboseLevel>1)
525 G4cout << "-> checking 2*j=" << j << G4endl;
526#endif
527 for (G4int l=lMin; l<=lMax; ++l)
528 {
529#ifdef G4VERBOSE
530 if (verboseLevel>1)
531 G4cout << " checking l=" << l << G4endl;
532#endif
533 if (l%2==0)
534 {
535 if (PParity == D1Parity*D2Parity) // check parity for this l
536 return l;
537 }
538 else
539 {
540 if (PParity == -1*D1Parity*D2Parity) // check parity for this l
541 return l;
542 }
543 }
544 }
545 }
546 else
547 {
548 G4Exception("G4VDecayChannel::GetAngularMomentum()",
549 "PART111", JustWarning,
550 "Sorry, can't handle 3 particle decays (up to now)");
551 return 0;
552 }
553 G4Exception ("G4VDecayChannel::GetAngularMomentum()",
554 "PART111", JustWarning,
555 "Can't find angular momentum for this decay");
556 return 0;
557}
@ JustWarning

References G4VDecayChannel::CheckAndFillDaughters(), G4cout, G4endl, G4Exception(), G4VDecayChannel::G4MT_daughters, G4VDecayChannel::G4MT_parent, G4ParticleDefinition::GetPDGiParity(), G4ParticleDefinition::GetPDGiSpin(), JustWarning, G4VDecayChannel::numberOfDaughters, and G4VDecayChannel::verboseLevel.

◆ GetBR()

G4double G4VDecayChannel::GetBR ( ) const
inlineinherited

◆ GetDaughter()

G4ParticleDefinition * G4VDecayChannel::GetDaughter ( G4int  anIndex)
inlineinherited

Definition at line 201 of file G4VDecayChannel.hh.

202{
203 // pointers to daughter particles are filled, if they are not set yet
205
206 // get the pointer to a daughter particle
207 if ( (anIndex>=0) && (anIndex<numberOfDaughters) )
208 {
209 return G4MT_daughters[anIndex];
210 }
211 else
212 {
213 if (verboseLevel>0)
214 G4cout << "G4VDecayChannel::GetDaughter index out of range "
215 << anIndex << G4endl;
216 return nullptr;
217 }
218}

References G4VDecayChannel::CheckAndFillDaughters(), G4cout, G4endl, G4VDecayChannel::G4MT_daughters, G4VDecayChannel::numberOfDaughters, and G4VDecayChannel::verboseLevel.

Referenced by G4IonTable::CreateIon(), G4KineticTrack::Decay(), G4KineticTrack::G4KineticTrack(), G4HtmlPPReporter::GeneratePropertyTable(), G4TextPPReporter::GeneratePropertyTable(), G4NuclearDecay::GetDaughterNucleus(), G4MTRunManagerKernel::SetUpDecayChannels(), and G4TaskRunManagerKernel::SetUpDecayChannels().

◆ GetDaughterMass()

G4double G4VDecayChannel::GetDaughterMass ( G4int  anIndex) const
inlineinherited

Definition at line 239 of file G4VDecayChannel.hh.

240{
241 if ( (anIndex>=0) && (anIndex<numberOfDaughters) )
242 {
243 return G4MT_daughters_mass[anIndex];
244 }
245 else
246 {
247 if (verboseLevel>0)
248 {
249 G4cout << "G4VDecayChannel::GetDaughterMass ";
250 G4cout << "index out of range " << anIndex << G4endl;
251 }
252 return 0.0;
253 }
254}

References G4cout, G4endl, G4VDecayChannel::G4MT_daughters_mass, G4VDecayChannel::numberOfDaughters, and G4VDecayChannel::verboseLevel.

◆ GetDaughterName()

const G4String & G4VDecayChannel::GetDaughterName ( G4int  anIndex) const
inlineinherited

◆ GetKinematicsName()

const G4String & G4VDecayChannel::GetKinematicsName ( ) const
inlineinherited

◆ GetNoName()

const G4String & G4VDecayChannel::GetNoName ( ) const
privateinherited

Definition at line 579 of file G4VDecayChannel.cc.

580{
581 return noName;
582}
static const G4String noName

References G4VDecayChannel::noName.

Referenced by G4VDecayChannel::GetDaughterName().

◆ GetNumberOfDaughters()

G4int G4VDecayChannel::GetNumberOfDaughters ( ) const
inlineinherited

◆ GetParent()

G4ParticleDefinition * G4VDecayChannel::GetParent ( )
inlineinherited

Definition at line 257 of file G4VDecayChannel.hh.

258{
259 // the pointer to the parent particle is filled, if it is not set yet
261 // get the pointer to the parent particle
262 return G4MT_parent;
263}

References G4VDecayChannel::CheckAndFillParent(), and G4VDecayChannel::G4MT_parent.

Referenced by G4DecayTable::Insert().

◆ GetParentMass()

G4double G4VDecayChannel::GetParentMass ( ) const
inlineinherited

Definition at line 272 of file G4VDecayChannel.hh.

273{
274 return G4MT_parent_mass;
275}

References G4VDecayChannel::G4MT_parent_mass.

◆ GetParentName()

const G4String & G4VDecayChannel::GetParentName ( ) const
inlineinherited

◆ GetPolarization()

const G4ThreeVector & G4VDecayChannel::GetPolarization ( ) const
inlineinherited

Definition at line 331 of file G4VDecayChannel.hh.

332{
333 return parent_polarization;
334}
G4ThreeVector parent_polarization

References G4VDecayChannel::parent_polarization.

◆ GetRangeMass()

G4double G4VDecayChannel::GetRangeMass ( ) const
inlineinherited

Definition at line 316 of file G4VDecayChannel.hh.

317{
318 return rangeMass;
319}

References G4VDecayChannel::rangeMass.

◆ GetVerboseLevel()

G4int G4VDecayChannel::GetVerboseLevel ( ) const
inlineinherited

◆ IsOKWithParentMass()

G4bool G4PhaseSpaceDecayChannel::IsOKWithParentMass ( G4double  parentMass)
virtual

Reimplemented from G4VDecayChannel.

Definition at line 853 of file G4PhaseSpaceDecayChannel.cc.

854{
856 return G4VDecayChannel::IsOKWithParentMass(parentMass);
857
860
861 G4double sumOfDaughterMassMin=0.0;
862 for (G4int index=0; index<numberOfDaughters; ++index)
863 {
864 sumOfDaughterMassMin += givenDaughterMasses[index];
865 }
866 return (parentMass >= sumOfDaughterMassMin);
867}
G4double givenDaughterMasses[MAX_N_DAUGHTERS]
virtual G4bool IsOKWithParentMass(G4double parentMass)

References G4VDecayChannel::CheckAndFillDaughters(), G4VDecayChannel::CheckAndFillParent(), givenDaughterMasses, G4VDecayChannel::IsOKWithParentMass(), G4VDecayChannel::numberOfDaughters, and useGivenDaughterMass.

◆ ManyBodyDecayIt()

G4DecayProducts * G4PhaseSpaceDecayChannel::ManyBodyDecayIt ( )
private

Definition at line 530 of file G4PhaseSpaceDecayChannel.cc.

531{
532 // Algorithm of this code originally written in FORTRAN by M.Asai
533 // **************************************************************
534 // NBODY - N-body phase space Monte-Carlo generator
535 // Makoto Asai - Hiroshima Institute of Technology
536 // Revised release : 19/Apr/1995
537
538 G4int index, index2;
539
540#ifdef G4VERBOSE
541 if (GetVerboseLevel()>1)
542 G4cout << "G4PhaseSpaceDecayChannel::ManyBodyDecayIt()" << G4endl;
543#endif
544
545 // parent mass
546 G4double parentmass = current_parent_mass.Get();
547
548 // parent particle
549 G4ThreeVector dummy;
550 G4DynamicParticle* parentparticle
551 = new G4DynamicParticle( G4MT_parent, dummy, 0.0, parentmass);
552
553 // create G4Decayproducts
554 G4DecayProducts* products = new G4DecayProducts(*parentparticle);
555 delete parentparticle;
556
557 // daughters'mass
558 G4double* daughtermass = new G4double[numberOfDaughters];
559
560 G4double sumofdaughtermass = 0.0;
561 for (index=0; index<numberOfDaughters; ++index)
562 {
564 {
565 daughtermass[index] = G4MT_daughters_mass[index];
566 }
567 else
568 {
569 daughtermass[index] = givenDaughterMasses[index];
570 }
571 sumofdaughtermass += daughtermass[index];
572 }
573 if (sumofdaughtermass >parentmass)
574 {
575#ifdef G4VERBOSE
576 if (GetVerboseLevel()>0)
577 {
578 G4cout << "G4PhaseSpaceDecayChannel::ManyBodyDecayIt()" << G4endl
579 << "Sum of daughter mass is larger than parent mass!" << G4endl;
580 G4cout << "Parent :" << G4MT_parent->GetParticleName()
581 << " " << current_parent_mass.Get()/GeV << G4endl;
582 G4cout << "Daughter 1 :" << G4MT_daughters[0]->GetParticleName()
583 << " " << daughtermass[0]/GeV << G4endl;
584 G4cout << "Daughter 2:" << G4MT_daughters[1]->GetParticleName()
585 << " " << daughtermass[1]/GeV << G4endl;
586 G4cout << "Daughter 3:" << G4MT_daughters[2]->GetParticleName()
587 << " " << daughtermass[2]/GeV << G4endl;
588 G4cout << "Daughter 4:" << G4MT_daughters[3]->GetParticleName()
589 << " " << daughtermass[3]/GeV << G4endl;
590 G4cout << "Daughter 5:" << G4MT_daughters[4]->GetParticleName()
591 << " " << daughtermass[4]/GeV << G4endl;
592 }
593#endif
594 G4Exception("G4PhaseSpaceDecayChannel::ManyBodyDecayIt()",
595 "PART112", JustWarning,
596 "Can not create decay products: sum of daughter mass \
597 is larger than parent mass!");
598 delete [] daughtermass;
599 return products;
600 }
601
602 // Calculate daughter momentum
603 G4double* daughtermomentum = new G4double[numberOfDaughters];
604 G4ThreeVector direction;
605 G4DynamicParticle** daughterparticle;
607 G4double tmas;
608 G4double weight = 1.0;
609 G4int numberOfTry = 0;
610 const std::size_t MAX_LOOP=10000;
611
612 for (std::size_t loop_counter=0; loop_counter<MAX_LOOP; ++loop_counter)
613 {
614 // Generate rundom number in descending order
615 G4double temp;
617 rd[0] = 1.0;
618 for( index=1; index<numberOfDaughters-1; ++index )
619 {
620 rd[index] = G4UniformRand();
621 }
622 rd[ numberOfDaughters -1] = 0.0;
623 for( index=1; index<numberOfDaughters-1; ++index)
624 {
625 for( index2=index+1; index2<numberOfDaughters; ++index2)
626 {
627 if (rd[index] < rd[index2])
628 {
629 temp = rd[index];
630 rd[index] = rd[index2];
631 rd[index2] = temp;
632 }
633 }
634 }
635
636 // calculate virtual mass
637 tmas = parentmass - sumofdaughtermass;
638 temp = sumofdaughtermass;
639 for( index=0; index<numberOfDaughters; ++index )
640 {
641 sm[index] = rd[index]*tmas + temp;
642 temp -= daughtermass[index];
643 if (GetVerboseLevel()>1)
644 {
645 G4cout << index << " rundom number:" << rd[index];
646 G4cout << " virtual mass:" << sm[index]/GeV << "[GeV/c/c]" << G4endl;
647 }
648 }
649 delete [] rd;
650
651 // Calculate daughter momentum
652 weight = 1.0;
653 G4bool smOK=true;
654 for( index=0; index<numberOfDaughters-1 && smOK; ++index )
655 {
656 smOK = (sm[index]-daughtermass[index]- sm[index+1] >=0.);
657 }
658 if (!smOK) continue;
659
660 index = numberOfDaughters-1;
661 daughtermomentum[index]= Pmx(sm[index-1],daughtermass[index-1],sm[index]);
662#ifdef G4VERBOSE
663 if (GetVerboseLevel()>1)
664 {
665 G4cout << " daughter " << index << ":" << *daughters_name[index];
666 G4cout << " momentum:" << daughtermomentum[index]/GeV << "[GeV/c]"
667 << G4endl;
668 }
669#endif
670 for( index=numberOfDaughters-2; index>=0; --index)
671 {
672 // calculate
673 daughtermomentum[index]= Pmx( sm[index],daughtermass[index],sm[index +1]);
674 if(daughtermomentum[index] < 0.0)
675 {
676 // !!! illegal momentum !!!
677#ifdef G4VERBOSE
678 if (GetVerboseLevel()>0)
679 {
680 G4cout << "G4PhaseSpaceDecayChannel::ManyBodyDecayIt()" << G4endl;
681 G4cout << " Cannot calculate daughter momentum " << G4endl;
682 G4cout << " parent:" << *parent_name;
683 G4cout << " mass:" << parentmass/GeV << "[GeV/c/c]" << G4endl;
684 G4cout << " daughter " << index << ":" << *daughters_name[index];
685 G4cout << " mass:" << daughtermass[index]/GeV << "[GeV/c/c]";
686 G4cout << " mass:" << daughtermomentum[index]/GeV << "[GeV/c]"
687 << G4endl;
689 {
690 G4cout << "Daughter Mass is given." << G4endl;
691 }
692 }
693#endif
694 delete [] sm;
695 delete [] daughtermass;
696 delete [] daughtermomentum;
697 delete products;
698 G4Exception("G4PhaseSpaceDecayChannel::ManyBodyDecayIt()",
699 "PART112", JustWarning,
700 "Can not create decay products: sum of daughter mass \
701 is larger than parent mass");
702 return nullptr; // Error detection
703
704 }
705 else
706 {
707 // calculate weight of this events
708 weight *= daughtermomentum[index]/sm[index];
709#ifdef G4VERBOSE
710 if (GetVerboseLevel()>1)
711 {
712 G4cout << " daughter " << index << ":" << *daughters_name[index];
713 G4cout << " momentum:" << daughtermomentum[index]/GeV << "[GeV/c]"
714 << G4endl;
715 }
716#endif
717 }
718 }
719
720#ifdef G4VERBOSE
721 if (GetVerboseLevel()>1)
722 {
723 G4cout << " weight: " << weight << G4endl;
724 }
725#endif
726
727 // exit if number of Try exceeds 100
728 if (++numberOfTry > 100)
729 {
730#ifdef G4VERBOSE
731 if (GetVerboseLevel()>0)
732 {
733 G4cout << "G4PhaseSpaceDecayChannel::ManyBodyDecayIt()" << G4endl;
734 G4cout << "Cannot determine Decay Kinematics " << G4endl;
735 G4cout << "parent : " << *parent_name << G4endl;
736 G4cout << "daughters : ";
737 for(index=0; index<numberOfDaughters; ++index)
738 {
739 G4cout << *daughters_name[index] << " , ";
740 }
741 G4cout << G4endl;
742 }
743#endif
744 G4Exception("G4PhaseSpaceDecayChannel::ManyBodyDecayIt()",
745 "PART113", JustWarning,
746 "Cannot decay: Decay Kinematics cannot be calculated");
747
748 delete [] sm;
749 delete [] daughtermass;
750 delete [] daughtermomentum;
751 delete products;
752 return nullptr; // Error detection
753 }
754 if ( weight < G4UniformRand()) break;
755 }
756
757#ifdef G4VERBOSE
758 if (GetVerboseLevel()>1)
759 {
760 G4cout << "Start calculation of daughters momentum vector " << G4endl;
761 }
762#endif
763
764 G4double costheta, sintheta, phi;
766 daughterparticle = new G4DynamicParticle*[numberOfDaughters];
767
768 index = numberOfDaughters -2;
769 costheta = 2.*G4UniformRand()-1.0;
770 sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
771 phi = twopi*G4UniformRand()*rad;
772 direction.setZ(costheta);
773 direction.setY(sintheta*std::sin(phi));
774 direction.setX(sintheta*std::cos(phi));
775 daughterparticle[index] = new G4DynamicParticle( G4MT_daughters[index],
776 direction*daughtermomentum[index] );
777 daughterparticle[index+1] = new G4DynamicParticle( G4MT_daughters[index+1],
778 direction*(-1.0*daughtermomentum[index]) );
779
780 for (index = numberOfDaughters-3; index >= 0; --index)
781 {
782 // calculate momentum direction
783 costheta = 2.*G4UniformRand()-1.0;
784 sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
785 phi = twopi*G4UniformRand()*rad;
786 direction.setZ(costheta);
787 direction.setY(sintheta*std::sin(phi));
788 direction.setX(sintheta*std::cos(phi));
789
790 // boost already created particles
791 beta = daughtermomentum[index];
792 beta /= std::sqrt( daughtermomentum[index]*daughtermomentum[index]
793 + sm[index+1]*sm[index+1] );
794 for (index2 = index+1; index2<numberOfDaughters; ++index2)
795 {
797 // make G4LorentzVector for secondaries
798 p4 = daughterparticle[index2]->Get4Momentum();
799
800 // boost secondaries to new frame
801 p4.boost( direction.x()*beta, direction.y()*beta, direction.z()*beta);
802
803 // change energy/momentum
804 daughterparticle[index2]->Set4Momentum(p4);
805 }
806 // create daughter G4DynamicParticle
807 daughterparticle[index] = new G4DynamicParticle( G4MT_daughters[index],
808 direction*(-1.0*daughtermomentum[index]));
809 }
810
811 // add daughters to G4Decayproducts
812 for (index = 0; index<numberOfDaughters; ++index)
813 {
814 products->PushProducts(daughterparticle[index]);
815 }
816
817#ifdef G4VERBOSE
818 if (GetVerboseLevel()>1)
819 {
820 G4cout << "G4PhaseSpaceDecayChannel::ManyBodyDecayIt() -";
821 G4cout << " create decay products in rest frame " << G4endl;
822 products->DumpInfo();
823 }
824#endif
825
826 delete [] daughterparticle;
827 delete [] daughtermomentum;
828 delete [] daughtermass;
829 delete [] sm;
830
831 return products;
832}
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double rad
Definition: G4SIunits.hh:129
bool G4bool
Definition: G4Types.hh:86
double z() const
double x() const
void setY(double)
double y() const
void setZ(double)
void setX(double)
HepLorentzVector & boost(double, double, double)
value_type & Get() const
Definition: G4Cache.hh:315
void DumpInfo() const
G4int PushProducts(G4DynamicParticle *aParticle)
G4LorentzVector Get4Momentum() const
void Set4Momentum(const G4LorentzVector &momentum)
static G4double Pmx(G4double e, G4double p1, G4double p2)

References anonymous_namespace{G4PionRadiativeDecayChannel.cc}::beta, CLHEP::HepLorentzVector::boost(), current_parent_mass, G4VDecayChannel::daughters_name, G4DecayProducts::DumpInfo(), G4cout, G4endl, G4Exception(), G4VDecayChannel::G4MT_daughters, G4VDecayChannel::G4MT_daughters_mass, G4VDecayChannel::G4MT_parent, G4UniformRand, G4Cache< VALTYPE >::Get(), G4DynamicParticle::Get4Momentum(), G4ParticleDefinition::GetParticleName(), G4VDecayChannel::GetVerboseLevel(), GeV, givenDaughterMasses, JustWarning, G4VDecayChannel::numberOfDaughters, G4VDecayChannel::parent_name, Pmx(), G4DecayProducts::PushProducts(), rad, G4DynamicParticle::Set4Momentum(), CLHEP::Hep3Vector::setX(), CLHEP::Hep3Vector::setY(), CLHEP::Hep3Vector::setZ(), G4InuclParticleNames::sm, twopi, useGivenDaughterMass, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by DecayIt().

◆ OneBodyDecayIt()

G4DecayProducts * G4PhaseSpaceDecayChannel::OneBodyDecayIt ( )
private

Definition at line 119 of file G4PhaseSpaceDecayChannel.cc.

120{
121#ifdef G4VERBOSE
122 if (GetVerboseLevel()>1)
123 G4cout << "G4PhaseSpaceDecayChannel::OneBodyDecayIt()" << G4endl;
124#endif
125 // parent mass
126 G4double parentmass = current_parent_mass.Get();
127
128 // create parent G4DynamicParticle at rest
129 G4ThreeVector dummy;
130 G4DynamicParticle* parentparticle
131 = new G4DynamicParticle( G4MT_parent, dummy, 0.0, parentmass);
132 // create G4Decayproducts
133 G4DecayProducts *products = new G4DecayProducts(*parentparticle);
134 delete parentparticle;
135
136 // create daughter G4DynamicParticle at rest
137 G4DynamicParticle* daughterparticle
138 = new G4DynamicParticle( G4MT_daughters[0], dummy, 0.0);
139 if (useGivenDaughterMass) daughterparticle->SetMass(givenDaughterMasses[0]);
140 products->PushProducts(daughterparticle);
141
142#ifdef G4VERBOSE
143 if (GetVerboseLevel()>1)
144 {
145 G4cout << "G4PhaseSpaceDecayChannel::OneBodyDecayIt() -";
146 G4cout << " create decay products in rest frame " << G4endl;
147 products->DumpInfo();
148 }
149#endif
150 return products;
151}
void SetMass(G4double mass)

References current_parent_mass, G4DecayProducts::DumpInfo(), G4cout, G4endl, G4VDecayChannel::G4MT_daughters, G4VDecayChannel::G4MT_parent, G4Cache< VALTYPE >::Get(), G4VDecayChannel::GetVerboseLevel(), givenDaughterMasses, G4DecayProducts::PushProducts(), G4DynamicParticle::SetMass(), and useGivenDaughterMass.

Referenced by DecayIt().

◆ operator!=()

G4bool G4VDecayChannel::operator!= ( const G4VDecayChannel r) const
inlineinherited

Definition at line 69 of file G4VDecayChannel.hh.

69{ return (this != &r); }

◆ operator<()

G4bool G4VDecayChannel::operator< ( const G4VDecayChannel right) const
inlineinherited

Definition at line 195 of file G4VDecayChannel.hh.

196{
197 return (this->rbranch < right.rbranch);
198}

References G4VDecayChannel::rbranch.

◆ operator==()

G4bool G4VDecayChannel::operator== ( const G4VDecayChannel r) const
inlineinherited

Definition at line 68 of file G4VDecayChannel.hh.

68{ return (this == &r); }

◆ Pmx()

G4double G4PhaseSpaceDecayChannel::Pmx ( G4double  e,
G4double  p1,
G4double  p2 
)
static

Definition at line 870 of file G4PhaseSpaceDecayChannel.cc.

871{
872 // calcurate momentum of daughter particles in two-body decay
873 G4double ppp = (e+p1+p2)*(e+p1-p2)*(e-p1+p2)*(e-p1-p2)/(4.0*e*e);
874 if (ppp>0) return std::sqrt(ppp);
875 else return -1.;
876}

Referenced by G4DalitzDecayChannel::DecayIt(), ManyBodyDecayIt(), and TwoBodyDecayIt().

◆ SampleDaughterMasses()

G4bool G4PhaseSpaceDecayChannel::SampleDaughterMasses ( )

Definition at line 846 of file G4PhaseSpaceDecayChannel.cc.

847{
848 useGivenDaughterMass = false;
849 return useGivenDaughterMass;
850}

References useGivenDaughterMass.

◆ SetBR()

void G4VDecayChannel::SetBR ( G4double  value)
inherited

◆ SetDaughter() [1/2]

void G4VDecayChannel::SetDaughter ( G4int  anIndex,
const G4ParticleDefinition particle_type 
)
inherited

◆ SetDaughter() [2/2]

void G4VDecayChannel::SetDaughter ( G4int  anIndex,
const G4String particle_name 
)
inherited

Definition at line 234 of file G4VDecayChannel.cc.

235{
236 // check numberOfDaughters is positive
237 if (numberOfDaughters<=0)
238 {
239#ifdef G4VERBOSE
240 if (verboseLevel>0)
241 {
242 G4cout << "G4VDecayChannel::SetDaughter() - "
243 << "Number of daughters is not defined" << G4endl;
244 }
245#endif
246 return;
247 }
248
249 // An analysis of this code, shows that this method is called
250 // only in the constructor of derived classes.
251 // The general idea of this method is probably to support
252 // the possibility to re-define daughters on the fly, however
253 // this design is extremely problematic for MT mode, we thus
254 // require (as practically happens) that the method is called only
255 // at construction, i.e. when G4MT_daughters == 0
256 // moreover this method can be called only after SetNumberOfDaugthers()
257 // has been called (see previous if), in such a case daughters_name != nullptr
258 //
259 if ( daughters_name == nullptr )
260 {
261 G4Exception("G4VDecayChannel::SetDaughter()", "PART112", FatalException,
262 "Trying to add a daughter without specifying number of secondaries!");
263 return;
264 }
265 if ( G4MT_daughters != nullptr )
266 {
267 G4Exception("G4VDecayChannel::SetDaughter()", "PART111", FatalException,
268 "Trying to modify a daughter of a decay channel, \
269 but decay channel already has daughters.");
270 return;
271 }
272
273 // check an index
274 if ( (anIndex<0) || (anIndex>=numberOfDaughters) )
275 {
276#ifdef G4VERBOSE
277 if (verboseLevel>0)
278 {
279 G4cout << "G4VDecayChannel::SetDaughter() - "
280 << "index out of range " << anIndex << G4endl;
281 }
282#endif
283 }
284 else
285 {
286 // fill the name
287 daughters_name[anIndex] = new G4String(particle_name);
288#ifdef G4VERBOSE
289 if (verboseLevel>1)
290 {
291 G4cout << "G4VDecayChannel::SetDaughter[" << anIndex <<"] :";
292 G4cout << daughters_name[anIndex] << ":"
293 << *daughters_name[anIndex]<<G4endl;
294 }
295#endif
296 }
297}

References G4VDecayChannel::daughters_name, FatalException, G4cout, G4endl, G4Exception(), G4VDecayChannel::G4MT_daughters, G4VDecayChannel::numberOfDaughters, and G4VDecayChannel::verboseLevel.

◆ SetDaughterMasses()

G4bool G4PhaseSpaceDecayChannel::SetDaughterMasses ( G4double  masses[])

Definition at line 835 of file G4PhaseSpaceDecayChannel.cc.

836{
837 for (G4int idx=0; idx<numberOfDaughters; ++idx)
838 {
839 givenDaughterMasses[idx] = masses[idx];
840 }
842 return useGivenDaughterMass;
843}

References givenDaughterMasses, G4VDecayChannel::numberOfDaughters, and useGivenDaughterMass.

◆ SetNumberOfDaughters()

void G4VDecayChannel::SetNumberOfDaughters ( G4int  value)
inherited

◆ SetParent() [1/2]

void G4VDecayChannel::SetParent ( const G4ParticleDefinition particle_type)
inherited

◆ SetParent() [2/2]

void G4VDecayChannel::SetParent ( const G4String particle_name)
inlineinherited

Definition at line 278 of file G4VDecayChannel.hh.

279{
280 if (parent_name != nullptr) delete parent_name;
281 parent_name = new G4String(particle_name);
282 G4MT_parent = nullptr;
283}

References G4VDecayChannel::G4MT_parent, and G4VDecayChannel::parent_name.

◆ SetPolarization()

void G4VDecayChannel::SetPolarization ( const G4ThreeVector polar)
inlineinherited

◆ SetRangeMass()

void G4VDecayChannel::SetRangeMass ( G4double  val)
inlineinherited

Definition at line 322 of file G4VDecayChannel.hh.

322{ if(val>=0.) rangeMass=val; }

References G4VDecayChannel::rangeMass.

◆ SetVerboseLevel()

void G4VDecayChannel::SetVerboseLevel ( G4int  value)
inlineinherited

Definition at line 304 of file G4VDecayChannel.hh.

305{
306 verboseLevel = value;
307}

References G4VDecayChannel::verboseLevel.

Referenced by G4Decay::DecayIt(), and G4MuonicAtomDecay::DecayIt().

◆ ThreeBodyDecayIt()

G4DecayProducts * G4PhaseSpaceDecayChannel::ThreeBodyDecayIt ( )
private

Definition at line 292 of file G4PhaseSpaceDecayChannel.cc.

293{
294 // Algorithm of this code originally written in GDECA3 of GEANT3
295
296#ifdef G4VERBOSE
297 if (GetVerboseLevel()>1)
298 G4cout << "G4PhaseSpaceDecayChannel::ThreeBodyDecayIt()" << G4endl;
299#endif
300 // parent mass
301 G4double parentmass = current_parent_mass.Get();
302 // daughters'mass
303 G4double daughtermass[3], daughterwidth[3];
304 G4double sumofdaughtermass = 0.0;
305 G4double sumofdaughterwidthsq = 0.0;
306 G4bool withWidth = false;
307 for (G4int index=0; index<3; ++index)
308 {
309 daughtermass[index] = G4MT_daughters_mass[index];
310 sumofdaughtermass += daughtermass[index];
311 daughterwidth[index] = G4MT_daughters_width[index];
312 sumofdaughterwidthsq += daughterwidth[index]*daughterwidth[index];
313 withWidth = withWidth ||(daughterwidth[index]>1.0e-3*daughtermass[index]);
314 }
315
316 // create parent G4DynamicParticle at rest
317 G4ThreeVector dummy;
318 G4DynamicParticle* parentparticle
319 = new G4DynamicParticle( G4MT_parent, dummy, 0.0, parentmass);
320
321 // create G4Decayproducts
322 G4DecayProducts* products = new G4DecayProducts(*parentparticle);
323 delete parentparticle;
324
326 {
327 if (withWidth)
328 {
329 G4double maxDev = (parentmass - sumofdaughtermass )
330 / std::sqrt(sumofdaughterwidthsq);
331 if (maxDev <= -1.0*rangeMass )
332 {
333#ifdef G4VERBOSE
334 if (GetVerboseLevel()>0)
335 {
336 G4cout << "G4PhaseSpaceDecayChannel::ThreeBodyDecayIt()" << G4endl
337 << "Sum of daughter mass is larger than parent mass!"
338 << G4endl;
339 G4cout << "Parent :" << G4MT_parent->GetParticleName()
340 << " " << current_parent_mass.Get()/GeV << G4endl;
341 G4cout << "Daughter 1 :" << G4MT_daughters[0]->GetParticleName()
342 << " " << daughtermass[0]/GeV << G4endl;
343 G4cout << "Daughter 2:" << G4MT_daughters[1]->GetParticleName()
344 << " " << daughtermass[1]/GeV << G4endl;
345 G4cout << "Daughter 3:" << G4MT_daughters[2]->GetParticleName()
346 << " " << daughtermass[2]/GeV << G4endl;
347 }
348#endif
349 G4Exception("G4PhaseSpaceDecayChannel::ThreeBodyDecayIt()",
350 "PART112", JustWarning,
351 "Cannot create decay products: sum of daughter mass \
352 is larger than parent mass!");
353 return products;
354 }
355 G4double dm1=daughtermass[0];
356 if (daughterwidth[0] > 0.)
357 dm1= DynamicalMass(daughtermass[0],daughterwidth[0], maxDev);
358 G4double dm2= daughtermass[1];
359 if (daughterwidth[1] > 0.)
360 dm2= DynamicalMass(daughtermass[1],daughterwidth[1], maxDev);
361 G4double dm3= daughtermass[2];
362 if (daughterwidth[2] > 0.)
363 dm3= DynamicalMass(daughtermass[2],daughterwidth[2], maxDev);
364 while (dm1+dm2+dm3>parentmass) // Loop checking, 09.08.2015, K.Kurashige
365 {
366 dm1= DynamicalMass(daughtermass[0],daughterwidth[0], maxDev);
367 dm2= DynamicalMass(daughtermass[1],daughterwidth[1], maxDev);
368 dm3= DynamicalMass(daughtermass[2],daughterwidth[2], maxDev);
369 }
370 daughtermass[0] = dm1;
371 daughtermass[1] = dm2;
372 daughtermass[2] = dm3;
373 sumofdaughtermass = dm1+dm2+dm3;
374 }
375 }
376 else
377 {
378 // use given daughter mass;
379 daughtermass[0] = givenDaughterMasses[0];
380 daughtermass[1] = givenDaughterMasses[1];
381 daughtermass[2] = givenDaughterMasses[2];
382 sumofdaughtermass = daughtermass[0] + daughtermass[1] + daughtermass[2];
383 }
384
385 if (sumofdaughtermass >parentmass)
386 {
387#ifdef G4VERBOSE
388 if (GetVerboseLevel()>0)
389 {
390 G4cout << "G4PhaseSpaceDecayChannel::ThreeBodyDecayIt()" << G4endl
391 << "Sum of daughter mass is larger than parent mass!" << G4endl;
392 G4cout << "Parent :" << G4MT_parent->GetParticleName()
393 << " " << current_parent_mass.Get()/GeV << G4endl;
394 G4cout << "Daughter 1 :" << G4MT_daughters[0]->GetParticleName()
395 << " " << daughtermass[0]/GeV << G4endl;
396 G4cout << "Daughter 2:" << G4MT_daughters[1]->GetParticleName()
397 << " " << daughtermass[1]/GeV << G4endl;
398 G4cout << "Daughter 3:" << G4MT_daughters[2]->GetParticleName()
399 << " " << daughtermass[2]/GeV << G4endl;
400 }
402 {
403 G4cout << "Daughter Mass is given." << G4endl;
404 }
405#endif
406 G4Exception("G4PhaseSpaceDecayChannel::ThreeBodyDecayIt",
407 "PART112", JustWarning,
408 "Can not create decay products: sum of daughter mass \
409 is larger than parent mass!");
410 return products;
411 }
412
413 // calculate daughter momentum
414 // Generate two
415 G4double rd1, rd2, rd;
416 G4double daughtermomentum[3];
417 G4double momentummax=0.0, momentumsum = 0.0;
419 const std::size_t MAX_LOOP=10000;
420
421 for (std::size_t loop_counter=0; loop_counter<MAX_LOOP; ++loop_counter)
422 {
423 rd1 = G4UniformRand();
424 rd2 = G4UniformRand();
425 if (rd2 > rd1)
426 {
427 rd = rd1;
428 rd1 = rd2;
429 rd2 = rd;
430 }
431 momentummax = 0.0;
432 momentumsum = 0.0;
433 // daughter 0
434 energy = rd2*(parentmass - sumofdaughtermass);
435 daughtermomentum[0] = std::sqrt(energy*energy
436 + 2.0*energy* daughtermass[0]);
437 if ( daughtermomentum[0] > momentummax )
438 momentummax = daughtermomentum[0];
439 momentumsum += daughtermomentum[0];
440 // daughter 1
441 energy = (1.-rd1)*(parentmass - sumofdaughtermass);
442 daughtermomentum[1] = std::sqrt(energy*energy
443 + 2.0*energy*daughtermass[1]);
444 if ( daughtermomentum[1] >momentummax )momentummax = daughtermomentum[1];
445 momentumsum += daughtermomentum[1];
446 // daughter 2
447 energy = (rd1-rd2)*(parentmass - sumofdaughtermass);
448 daughtermomentum[2] = std::sqrt(energy*energy
449 + 2.0*energy* daughtermass[2]);
450 if ( daughtermomentum[2] > momentummax )
451 momentummax = daughtermomentum[2];
452 momentumsum += daughtermomentum[2];
453 if (momentummax <= momentumsum - momentummax ) break;
454 }
455
456 // output message
457#ifdef G4VERBOSE
458 if (GetVerboseLevel()>1)
459 {
460 G4cout << " daughter 0:" << daughtermomentum[0]/GeV << "[GeV/c]"
461 << G4endl;
462 G4cout << " daughter 1:" << daughtermomentum[1]/GeV << "[GeV/c]"
463 << G4endl;
464 G4cout << " daughter 2:" << daughtermomentum[2]/GeV << "[GeV/c]"
465 << G4endl;
466 G4cout << " momentum sum:" << momentumsum/GeV << "[GeV/c]"
467 << G4endl;
468 }
469#endif
470
471 // create daughter G4DynamicParticle
472 G4double costheta, sintheta, phi, sinphi, cosphi;
473 G4double costhetan, sinthetan, phin, sinphin, cosphin;
474 costheta = 2.*G4UniformRand()-1.0;
475 sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
476 phi = twopi*G4UniformRand()*rad;
477 sinphi = std::sin(phi);
478 cosphi = std::cos(phi);
479
480 G4ThreeVector direction0(sintheta*cosphi,sintheta*sinphi,costheta);
481 G4double Ekin = std::sqrt(daughtermomentum[0]*daughtermomentum[0]
482 + daughtermass[0]*daughtermass[0]) - daughtermass[0];
483 G4DynamicParticle* daughterparticle
484 = new G4DynamicParticle(G4MT_daughters[0],direction0,Ekin,daughtermass[0]);
485 products->PushProducts(daughterparticle);
486
487 costhetan = (daughtermomentum[1]*daughtermomentum[1]
488 - daughtermomentum[2]*daughtermomentum[2]
489 - daughtermomentum[0]*daughtermomentum[0])
490 / (2.0*daughtermomentum[2]*daughtermomentum[0]);
491 sinthetan = std::sqrt((1.0-costhetan)*(1.0+costhetan));
492 phin = twopi*G4UniformRand()*rad;
493 sinphin = std::sin(phin);
494 cosphin = std::cos(phin);
495 G4ThreeVector direction2;
496 direction2.setX( sinthetan*cosphin*costheta*cosphi
497 - sinthetan*sinphin*sinphi + costhetan*sintheta*cosphi);
498 direction2.setY( sinthetan*cosphin*costheta*sinphi
499 + sinthetan*sinphin*cosphi + costhetan*sintheta*sinphi);
500 direction2.setZ( -sinthetan*cosphin*sintheta + costhetan*costheta);
501 G4ThreeVector pmom = daughtermomentum[2]*direction2/direction2.mag();
502 Ekin = std::sqrt(pmom.mag2() + daughtermass[2]*daughtermass[2])
503 - daughtermass[2];
504 daughterparticle = new G4DynamicParticle( G4MT_daughters[2],
505 pmom/pmom.mag(),
506 Ekin, daughtermass[2] );
507 products->PushProducts(daughterparticle);
508
509 pmom = (direction0*daughtermomentum[0]
510 + direction2*(daughtermomentum[2]/direction2.mag()))*(-1.0);
511 Ekin = std::sqrt(pmom.mag2() + daughtermass[1]*daughtermass[1])
512 - daughtermass[1];
513 daughterparticle = new G4DynamicParticle( G4MT_daughters[1],
514 pmom/pmom.mag(),
515 Ekin, daughtermass[1] );
516 products->PushProducts(daughterparticle);
517
518#ifdef G4VERBOSE
519 if (GetVerboseLevel()>1)
520 {
521 G4cout << "G4PhaseSpaceDecayChannel::ThreeBodyDecayIt -";
522 G4cout << " create decay products in rest frame " << G4endl;
523 products->DumpInfo();
524 }
525#endif
526 return products;
527}
double mag2() const
double mag() const
G4double DynamicalMass(G4double massPDG, G4double width, G4double maxDev=1.0) const
G4double energy(const ThreeVector &p, const G4double m)

References current_parent_mass, G4DecayProducts::DumpInfo(), G4VDecayChannel::DynamicalMass(), G4INCL::KinematicsUtils::energy(), G4cout, G4endl, G4Exception(), G4VDecayChannel::G4MT_daughters, G4VDecayChannel::G4MT_daughters_mass, G4VDecayChannel::G4MT_daughters_width, G4VDecayChannel::G4MT_parent, G4UniformRand, G4Cache< VALTYPE >::Get(), G4ParticleDefinition::GetParticleName(), G4VDecayChannel::GetVerboseLevel(), GeV, givenDaughterMasses, JustWarning, CLHEP::Hep3Vector::mag(), CLHEP::Hep3Vector::mag2(), G4DecayProducts::PushProducts(), rad, G4VDecayChannel::rangeMass, CLHEP::Hep3Vector::setX(), CLHEP::Hep3Vector::setY(), CLHEP::Hep3Vector::setZ(), twopi, and useGivenDaughterMass.

Referenced by DecayIt().

◆ TwoBodyDecayIt()

G4DecayProducts * G4PhaseSpaceDecayChannel::TwoBodyDecayIt ( )
private

Definition at line 154 of file G4PhaseSpaceDecayChannel.cc.

155{
156#ifdef G4VERBOSE
157 if (GetVerboseLevel()>1)
158 G4cout << "G4PhaseSpaceDecayChannel::TwoBodyDecayIt()" << G4endl;
159#endif
160 // parent mass
161 G4double parentmass = current_parent_mass.Get();
162
163 // daughters'mass, width
164 G4double daughtermass[2], daughterwidth[2];
165 daughtermass[0] = G4MT_daughters_mass[0];
166 daughtermass[1] = G4MT_daughters_mass[1];
167 daughterwidth[0] = G4MT_daughters_width[0];
168 daughterwidth[1] = G4MT_daughters_width[1];
169
170 // create parent G4DynamicParticle at rest
171 G4ThreeVector dummy;
172 G4DynamicParticle* parentparticle
173 = new G4DynamicParticle( G4MT_parent, dummy, 0.0, parentmass);
174 // create G4Decayproducts
175 G4DecayProducts* products = new G4DecayProducts(*parentparticle);
176 delete parentparticle;
177
179 {
180 G4bool withWidth = (daughterwidth[0]>1.0e-3*daughtermass[0])
181 || (daughterwidth[1]>1.0e-3*daughtermass[1]);
182 if (withWidth)
183 {
184 G4double sumofdaughterwidthsq = daughterwidth[0]*daughterwidth[0]
185 + daughterwidth[1]*daughterwidth[1];
186 // check parent mass and daughter mass
187 G4double maxDev = (parentmass - daughtermass[0] - daughtermass[1] )
188 / std::sqrt(sumofdaughterwidthsq);
189 if (maxDev <= -1.0*rangeMass )
190 {
191#ifdef G4VERBOSE
192 if (GetVerboseLevel()>0)
193 {
194 G4cout << "G4PhaseSpaceDecayChannel::TwoBodyDecayIt()" << G4endl
195 << "Sum of daughter mass is larger than parent mass!"
196 << G4endl;
197 G4cout << "Parent :" << G4MT_parent->GetParticleName()
198 << " " << current_parent_mass.Get()/GeV << G4endl;
199 G4cout << "Daughter 1 :" << G4MT_daughters[0]->GetParticleName()
200 << " " << daughtermass[0]/GeV << G4endl;
201 G4cout << "Daughter 2:" << G4MT_daughters[1]->GetParticleName()
202 << " " << daughtermass[1]/GeV << G4endl;
203 }
204#endif
205 G4Exception("G4PhaseSpaceDecayChannel::TwoBodyDecayIt()",
206 "PART112", JustWarning,
207 "Cannot create decay products: sum of daughter mass is \
208 larger than parent mass!");
209 return products;
210 }
211 G4double dm1=daughtermass[0];
212 if (daughterwidth[0] > 0.)
213 dm1 = DynamicalMass(daughtermass[0],daughterwidth[0], maxDev);
214 G4double dm2= daughtermass[1];
215 if (daughterwidth[1] > 0.)
216 dm2 = DynamicalMass(daughtermass[1],daughterwidth[1], maxDev);
217 while (dm1+dm2>parentmass) // Loop checking, 09.08.2015, K.Kurashige
218 {
219 dm1 = DynamicalMass(daughtermass[0],daughterwidth[0], maxDev);
220 dm2 = DynamicalMass(daughtermass[1],daughterwidth[1], maxDev);
221 }
222 daughtermass[0] = dm1;
223 daughtermass[1] = dm2;
224 }
225 }
226 else
227 {
228 // use given daughter mass;
229 daughtermass[0] = givenDaughterMasses[0];
230 daughtermass[1] = givenDaughterMasses[1];
231 }
232 if (parentmass < daughtermass[0] + daughtermass[1] )
233 {
234#ifdef G4VERBOSE
235 if (GetVerboseLevel()>0)
236 {
237 G4cout << "G4PhaseSpaceDecayChannel::TwoBodyDecayIt()" << G4endl
238 << "Sum of daughter mass is larger than parent mass!" << G4endl;
239 G4cout << "Parent :" << G4MT_parent->GetParticleName()
240 << " " << current_parent_mass.Get()/GeV << G4endl;
241 G4cout << "Daughter 1 :" << G4MT_daughters[0]->GetParticleName()
242 << " " << daughtermass[0]/GeV << G4endl;
243 G4cout << "Daughter 2:" << G4MT_daughters[1]->GetParticleName()
244 << " " << daughtermass[1]/GeV << G4endl;
246 {
247 G4cout << "Daughter Mass is given." << G4endl;
248 }
249 }
250#endif
251 G4Exception("G4PhaseSpaceDecayChannel::TwoBodyDecayIt()",
252 "PART112", JustWarning,
253 "Cannot create decay products: sum of daughter mass is \
254 larger than parent mass!");
255 return products;
256 }
257
258 // calculate daughter momentum
259 G4double daughtermomentum = Pmx(parentmass,daughtermass[0],daughtermass[1]);
260
261 G4double costheta = 2.*G4UniformRand()-1.0;
262 G4double sintheta = std::sqrt((1.0 - costheta)*(1.0 + costheta));
264 G4ThreeVector direction(sintheta*std::cos(phi),
265 sintheta*std::sin(phi), costheta);
266
267 // create daughter G4DynamicParticle
268 G4double Ekin = std::sqrt(daughtermomentum*daughtermomentum
269 + daughtermass[0]*daughtermass[0]) - daughtermass[0];
270 G4DynamicParticle* daughterparticle
271 = new G4DynamicParticle(G4MT_daughters[0],direction,Ekin,daughtermass[0]);
272 products->PushProducts(daughterparticle);
273 Ekin = std::sqrt(daughtermomentum*daughtermomentum
274 + daughtermass[1]*daughtermass[1]) - daughtermass[1];
275 daughterparticle
276 = new G4DynamicParticle(G4MT_daughters[1], -1.0*direction,
277 Ekin, daughtermass[1]);
278 products->PushProducts(daughterparticle);
279
280#ifdef G4VERBOSE
281 if (GetVerboseLevel()>1)
282 {
283 G4cout << "G4PhaseSpaceDecayChannel::TwoBodyDecayIt() -";
284 G4cout << " Create decay products in rest frame " << G4endl;
285 products->DumpInfo();
286 }
287#endif
288 return products;
289}

References current_parent_mass, G4DecayProducts::DumpInfo(), G4VDecayChannel::DynamicalMass(), G4cout, G4endl, G4Exception(), G4VDecayChannel::G4MT_daughters, G4VDecayChannel::G4MT_daughters_mass, G4VDecayChannel::G4MT_daughters_width, G4VDecayChannel::G4MT_parent, G4UniformRand, G4Cache< VALTYPE >::Get(), G4ParticleDefinition::GetParticleName(), G4VDecayChannel::GetVerboseLevel(), GeV, givenDaughterMasses, JustWarning, Pmx(), G4DecayProducts::PushProducts(), rad, G4VDecayChannel::rangeMass, twopi, and useGivenDaughterMass.

Referenced by DecayIt().

Field Documentation

◆ current_parent_mass

G4Cache<G4double> G4PhaseSpaceDecayChannel::current_parent_mass
private

◆ daughters_name

G4String** G4VDecayChannel::daughters_name = nullptr
protectedinherited

◆ daughtersMutex

G4Mutex G4VDecayChannel::daughtersMutex
protectedinherited

◆ G4MT_daughters

G4ParticleDefinition** G4VDecayChannel::G4MT_daughters = nullptr
protectedinherited

◆ G4MT_daughters_mass

G4double* G4VDecayChannel::G4MT_daughters_mass = nullptr
protectedinherited

◆ G4MT_daughters_width

G4double* G4VDecayChannel::G4MT_daughters_width = nullptr
protectedinherited

◆ G4MT_parent

G4ParticleDefinition* G4VDecayChannel::G4MT_parent = nullptr
protectedinherited

◆ G4MT_parent_mass

G4double G4VDecayChannel::G4MT_parent_mass = 0.0
protectedinherited

◆ givenDaughterMasses

G4double G4PhaseSpaceDecayChannel::givenDaughterMasses[MAX_N_DAUGHTERS]
private

◆ kinematics_name

G4String G4VDecayChannel::kinematics_name = ""
protectedinherited

◆ noName

const G4String G4VDecayChannel::noName = " "
staticprotectedinherited

Definition at line 170 of file G4VDecayChannel.hh.

Referenced by G4VDecayChannel::GetNoName().

◆ numberOfDaughters

G4int G4VDecayChannel::numberOfDaughters = 0
protectedinherited

◆ parent_name

G4String* G4VDecayChannel::parent_name = nullptr
protectedinherited

◆ parent_polarization

G4ThreeVector G4VDecayChannel::parent_polarization
protectedinherited

◆ parentMutex

G4Mutex G4VDecayChannel::parentMutex
protectedinherited

◆ particletable

G4ParticleTable* G4VDecayChannel::particletable = nullptr
protectedinherited

◆ rangeMass

G4double G4VDecayChannel::rangeMass = 2.5
protectedinherited

◆ rbranch

G4double G4VDecayChannel::rbranch = 0.0
protectedinherited

◆ useGivenDaughterMass

G4bool G4PhaseSpaceDecayChannel::useGivenDaughterMass = false
private

◆ verboseLevel

G4int G4VDecayChannel::verboseLevel = 1
protectedinherited

The documentation for this class was generated from the following files: