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

#include <G4KL3DecayChannel.hh>

Inheritance diagram for G4KL3DecayChannel:
G4VDecayChannel

Public Member Functions

virtual G4DecayProductsDecayIt (G4double)
 
void DumpInfo ()
 
 G4KL3DecayChannel (const G4String &theParentName, G4double theBR, const G4String &thePionName, const G4String &theLeptonName, const G4String &theNutrinoName)
 
G4int GetAngularMomentum ()
 
G4double GetBR () const
 
G4double GetDalitzParameterLambda () const
 
G4double GetDalitzParameterXi () 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 SetDalitzParameter (G4double aLambda, G4double aXi)
 
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 SetPolarization (const G4ThreeVector &)
 
void SetRangeMass (G4double val)
 
void SetVerboseLevel (G4int value)
 
virtual ~G4KL3DecayChannel ()
 

Protected Types

enum  { idPi =0 , idLepton =1 , idNutrino =2 }
 

Protected Member Functions

void CheckAndFillDaughters ()
 
void CheckAndFillParent ()
 
void ClearDaughtersName ()
 
G4double DalitzDensity (G4double parentmass, G4double Epi, G4double El, G4double Enu, G4double massPi, G4double massL, G4double massNu)
 
G4double DynamicalMass (G4double massPDG, G4double width, G4double maxDev=1.0) const
 
 G4KL3DecayChannel (const G4KL3DecayChannel &)
 
G4KL3DecayChanneloperator= (const G4KL3DecayChannel &)
 
void PhaseSpace (G4double Mparent, const G4double *Mdaughter, G4double *Edaughter, G4double *Pdaughter)
 

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 ()
 
 G4KL3DecayChannel ()
 
const G4StringGetNoName () const
 

Private Attributes

G4double pLambda = 0.0
 
G4double pXi0 = 0.0
 

Detailed Description

Definition at line 37 of file G4KL3DecayChannel.hh.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
protected
Enumerator
idPi 
idLepton 
idNutrino 

Definition at line 61 of file G4KL3DecayChannel.hh.

Constructor & Destructor Documentation

◆ G4KL3DecayChannel() [1/3]

G4KL3DecayChannel::G4KL3DecayChannel ( const G4String theParentName,
G4double  theBR,
const G4String thePionName,
const G4String theLeptonName,
const G4String theNutrinoName 
)

Definition at line 47 of file G4KL3DecayChannel.cc.

53 : G4VDecayChannel("KL3 Decay", theParentName, theBR, 3,
54 thePionName, theLeptonName, theNutrinoName)
55{
56 static const G4String K_plus("kaon+");
57 static const G4String K_minus("kaon-");
58 static const G4String K_L("kaon0L");
59 static const G4String Mu_plus("mu+");
60 static const G4String Mu_minus("mu-");
61 static const G4String E_plus("e+");
62 static const G4String E_minus("e-");
63
64 // check modes
65 if ( ((theParentName == K_plus)&&(theLeptonName == E_plus)) ||
66 ((theParentName == K_minus)&&(theLeptonName == E_minus)) )
67 {
68 // K+- (Ke3)
69 pLambda = 0.0286;
70 pXi0 = -0.35;
71 }
72 else if ( ((theParentName == K_plus)&&(theLeptonName == Mu_plus)) ||
73 ((theParentName == K_minus)&&(theLeptonName == Mu_minus)) )
74 {
75 // K+- (Kmu3)
76 pLambda = 0.033;
77 pXi0 = -0.35;
78 }
79 else if ( (theParentName == K_L) &&
80 ((theLeptonName == E_plus) || (theLeptonName == E_minus)) )
81 {
82 // K0L (Ke3)
83 pLambda = 0.0300;
84 pXi0 = -0.11;
85 }
86 else if ( (theParentName == K_L) &&
87 ((theLeptonName == Mu_plus) || (theLeptonName == Mu_minus)) )
88 {
89 // K0L (Kmu3)
90 pLambda = 0.034;
91 pXi0 = -0.11;
92 }
93 else
94 {
95#ifdef G4VERBOSE
96 if (GetVerboseLevel()>2)
97 {
98 G4cout << "G4KL3DecayChannel:: constructor :";
99 G4cout << "illegal arguments " << G4endl;;
100 DumpInfo();
101 }
102#endif
103 // set values for K0L (Ke3) temporarily
104 pLambda = 0.0300;
105 pXi0 = -0.11;
106 }
107}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4int GetVerboseLevel() const

