Geant4  10.02
G4VDecayChannel.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4VDecayChannel.cc 93024 2015-09-30 15:59:32Z gcosmo $
28 //
29 //
30 // ------------------------------------------------------------
31 // GEANT 4 class header file
32 //
33 // History: first implementation, based on object model of
34 // 27 July 1996 H.Kurashige
35 // 30 May 1997 H.Kurashige
36 // 23 Mar. 2000 H.Weber : add GetAngularMomentum
37 // ------------------------------------------------------------
38 
39 #include "G4ParticleDefinition.hh"
40 #include "G4SystemOfUnits.hh"
41 #include "G4ParticleTable.hh"
42 #include "G4DecayTable.hh"
43 #include "G4DecayProducts.hh"
44 #include "G4VDecayChannel.hh"
45 
61 
63 
65  :kinematics_name(""),
66  rbranch(0.0),
67  numberOfDaughters(0),
68  parent_name(0), daughters_name(0),
69  rangeMass(2.5),
70  parent_polarization(),
71  particletable(0),
72  verboseLevel(1)
73 {
75  G4MT_parent = 0;
76  G4MT_daughters = 0;
77  G4MT_parent_mass = 0.0;
80 
81  // set pointer to G4ParticleTable (static and singleton object)
83 }
84 
86  :kinematics_name(aName),
87  rbranch(0.0),
88  numberOfDaughters(0),
89  parent_name(0), daughters_name(0),
90  rangeMass(2.5),
91  parent_polarization(),
92  particletable(0),
93  verboseLevel(Verbose)
94 {
96  G4MT_parent = 0;
97  G4MT_daughters = 0;
98  G4MT_parent_mass = 0.0;
101 
102  // set pointer to G4ParticleTable (static and singleton object)
104 }
105 
107  const G4String& theParentName,
108  G4double theBR,
109  G4int theNumberOfDaughters,
110  const G4String& theDaughterName1,
111  const G4String& theDaughterName2,
112  const G4String& theDaughterName3,
113  const G4String& theDaughterName4 )
114  :kinematics_name(aName),
115  rbranch(theBR),
116  numberOfDaughters(theNumberOfDaughters),
117  parent_name(0), daughters_name(0),
118  rangeMass(1.0),
119  parent_polarization(),
120  particletable(0),
121  verboseLevel(1)
122 {
124  G4MT_parent = 0;
125  G4MT_daughters = 0;
126  G4MT_parent_mass = 0.0;
129 
130  // set pointer to G4ParticleTable (static and singleton object)
132 
133  // parent name
134  parent_name = new G4String(theParentName);
135 
136  // cleate array
138  for (G4int index=0;index<numberOfDaughters;index++) daughters_name[index]=0;
139 
140  // daughters' name
141  if (numberOfDaughters>0) daughters_name[0] = new G4String(theDaughterName1);
142  if (numberOfDaughters>1) daughters_name[1] = new G4String(theDaughterName2);
143  if (numberOfDaughters>2) daughters_name[2] = new G4String(theDaughterName3);
144  if (numberOfDaughters>3) daughters_name[3] = new G4String(theDaughterName4);
145 }
146 
148 {
150 
152  verboseLevel = right.verboseLevel;
153  rbranch = right.rbranch;
154  rangeMass = right.rangeMass;
155 
156  // copy parent name
157  parent_name = new G4String(*right.parent_name);
158  G4MT_parent = 0;
159  G4MT_parent_mass = 0.0;
160 
161  //create array
163 
164  daughters_name =0;
165  if ( numberOfDaughters >0 ) {
167  //copy daughters name
168  for (G4int index=0; index < numberOfDaughters; index++)
169  {
170  daughters_name[index] = new G4String(*right.daughters_name[index]);
171  }
172  }
173 
174  //
176  G4MT_daughters = 0;
178 
179  // particle table
181 
183 
184 }
185 
187 {
188  if (this != &right) {
190  verboseLevel = right.verboseLevel;
191  rbranch = right.rbranch;
192  rangeMass = right.rangeMass;
194  // copy parent name
195  parent_name = new G4String(*right.parent_name);
196 
197  // clear daughters_name array
199 
200  // recreate array
202  if ( numberOfDaughters >0 ) {
205  //copy daughters name
206  for (G4int index=0; index < numberOfDaughters; index++) {
207  daughters_name[index] = new G4String(*right.daughters_name[index]);
208  }
209  }
210  }
211 
212  //
213  G4MT_parent = 0;
214  G4MT_daughters = 0;
215  G4MT_parent_mass = 0.0;
218 
219  // particle table
221 
222  return *this;
223 }
224 
226 {
228  if (parent_name != 0) delete parent_name;
229  parent_name = 0;
230  if (G4MT_daughters_mass != 0) delete [] G4MT_daughters_mass;
232  if (G4MT_daughters_width != 0) delete [] G4MT_daughters_width;
234 }
235 
240 
242 {
243  if ( daughters_name != 0) {
244  if (numberOfDaughters>0) {
245 #ifdef G4VERBOSE
246  if (verboseLevel>1) {
247  G4cout << "G4VDecayChannel::ClearDaughtersName "
248  << " for " << *parent_name << G4endl;
249  }
250 #endif
251  for (G4int index=0; index < numberOfDaughters; index++) {
252  if (daughters_name[index] != 0) delete daughters_name[index];
253  }
254  }
255  delete [] daughters_name;
256  daughters_name = 0;
257  }
258  //
259  if (G4MT_daughters != 0) delete [] G4MT_daughters;
260  if (G4MT_daughters_mass != 0) delete [] G4MT_daughters_mass;
261  if (G4MT_daughters_width != 0) delete [] G4MT_daughters_width;
263  G4MT_daughters = 0;
266 
267  numberOfDaughters = 0;
268 }
269 
271 {
272  if (size >0) {
273  // remove old contents
275  // cleate array
276  daughters_name = new G4String*[size];
277  for (G4int index=0;index<size;index++) daughters_name[index]=0;
278  numberOfDaughters = size;
279  }
280 }
281 
283  const G4String &particle_name)
284 {
285  // check numberOfDaughters is positive
286  if (numberOfDaughters<=0) {
287 #ifdef G4VERBOSE
288  if (verboseLevel>0) {
289  G4cout << "G4VDecayChannel::SetDaughter: "
290  << "Number of daughters is not defined" << G4endl;
291  }
292 #endif
293  return;
294  }
295 
296  // check existence of daughters_name array
297  if (daughters_name == 0) {
298  // cleate array
300  for (G4int index=0;index<numberOfDaughters;index++) {
301  daughters_name[index]=0;
302  }
303  }
304 
305  // check an index
306  if ( (anIndex<0) || (anIndex>=numberOfDaughters) ) {
307 #ifdef G4VERBOSE
308  if (verboseLevel>0) {
309  G4cout << "G4VDecayChannel::SetDaughter"
310  << "index out of range " << anIndex << G4endl;
311  }
312 #endif
313  } else {
314  // delete the old name if it exists
315  if (daughters_name[anIndex]!=0) delete daughters_name[anIndex];
316  // fill the name
317  daughters_name[anIndex] = new G4String(particle_name);
318  // refill the array of daughters[] if it exists
319  if (G4MT_daughters != 0) FillDaughters();
320 #ifdef G4VERBOSE
321  if (verboseLevel>1) {
322  G4cout << "G4VDecayChannel::SetDaughter[" << anIndex <<"] :";
323  G4cout << daughters_name[anIndex] << ":" << *daughters_name[anIndex]<<G4endl;
324  }
325 #endif
326  }
327 }
328 
329 void G4VDecayChannel::SetDaughter(G4int anIndex, const G4ParticleDefinition * parent_type)
330 {
331  if (parent_type != 0) SetDaughter(anIndex, parent_type->GetParticleName());
332 }
333 
335 {
336  G4int index;
337 
338 #ifdef G4VERBOSE
339  if (verboseLevel>1) G4cout << "G4VDecayChannel::FillDaughters()" <<G4endl;
340 #endif
341  if (G4MT_daughters != 0) {
342  delete [] G4MT_daughters;
343  G4MT_daughters = 0;
344  }
345 
346  // parent mass
347  if (G4MT_parent == 0) FillParent();
348  G4double parentmass = G4MT_parent->GetPDGMass();
349 
350  //
351  G4double sumofdaughtermass = 0.0;
352  G4double sumofdaughterwidthsq = 0.0;
353 
354  if ((numberOfDaughters <=0) || (daughters_name == 0) ){
355 #ifdef G4VERBOSE
356  if (verboseLevel>0) {
357  G4cout << "G4VDecayChannel::FillDaughters "
358  << "[ " << G4MT_parent->GetParticleName() << " ]"
359  << "numberOfDaughters is not defined yet";
360  }
361 #endif
362  G4MT_daughters = 0;
363  G4Exception("G4VDecayChannel::FillDaughters",
364  "PART011", FatalException,
365  "Can not fill daughters: numberOfDaughters is not defined yet");
366  }
367 
368  //create and set the array of pointers to daughter particles
370  if (G4MT_daughters_mass != 0) delete [] G4MT_daughters_mass;
371  if (G4MT_daughters_width != 0) delete [] G4MT_daughters_width;
374  // loop over all daughters
375  for (index=0; index < numberOfDaughters; index++) {
376  if (daughters_name[index] == 0) {
377  // daughter name is not defined
378 #ifdef G4VERBOSE
379  if (verboseLevel>0) {
380  G4cout << "G4VDecayChannel::FillDaughters "
381  << "[ " << G4MT_parent->GetParticleName() << " ]"
382  << index << "-th daughter is not defined yet" << G4endl;
383  }
384 #endif
385  G4MT_daughters[index] = 0;
386  G4Exception("G4VDecayChannel::FillDaughters",
387  "PART011", FatalException,
388  "Can not fill daughters: name of a daughter is not defined yet");
389  }
390  //search daughter particles in the particle table
392  if (G4MT_daughters[index] == 0) {
393  // can not find the daughter particle
394 #ifdef G4VERBOSE
395  if (verboseLevel>0) {
396  G4cout << "G4VDecayChannel::FillDaughters "
397  << "[ " << G4MT_parent->GetParticleName() << " ]"
398  << index << ":" << *daughters_name[index]
399  << " is not defined !!" << G4endl;
400  G4cout << " The BR of this decay mode is set to zero " << G4endl;
401  }
402 #endif
403  SetBR(0.0);
404  return;
405  }
406 #ifdef G4VERBOSE
407  if (verboseLevel>1) {
408  G4cout << index << ":" << *daughters_name[index];
409  G4cout << ":" << G4MT_daughters[index] << G4endl;
410  }
411 #endif
412  G4MT_daughters_mass[index] = G4MT_daughters[index]->GetPDGMass();
413  G4double d_width = G4MT_daughters[index]->GetPDGWidth();
414  G4MT_daughters_width[index] = d_width;
415  sumofdaughtermass += G4MT_daughters[index]->GetPDGMass();
416  sumofdaughterwidthsq += d_width*d_width;
417  } // end loop over all daughters
418 
419  // check sum of daghter mass
420  G4double widthMass = std::sqrt(G4MT_parent->GetPDGWidth()*G4MT_parent->GetPDGWidth()+sumofdaughterwidthsq);
421  if ( (G4MT_parent->GetParticleType() != "nucleus") &&
422  (sumofdaughtermass > parentmass + rangeMass*widthMass) ){
423  // !!! illegal mass !!!
424 #ifdef G4VERBOSE
425  if (GetVerboseLevel()>0) {
426  G4cout << "G4VDecayChannel::FillDaughters "
427  << "[ " << G4MT_parent->GetParticleName() << " ]"
428  << " Energy/Momentum conserevation breaks " <<G4endl;
429  if (GetVerboseLevel()>1) {
430  G4cout << " parent:" << *parent_name
431  << " mass:" << parentmass/GeV << "[GeV/c/c]" <<G4endl;
432  for (index=0; index < numberOfDaughters; index++){
433  G4cout << " daughter " << index << ":" << *daughters_name[index]
434  << " mass:" << G4MT_daughters[index]->GetPDGMass()/GeV
435  << "[GeV/c/c]" <<G4endl;
436  }
437  }
438  }
439 #endif
440  }
441 }
442 
443 
445 {
446  if (parent_name == 0) {
447  // parent name is not defined
448 #ifdef G4VERBOSE
449  if (verboseLevel>0) {
450  G4cout << "G4VDecayChannel::FillParent "
451  << ": parent name is not defined !!" << G4endl;
452  }
453 #endif
454  G4MT_parent = 0;
455  G4Exception("G4VDecayChannel::FillParent()",
456  "PART012", FatalException,
457  "Can not fill parent: parent name is not defined yet");
458  return;
459  }
460  // search parent particle in the particle table
462  if (G4MT_parent == 0) {
463  // parent particle does not exist
464 #ifdef G4VERBOSE
465  if (verboseLevel>0) {
466  G4cout << "G4VDecayChannel::FillParent "
467  << *parent_name << " does not exist !!" << G4endl;
468  }
469 #endif
470  G4Exception("G4VDecayChannel::FillParent()",
471  "PART012", FatalException,
472  "Can not fill parent: parent does not exist");
473  return;
474  }
476 }
477 
479 {
480  if (parent_type != 0) SetParent(parent_type->GetParticleName());
481 }
482 
484 {
485  // determine angular momentum
486 
487  // fill pointers to daughter particles if not yet set
488  if (G4MT_daughters == 0) FillDaughters();
489 
490  const G4int PiSpin = G4MT_parent->GetPDGiSpin();
491  const G4int PParity = G4MT_parent->GetPDGiParity();
492  if (2==numberOfDaughters) { // up to now we can only handle two particle decays
493  const G4int D1iSpin = G4MT_daughters[0]->GetPDGiSpin();
494  const G4int D1Parity = G4MT_daughters[0]->GetPDGiParity();
495  const G4int D2iSpin = G4MT_daughters[1]->GetPDGiSpin();
496  const G4int D2Parity = G4MT_daughters[1]->GetPDGiParity();
497  const G4int MiniSpin = std::abs (D1iSpin - D2iSpin);
498  const G4int MaxiSpin = D1iSpin + D2iSpin;
499  const G4int lMax = (PiSpin+D1iSpin+D2iSpin)/2; // l is allways int
500  G4int lMin;
501 #ifdef G4VERBOSE
502  if (verboseLevel>1) {
503  G4cout << "iSpin: " << PiSpin << " -> " << D1iSpin << " + " << D2iSpin << G4endl;
504  G4cout << "2*jmin, 2*jmax, lmax " << MiniSpin << " " << MaxiSpin << " " << lMax << G4endl;
505  }
506 #endif
507  for (G4int j=MiniSpin; j<=MaxiSpin; j+=2){ // loop over all possible spin couplings
508  lMin = std::abs(PiSpin-j)/2;
509 #ifdef G4VERBOSE
510  if (verboseLevel>1)
511  G4cout << "-> checking 2*j=" << j << G4endl;
512 #endif
513  for (G4int l=lMin; l<=lMax; l++) {
514 #ifdef G4VERBOSE
515  if (verboseLevel>1)
516  G4cout << " checking l=" << l << G4endl;
517 #endif
518  if (l%2==0) {
519  if (PParity == D1Parity*D2Parity) { // check parity for this l
520  return l;
521  }
522  } else {
523  if (PParity == -1*D1Parity*D2Parity) { // check parity for this l
524  return l;
525  }
526  }
527  }
528  }
529  } else {
530  G4Exception("G4VDecayChannel::GetAngularMomentum",
531  "PART111", JustWarning,
532  "Sorry, can't handle 3 particle decays (up to now)");
533  return 0;
534  }
535  G4Exception ("G4VDecayChannel::GetAngularMomentum",
536  "PART111", JustWarning,
537  "Can't find angular momentum for this decay");
538  return 0;
539 }
540 
542 {
543  G4cout << " BR: " << rbranch << " [" << kinematics_name << "]";
544  G4cout << " : " ;
545  for (G4int index=0; index < numberOfDaughters; index++){
546  if(daughters_name[index] != 0) {
547  G4cout << " " << *(daughters_name[index]);
548  } else {
549  G4cout << " not defined ";
550  }
551  }
552  G4cout << G4endl;
553 }
554 
556 {
557  return noName;
558 }
559 
560 #include "Randomize.hh"
562 {
563  if (maxDev >rangeMass) maxDev = rangeMass;
564  if (maxDev <=-1.*rangeMass) return massPDG; // can not calculate
565 
566  G4double x = G4UniformRand()*(maxDev+rangeMass) - rangeMass;
567  G4double y = G4UniformRand();
568  const size_t MAX_LOOP=10000;
569  for (size_t loop_counter=0; loop_counter <MAX_LOOP; ++loop_counter){
570  if ( y * (width*width*x*x + massPDG*massPDG*width*width) <= massPDG*massPDG*width*width ) break;
571  x = G4UniformRand()*(maxDev+rangeMass) - rangeMass;
572  y = G4UniformRand();
573  }
574  G4double mass = massPDG + x*width;
575  return mass;
576 }
577 
579 {
580  G4double sumOfDaughterMassMin=0.0;
581  if (G4MT_parent == 0) FillParent();
582  if (G4MT_daughters == 0) FillDaughters();
583 
584  for (G4int index=0; index < numberOfDaughters; index++) {
585  sumOfDaughterMassMin +=
587  }
588  return (parentMass > sumOfDaughterMassMin);
589 }
const G4String & GetNoName() const
static const G4String noName
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
void SetBR(G4double value)
G4double * G4MT_daughters_mass
G4ParticleDefinition * G4MT_parent
G4ParticleDefinition ** G4MT_daughters
#define width
G4double G4MT_parent_mass
virtual G4bool IsOKWithParentMass(G4double parentMass)
G4ParticleTable * particletable
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4int GetAngularMomentum()
G4double GetPDGWidth() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
G4String kinematics_name
void SetNumberOfDaughters(G4int value)
G4String * parent_name
static const double GeV
Definition: G4SIunits.hh:214
const G4String & GetParticleType() const
G4VDecayChannel & operator=(const G4VDecayChannel &)
G4ThreeVector parent_polarization
virtual ~G4VDecayChannel()
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4int GetVerboseLevel() const
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
G4double * G4MT_daughters_width
const G4double x[NPOINTSGL]
void SetParent(const G4ParticleDefinition *particle_type)
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
#define G4endl
Definition: G4ios.hh:61
G4double DynamicalMass(G4double massPDG, G4double width, G4double maxDev=+1.) const
double G4double
Definition: G4Types.hh:76
G4String ** daughters_name