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

#include <G4GeneralPhaseSpaceDecay.hh>

Inheritance diagram for G4GeneralPhaseSpaceDecay:
G4VDecayChannel

Public Member Functions

virtual G4DecayProductsDecayIt (G4double mass=0.0)
 
void DumpInfo ()
 
 G4GeneralPhaseSpaceDecay (const G4String &theParentName, G4double theBR, G4int theNumberOfDaughters, const G4String &theDaughterName1, const G4String &theDaughterName2="", const G4String &theDaughterName3="")
 
 G4GeneralPhaseSpaceDecay (const G4String &theParentName, G4double theParentMass, G4double theBR, G4int theNumberOfDaughters, const G4String &theDaughterName1, const G4String &theDaughterName2, const G4String &theDaughterName3, const G4double *masses)
 
 G4GeneralPhaseSpaceDecay (const G4String &theParentName, G4double theParentMass, G4double theBR, G4int theNumberOfDaughters, const G4String &theDaughterName1, const G4String &theDaughterName2, const G4String &theDaughterName3, const G4String &theDaughterName4, const G4double *masses)
 
 G4GeneralPhaseSpaceDecay (const G4String &theParentName, G4double theParentMass, G4double theBR, G4int theNumberOfDaughters, const G4String &theDaughterName1, const G4String &theDaughterName2="", const G4String &theDaughterName3="")
 
 G4GeneralPhaseSpaceDecay (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
 
virtual G4bool IsOKWithParentMass (G4double parentMass)
 
G4bool operator!= (const G4VDecayChannel &r) const
 
G4bool operator< (const G4VDecayChannel &right) const
 
G4bool operator== (const G4VDecayChannel &r) const
 
void SetBR (G4double value)
 
void SetDaughter (G4int anIndex, const G4ParticleDefinition *particle_type)
 
void SetDaughter (G4int anIndex, const G4String &particle_name)
 
void SetNumberOfDaughters (G4int value)
 
void SetParent (const G4ParticleDefinition *particle_type)
 
void SetParent (const G4String &particle_name)
 
void SetParentMass (const G4double aParentMass)
 
void SetPolarization (const G4ThreeVector &)
 
void SetRangeMass (G4double val)
 
void SetVerboseLevel (G4int value)
 
virtual ~G4GeneralPhaseSpaceDecay ()
 

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
 
G4DecayProductsManyBodyDecayIt ()
 
G4DecayProductsOneBodyDecayIt ()
 
G4DecayProductsThreeBodyDecayIt ()
 
G4DecayProductsTwoBodyDecayIt ()
 

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
 

Private Attributes

G4double parentmass
 
const G4doubletheDaughterMasses
 

Detailed Description

Definition at line 44 of file G4GeneralPhaseSpaceDecay.hh.

Constructor & Destructor Documentation

◆ G4GeneralPhaseSpaceDecay() [1/5]

G4GeneralPhaseSpaceDecay::G4GeneralPhaseSpaceDecay ( G4int  Verbose = 1)

Definition at line 48 of file G4GeneralPhaseSpaceDecay.cc.

48 :
49 G4VDecayChannel("Phase Space", Verbose),
51{
52 if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay:: constructor " << G4endl;
53}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4int GetVerboseLevel() const

References G4cout, G4endl, and G4VDecayChannel::GetVerboseLevel().

◆ G4GeneralPhaseSpaceDecay() [2/5]

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

Definition at line 55 of file G4GeneralPhaseSpaceDecay.cc.

60 :
61 G4VDecayChannel("Phase Space",
62 theParentName,theBR,
63 theNumberOfDaughters,
64 theDaughterName1,
65 theDaughterName2,
66 theDaughterName3),
68{
69 if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay:: constructor " << G4endl;
70
71 // Set the parent particle (resonance) mass to the (default) PDG vale
72 if (G4MT_parent != NULL)
73 {
75 } else {
76 parentmass=0.;
77 }
78
79}
G4ParticleDefinition * G4MT_parent

References G4cout, G4endl, G4VDecayChannel::G4MT_parent, G4ParticleDefinition::GetPDGMass(), G4VDecayChannel::GetVerboseLevel(), and parentmass.

◆ G4GeneralPhaseSpaceDecay() [3/5]

G4GeneralPhaseSpaceDecay::G4GeneralPhaseSpaceDecay ( const G4String theParentName,
G4double  theParentMass,
G4double  theBR,
G4int  theNumberOfDaughters,
const G4String theDaughterName1,
const G4String theDaughterName2 = "",
const G4String theDaughterName3 = "" 
)

Definition at line 81 of file G4GeneralPhaseSpaceDecay.cc.

87 :
88 G4VDecayChannel("Phase Space",
89 theParentName,theBR,
90 theNumberOfDaughters,
91 theDaughterName1,
92 theDaughterName2,
93 theDaughterName3),
94 parentmass(theParentMass),
96{
97 if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay:: constructor " << G4endl;
98}

References G4cout, G4endl, and G4VDecayChannel::GetVerboseLevel().

◆ G4GeneralPhaseSpaceDecay() [4/5]

G4GeneralPhaseSpaceDecay::G4GeneralPhaseSpaceDecay ( const G4String theParentName,
G4double  theParentMass,
G4double  theBR,
G4int  theNumberOfDaughters,
const G4String theDaughterName1,
const G4String theDaughterName2,
const G4String theDaughterName3,
const G4double masses 
)

Definition at line 100 of file G4GeneralPhaseSpaceDecay.cc.

107 :
108 G4VDecayChannel("Phase Space",
109 theParentName,theBR,
110 theNumberOfDaughters,
111 theDaughterName1,
112 theDaughterName2,
113 theDaughterName3),
114 parentmass(theParentMass),
115 theDaughterMasses(masses)
116{
117 if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay:: constructor " << G4endl;
118}

References G4cout, G4endl, and G4VDecayChannel::GetVerboseLevel().

◆ G4GeneralPhaseSpaceDecay() [5/5]

G4GeneralPhaseSpaceDecay::G4GeneralPhaseSpaceDecay ( const G4String theParentName,
G4double  theParentMass,
G4double  theBR,
G4int  theNumberOfDaughters,
const G4String theDaughterName1,
const G4String theDaughterName2,
const G4String theDaughterName3,
const G4String theDaughterName4,
const G4double masses 
)

Definition at line 120 of file G4GeneralPhaseSpaceDecay.cc.

128 :
129 G4VDecayChannel("Phase Space",
130 theParentName,theBR,
131 theNumberOfDaughters,
132 theDaughterName1,
133 theDaughterName2,
134 theDaughterName3,
135 theDaughterName4),
136 parentmass(theParentMass),
137 theDaughterMasses(masses)
138{
139 if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay:: constructor " << G4endl;
140}

References G4cout, G4endl, and G4VDecayChannel::GetVerboseLevel().

◆ ~G4GeneralPhaseSpaceDecay()

G4GeneralPhaseSpaceDecay::~G4GeneralPhaseSpaceDecay ( )
virtual

Definition at line 142 of file G4GeneralPhaseSpaceDecay.cc.

143{
144}

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
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 * G4GeneralPhaseSpaceDecay::DecayIt ( G4double  mass = 0.0)
virtual

Implements G4VDecayChannel.

Definition at line 146 of file G4GeneralPhaseSpaceDecay.cc.

147{
148 if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay::DecayIt ";
149 G4DecayProducts * products = NULL;
150
153
154 switch (numberOfDaughters){
155 case 0:
156 if (GetVerboseLevel()>0) {
157 G4cout << "G4GeneralPhaseSpaceDecay::DecayIt ";
158 G4cout << " daughters not defined " <<G4endl;
159 }
160 break;
161 case 1:
162 products = OneBodyDecayIt();
163 break;
164 case 2:
165 products = TwoBodyDecayIt();
166 break;
167 case 3:
168 products = ThreeBodyDecayIt();
169 break;
170 default:
171 products = ManyBodyDecayIt();
172 break;
173 }
174 if ((products == NULL) && (GetVerboseLevel()>0)) {
175 G4cout << "G4GeneralPhaseSpaceDecay::DecayIt ";
176 G4cout << *parent_name << " can not decay " << G4endl;
177 DumpInfo();
178 }
179 return products;
180}
void CheckAndFillDaughters()

References G4VDecayChannel::CheckAndFillDaughters(), G4VDecayChannel::CheckAndFillParent(), G4VDecayChannel::DumpInfo(), G4cout, G4endl, G4VDecayChannel::GetVerboseLevel(), ManyBodyDecayIt(), G4VDecayChannel::numberOfDaughters, OneBodyDecayIt(), G4VDecayChannel::parent_name, ThreeBodyDecayIt(), and TwoBodyDecayIt().

Referenced by G4KineticTrack::Decay().

◆ 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 DecayIt(), G4PhaseSpaceDecayChannel::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 G4PhaseSpaceDecayChannel::ThreeBodyDecayIt(), and G4PhaseSpaceDecayChannel::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}
G4double G4MT_parent_mass

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 G4GeneralPhaseSpaceDecay::GetParentMass ( ) const
inline