References G4VDecayChannel::DumpInfo(), G4cout, G4endl, G4VDecayChannel::GetVerboseLevel(), pLambda, and pXi0.

◆ ~G4KL3DecayChannel()

G4KL3DecayChannel::~G4KL3DecayChannel ( )
virtual

Definition at line 109 of file G4KL3DecayChannel.cc.

110{
111}

◆ G4KL3DecayChannel() [2/3]

G4KL3DecayChannel::G4KL3DecayChannel ( const G4KL3DecayChannel right)
protected

Definition at line 113 of file G4KL3DecayChannel.cc.

114 : G4VDecayChannel(right),
115 pLambda(right.pLambda),
116 pXi0(right.pXi0)
117{
118}

◆ G4KL3DecayChannel() [3/3]

G4KL3DecayChannel::G4KL3DecayChannel ( )
private

Definition at line 41 of file G4KL3DecayChannel.cc.

43{
44}

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=(), operator=(), G4MuonDecayChannel::operator=(), G4MuonDecayChannelWithSpin::operator=(), G4MuonRadiativeDecayChannelWithSpin::operator=(), G4NeutronBetaDecayChannel::operator=(), G4PionRadiativeDecayChannel::operator=(), G4TauLeptonicDecayChannel::operator=(), G4VDecayChannel::operator=(), G4VDecayChannel::SetNumberOfDaughters(), and G4VDecayChannel::~G4VDecayChannel().

◆ DalitzDensity()

G4double G4KL3DecayChannel::DalitzDensity ( G4double  parentmass,
G4double  Epi,
G4double  El,
G4double  Enu,
G4double  massPi,
G4double  massL,
G4double  massNu 
)
protected

Definition at line 333 of file G4KL3DecayChannel.cc.

336{
337 // KL3 decay - Dalitz Plot Density, see Chounet et al Phys. Rep. 4, 201
338 // Arguments
339 // Epi: kinetic enregy of pion
340 // El: kinetic enregy of lepton (e or mu)
341 // Enu: kinetic energy of nutrino
342 // Constants
343 // pLambda : linear energy dependence of f+
344 // pXi0 : = f+(0)/f-
345 // pNorm : normalization factor
346 // Variables
347 // Epi: total energy of pion
348 // El: total energy of lepton (e or mu)
349 // Enu: total energy of nutrino
350
351 // calculate total energy
352 Epi = Epi + massPi;
353 El = El + massL;
354 Enu = Enu + massNu;
355
356 G4double Epi_max = (massK*massK+massPi*massPi-massL*massL)/2.0/massK;
357 G4double E = Epi_max - Epi;
358 G4double q2 = massK*massK + massPi*massPi - 2.0*massK*Epi;
359
360 G4double F = 1.0 + pLambda*q2/massPi/massPi;
361 G4double Fmax = 1.0;
362 if (pLambda >0.0) Fmax = (1.0 + pLambda*(massK*massK/massPi/massPi+1.0));
363
364 G4double Xi = pXi0*(1.0 + pLambda*q2/massPi/massPi);
365
366 G4double coeffA = massK*(2.0*El*Enu-massK*E)+massL*massL*(E/4.0-Enu);
367 G4double coeffB = massL*massL*(Enu-E/2.0);
368 G4double coeffC = massL*massL*E/4.0;
369
370 G4double RhoMax = (Fmax*Fmax)*(massK*massK*massK/8.0);
371
372 G4double Rho = (F*F)*(coeffA + coeffB*Xi + coeffC*Xi*Xi);
373
374#ifdef G4VERBOSE
375 if (GetVerboseLevel()>2)
376 {
377 G4cout << "G4KL3DecayChannel::DalitzDensity " <<G4endl;
378 G4cout << " Pi[" << massPi/GeV <<"GeV/c/c] :" << Epi/GeV << "GeV" <<G4endl;
379 G4cout << " L[" << massL/GeV <<"GeV/c/c] :" << El/GeV << "GeV" <<G4endl;
380 G4cout << " Nu[" << massNu/GeV <<"GeV/c/c] :" << Enu/GeV << "GeV" <<G4endl;
381 G4cout << " F :" << F << " Fmax :" << Fmax << " Xi :" << Xi << G4endl;
382 G4cout << " A :" << coeffA << " B :" << coeffB << " C :"<< coeffC
383 << G4endl;
384 G4cout << " Rho :" << Rho << " RhoMax :" << RhoMax << G4endl;
385 }
386#endif
387 return (Rho/RhoMax);
388}
static constexpr double GeV
Definition: G4SIunits.hh:203
double G4double
Definition: G4Types.hh:83

