#include <G4VDecayChannel.hh>
Inheritance diagram for G4VDecayChannel:
Definition at line 49 of file G4VDecayChannel.hh.
Definition at line 63 of file G4VDecayChannel.cc.
References G4ParticleTable::GetParticleTable(), and particletable.
00064 :kinematics_name(aName), 00065 rbranch(0.0), 00066 numberOfDaughters(0), 00067 parent_name(0), daughters_name(0), 00068 particletable(0), 00069 parent(0), daughters(0), 00070 parent_mass(0.0), daughters_mass(0), 00071 verboseLevel(Verbose) 00072 { 00073 // set pointer to G4ParticleTable (static and singleton object) 00074 particletable = G4ParticleTable::GetParticleTable(); 00075 }
G4VDecayChannel::G4VDecayChannel | ( | const G4String & | aName, | |
const G4String & | theParentName, | |||
G4double | theBR, | |||
G4int | theNumberOfDaughters, | |||
const G4String & | theDaughterName1, | |||
const G4String & | theDaughterName2 = "" , |
|||
const G4String & | theDaughterName3 = "" , |
|||
const G4String & | theDaughterName4 = "" | |||
) |
Definition at line 77 of file G4VDecayChannel.cc.
References daughters_name, G4ParticleTable::GetParticleTable(), numberOfDaughters, parent_name, and particletable.
00085 :kinematics_name(aName), 00086 rbranch(theBR), 00087 numberOfDaughters(theNumberOfDaughters), 00088 parent_name(0), daughters_name(0), 00089 particletable(0), 00090 parent(0), daughters(0), 00091 parent_mass(0.0), daughters_mass(0), 00092 verboseLevel(1) 00093 { 00094 // set pointer to G4ParticleTable (static and singleton object) 00095 particletable = G4ParticleTable::GetParticleTable(); 00096 00097 // parent name 00098 parent_name = new G4String(theParentName); 00099 00100 // cleate array 00101 daughters_name = new G4String*[numberOfDaughters]; 00102 for (G4int index=0;index<numberOfDaughters;index++) daughters_name[index]=0; 00103 00104 // daughters' name 00105 if (numberOfDaughters>0) daughters_name[0] = new G4String(theDaughterName1); 00106 if (numberOfDaughters>1) daughters_name[1] = new G4String(theDaughterName2); 00107 if (numberOfDaughters>2) daughters_name[2] = new G4String(theDaughterName3); 00108 if (numberOfDaughters>3) daughters_name[3] = new G4String(theDaughterName4); 00109 }
G4VDecayChannel::~G4VDecayChannel | ( | ) | [virtual] |
Definition at line 183 of file G4VDecayChannel.cc.
References ClearDaughtersName(), daughters_mass, and parent_name.
00184 { 00185 ClearDaughtersName(); 00186 if (parent_name != 0) delete parent_name; 00187 parent_name = 0; 00188 if (daughters_mass != 0) delete [] daughters_mass; 00189 daughters_mass =0; 00190 }
G4VDecayChannel::G4VDecayChannel | ( | ) | [protected] |
Definition at line 48 of file G4VDecayChannel.cc.
References G4ParticleTable::GetParticleTable(), and particletable.
00049 :kinematics_name(""), 00050 rbranch(0.0), 00051 numberOfDaughters(0), 00052 parent_name(0), daughters_name(0), 00053 particletable(0), 00054 parent(0), daughters(0), 00055 parent_mass(0.0), daughters_mass(0), 00056 verboseLevel(1) 00057 { 00058 // set pointer to G4ParticleTable (static and singleton object) 00059 particletable = G4ParticleTable::GetParticleTable(); 00060 }
G4VDecayChannel::G4VDecayChannel | ( | const G4VDecayChannel & | ) | [protected] |
Definition at line 113 of file G4VDecayChannel.cc.
References daughters, daughters_mass, daughters_name, G4ParticleTable::GetParticleTable(), kinematics_name, numberOfDaughters, parent, parent_mass, parent_name, particletable, rbranch, and verboseLevel.
00114 { 00115 kinematics_name = right.kinematics_name; 00116 verboseLevel = right.verboseLevel; 00117 rbranch = right.rbranch; 00118 00119 // copy parent name 00120 parent_name = new G4String(*right.parent_name); 00121 parent = 0; 00122 parent_mass = 0.0; 00123 00124 //create array 00125 numberOfDaughters = right.numberOfDaughters; 00126 00127 daughters_name =0; 00128 if ( numberOfDaughters >0 ) { 00129 daughters_name = new G4String*[numberOfDaughters]; 00130 //copy daughters name 00131 for (G4int index=0; index < numberOfDaughters; index++) 00132 { 00133 daughters_name[index] = new G4String(*right.daughters_name[index]); 00134 } 00135 } 00136 00137 // 00138 daughters_mass = 0; 00139 daughters = 0; 00140 00141 // particle table 00142 particletable = G4ParticleTable::GetParticleTable(); 00143 }
void G4VDecayChannel::ClearDaughtersName | ( | ) | [protected] |
Definition at line 192 of file G4VDecayChannel.cc.
References daughters, daughters_mass, daughters_name, G4cerr, G4endl, numberOfDaughters, parent_name, and verboseLevel.
Referenced by operator=(), G4TauLeptonicDecayChannel::operator=(), G4PionRadiativeDecayChannel::operator=(), G4NeutronBetaDecayChannel::operator=(), G4MuonRadiativeDecayChannelWithSpin::operator=(), G4MuonDecayChannelWithSpin::operator=(), G4MuonDecayChannel::operator=(), G4KL3DecayChannel::operator=(), G4DalitzDecayChannel::operator=(), SetNumberOfDaughters(), and ~G4VDecayChannel().
00193 { 00194 if ( daughters_name != 0) { 00195 if (numberOfDaughters>0) { 00196 #ifdef G4VERBOSE 00197 if (verboseLevel>1) { 00198 G4cerr << "G4VDecayChannel::ClearDaughtersName " 00199 << " for " << *parent_name << G4endl; 00200 } 00201 #endif 00202 for (G4int index=0; index < numberOfDaughters; index++) { 00203 if (daughters_name[index] != 0) delete daughters_name[index]; 00204 } 00205 } 00206 delete [] daughters_name; 00207 daughters_name = 0; 00208 } 00209 // 00210 if (daughters != 0) delete [] daughters; 00211 if (daughters_mass != 0) delete [] daughters_mass; 00212 daughters = 0; 00213 daughters_mass = 0; 00214 00215 numberOfDaughters = 0; 00216 }
virtual G4DecayProducts* G4VDecayChannel::DecayIt | ( | G4double | parentMass = -1.0 |
) | [pure virtual] |
Implemented in G4DalitzDecayChannel, G4KL3DecayChannel, G4MuonDecayChannel, G4MuonDecayChannelWithSpin, G4MuonRadiativeDecayChannelWithSpin, G4NeutronBetaDecayChannel, G4PhaseSpaceDecayChannel, G4PionRadiativeDecayChannel, G4TauLeptonicDecayChannel, G4NuclearDecayChannel, and G4GeneralPhaseSpaceDecay.
Referenced by G4Decay::DecayIt(), G4IntraNucleiCascader::decayTrappedParticle(), and G4RadioactiveDecay::DoDecay().
void G4VDecayChannel::DumpInfo | ( | ) |
Definition at line 480 of file G4VDecayChannel.cc.
References daughters_name, G4cout, G4endl, kinematics_name, numberOfDaughters, and rbranch.
Referenced by G4PhaseSpaceDecayChannel::DecayIt(), G4NuclearDecayChannel::DecayIt(), G4GeneralPhaseSpaceDecay::DecayIt(), and G4KL3DecayChannel::G4KL3DecayChannel().
00481 { 00482 G4cout << " BR: " << rbranch << " [" << kinematics_name << "]"; 00483 G4cout << " : " ; 00484 for (G4int index=0; index < numberOfDaughters; index++) 00485 { 00486 if(daughters_name[index] != 0) { 00487 G4cout << " " << *(daughters_name[index]); 00488 } else { 00489 G4cout << " not defined "; 00490 } 00491 } 00492 G4cout << G4endl; 00493 }
void G4VDecayChannel::FillDaughters | ( | ) | [protected] |
Definition at line 282 of file G4VDecayChannel.cc.
References daughters, daughters_mass, daughters_name, FatalException, FillParent(), G4ParticleTable::FindParticle(), G4cerr, G4cout, G4endl, G4Exception(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetParticleType(), G4ParticleDefinition::GetPDGMass(), G4ParticleDefinition::GetPDGWidth(), GetVerboseLevel(), numberOfDaughters, parent, parent_name, particletable, SetBR(), and verboseLevel.
Referenced by G4TauLeptonicDecayChannel::DecayIt(), G4PionRadiativeDecayChannel::DecayIt(), G4PhaseSpaceDecayChannel::DecayIt(), G4NuclearDecayChannel::DecayIt(), G4NeutronBetaDecayChannel::DecayIt(), G4MuonRadiativeDecayChannelWithSpin::DecayIt(), G4MuonDecayChannelWithSpin::DecayIt(), G4MuonDecayChannel::DecayIt(), G4KL3DecayChannel::DecayIt(), G4GeneralPhaseSpaceDecay::DecayIt(), G4DalitzDecayChannel::DecayIt(), GetAngularMomentum(), GetDaughter(), and SetDaughter().
00283 { 00284 G4int index; 00285 00286 #ifdef G4VERBOSE 00287 if (verboseLevel>1) G4cout << "G4VDecayChannel::FillDaughters()" <<G4endl; 00288 #endif 00289 if (daughters != 0) { 00290 delete [] daughters; 00291 daughters = 0; 00292 } 00293 00294 // parent mass 00295 if (parent == 0) FillParent(); 00296 G4double parentmass = parent->GetPDGMass(); 00297 00298 // 00299 G4double sumofdaughtermass = 0.0; 00300 if ((numberOfDaughters <=0) || (daughters_name == 0) ){ 00301 #ifdef G4VERBOSE 00302 if (verboseLevel>0) { 00303 G4cerr << "G4VDecayChannel::FillDaughters " 00304 << "[ " << parent->GetParticleName() << " ]" 00305 << "numberOfDaughters is not defined yet"; 00306 } 00307 #endif 00308 daughters = 0; 00309 G4Exception("G4VDecayChannel::FillDaughters", 00310 "PART011", FatalException, 00311 "Can not fill daughters: numberOfDaughters is not defined yet"); 00312 } 00313 00314 //create and set the array of pointers to daughter particles 00315 daughters = new G4ParticleDefinition*[numberOfDaughters]; 00316 if (daughters_mass != 0) delete [] daughters_mass; 00317 daughters_mass = new G4double[numberOfDaughters]; 00318 // loop over all daughters 00319 for (index=0; index < numberOfDaughters; index++) { 00320 if (daughters_name[index] == 0) { 00321 // daughter name is not defined 00322 #ifdef G4VERBOSE 00323 if (verboseLevel>0) { 00324 G4cerr << "G4VDecayChannel::FillDaughters " 00325 << "[ " << parent->GetParticleName() << " ]" 00326 << index << "-th daughter is not defined yet" << G4endl; 00327 } 00328 #endif 00329 daughters[index] = 0; 00330 G4Exception("G4VDecayChannel::FillDaughters", 00331 "PART011", FatalException, 00332 "Can not fill daughters: name of a daughter is not defined yet"); 00333 } 00334 //search daughter particles in the particle table 00335 daughters[index] = particletable->FindParticle(*daughters_name[index]); 00336 if (daughters[index] == 0) { 00337 // can not find the daughter particle 00338 #ifdef G4VERBOSE 00339 if (verboseLevel>0) { 00340 G4cerr << "G4VDecayChannel::FillDaughters " 00341 << "[ " << parent->GetParticleName() << " ]" 00342 << index << ":" << *daughters_name[index] 00343 << " is not defined !!" << G4endl; 00344 G4cerr << " The BR of this decay mode is set to zero " << G4endl; 00345 } 00346 #endif 00347 SetBR(0.0); 00348 return; 00349 } 00350 #ifdef G4VERBOSE 00351 if (verboseLevel>1) { 00352 G4cout << index << ":" << *daughters_name[index]; 00353 G4cout << ":" << daughters[index] << G4endl; 00354 } 00355 #endif 00356 daughters_mass[index] = daughters[index]->GetPDGMass(); 00357 sumofdaughtermass += daughters[index]->GetPDGMass(); 00358 } // end loop over all daughters 00359 00360 // check sum of daghter mass 00361 G4double widthMass = parent->GetPDGWidth(); 00362 if ( (parent->GetParticleType() != "nucleus") && 00363 (sumofdaughtermass > parentmass + 5*widthMass) ){ 00364 // !!! illegal mass !!! 00365 #ifdef G4VERBOSE 00366 if (GetVerboseLevel()>0) { 00367 G4cerr << "G4VDecayChannel::FillDaughters " 00368 << "[ " << parent->GetParticleName() << " ]" 00369 << " Energy/Momentum conserevation breaks " <<G4endl; 00370 if (GetVerboseLevel()>1) { 00371 G4cerr << " parent:" << *parent_name 00372 << " mass:" << parentmass/GeV << "[GeV/c/c]" <<G4endl; 00373 for (index=0; index < numberOfDaughters; index++){ 00374 G4cerr << " daughter " << index << ":" << *daughters_name[index] 00375 << " mass:" << daughters[index]->GetPDGMass()/GeV 00376 << "[GeV/c/c]" <<G4endl; 00377 } 00378 } 00379 } 00380 #endif 00381 } 00382 }
void G4VDecayChannel::FillParent | ( | ) | [protected] |
Definition at line 385 of file G4VDecayChannel.cc.
References FatalException, G4ParticleTable::FindParticle(), G4cerr, G4endl, G4Exception(), G4ParticleDefinition::GetPDGMass(), parent, parent_mass, parent_name, particletable, and verboseLevel.
Referenced by G4TauLeptonicDecayChannel::DecayIt(), G4PionRadiativeDecayChannel::DecayIt(), G4PhaseSpaceDecayChannel::DecayIt(), G4NuclearDecayChannel::DecayIt(), G4NeutronBetaDecayChannel::DecayIt(), G4MuonRadiativeDecayChannelWithSpin::DecayIt(), G4MuonDecayChannelWithSpin::DecayIt(), G4MuonDecayChannel::DecayIt(), G4KL3DecayChannel::DecayIt(), G4GeneralPhaseSpaceDecay::DecayIt(), G4DalitzDecayChannel::DecayIt(), FillDaughters(), G4NuclearDecayChannel::G4NuclearDecayChannel(), and GetParent().
00386 { 00387 if (parent_name == 0) { 00388 // parent name is not defined 00389 #ifdef G4VERBOSE 00390 if (verboseLevel>0) { 00391 G4cerr << "G4VDecayChannel::FillParent " 00392 << ": parent name is not defined !!" << G4endl; 00393 } 00394 #endif 00395 parent = 0; 00396 G4Exception("G4VDecayChannel::FillParent()", 00397 "PART012", FatalException, 00398 "Can not fill parent: parent name is not defined yet"); 00399 } 00400 // search parent particle in the particle table 00401 parent = particletable->FindParticle(*parent_name); 00402 if (parent == 0) { 00403 // parent particle does not exist 00404 #ifdef G4VERBOSE 00405 if (verboseLevel>0) { 00406 G4cerr << "G4VDecayChannel::FillParent " 00407 << *parent_name << " does not exist !!" << G4endl; 00408 } 00409 #endif 00410 G4Exception("G4VDecayChannel::FillParent()", 00411 "PART012", FatalException, 00412 "Can not fill parent: parent does not exist"); 00413 } 00414 parent_mass = parent->GetPDGMass(); 00415 }
G4int G4VDecayChannel::GetAngularMomentum | ( | ) |
Definition at line 422 of file G4VDecayChannel.cc.
References daughters, FillDaughters(), G4cout, G4endl, G4Exception(), G4ParticleDefinition::GetPDGiParity(), G4ParticleDefinition::GetPDGiSpin(), JustWarning, numberOfDaughters, parent, and verboseLevel.
00423 { 00424 // determine angular momentum 00425 00426 // fill pointers to daughter particles if not yet set 00427 if (daughters == 0) FillDaughters(); 00428 00429 const G4int PiSpin = parent->GetPDGiSpin(); 00430 const G4int PParity = parent->GetPDGiParity(); 00431 if (2==numberOfDaughters) { // up to now we can only handle two particle decays 00432 const G4int D1iSpin = daughters[0]->GetPDGiSpin(); 00433 const G4int D1Parity = daughters[0]->GetPDGiParity(); 00434 const G4int D2iSpin = daughters[1]->GetPDGiSpin(); 00435 const G4int D2Parity = daughters[1]->GetPDGiParity(); 00436 const G4int MiniSpin = std::abs (D1iSpin - D2iSpin); 00437 const G4int MaxiSpin = D1iSpin + D2iSpin; 00438 const G4int lMax = (PiSpin+D1iSpin+D2iSpin)/2; // l is allways int 00439 G4int lMin; 00440 #ifdef G4VERBOSE 00441 if (verboseLevel>1) { 00442 G4cout << "iSpin: " << PiSpin << " -> " << D1iSpin << " + " << D2iSpin << G4endl; 00443 G4cout << "2*jmin, 2*jmax, lmax " << MiniSpin << " " << MaxiSpin << " " << lMax << G4endl; 00444 } 00445 #endif 00446 for (G4int j=MiniSpin; j<=MaxiSpin; j+=2){ // loop over all possible spin couplings 00447 lMin = std::abs(PiSpin-j)/2; 00448 #ifdef G4VERBOSE 00449 if (verboseLevel>1) 00450 G4cout << "-> checking 2*j=" << j << G4endl; 00451 #endif 00452 for (G4int l=lMin; l<=lMax; l++) { 00453 #ifdef G4VERBOSE 00454 if (verboseLevel>1) 00455 G4cout << " checking l=" << l << G4endl; 00456 #endif 00457 if (l%2==0) { 00458 if (PParity == D1Parity*D2Parity) { // check parity for this l 00459 return l; 00460 } 00461 } else { 00462 if (PParity == -1*D1Parity*D2Parity) { // check parity for this l 00463 return l; 00464 } 00465 } 00466 } 00467 } 00468 } else { 00469 G4Exception("G4VDecayChannel::GetAngularMomentum", 00470 "PART111", JustWarning, 00471 "Sorry, can't handle 3 particle decays (up to now)"); 00472 return 0; 00473 } 00474 G4Exception ("G4VDecayChannel::GetAngularMomentum", 00475 "PART111", JustWarning, 00476 "Can't find angular momentum for this decay"); 00477 return 0; 00478 }
G4double G4VDecayChannel::GetBR | ( | ) | const [inline] |
Definition at line 270 of file G4VDecayChannel.hh.
References rbranch.
Referenced by G4RadioactiveDecay::AddDecayRateTable(), G4TextPPReporter::GeneratePropertyTable(), G4DecayTableMessenger::GetCurrentValue(), G4DecayTable::Insert(), and G4RadioactiveDecay::LoadDecayTable().
00270 { return rbranch; }
G4ParticleDefinition * G4VDecayChannel::GetDaughter | ( | G4int | anIndex | ) | [inline] |
Definition at line 184 of file G4VDecayChannel.hh.
References daughters, FillDaughters(), G4cout, G4endl, numberOfDaughters, and verboseLevel.
Referenced by G4KineticTrack::Decay().
00185 { 00186 //pointers to daughter particles are filled, if they are not set yet 00187 if (daughters == 0) FillDaughters(); 00188 00189 //get the pointer to a daughter particle 00190 if ( (anIndex>=0) && (anIndex<numberOfDaughters) ) { 00191 return daughters[anIndex]; 00192 } else { 00193 if (verboseLevel>0) 00194 G4cout << "G4VDecayChannel::GetDaughter index out of range "<<anIndex<<G4endl; 00195 return 0; 00196 } 00197 }
Definition at line 214 of file G4VDecayChannel.hh.
References daughters_mass, G4cout, G4endl, numberOfDaughters, and verboseLevel.
00215 { 00216 if ( (anIndex>=0) && (anIndex<numberOfDaughters) ) { 00217 return daughters_mass[anIndex]; 00218 } else { 00219 if (verboseLevel>0){ 00220 G4cout << "G4VDecayChannel::GetDaughterMass "; 00221 G4cout << "index out of range " << anIndex << G4endl; 00222 } 00223 return 0.0; 00224 } 00225 }
Definition at line 200 of file G4VDecayChannel.hh.
References daughters_name, G4cout, G4endl, numberOfDaughters, and verboseLevel.
Referenced by G4KineticTrack::Decay().
00201 { 00202 if ( (anIndex>=0) && (anIndex<numberOfDaughters) ) { 00203 return *daughters_name[anIndex]; 00204 } else { 00205 if (verboseLevel>0){ 00206 G4cout << "G4VDecayChannel::GetDaughterName "; 00207 G4cout << "index out of range " << anIndex << G4endl; 00208 } 00209 return GetNoName(); 00210 } 00211 }
const G4String & G4VDecayChannel::GetKinematicsName | ( | ) | const [inline] |
Definition at line 264 of file G4VDecayChannel.hh.
References kinematics_name.
00264 { return kinematics_name; }
G4int G4VDecayChannel::GetNumberOfDaughters | ( | ) | const [inline] |
Definition at line 258 of file G4VDecayChannel.hh.
References numberOfDaughters.
Referenced by G4KineticTrack::Decay(), and G4TextPPRetriever::ModifyPropertyTable().
00259 { 00260 return numberOfDaughters; 00261 }
G4ParticleDefinition * G4VDecayChannel::GetParent | ( | ) | [inline] |
Definition at line 228 of file G4VDecayChannel.hh.
References FillParent(), and parent.
Referenced by G4DecayTable::Insert().
00229 { 00230 //the pointer to the parent particle is filled, if it is not set yet 00231 if (parent == 0) FillParent(); 00232 //get the pointer to the parent particle 00233 return parent; 00234 }
G4double G4VDecayChannel::GetParentMass | ( | ) | const [inline] |
Reimplemented in G4GeneralPhaseSpaceDecay.
Definition at line 243 of file G4VDecayChannel.hh.
References parent_mass.
00244 { 00245 return parent_mass; 00246 }
const G4String & G4VDecayChannel::GetParentName | ( | ) | const [inline] |
Definition at line 237 of file G4VDecayChannel.hh.
References parent_name.
Referenced by G4KineticTrack::Decay(), and G4MuonRadiativeDecayChannelWithSpin::DecayIt().
00238 { 00239 return *parent_name; 00240 }
G4int G4VDecayChannel::GetVerboseLevel | ( | ) | const [inline] |
Definition at line 276 of file G4VDecayChannel.hh.
References verboseLevel.
Referenced by G4KL3DecayChannel::DalitzDensity(), G4TauLeptonicDecayChannel::DecayIt(), G4PionRadiativeDecayChannel::DecayIt(), G4PhaseSpaceDecayChannel::DecayIt(), G4NuclearDecayChannel::DecayIt(), G4NeutronBetaDecayChannel::DecayIt(), G4MuonRadiativeDecayChannelWithSpin::DecayIt(), G4MuonDecayChannelWithSpin::DecayIt(), G4MuonDecayChannel::DecayIt(), G4KL3DecayChannel::DecayIt(), G4GeneralPhaseSpaceDecay::DecayIt(), G4Decay::DecayIt(), G4DalitzDecayChannel::DecayIt(), FillDaughters(), G4AlphaDecayChannel::G4AlphaDecayChannel(), G4BetaMinusDecayChannel::G4BetaMinusDecayChannel(), G4BetaPlusDecayChannel::G4BetaPlusDecayChannel(), G4GeneralPhaseSpaceDecay::G4GeneralPhaseSpaceDecay(), G4ITDecayChannel::G4ITDecayChannel(), G4KL3DecayChannel::G4KL3DecayChannel(), G4KshellECDecayChannel::G4KshellECDecayChannel(), G4LshellECDecayChannel::G4LshellECDecayChannel(), G4MshellECDecayChannel::G4MshellECDecayChannel(), G4MuonDecayChannel::G4MuonDecayChannel(), G4MuonRadiativeDecayChannelWithSpin::G4MuonRadiativeDecayChannelWithSpin(), G4NeutronBetaDecayChannel::G4NeutronBetaDecayChannel(), G4NuclearDecayChannel::G4NuclearDecayChannel(), G4PionRadiativeDecayChannel::G4PionRadiativeDecayChannel(), G4TauLeptonicDecayChannel::G4TauLeptonicDecayChannel(), G4GeneralPhaseSpaceDecay::ManyBodyDecayIt(), G4GeneralPhaseSpaceDecay::OneBodyDecayIt(), G4KL3DecayChannel::PhaseSpace(), G4GeneralPhaseSpaceDecay::ThreeBodyDecayIt(), and G4GeneralPhaseSpaceDecay::TwoBodyDecayIt().
00276 { return verboseLevel; }
G4int G4VDecayChannel::operator!= | ( | const G4VDecayChannel & | right | ) | const [inline] |
G4int G4VDecayChannel::operator< | ( | const G4VDecayChannel & | right | ) | const [inline] |
Definition at line 178 of file G4VDecayChannel.hh.
References rbranch.
00179 { 00180 return (this->rbranch < right.rbranch); 00181 }
G4VDecayChannel & G4VDecayChannel::operator= | ( | const G4VDecayChannel & | ) | [protected] |
Definition at line 145 of file G4VDecayChannel.cc.
References ClearDaughtersName(), daughters, daughters_mass, daughters_name, G4ParticleTable::GetParticleTable(), kinematics_name, numberOfDaughters, parent, parent_mass, parent_name, particletable, rbranch, and verboseLevel.
00146 { 00147 if (this != &right) { 00148 kinematics_name = right.kinematics_name; 00149 verboseLevel = right.verboseLevel; 00150 rbranch = right.rbranch; 00151 00152 // copy parent name 00153 parent_name = new G4String(*right.parent_name); 00154 00155 // clear daughters_name array 00156 ClearDaughtersName(); 00157 00158 // recreate array 00159 numberOfDaughters = right.numberOfDaughters; 00160 if ( numberOfDaughters >0 ) { 00161 if (daughters_name !=0) ClearDaughtersName(); 00162 daughters_name = new G4String*[numberOfDaughters]; 00163 //copy daughters name 00164 for (G4int index=0; index < numberOfDaughters; index++) { 00165 daughters_name[index] = new G4String(*right.daughters_name[index]); 00166 } 00167 } 00168 } 00169 00170 // 00171 parent = 0; 00172 daughters = 0; 00173 parent_mass = 0.0; 00174 daughters_mass = 0; 00175 00176 // particle table 00177 particletable = G4ParticleTable::GetParticleTable(); 00178 00179 return *this; 00180 }
G4int G4VDecayChannel::operator== | ( | const G4VDecayChannel & | right | ) | const [inline] |
void G4VDecayChannel::SetBR | ( | G4double | value | ) | [inline] |
Definition at line 267 of file G4VDecayChannel.hh.
References rbranch.
Referenced by G4TauPlus::Definition(), G4TauMinus::Definition(), FillDaughters(), G4DalitzDecayChannel::G4DalitzDecayChannel(), G4MuonDecayChannel::G4MuonDecayChannel(), G4MuonRadiativeDecayChannelWithSpin::G4MuonRadiativeDecayChannelWithSpin(), G4NeutronBetaDecayChannel::G4NeutronBetaDecayChannel(), G4NuclearDecayChannel::G4NuclearDecayChannel(), G4PionRadiativeDecayChannel::G4PionRadiativeDecayChannel(), G4TauLeptonicDecayChannel::G4TauLeptonicDecayChannel(), G4RadioactiveDecay::LoadDecayTable(), G4TextPPRetriever::ModifyPropertyTable(), and G4DecayTableMessenger::SetNewValue().
00267 { rbranch = value; }
Definition at line 230 of file G4VDecayChannel.cc.
References daughters, daughters_name, FillDaughters(), G4cerr, G4cout, G4endl, numberOfDaughters, and verboseLevel.
00232 { 00233 // check numberOfDaughters is positive 00234 if (numberOfDaughters<=0) { 00235 #ifdef G4VERBOSE 00236 if (verboseLevel>0) { 00237 G4cerr << "G4VDecayChannel::SetDaughter: " 00238 << "Number of daughters is not defined" << G4endl; 00239 } 00240 #endif 00241 return; 00242 } 00243 00244 // check existence of daughters_name array 00245 if (daughters_name == 0) { 00246 // cleate array 00247 daughters_name = new G4String*[numberOfDaughters]; 00248 for (G4int index=0;index<numberOfDaughters;index++) { 00249 daughters_name[index]=0; 00250 } 00251 } 00252 00253 // check an index 00254 if ( (anIndex<0) || (anIndex>=numberOfDaughters) ) { 00255 #ifdef G4VERBOSE 00256 if (verboseLevel>0) { 00257 G4cerr << "G4VDecayChannel::SetDaughter" 00258 << "index out of range " << anIndex << G4endl; 00259 } 00260 #endif 00261 } else { 00262 // delete the old name if it exists 00263 if (daughters_name[anIndex]!=0) delete daughters_name[anIndex]; 00264 // fill the name 00265 daughters_name[anIndex] = new G4String(particle_name); 00266 // refill the array of daughters[] if it exists 00267 if (daughters != 0) FillDaughters(); 00268 #ifdef G4VERBOSE 00269 if (verboseLevel>1) { 00270 G4cout << "G4VDecayChannel::SetDaughter[" << anIndex <<"] :"; 00271 G4cout << daughters_name[anIndex] << ":" << *daughters_name[anIndex]<<G4endl; 00272 } 00273 #endif 00274 } 00275 }
void G4VDecayChannel::SetDaughter | ( | G4int | anIndex, | |
const G4ParticleDefinition * | particle_type | |||
) |
Definition at line 277 of file G4VDecayChannel.cc.
References G4ParticleDefinition::GetParticleName().
Referenced by G4TauPlus::Definition(), G4TauMinus::Definition(), G4DalitzDecayChannel::G4DalitzDecayChannel(), G4MuonDecayChannel::G4MuonDecayChannel(), G4MuonRadiativeDecayChannelWithSpin::G4MuonRadiativeDecayChannelWithSpin(), G4NeutronBetaDecayChannel::G4NeutronBetaDecayChannel(), G4NuclearDecayChannel::G4NuclearDecayChannel(), G4PionRadiativeDecayChannel::G4PionRadiativeDecayChannel(), and G4TauLeptonicDecayChannel::G4TauLeptonicDecayChannel().
00278 { 00279 if (parent_type != 0) SetDaughter(anIndex, parent_type->GetParticleName()); 00280 }
void G4VDecayChannel::SetNumberOfDaughters | ( | G4int | value | ) |
Definition at line 218 of file G4VDecayChannel.cc.
References ClearDaughtersName(), daughters_name, and numberOfDaughters.
Referenced by G4TauPlus::Definition(), G4TauMinus::Definition(), G4DalitzDecayChannel::G4DalitzDecayChannel(), G4MuonDecayChannel::G4MuonDecayChannel(), G4MuonRadiativeDecayChannelWithSpin::G4MuonRadiativeDecayChannelWithSpin(), G4NeutronBetaDecayChannel::G4NeutronBetaDecayChannel(), G4NuclearDecayChannel::G4NuclearDecayChannel(), G4PionRadiativeDecayChannel::G4PionRadiativeDecayChannel(), and G4TauLeptonicDecayChannel::G4TauLeptonicDecayChannel().
00219 { 00220 if (size >0) { 00221 // remove old contents 00222 ClearDaughtersName(); 00223 // cleate array 00224 daughters_name = new G4String*[size]; 00225 for (G4int index=0;index<size;index++) daughters_name[index]=0; 00226 numberOfDaughters = size; 00227 } 00228 }
void G4VDecayChannel::SetParent | ( | const G4String & | particle_name | ) | [inline] |
Definition at line 250 of file G4VDecayChannel.hh.
References parent, and parent_name.
00251 { 00252 if (parent_name != 0) delete parent_name; 00253 parent_name = new G4String(particle_name); 00254 parent = 0; 00255 }
void G4VDecayChannel::SetParent | ( | const G4ParticleDefinition * | particle_type | ) |
Definition at line 417 of file G4VDecayChannel.cc.
References G4ParticleDefinition::GetParticleName().
Referenced by G4TauPlus::Definition(), G4TauMinus::Definition(), G4DalitzDecayChannel::G4DalitzDecayChannel(), G4MuonDecayChannel::G4MuonDecayChannel(), G4MuonRadiativeDecayChannelWithSpin::G4MuonRadiativeDecayChannelWithSpin(), G4NeutronBetaDecayChannel::G4NeutronBetaDecayChannel(), G4NuclearDecayChannel::G4NuclearDecayChannel(), G4PionRadiativeDecayChannel::G4PionRadiativeDecayChannel(), and G4TauLeptonicDecayChannel::G4TauLeptonicDecayChannel().
00418 { 00419 if (parent_type != 0) SetParent(parent_type->GetParticleName()); 00420 }
void G4VDecayChannel::SetVerboseLevel | ( | G4int | value | ) | [inline] |
Definition at line 273 of file G4VDecayChannel.hh.
References verboseLevel.
Referenced by G4Decay::DecayIt().
00273 { verboseLevel = value; }
G4ParticleDefinition** G4VDecayChannel::daughters [protected] |
Definition at line 147 of file G4VDecayChannel.hh.
Referenced by ClearDaughtersName(), G4TauLeptonicDecayChannel::DecayIt(), G4PionRadiativeDecayChannel::DecayIt(), G4PhaseSpaceDecayChannel::DecayIt(), G4NuclearDecayChannel::DecayIt(), G4NeutronBetaDecayChannel::DecayIt(), G4MuonRadiativeDecayChannelWithSpin::DecayIt(), G4MuonDecayChannelWithSpin::DecayIt(), G4MuonDecayChannel::DecayIt(), G4KL3DecayChannel::DecayIt(), G4GeneralPhaseSpaceDecay::DecayIt(), G4DalitzDecayChannel::DecayIt(), FillDaughters(), G4VDecayChannel(), GetAngularMomentum(), GetDaughter(), G4GeneralPhaseSpaceDecay::ManyBodyDecayIt(), G4GeneralPhaseSpaceDecay::OneBodyDecayIt(), operator=(), SetDaughter(), G4GeneralPhaseSpaceDecay::ThreeBodyDecayIt(), and G4GeneralPhaseSpaceDecay::TwoBodyDecayIt().
G4double* G4VDecayChannel::daughters_mass [protected] |
Definition at line 151 of file G4VDecayChannel.hh.
Referenced by ClearDaughtersName(), FillDaughters(), G4VDecayChannel(), GetDaughterMass(), operator=(), and ~G4VDecayChannel().
G4String** G4VDecayChannel::daughters_name [protected] |
Definition at line 135 of file G4VDecayChannel.hh.
Referenced by ClearDaughtersName(), G4KL3DecayChannel::DecayIt(), DumpInfo(), FillDaughters(), G4VDecayChannel(), GetDaughterName(), G4GeneralPhaseSpaceDecay::ManyBodyDecayIt(), operator=(), G4TauLeptonicDecayChannel::operator=(), G4PionRadiativeDecayChannel::operator=(), G4NeutronBetaDecayChannel::operator=(), G4MuonRadiativeDecayChannelWithSpin::operator=(), G4MuonDecayChannelWithSpin::operator=(), G4MuonDecayChannel::operator=(), G4KL3DecayChannel::operator=(), G4DalitzDecayChannel::operator=(), SetDaughter(), and SetNumberOfDaughters().
G4String G4VDecayChannel::kinematics_name [protected] |
Definition at line 127 of file G4VDecayChannel.hh.
Referenced by DumpInfo(), G4VDecayChannel(), GetKinematicsName(), operator=(), G4TauLeptonicDecayChannel::operator=(), G4PionRadiativeDecayChannel::operator=(), G4NeutronBetaDecayChannel::operator=(), G4MuonRadiativeDecayChannelWithSpin::operator=(), G4MuonDecayChannelWithSpin::operator=(), G4MuonDecayChannel::operator=(), G4KL3DecayChannel::operator=(), and G4DalitzDecayChannel::operator=().
const G4String G4VDecayChannel::noName = " " [static, protected] |
Definition at line 174 of file G4VDecayChannel.hh.
G4int G4VDecayChannel::numberOfDaughters [protected] |
Definition at line 131 of file G4VDecayChannel.hh.
Referenced by ClearDaughtersName(), G4PhaseSpaceDecayChannel::DecayIt(), G4NuclearDecayChannel::DecayIt(), G4GeneralPhaseSpaceDecay::DecayIt(), DumpInfo(), FillDaughters(), G4VDecayChannel(), GetAngularMomentum(), GetDaughter(), GetDaughterMass(), GetDaughterName(), GetNumberOfDaughters(), G4GeneralPhaseSpaceDecay::ManyBodyDecayIt(), operator=(), G4TauLeptonicDecayChannel::operator=(), G4PionRadiativeDecayChannel::operator=(), G4NeutronBetaDecayChannel::operator=(), G4MuonRadiativeDecayChannelWithSpin::operator=(), G4MuonDecayChannelWithSpin::operator=(), G4MuonDecayChannel::operator=(), G4KL3DecayChannel::operator=(), G4DalitzDecayChannel::operator=(), SetDaughter(), and SetNumberOfDaughters().
G4ParticleDefinition* G4VDecayChannel::parent [protected] |
Definition at line 146 of file G4VDecayChannel.hh.
Referenced by G4TauLeptonicDecayChannel::DecayIt(), G4PionRadiativeDecayChannel::DecayIt(), G4PhaseSpaceDecayChannel::DecayIt(), G4NuclearDecayChannel::DecayIt(), G4NeutronBetaDecayChannel::DecayIt(), G4MuonRadiativeDecayChannelWithSpin::DecayIt(), G4MuonDecayChannelWithSpin::DecayIt(), G4MuonDecayChannel::DecayIt(), G4KL3DecayChannel::DecayIt(), G4GeneralPhaseSpaceDecay::DecayIt(), G4DalitzDecayChannel::DecayIt(), FillDaughters(), FillParent(), G4GeneralPhaseSpaceDecay::G4GeneralPhaseSpaceDecay(), G4VDecayChannel(), GetAngularMomentum(), GetParent(), G4GeneralPhaseSpaceDecay::ManyBodyDecayIt(), G4GeneralPhaseSpaceDecay::OneBodyDecayIt(), operator=(), SetParent(), G4GeneralPhaseSpaceDecay::ThreeBodyDecayIt(), and G4GeneralPhaseSpaceDecay::TwoBodyDecayIt().
G4double G4VDecayChannel::parent_mass [protected] |
Definition at line 150 of file G4VDecayChannel.hh.
Referenced by G4PhaseSpaceDecayChannel::DecayIt(), FillParent(), G4NuclearDecayChannel::G4NuclearDecayChannel(), G4VDecayChannel(), GetParentMass(), and operator=().
G4String* G4VDecayChannel::parent_name [protected] |
Definition at line 133 of file G4VDecayChannel.hh.
Referenced by ClearDaughtersName(), G4PhaseSpaceDecayChannel::DecayIt(), G4NuclearDecayChannel::DecayIt(), G4GeneralPhaseSpaceDecay::DecayIt(), FillDaughters(), FillParent(), G4VDecayChannel(), GetParentName(), G4GeneralPhaseSpaceDecay::ManyBodyDecayIt(), operator=(), G4TauLeptonicDecayChannel::operator=(), G4PionRadiativeDecayChannel::operator=(), G4NeutronBetaDecayChannel::operator=(), G4MuonRadiativeDecayChannelWithSpin::operator=(), G4MuonDecayChannelWithSpin::operator=(), G4MuonDecayChannel::operator=(), G4KL3DecayChannel::operator=(), G4DalitzDecayChannel::operator=(), SetParent(), and ~G4VDecayChannel().
G4ParticleTable* G4VDecayChannel::particletable [protected] |
Definition at line 143 of file G4VDecayChannel.hh.
Referenced by FillDaughters(), FillParent(), G4VDecayChannel(), and operator=().
G4double G4VDecayChannel::rbranch [protected] |
Definition at line 129 of file G4VDecayChannel.hh.
Referenced by DumpInfo(), G4VDecayChannel(), GetBR(), operator<(), operator=(), G4TauLeptonicDecayChannel::operator=(), G4PionRadiativeDecayChannel::operator=(), G4NeutronBetaDecayChannel::operator=(), G4MuonRadiativeDecayChannelWithSpin::operator=(), G4MuonDecayChannelWithSpin::operator=(), G4MuonDecayChannel::operator=(), G4KL3DecayChannel::operator=(), G4DalitzDecayChannel::operator=(), and SetBR().
G4int G4VDecayChannel::verboseLevel [protected] |
Definition at line 169 of file G4VDecayChannel.hh.
Referenced by ClearDaughtersName(), FillDaughters(), FillParent(), G4VDecayChannel(), GetAngularMomentum(), GetDaughter(), GetDaughterMass(), GetDaughterName(), GetVerboseLevel(), operator=(), G4TauLeptonicDecayChannel::operator=(), G4PionRadiativeDecayChannel::operator=(), G4NeutronBetaDecayChannel::operator=(), G4MuonRadiativeDecayChannelWithSpin::operator=(), G4MuonDecayChannelWithSpin::operator=(), G4MuonDecayChannel::operator=(), G4KL3DecayChannel::operator=(), G4DalitzDecayChannel::operator=(), SetDaughter(), and SetVerboseLevel().