Definition at line 107 of file G4GeneralPhaseSpaceDecay.hh.

108{
109 return parentmass;
110}

References parentmass.

◆ 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 G4VDecayChannel::IsOKWithParentMass ( G4double  parentMass)
virtualinherited

Reimplemented in G4PhaseSpaceDecayChannel.

Definition at line 607 of file G4VDecayChannel.cc.

608{
609 G4double sumOfDaughterMassMin = 0.0;
612 // skip one body decay
613 if (numberOfDaughters==1) return true;
614
615 for (G4int index=0; index<numberOfDaughters; ++index)
616 {
617 sumOfDaughterMassMin +=
619 }
620 return (parentMass >= sumOfDaughterMassMin);
621}

References G4VDecayChannel::CheckAndFillDaughters(), G4VDecayChannel::CheckAndFillParent(), G4VDecayChannel::G4MT_daughters_mass, G4VDecayChannel::G4MT_daughters_width, G4VDecayChannel::numberOfDaughters, and G4VDecayChannel::rangeMass.

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

◆ ManyBodyDecayIt()

G4DecayProducts * G4GeneralPhaseSpaceDecay::ManyBodyDecayIt ( )
protected

Definition at line 381 of file G4GeneralPhaseSpaceDecay.cc.

391{
392 //return value
393 G4DecayProducts *products;
394
395 if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay::ManyBodyDecayIt()"<<G4endl;
396
397 //daughters'mass
398 G4double *daughtermass = new G4double[numberOfDaughters];
399 G4double sumofdaughtermass = 0.0;
400 for (G4int index=0; index<numberOfDaughters; index++){
401 daughtermass[index] = G4MT_daughters[index]->GetPDGMass();
402 sumofdaughtermass += daughtermass[index];
403 }
404
405 //Calculate daughter momentum
406 G4double *daughtermomentum = new G4double[numberOfDaughters];
407 G4ParticleMomentum direction;
408 G4DynamicParticle **daughterparticle;
410 G4double tmas;
411 G4double weight = 1.0;
412 G4int numberOfTry = 0;
413 G4int index1;
414
415 do {
416 //Generate rundom number in descending order
417 G4double temp;
419 rd[0] = 1.0;
420 for(index1 =1; index1 < numberOfDaughters -1; index1++)
421 rd[index1] = G4UniformRand();
422 rd[ numberOfDaughters -1] = 0.0;
423 for(index1 =1; index1 < numberOfDaughters -1; index1++) {
424 for(G4int index2 = index1+1; index2 < numberOfDaughters; index2++) {
425 if (rd[index1] < rd[index2]){
426 temp = rd[index1];
427 rd[index1] = rd[index2];
428 rd[index2] = temp;
429 }
430 }
431 }
432
433 //calcurate virtual mass
434 tmas = parentmass - sumofdaughtermass;
435 temp = sumofdaughtermass;
436 for(index1 =0; index1 < numberOfDaughters; index1++) {
437 sm[index1] = rd[index1]*tmas + temp;
438 temp -= daughtermass[index1];
439 if (GetVerboseLevel()>1) {
440 G4cout << index1 << " rundom number:" << rd[index1];
441 G4cout << " virtual mass:" << sm[index1]/GeV << "[GeV/c/c]" <<G4endl;
442 }
443 }
444 delete [] rd;
445
446 //Calculate daughter momentum
447 weight = 1.0;
448 index1 =numberOfDaughters-1;
449 daughtermomentum[index1]= Pmx( sm[index1-1],daughtermass[index1-1],sm[index1]);
450 if (GetVerboseLevel()>1) {
451 G4cout << " daughter " << index1 << ":" << *daughters_name[index1];
452 G4cout << " momentum:" << daughtermomentum[index1]/GeV << "[GeV/c]" <<G4endl;
453 }
454 for(index1 =numberOfDaughters-2; index1>=0; index1--) {
455 // calculate
456 daughtermomentum[index1]= Pmx( sm[index1],daughtermass[index1], sm[index1 +1]);
457 if(daughtermomentum[index1] < 0.0) {
458 // !!! illegal momentum !!!
459 if (GetVerboseLevel()>0) {
460 G4cout << "G4GeneralPhaseSpaceDecay::ManyBodyDecayIt ";
461 G4cout << " can not calculate daughter momentum " <<G4endl;
462 G4cout << " parent:" << *parent_name;
463 G4cout << " mass:" << parentmass/GeV << "[GeV/c/c]" <<G4endl;
464 G4cout << " daughter " << index1 << ":" << *daughters_name[index1];
465 G4cout << " mass:" << daughtermass[index1]/GeV << "[GeV/c/c]" ;
466 G4cout << " mass:" << daughtermomentum[index1]/GeV << "[GeV/c]" <<G4endl;
467 }
468 delete [] sm;
469 delete [] daughtermass;
470 delete [] daughtermomentum;
471 return NULL; // Error detection
472
473 } else {
474 // calculate weight of this events
475 weight *= daughtermomentum[index1]/sm[index1];
476 if (GetVerboseLevel()>1) {
477 G4cout << " daughter " << index1 << ":" << *daughters_name[index1];
478 G4cout << " momentum:" << daughtermomentum[index1]/GeV << "[GeV/c]" <<G4endl;
479 }
480 }
481 }
482 if (GetVerboseLevel()>1) {
483 G4cout << " weight: " << weight <<G4endl;
484 }
485
486 // exit if number of Try exceeds 100
487 if (numberOfTry++ >100) {
488 if (GetVerboseLevel()>0) {
489 G4cout << "G4GeneralPhaseSpaceDecay::ManyBodyDecayIt: ";
490 G4cout << " can not determine Decay Kinematics " << G4endl;
491 }
492 delete [] sm;
493 delete [] daughtermass;
494 delete [] daughtermomentum;
495 return NULL; // Error detection
496 }
497 } while ( weight > G4UniformRand()); /* Loop checking, 02.11.2015, A.Ribon */
498 if (GetVerboseLevel()>1) {
499 G4cout << "Start calculation of daughters momentum vector "<<G4endl;
500 }
501
502 G4double costheta, sintheta, phi;
504 daughterparticle = new G4DynamicParticle*[numberOfDaughters];
505
506 index1 = numberOfDaughters -2;
507 costheta = 2.*G4UniformRand()-1.0;
508 sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
509 phi = twopi*G4UniformRand()*rad;
510 direction.setZ(costheta);
511 direction.setY(sintheta*std::sin(phi));
512 direction.setX(sintheta*std::cos(phi));
513 daughterparticle[index1] = new G4DynamicParticle( G4MT_daughters[index1], direction*daughtermomentum[index1] );
514 daughterparticle[index1+1] = new G4DynamicParticle( G4MT_daughters[index1+1], direction*(-1.0*daughtermomentum[index1]) );
515
516 for (index1 = numberOfDaughters -3; index1 >= 0; index1--) {
517 //calculate momentum direction
518 costheta = 2.*G4UniformRand()-1.0;
519 sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
520 phi = twopi*G4UniformRand()*rad;
521 direction.setZ(costheta);
522 direction.setY(sintheta*std::sin(phi));
523 direction.setX(sintheta*std::cos(phi));
524
525 // boost already created particles
526 beta = daughtermomentum[index1];
527 beta /= std::sqrt( daughtermomentum[index1]*daughtermomentum[index1] + sm[index1+1]*sm[index1+1] );
528 for (G4int index2 = index1+1; index2<numberOfDaughters; index2++) {
530 // make G4LorentzVector for secondaries
531 p4 = daughterparticle[index2]->Get4Momentum();
532
533 // boost secondaries to new frame
534 p4.boost( direction.x()*beta, direction.y()*beta, direction.z()*beta);
535
536 // change energy/momentum
537 daughterparticle[index2]->Set4Momentum(p4);
538 }
539 //create daughter G4DynamicParticle
540 daughterparticle[index1]= new G4DynamicParticle( G4MT_daughters[index1], direction*(-1.0*daughtermomentum[index1]));
541 }
542
543 //create G4Decayproducts
544 G4DynamicParticle *parentparticle;
545 direction.setX(1.0); direction.setY(0.0); direction.setZ(0.0);
546 parentparticle = new G4DynamicParticle( G4MT_parent, direction, 0.0);
547 products = new G4DecayProducts(*parentparticle);
548 delete parentparticle;
549 for (index1 = 0; index1<numberOfDaughters; index1++) {
550 products->PushProducts(daughterparticle[index1]);
551 }
552 if (GetVerboseLevel()>1) {
553 G4cout << "G4GeneralPhaseSpaceDecay::ManyBodyDecayIt ";
554 G4cout << " create decay products in rest frame " << G4endl;
555 products->DumpInfo();
556 }
557
558 delete [] daughterparticle;
559 delete [] daughtermomentum;
560 delete [] daughtermass;
561 delete [] sm;
562
563 return products;
564}
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double rad
Definition: G4SIunits.hh:129
double z() const
double x() const
void setY(double)
double y() const
void setZ(double)
void setX(double)
HepLorentzVector & boost(double, double, double)
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(), G4VDecayChannel::daughters_name, G4DecayProducts::DumpInfo(), G4cout, G4endl, G4VDecayChannel::G4MT_daughters, G4VDecayChannel::G4MT_parent, G4UniformRand, G4DynamicParticle::Get4Momentum(), G4ParticleDefinition::GetPDGMass(), G4VDecayChannel::GetVerboseLevel(), GeV, G4VDecayChannel::numberOfDaughters, G4VDecayChannel::parent_name, parentmass, Pmx(), G4DecayProducts::PushProducts(), rad, G4DynamicParticle::Set4Momentum(), CLHEP::Hep3Vector::setX(), CLHEP::Hep3Vector::setY(), CLHEP::Hep3Vector::setZ(), G4InuclParticleNames::sm, twopi, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by DecayIt().