References G4cout, G4endl, G4VDecayChannel::GetVerboseLevel(), GeV, pLambda, and pXi0.

Referenced by DecayIt().

◆ DecayIt()

G4DecayProducts * G4KL3DecayChannel::DecayIt ( G4double  )
virtual

Implements G4VDecayChannel.

Definition at line 152 of file G4KL3DecayChannel.cc.

153{
154 // this version neglects muon polarization
155 // assumes the pure V-A coupling
156 // gives incorrect energy spectrum for Neutrinos
157#ifdef G4VERBOSE
158 if (GetVerboseLevel()>1) G4cout << "G4KL3DecayChannel::DecayIt " << G4endl;
159#endif
160
161 // fill parent particle and its mass
164
165 // fill daughter particles and their mass
167 G4double daughterM[3];
168 daughterM[idPi] = G4MT_daughters[idPi]->GetPDGMass();
171
172 // determine momentum/energy of daughters according to DalitzDensity
173 G4double daughterP[3], daughterE[3];
174 G4double w;
175 G4double r;
176 const size_t MAX_LOOP = 10000;
177 for (std::size_t loop_counter=0; loop_counter<MAX_LOOP; ++loop_counter)
178 {
179 r = G4UniformRand();
180 PhaseSpace(massK, &daughterM[0], &daughterE[0], &daughterP[0]);
181 w = DalitzDensity(massK, daughterE[idPi], daughterE[idLepton],
182 daughterE[idNutrino], daughterM[idPi],
183 daughterM[idLepton], daughterM[idNutrino]);
184 if ( r <= w) break;
185 }
186
187 // output message
188#ifdef G4VERBOSE
189 if (GetVerboseLevel()>1)
190 {
191 G4cout << *daughters_name[0] << ":" << daughterP[0]/GeV << "[GeV/c]"
192 << G4endl;
193 G4cout << *daughters_name[1] << ":" << daughterP[1]/GeV << "[GeV/c]"
194 << G4endl;
195 G4cout << *daughters_name[2] << ":" << daughterP[2]/GeV << "[GeV/c]"
196 << G4endl;
197 }
198#endif
199
200 // create parent G4DynamicParticle at rest
201 G4ThreeVector* direction = new G4ThreeVector(1.0,0.0,0.0);
202 G4DynamicParticle* parentparticle
203 = new G4DynamicParticle( G4MT_parent, *direction, 0.0);
204 delete direction;
205
206 // create G4Decayproducts
207 G4DecayProducts *products
208 = new G4DecayProducts(*parentparticle);
209 delete parentparticle;
210
211 // create daughter G4DynamicParticle
212 G4double costheta, sintheta, phi, sinphi, cosphi;
213 G4double costhetan, sinthetan, phin, sinphin, cosphin;
214
215 // pion
216 costheta = 2.*G4UniformRand()-1.0;
217 sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
218 phi = twopi*G4UniformRand()*rad;
219 sinphi = std::sin(phi);
220 cosphi = std::cos(phi);
221 direction = new G4ThreeVector(sintheta*cosphi,sintheta*sinphi,costheta);
222 G4ThreeVector momentum0 = (*direction)*daughterP[0];
223 G4DynamicParticle * daughterparticle
224 = new G4DynamicParticle(G4MT_daughters[0], momentum0);
225 products->PushProducts(daughterparticle);
226
227 // neutrino
228 costhetan = (daughterP[1]*daughterP[1]
229 - daughterP[2]*daughterP[2]
230 - daughterP[0]*daughterP[0])/(2.0*daughterP[2]*daughterP[0]);
231 sinthetan = std::sqrt((1.0-costhetan)*(1.0+costhetan));
232 phin = twopi*G4UniformRand()*rad;
233 sinphin = std::sin(phin);
234 cosphin = std::cos(phin);
235 direction->setX( sinthetan*cosphin*costheta*cosphi
236 - sinthetan*sinphin*sinphi + costhetan*sintheta*cosphi);
237 direction->setY( sinthetan*cosphin*costheta*sinphi
238 + sinthetan*sinphin*cosphi + costhetan*sintheta*sinphi);
239 direction->setZ( -sinthetan*cosphin*sintheta + costhetan*costheta);
240
241 G4ThreeVector momentum2 = (*direction)*daughterP[2];
242 daughterparticle = new G4DynamicParticle( G4MT_daughters[2], momentum2);
243 products->PushProducts(daughterparticle);
244
245 // lepton
246 G4ThreeVector momentum1 = (momentum0 + momentum2) * (-1.0);
247 daughterparticle = new G4DynamicParticle( G4MT_daughters[1], momentum1);
248 products->PushProducts(daughterparticle);
249
250#ifdef G4VERBOSE
251 if (GetVerboseLevel()>1)
252 {
253 G4cout << "G4KL3DecayChannel::DecayIt ";
254 G4cout << " create decay products in rest frame " <<G4endl;
255 G4cout << " decay products address=" << products << G4endl;
256 products->DumpInfo();
257 }
258#endif
259 delete direction;
260 return products;
261}
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double rad
Definition: G4SIunits.hh:129
CLHEP::Hep3Vector G4ThreeVector
#define G4UniformRand()
Definition: Randomize.hh:52
void setY(double)
void setZ(double)
void setX(double)
void DumpInfo() const
G4int PushProducts(G4DynamicParticle *aParticle)
void PhaseSpace(G4double Mparent, const G4double *Mdaughter, G4double *Edaughter, G4double *Pdaughter)
G4double DalitzDensity(G4double parentmass, G4double Epi, G4double El, G4double Enu, G4double massPi, G4double massL, G4double massNu)
void CheckAndFillDaughters()

