Geant4_10
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 77085 2013-11-21 10:37:09Z 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  particletable(0),
70  verboseLevel(1)
71 {
73  G4MT_parent = 0;
74  G4MT_daughters = 0;
75  G4MT_parent_mass = 0.0;
77 
78  // set pointer to G4ParticleTable (static and singleton object)
80 }
81 
83  :kinematics_name(aName),
84  rbranch(0.0),
85  numberOfDaughters(0),
86  parent_name(0), daughters_name(0),
87  particletable(0),
88  verboseLevel(Verbose)
89 {
91  G4MT_parent = 0;
92  G4MT_daughters = 0;
93  G4MT_parent_mass = 0.0;
95 
96  // set pointer to G4ParticleTable (static and singleton object)
98 }
99 
101  const G4String& theParentName,
102  G4double theBR,
103  G4int theNumberOfDaughters,
104  const G4String& theDaughterName1,
105  const G4String& theDaughterName2,
106  const G4String& theDaughterName3,
107  const G4String& theDaughterName4 )
108  :kinematics_name(aName),
109  rbranch(theBR),
110  numberOfDaughters(theNumberOfDaughters),
111  parent_name(0), daughters_name(0),
112  particletable(0),
113  verboseLevel(1)
114 {
116  G4MT_parent = 0;
117  G4MT_daughters = 0;
118  G4MT_parent_mass = 0.0;
120 
121  // set pointer to G4ParticleTable (static and singleton object)
123 
124  // parent name
125  parent_name = new G4String(theParentName);
126 
127  // cleate array
130 
131  // daughters' name
132  if (numberOfDaughters>0) daughters_name[0] = new G4String(theDaughterName1);
133  if (numberOfDaughters>1) daughters_name[1] = new G4String(theDaughterName2);
134  if (numberOfDaughters>2) daughters_name[2] = new G4String(theDaughterName3);
135  if (numberOfDaughters>3) daughters_name[3] = new G4String(theDaughterName4);
136 }
137 
139 {
141 
143  verboseLevel = right.verboseLevel;
144  rbranch = right.rbranch;
145 
146  // copy parent name
147  parent_name = new G4String(*right.parent_name);
148  G4MT_parent = 0;
149  G4MT_parent_mass = 0.0;
150 
151  //create array
153 
154  daughters_name =0;
155  if ( numberOfDaughters >0 ) {
157  //copy daughters name
158  for (G4int index=0; index < numberOfDaughters; index++)
159  {
161  }
162  }
163 
164  //
166  G4MT_daughters = 0;
167 
168  // particle table
170 }
171 
173 {
174  if (this != &right) {
176  verboseLevel = right.verboseLevel;
177  rbranch = right.rbranch;
178 
179  // copy parent name
180  parent_name = new G4String(*right.parent_name);
181 
182  // clear daughters_name array
184 
185  // recreate array
187  if ( numberOfDaughters >0 ) {
190  //copy daughters name
191  for (G4int index=0; index < numberOfDaughters; index++) {
193  }
194  }
195  }
196 
197  //
198  G4MT_parent = 0;
199  G4MT_daughters = 0;
200  G4MT_parent_mass = 0.0;
202 
203  // particle table
205 
206  return *this;
207 }
208 
210 {
212  if (parent_name != 0) delete parent_name;
213  parent_name = 0;
214  if (G4MT_daughters_mass != 0) delete [] G4MT_daughters_mass;
216 }
217 
222 
224 {
225  if ( daughters_name != 0) {
226  if (numberOfDaughters>0) {
227 #ifdef G4VERBOSE
228  if (verboseLevel>1) {
229  G4cout << "G4VDecayChannel::ClearDaughtersName "
230  << " for " << *parent_name << G4endl;
231  }
232 #endif
233  for (G4int index=0; index < numberOfDaughters; index++) {
234  if (daughters_name[index] != 0) delete daughters_name[index];
235  }
236  }
237  delete [] daughters_name;
238  daughters_name = 0;
239  }
240  //
241  if (G4MT_daughters != 0) delete [] G4MT_daughters;
242  if (G4MT_daughters_mass != 0) delete [] G4MT_daughters_mass;
243  G4MT_daughters = 0;
245 
246  numberOfDaughters = 0;
247 }
248 
250 {
251  if (size >0) {
252  // remove old contents
254  // cleate array
255  daughters_name = new G4String*[size];
256  for (G4int index=0;index<size;index++) daughters_name[index]=0;
257  numberOfDaughters = size;
258  }
259 }
260 
262  const G4String &particle_name)
263 {
264  // check numberOfDaughters is positive
265  if (numberOfDaughters<=0) {
266 #ifdef G4VERBOSE
267  if (verboseLevel>0) {
268  G4cout << "G4VDecayChannel::SetDaughter: "
269  << "Number of daughters is not defined" << G4endl;
270  }
271 #endif
272  return;
273  }
274 
275  // check existence of daughters_name array
276  if (daughters_name == 0) {
277  // cleate array
281  }
282  }
283 
284  // check an index
285  if ( (anIndex<0) || (anIndex>=numberOfDaughters) ) {
286 #ifdef G4VERBOSE
287  if (verboseLevel>0) {
288  G4cout << "G4VDecayChannel::SetDaughter"
289  << "index out of range " << anIndex << G4endl;
290  }
291 #endif
292  } else {
293  // delete the old name if it exists
294  if (daughters_name[anIndex]!=0) delete daughters_name[anIndex];
295  // fill the name
296  daughters_name[anIndex] = new G4String(particle_name);
297  // refill the array of daughters[] if it exists
298  if (G4MT_daughters != 0) FillDaughters();
299 #ifdef G4VERBOSE
300  if (verboseLevel>1) {
301  G4cout << "G4VDecayChannel::SetDaughter[" << anIndex <<"] :";
302  G4cout << daughters_name[anIndex] << ":" << *daughters_name[anIndex]<<G4endl;
303  }
304 #endif
305  }
306 }
307 
308 void G4VDecayChannel::SetDaughter(G4int anIndex, const G4ParticleDefinition * parent_type)
309 {
310  if (parent_type != 0) SetDaughter(anIndex, parent_type->GetParticleName());
311 }
312 
314 {
315  G4int index;
316 
317 #ifdef G4VERBOSE
318  if (verboseLevel>1) G4cout << "G4VDecayChannel::FillDaughters()" <<G4endl;
319 #endif
320  if (G4MT_daughters != 0) {
321  delete [] G4MT_daughters;
322  G4MT_daughters = 0;
323  }
324 
325  // parent mass
326  if (G4MT_parent == 0) FillParent();
327  G4double parentmass = G4MT_parent->GetPDGMass();
328 
329  //
330  G4double sumofdaughtermass = 0.0;
331  if ((numberOfDaughters <=0) || (daughters_name == 0) ){
332 #ifdef G4VERBOSE
333  if (verboseLevel>0) {
334  G4cout << "G4VDecayChannel::FillDaughters "
335  << "[ " << G4MT_parent->GetParticleName() << " ]"
336  << "numberOfDaughters is not defined yet";
337  }
338 #endif
339  G4MT_daughters = 0;
340  G4Exception("G4VDecayChannel::FillDaughters",
341  "PART011", FatalException,
342  "Can not fill daughters: numberOfDaughters is not defined yet");
343  }
344 
345  //create and set the array of pointers to daughter particles
347  if (G4MT_daughters_mass != 0) delete [] G4MT_daughters_mass;
349  // loop over all daughters
350  for (index=0; index < numberOfDaughters; index++) {
351  if (daughters_name[index] == 0) {
352  // daughter name is not defined
353 #ifdef G4VERBOSE
354  if (verboseLevel>0) {
355  G4cout << "G4VDecayChannel::FillDaughters "
356  << "[ " << G4MT_parent->GetParticleName() << " ]"
357  << index << "-th daughter is not defined yet" << G4endl;
358  }
359 #endif
360  G4MT_daughters[index] = 0;
361  G4Exception("G4VDecayChannel::FillDaughters",
362  "PART011", FatalException,
363  "Can not fill daughters: name of a daughter is not defined yet");
364  }
365  //search daughter particles in the particle table
367  if (G4MT_daughters[index] == 0) {
368  // can not find the daughter particle
369 #ifdef G4VERBOSE
370  if (verboseLevel>0) {
371  G4cout << "G4VDecayChannel::FillDaughters "
372  << "[ " << G4MT_parent->GetParticleName() << " ]"
373  << index << ":" << *daughters_name[index]
374  << " is not defined !!" << G4endl;
375  G4cout << " The BR of this decay mode is set to zero " << G4endl;
376  }
377 #endif
378  SetBR(0.0);
379  return;
380  }
381 #ifdef G4VERBOSE
382  if (verboseLevel>1) {
383  G4cout << index << ":" << *daughters_name[index];
384  G4cout << ":" << G4MT_daughters[index] << G4endl;
385  }
386 #endif
388  sumofdaughtermass += G4MT_daughters[index]->GetPDGMass();
389  } // end loop over all daughters
390 
391  // check sum of daghter mass
392  G4double widthMass = G4MT_parent->GetPDGWidth();
393  if ( (G4MT_parent->GetParticleType() != "nucleus") &&
394  (sumofdaughtermass > parentmass + 5*widthMass) ){
395  // !!! illegal mass !!!
396 #ifdef G4VERBOSE
397  if (GetVerboseLevel()>0) {
398  G4cout << "G4VDecayChannel::FillDaughters "
399  << "[ " << G4MT_parent->GetParticleName() << " ]"
400  << " Energy/Momentum conserevation breaks " <<G4endl;
401  if (GetVerboseLevel()>1) {
402  G4cout << " parent:" << *parent_name
403  << " mass:" << parentmass/GeV << "[GeV/c/c]" <<G4endl;
404  for (index=0; index < numberOfDaughters; index++){
405  G4cout << " daughter " << index << ":" << *daughters_name[index]
406  << " mass:" << G4MT_daughters[index]->GetPDGMass()/GeV
407  << "[GeV/c/c]" <<G4endl;
408  }
409  }
410  }
411 #endif
412  }
413 }
414 
415 
417 {
418  if (parent_name == 0) {
419  // parent name is not defined
420 #ifdef G4VERBOSE
421  if (verboseLevel>0) {
422  G4cout << "G4VDecayChannel::FillParent "
423  << ": parent name is not defined !!" << G4endl;
424  }
425 #endif
426  G4MT_parent = 0;
427  G4Exception("G4VDecayChannel::FillParent()",
428  "PART012", FatalException,
429  "Can not fill parent: parent name is not defined yet");
430  return;
431  }
432  // search parent particle in the particle table
434  if (G4MT_parent == 0) {
435  // parent particle does not exist
436 #ifdef G4VERBOSE
437  if (verboseLevel>0) {
438  G4cout << "G4VDecayChannel::FillParent "
439  << *parent_name << " does not exist !!" << G4endl;
440  }
441 #endif
442  G4Exception("G4VDecayChannel::FillParent()",
443  "PART012", FatalException,
444  "Can not fill parent: parent does not exist");
445  return;
446  }
448 }
449 
451 {
452  if (parent_type != 0) SetParent(parent_type->GetParticleName());
453 }
454 
456 {
457  // determine angular momentum
458 
459  // fill pointers to daughter particles if not yet set
460  if (G4MT_daughters == 0) FillDaughters();
461 
462  const G4int PiSpin = G4MT_parent->GetPDGiSpin();
463  const G4int PParity = G4MT_parent->GetPDGiParity();
464  if (2==numberOfDaughters) { // up to now we can only handle two particle decays
465  const G4int D1iSpin = G4MT_daughters[0]->GetPDGiSpin();
466  const G4int D1Parity = G4MT_daughters[0]->GetPDGiParity();
467  const G4int D2iSpin = G4MT_daughters[1]->GetPDGiSpin();
468  const G4int D2Parity = G4MT_daughters[1]->GetPDGiParity();
469  const G4int MiniSpin = std::abs (D1iSpin - D2iSpin);
470  const G4int MaxiSpin = D1iSpin + D2iSpin;
471  const G4int lMax = (PiSpin+D1iSpin+D2iSpin)/2; // l is allways int
472  G4int lMin;
473 #ifdef G4VERBOSE
474  if (verboseLevel>1) {
475  G4cout << "iSpin: " << PiSpin << " -> " << D1iSpin << " + " << D2iSpin << G4endl;
476  G4cout << "2*jmin, 2*jmax, lmax " << MiniSpin << " " << MaxiSpin << " " << lMax << G4endl;
477  }
478 #endif
479  for (G4int j=MiniSpin; j<=MaxiSpin; j+=2){ // loop over all possible spin couplings
480  lMin = std::abs(PiSpin-j)/2;
481 #ifdef G4VERBOSE
482  if (verboseLevel>1)
483  G4cout << "-> checking 2*j=" << j << G4endl;
484 #endif
485  for (G4int l=lMin; l<=lMax; l++) {
486 #ifdef G4VERBOSE
487  if (verboseLevel>1)
488  G4cout << " checking l=" << l << G4endl;
489 #endif
490  if (l%2==0) {
491  if (PParity == D1Parity*D2Parity) { // check parity for this l
492  return l;
493  }
494  } else {
495  if (PParity == -1*D1Parity*D2Parity) { // check parity for this l
496  return l;
497  }
498  }
499  }
500  }
501  } else {
502  G4Exception("G4VDecayChannel::GetAngularMomentum",
503  "PART111", JustWarning,
504  "Sorry, can't handle 3 particle decays (up to now)");
505  return 0;
506  }
507  G4Exception ("G4VDecayChannel::GetAngularMomentum",
508  "PART111", JustWarning,
509  "Can't find angular momentum for this decay");
510  return 0;
511 }
512 
514 {
515  G4cout << " BR: " << rbranch << " [" << kinematics_name << "]";
516  G4cout << " : " ;
517  for (G4int index=0; index < numberOfDaughters; index++)
518  {
519  if(daughters_name[index] != 0) {
520  G4cout << " " << *(daughters_name[index]);
521  } else {
522  G4cout << " not defined ";
523  }
524  }
525  G4cout << G4endl;
526 }
527 
528 const G4String& G4VDecayChannel::GetNoName() const
529 {
530  return noName;
531 }
static const G4String noName
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
void SetBR(G4double value)
G4double * G4MT_daughters_mass
Int_t index
Definition: macro.C:9
G4ParticleDefinition * G4MT_parent
G4ParticleDefinition ** G4MT_daughters
G4double G4MT_parent_mass
G4ParticleTable * particletable
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4int GetAngularMomentum()
G4double GetPDGWidth() const
G4GLOB_DLL std::ostream G4cout
G4String kinematics_name
void SetNumberOfDaughters(G4int value)
G4String * parent_name
const G4String & GetParticleType() const
G4VDecayChannel & operator=(const G4VDecayChannel &)
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()
void SetParent(const G4ParticleDefinition *particle_type)
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4String ** daughters_name