◆ OneBodyDecayIt()

G4DecayProducts * G4GeneralPhaseSpaceDecay::OneBodyDecayIt ( )
protected

Definition at line 182 of file G4GeneralPhaseSpaceDecay.cc.

183{
184 if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay::OneBodyDecayIt()"<<G4endl;
185
186// G4double daughtermass = daughters[0]->GetPDGMass();
187
188 //create parent G4DynamicParticle at rest
189 G4ParticleMomentum dummy;
190 G4DynamicParticle * parentparticle = new G4DynamicParticle(G4MT_parent, dummy, 0.0);
191
192 //create G4Decayproducts
193 G4DecayProducts *products = new G4DecayProducts(*parentparticle);
194 delete parentparticle;
195
196 //create daughter G4DynamicParticle at rest
197 G4DynamicParticle * daughterparticle = new G4DynamicParticle(G4MT_daughters[0], dummy, 0.0);
198 products->PushProducts(daughterparticle);
199
200 if (GetVerboseLevel()>1)
201 {
202 G4cout << "G4GeneralPhaseSpaceDecay::OneBodyDecayIt ";
203 G4cout << " create decay products in rest frame " <<G4endl;
204 products->DumpInfo();
205 }
206 return products;
207}

References G4DecayProducts::DumpInfo(), G4cout, G4endl, G4VDecayChannel::G4MT_daughters, G4VDecayChannel::G4MT_parent, G4VDecayChannel::GetVerboseLevel(), and G4DecayProducts::PushProducts().

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 G4GeneralPhaseSpaceDecay::Pmx ( G4double  e,
G4double  p1,
G4double  p2 
)
inlinestatic