References G4VDecayChannel::CheckAndFillDaughters(), G4VDecayChannel::CheckAndFillParent(), DalitzDensity(), G4VDecayChannel::daughters_name, G4DecayProducts::DumpInfo(), G4cout, G4endl, G4VDecayChannel::G4MT_daughters, G4VDecayChannel::G4MT_parent, G4UniformRand, G4ParticleDefinition::GetPDGMass(), G4VDecayChannel::GetVerboseLevel(), GeV, idLepton, idNutrino, idPi, PhaseSpace(), G4DecayProducts::PushProducts(), rad, CLHEP::Hep3Vector::setX(), CLHEP::Hep3Vector::setY(), CLHEP::Hep3Vector::setZ(), and twopi.

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

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

◆ GetDalitzParameterLambda()

G4double G4KL3DecayChannel::GetDalitzParameterLambda ( ) const
inline

Definition at line 107 of file G4KL3DecayChannel.hh.

108{
109 return pLambda;
110}

References pLambda.

◆ GetDalitzParameterXi()

G4double G4KL3DecayChannel::GetDalitzParameterXi ( ) const
inline

Definition at line 113 of file G4KL3DecayChannel.hh.

114{
115 return pXi0;
116}

References pXi0.

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

◆ 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=()

G4KL3DecayChannel & G4KL3DecayChannel::operator= ( const G4KL3DecayChannel right)
protected

Definition at line 120 of file G4KL3DecayChannel.cc.

121{
122 if (this != &right)
123 {
126 rbranch = right.rbranch;
127
128 // copy parent name
129 parent_name = new G4String(*right.parent_name);
130
131 // clear daughters_name array
133
134 // recreate array
136 if ( numberOfDaughters >0 )
137 {
140 //copy daughters name
141 for (G4int index=0; index<numberOfDaughters; ++index)
142 {
143 daughters_name[index] = new G4String(*right.daughters_name[index]);
144 }
145 }
146 pLambda = right.pLambda;
147 pXi0 = right.pXi0;
148 }
149 return *this;
150}

References G4VDecayChannel::ClearDaughtersName(), G4VDecayChannel::daughters_name, G4VDecayChannel::kinematics_name, G4VDecayChannel::numberOfDaughters, G4VDecayChannel::parent_name, pLambda, pXi0, G4VDecayChannel::rbranch, and G4VDecayChannel::verboseLevel.

◆ operator==()

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