Definition at line 120 of file G4GeneralPhaseSpaceDecay.hh.

121{
122 // calculate momentum of daughter particles in two-body decay
123 if (e-p1-p2 < 0 )
124 {
125 throw G4HadronicException(__FILE__, __LINE__, "G4GeneralPhaseSpaceDecay::Pmx energy in cms < mass1+mass2");
126 }
127 G4double ppp = (e+p1+p2)*(e+p1-p2)*(e-p1+p2)*(e-p1-p2)/(4.0*e*e);
128 if (ppp>0) return std::sqrt(ppp);
129 else return -1.;
130}

Referenced by ManyBodyDecayIt(), and TwoBodyDecayIt().

◆ 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.

◆ 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.

◆ SetParentMass()

void G4GeneralPhaseSpaceDecay::SetParentMass ( const G4double  aParentMass)
inline

Definition at line 112 of file G4GeneralPhaseSpaceDecay.hh.

113{
114 parentmass = aParentMass;
115}

References parentmass.

◆ 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 * G4GeneralPhaseSpaceDecay::ThreeBodyDecayIt ( )
protected

Definition at line 259 of file G4GeneralPhaseSpaceDecay.cc.

261{
262 if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay::ThreeBodyDecayIt()"<<G4endl;
263
264 //daughters'mass
265 G4double daughtermass[3];
266 G4double sumofdaughtermass = 0.0;
267 for (G4int index=0; index<3; index++)
268 {
269 if ( theDaughterMasses )
270 {
271 daughtermass[index]= *(theDaughterMasses+index);
272 } else {
273 daughtermass[index] = G4MT_daughters[index]->GetPDGMass();
274 }
275 sumofdaughtermass += daughtermass[index];
276 }
277
278 //create parent G4DynamicParticle at rest
279 G4ParticleMomentum dummy;
280 G4DynamicParticle * parentparticle = new G4DynamicParticle( G4MT_parent, dummy, 0.0);
281
282 //create G4Decayproducts
283 G4DecayProducts *products = new G4DecayProducts(*parentparticle);
284 delete parentparticle;
285
286 //calculate daughter momentum
287 // Generate two
288 G4double rd1, rd2, rd;
289 G4double daughtermomentum[3];
290 G4double momentummax=0.0, momentumsum = 0.0;
292 const G4int maxNumberOfLoops = 10000;
293 G4int loopCounter = 0;
294
295 do
296 {
297 rd1 = G4UniformRand();
298 rd2 = G4UniformRand();
299 if (rd2 > rd1)
300 {
301 rd = rd1;
302 rd1 = rd2;
303 rd2 = rd;
304 }
305 momentummax = 0.0;
306 momentumsum = 0.0;
307 // daughter 0
308
309 energy = rd2*(parentmass - sumofdaughtermass);
310 daughtermomentum[0] = std::sqrt(energy*energy + 2.0*energy* daughtermass[0]);
311 if ( daughtermomentum[0] >momentummax )momentummax = daughtermomentum[0];
312 momentumsum += daughtermomentum[0];
313
314 // daughter 1
315 energy = (1.-rd1)*(parentmass - sumofdaughtermass);
316 daughtermomentum[1] = std::sqrt(energy*energy + 2.0*energy* daughtermass[1]);
317 if ( daughtermomentum[1] >momentummax )momentummax = daughtermomentum[1];
318 momentumsum += daughtermomentum[1];
319
320 // daughter 2
321 energy = (rd1-rd2)*(parentmass - sumofdaughtermass);
322 daughtermomentum[2] = std::sqrt(energy*energy + 2.0*energy* daughtermass[2]);
323 if ( daughtermomentum[2] >momentummax )momentummax = daughtermomentum[2];
324 momentumsum += daughtermomentum[2];
325 } while ( ( momentummax > momentumsum - momentummax ) && /* Loop checking, 02.11.2015, A.Ribon */
326 ++loopCounter < maxNumberOfLoops );
327 if ( loopCounter >= maxNumberOfLoops ) {
329 ed << " Failed sampling after maxNumberOfLoops attempts : forced exit" << G4endl;
330 G4Exception( " G4GeneralPhaseSpaceDecay::ThreeBodyDecayIt ", "HAD_PHASESPACE_001", FatalException, ed );
331 }
332
333 // output message
334 if (GetVerboseLevel()>1) {
335 G4cout << " daughter 0:" << daughtermomentum[0]/GeV << "[GeV/c]" <<G4endl;
336 G4cout << " daughter 1:" << daughtermomentum[1]/GeV << "[GeV/c]" <<G4endl;
337 G4cout << " daughter 2:" << daughtermomentum[2]/GeV << "[GeV/c]" <<G4endl;
338 G4cout << " momentum sum:" << momentumsum/GeV << "[GeV/c]" <<G4endl;
339 }
340
341 //create daughter G4DynamicParticle
342 G4double costheta, sintheta, phi, sinphi, cosphi;
343 G4double costhetan, sinthetan, phin, sinphin, cosphin;
344 costheta = 2.*G4UniformRand()-1.0;
345 sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
346 phi = twopi*G4UniformRand()*rad;
347 sinphi = std::sin(phi);
348 cosphi = std::cos(phi);
349 G4ParticleMomentum direction0(sintheta*cosphi,sintheta*sinphi,costheta);
350 G4double Etotal=std::sqrt( daughtermass[0]*daughtermass[0] + daughtermomentum[0]*daughtermomentum[0]);
351 G4DynamicParticle * daughterparticle
352 = new G4DynamicParticle( G4MT_daughters[0], Etotal, direction0*daughtermomentum[0]);
353 products->PushProducts(daughterparticle);
354
355 costhetan = (daughtermomentum[1]*daughtermomentum[1]-daughtermomentum[2]*daughtermomentum[2]-daughtermomentum[0]*daughtermomentum[0])/(2.0*daughtermomentum[2]*daughtermomentum[0]);
356 sinthetan = std::sqrt((1.0-costhetan)*(1.0+costhetan));
357 phin = twopi*G4UniformRand()*rad;
358 sinphin = std::sin(phin);
359 cosphin = std::cos(phin);
360 G4ParticleMomentum direction2;
361 direction2.setX( sinthetan*cosphin*costheta*cosphi - sinthetan*sinphin*sinphi + costhetan*sintheta*cosphi);
362 direction2.setY( sinthetan*cosphin*costheta*sinphi + sinthetan*sinphin*cosphi + costhetan*sintheta*sinphi);
363 direction2.setZ( -sinthetan*cosphin*sintheta + costhetan*costheta);
364 Etotal=std::sqrt( daughtermass[2]*daughtermass[2] + daughtermomentum[2]*daughtermomentum[2]/direction2.mag2());
365 daughterparticle = new G4DynamicParticle( G4MT_daughters[2],Etotal, direction2*(daughtermomentum[2]/direction2.mag()));
366 products->PushProducts(daughterparticle);
367 G4ThreeVector mom=(direction0*daughtermomentum[0] + direction2*(daughtermomentum[2]/direction2.mag()))*(-1.0);
368 Etotal= std::sqrt( daughtermass[1]*daughtermass[1] + mom.mag2() );
369 daughterparticle =
370 new G4DynamicParticle(G4MT_daughters[1], Etotal, mom);
371 products->PushProducts(daughterparticle);
372
373 if (GetVerboseLevel()>1) {
374 G4cout << "G4GeneralPhaseSpaceDecay::ThreeBodyDecayIt ";
375 G4cout << " create decay products in rest frame " <<G4endl;
376 products->DumpInfo();
377 }
378 return products;
379}
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
double mag2() const
double mag() const
G4double energy(const ThreeVector &p, const G4double m)

References G4DecayProducts::DumpInfo(), G4INCL::KinematicsUtils::energy(), FatalException, G4cout, G4endl, G4Exception(), G4VDecayChannel::G4MT_daughters, G4VDecayChannel::G4MT_parent, G4UniformRand, G4ParticleDefinition::GetPDGMass(), G4VDecayChannel::GetVerboseLevel(), GeV, CLHEP::Hep3Vector::mag(), CLHEP::Hep3Vector::mag2(), parentmass, G4DecayProducts::PushProducts(), rad, CLHEP::Hep3Vector::setX(), CLHEP::Hep3Vector::setY(), CLHEP::Hep3Vector::setZ(), theDaughterMasses, and twopi.

Referenced by DecayIt().

◆ TwoBodyDecayIt()

G4DecayProducts * G4GeneralPhaseSpaceDecay::TwoBodyDecayIt ( )
protected

Definition at line 209 of file G4GeneralPhaseSpaceDecay.cc.

210{
211 if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay::TwoBodyDecayIt()"<<G4endl;
212
213 //daughters'mass
214 G4double daughtermass[2];
215 G4double daughtermomentum;
216 if ( theDaughterMasses )
217 {
218 daughtermass[0]= *(theDaughterMasses);
219 daughtermass[1] = *(theDaughterMasses+1);
220 } else {
221 daughtermass[0] = G4MT_daughters[0]->GetPDGMass();
222 daughtermass[1] = G4MT_daughters[1]->GetPDGMass();
223 }
224
225// G4double sumofdaughtermass = daughtermass[0] + daughtermass[1];
226
227 //create parent G4DynamicParticle at rest
228 G4ParticleMomentum dummy;
229 G4DynamicParticle * parentparticle = new G4DynamicParticle( G4MT_parent, dummy, 0.0);
230
231 //create G4Decayproducts @@GF why dummy parentparticle?
232 G4DecayProducts *products = new G4DecayProducts(*parentparticle);
233 delete parentparticle;
234
235 //calculate daughter momentum
236 daughtermomentum = Pmx(parentmass,daughtermass[0],daughtermass[1]);
237 G4double costheta = 2.*G4UniformRand()-1.0;
238 G4double sintheta = std::sqrt((1.0 - costheta)*(1.0 + costheta));
240 G4ParticleMomentum direction(sintheta*std::cos(phi),sintheta*std::sin(phi),costheta);
241
242 //create daughter G4DynamicParticle
243 G4double Etotal= std::sqrt(daughtermass[0]*daughtermass[0] + daughtermomentum*daughtermomentum);
244 G4DynamicParticle * daughterparticle = new G4DynamicParticle( G4MT_daughters[0],Etotal, direction*daughtermomentum);
245 products->PushProducts(daughterparticle);
246 Etotal= std::sqrt(daughtermass[1]*daughtermass[1] + daughtermomentum*daughtermomentum);
247 daughterparticle = new G4DynamicParticle( G4MT_daughters[1],Etotal, direction*(-1.0*daughtermomentum));
248 products->PushProducts(daughterparticle);
249
250 if (GetVerboseLevel()>1)
251 {
252 G4cout << "G4GeneralPhaseSpaceDecay::TwoBodyDecayIt ";
253 G4cout << " create decay products in rest frame " <<G4endl;
254 products->DumpInfo();
255 }
256 return products;
257}

References G4DecayProducts::DumpInfo(), G4cout, G4endl, G4VDecayChannel::G4MT_daughters, G4VDecayChannel::G4MT_parent, G4UniformRand, G4ParticleDefinition::GetPDGMass(), G4VDecayChannel::GetVerboseLevel(), parentmass, Pmx(), G4DecayProducts::PushProducts(), rad, theDaughterMasses, and twopi.

Referenced by DecayIt().

Field Documentation

◆ 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

◆ 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

◆ parentmass

G4double G4GeneralPhaseSpaceDecay::parentmass
private

◆ 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

◆ theDaughterMasses

const G4double* G4GeneralPhaseSpaceDecay::theDaughterMasses
private

Definition at line 101 of file G4GeneralPhaseSpaceDecay.hh.

Referenced by ThreeBodyDecayIt(), and TwoBodyDecayIt().

◆ verboseLevel

G4int G4VDecayChannel::verboseLevel = 1
protectedinherited

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