Definition at line 68 of file G4VDecayChannel.hh.

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

◆ PhaseSpace()

void G4KL3DecayChannel::PhaseSpace ( G4double  Mparent,
const G4double Mdaughter,
G4double Edaughter,
G4double Pdaughter 
)
protected

Definition at line 263 of file G4KL3DecayChannel.cc.

267{
268 // Algorithm in this code was originally written in GDECA3 in GEANT3
269
270 // sum of daughters'mass
271 G4double sumofdaughtermass = 0.0;
272 G4int index;
273 const G4int N_DAUGHTER=3;
274
275 for (index=0; index<N_DAUGHTER; ++index)
276 {
277 sumofdaughtermass += M[index];
278 }
279
280 // calculate daughter momentum. Generate two
281 G4double rd1, rd2, rd;
282 G4double momentummax=0.0, momentumsum = 0.0;
284 const size_t MAX_LOOP=10000;
285 for (std::size_t loop_counter=0; loop_counter<MAX_LOOP; ++loop_counter)
286 {
287 rd1 = G4UniformRand();
288 rd2 = G4UniformRand();
289 if (rd2 > rd1)
290 {
291 rd = rd1;
292 rd1 = rd2;
293 rd2 = rd;
294 }
295 momentummax = 0.0;
296 momentumsum = 0.0;
297 // daughter 0
298 energy = rd2*(parentM - sumofdaughtermass);
299 P[0] = std::sqrt(energy*energy + 2.0*energy*M[0]);
300 E[0] = energy;
301 if ( P[0] >momentummax )momentummax = P[0];
302 momentumsum += P[0];
303 // daughter 1
304 energy = (1.-rd1)*(parentM - sumofdaughtermass);
305 P[1] = std::sqrt(energy*energy + 2.0*energy*M[1]);
306 E[1] = energy;
307 if ( P[1] >momentummax )momentummax = P[1];
308 momentumsum += P[1];
309 // daughter 2
310 energy = (rd1-rd2)*(parentM - sumofdaughtermass);
311 P[2] = std::sqrt(energy*energy + 2.0*energy*M[2]);
312 E[2] = energy;
313 if ( P[2] >momentummax )momentummax = P[2];
314 momentumsum += P[2];
315 if (momentummax <= momentumsum - momentummax ) break;
316 }
317#ifdef G4VERBOSE
318 if (GetVerboseLevel()>2)
319 {
320 G4cout << "G4KL3DecayChannel::PhaseSpace ";
321 G4cout << "Kon mass:" << parentM/GeV << "GeV/c/c" << G4endl;
322 for (index=0; index<3; ++index)
323 {
324 G4cout << index << " : " << M[index]/GeV << "GeV/c/c ";
325 G4cout << " : " << E[index]/GeV << "GeV ";
326 G4cout << " : " << P[index]/GeV << "GeV/c " << G4endl;
327 }
328 }
329#endif
330}
#define M(row, col)
G4double energy(const ThreeVector &p, const G4double m)
static double P[]

References G4INCL::KinematicsUtils::energy(), G4cout, G4endl, G4UniformRand, G4VDecayChannel::GetVerboseLevel(), GeV, M, and P.

Referenced by DecayIt().

◆ SetBR()

void G4VDecayChannel::SetBR ( G4double  value)
inherited

◆ SetDalitzParameter()

void G4KL3DecayChannel::SetDalitzParameter ( G4double  aLambda,
G4double  aXi 
)
inline

Definition at line 100 of file G4KL3DecayChannel.hh.

101{
102 pLambda = aLambda;
103 pXi0 = aXi;
104}

References pLambda, and pXi0.

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

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

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

◆ parentMutex

G4Mutex G4VDecayChannel::parentMutex
protectedinherited

◆ particletable

G4ParticleTable* G4VDecayChannel::particletable = nullptr
protectedinherited

◆ pLambda

G4double G4KL3DecayChannel::pLambda = 0.0
private

◆ pXi0

G4double G4KL3DecayChannel::pXi0 = 0.0
private

◆ rangeMass

G4double G4VDecayChannel::rangeMass = 2.5
protectedinherited

◆ rbranch

G4double G4VDecayChannel::rbranch = 0.0
protectedinherited

◆ verboseLevel

G4int G4VDecayChannel::verboseLevel = 1
protectedinherited

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