Geant4  10.00.p02
G4GeneralPhaseSpaceDecay.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: G4GeneralPhaseSpaceDecay.cc 67984 2013-03-13 10:44:01Z gcosmo $
28 // ----------------------------------------------------------------
29 // GEANT 4 class header file
30 //
31 // History: first implementation, A. Feliciello, 21st May 1998
32 //
33 // Note: this class is a generalization of the
34 // G4PhaseSpaceDecayChannel one
35 // ----------------------------------------------------------------
36 
37 #include "G4ParticleDefinition.hh"
38 #include "G4DecayProducts.hh"
39 #include "G4VDecayChannel.hh"
41 #include "G4PhysicalConstants.hh"
42 #include "G4SystemOfUnits.hh"
43 #include "Randomize.hh"
44 #include "G4LorentzVector.hh"
45 #include "G4LorentzRotation.hh"
46 #include "G4ios.hh"
47 
48 
50  G4VDecayChannel("Phase Space", Verbose),
51  parentmass(0.), theDaughterMasses(0)
52 {
53  if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay:: constructor " << G4endl;
54 }
55 
57  G4double theBR,
58  G4int theNumberOfDaughters,
59  const G4String& theDaughterName1,
60  const G4String& theDaughterName2,
61  const G4String& theDaughterName3) :
62  G4VDecayChannel("Phase Space",
63  theParentName,theBR,
64  theNumberOfDaughters,
65  theDaughterName1,
66  theDaughterName2,
67  theDaughterName3),
68  theDaughterMasses(0)
69 {
70  if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay:: constructor " << G4endl;
71 
72  // Set the parent particle (resonance) mass to the (default) PDG vale
73  if (G4MT_parent != NULL)
74  {
76  } else {
77  parentmass=0.;
78  }
79 
80 }
81 
83  G4double theParentMass,
84  G4double theBR,
85  G4int theNumberOfDaughters,
86  const G4String& theDaughterName1,
87  const G4String& theDaughterName2,
88  const G4String& theDaughterName3) :
89  G4VDecayChannel("Phase Space",
90  theParentName,theBR,
91  theNumberOfDaughters,
92  theDaughterName1,
93  theDaughterName2,
94  theDaughterName3),
95  parentmass(theParentMass),
96  theDaughterMasses(0)
97 {
98  if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay:: constructor " << G4endl;
99 }
100 
102  G4double theParentMass,
103  G4double theBR,
104  G4int theNumberOfDaughters,
105  const G4String& theDaughterName1,
106  const G4String& theDaughterName2,
107  const G4String& theDaughterName3,
108  const G4double *masses) :
109  G4VDecayChannel("Phase Space",
110  theParentName,theBR,
111  theNumberOfDaughters,
112  theDaughterName1,
113  theDaughterName2,
114  theDaughterName3),
115  parentmass(theParentMass),
116  theDaughterMasses(masses)
117 {
118  if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay:: constructor " << G4endl;
119 }
120 
122 {
123 }
124 
126 {
127  if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay::DecayIt ";
128  G4DecayProducts * products = NULL;
129 
130  if (G4MT_parent == NULL) FillParent();
131  if (G4MT_daughters == NULL) FillDaughters();
132 
133  switch (numberOfDaughters){
134  case 0:
135  if (GetVerboseLevel()>0) {
136  G4cout << "G4GeneralPhaseSpaceDecay::DecayIt ";
137  G4cout << " daughters not defined " <<G4endl;
138  }
139  break;
140  case 1:
141  products = OneBodyDecayIt();
142  break;
143  case 2:
144  products = TwoBodyDecayIt();
145  break;
146  case 3:
147  products = ThreeBodyDecayIt();
148  break;
149  default:
150  products = ManyBodyDecayIt();
151  break;
152  }
153  if ((products == NULL) && (GetVerboseLevel()>0)) {
154  G4cout << "G4GeneralPhaseSpaceDecay::DecayIt ";
155  G4cout << *parent_name << " can not decay " << G4endl;
156  DumpInfo();
157  }
158  return products;
159 }
160 
162 {
163  if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay::OneBodyDecayIt()"<<G4endl;
164 
165 // G4double daughtermass = daughters[0]->GetPDGMass();
166 
167  //create parent G4DynamicParticle at rest
168  G4ParticleMomentum dummy;
169  G4DynamicParticle * parentparticle = new G4DynamicParticle(G4MT_parent, dummy, 0.0);
170 
171  //create G4Decayproducts
172  G4DecayProducts *products = new G4DecayProducts(*parentparticle);
173  delete parentparticle;
174 
175  //create daughter G4DynamicParticle at rest
176  G4DynamicParticle * daughterparticle = new G4DynamicParticle(G4MT_daughters[0], dummy, 0.0);
177  products->PushProducts(daughterparticle);
178 
179  if (GetVerboseLevel()>1)
180  {
181  G4cout << "G4GeneralPhaseSpaceDecay::OneBodyDecayIt ";
182  G4cout << " create decay products in rest frame " <<G4endl;
183  products->DumpInfo();
184  }
185  return products;
186 }
187 
189 {
190  if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay::TwoBodyDecayIt()"<<G4endl;
191 
192  //daughters'mass
193  G4double daughtermass[2];
194  G4double daughtermomentum;
195  if ( theDaughterMasses )
196  {
197  daughtermass[0]= *(theDaughterMasses);
198  daughtermass[1] = *(theDaughterMasses+1);
199  } else {
200  daughtermass[0] = G4MT_daughters[0]->GetPDGMass();
201  daughtermass[1] = G4MT_daughters[1]->GetPDGMass();
202  }
203 
204 // G4double sumofdaughtermass = daughtermass[0] + daughtermass[1];
205 
206  //create parent G4DynamicParticle at rest
207  G4ParticleMomentum dummy;
208  G4DynamicParticle * parentparticle = new G4DynamicParticle( G4MT_parent, dummy, 0.0);
209 
210  //create G4Decayproducts @@GF why dummy parentparticle?
211  G4DecayProducts *products = new G4DecayProducts(*parentparticle);
212  delete parentparticle;
213 
214  //calculate daughter momentum
215  daughtermomentum = Pmx(parentmass,daughtermass[0],daughtermass[1]);
216  G4double costheta = 2.*G4UniformRand()-1.0;
217  G4double sintheta = std::sqrt((1.0 - costheta)*(1.0 + costheta));
218  G4double phi = twopi*G4UniformRand()*rad;
219  G4ParticleMomentum direction(sintheta*std::cos(phi),sintheta*std::sin(phi),costheta);
220 
221  //create daughter G4DynamicParticle
222  G4double Etotal= std::sqrt(daughtermass[0]*daughtermass[0] + daughtermomentum*daughtermomentum);
223  G4DynamicParticle * daughterparticle = new G4DynamicParticle( G4MT_daughters[0],Etotal, direction*daughtermomentum);
224  products->PushProducts(daughterparticle);
225  Etotal= std::sqrt(daughtermass[1]*daughtermass[1] + daughtermomentum*daughtermomentum);
226  daughterparticle = new G4DynamicParticle( G4MT_daughters[1],Etotal, direction*(-1.0*daughtermomentum));
227  products->PushProducts(daughterparticle);
228 
229  if (GetVerboseLevel()>1)
230  {
231  G4cout << "G4GeneralPhaseSpaceDecay::TwoBodyDecayIt ";
232  G4cout << " create decay products in rest frame " <<G4endl;
233  products->DumpInfo();
234  }
235  return products;
236 }
237 
239 // algorism of this code is originally written in GDECA3 of GEANT3
240 {
241  if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay::ThreeBodyDecayIt()"<<G4endl;
242 
243  //daughters'mass
244  G4double daughtermass[3];
245  G4double sumofdaughtermass = 0.0;
246  for (G4int index=0; index<3; index++)
247  {
248  if ( theDaughterMasses )
249  {
250  daughtermass[index]= *(theDaughterMasses+index);
251  } else {
252  daughtermass[index] = G4MT_daughters[index]->GetPDGMass();
253  }
254  sumofdaughtermass += daughtermass[index];
255  }
256 
257  //create parent G4DynamicParticle at rest
258  G4ParticleMomentum dummy;
259  G4DynamicParticle * parentparticle = new G4DynamicParticle( G4MT_parent, dummy, 0.0);
260 
261  //create G4Decayproducts
262  G4DecayProducts *products = new G4DecayProducts(*parentparticle);
263  delete parentparticle;
264 
265  //calculate daughter momentum
266  // Generate two
267  G4double rd1, rd2, rd;
268  G4double daughtermomentum[3];
269  G4double momentummax=0.0, momentumsum = 0.0;
271 
272  do
273  {
274  rd1 = G4UniformRand();
275  rd2 = G4UniformRand();
276  if (rd2 > rd1)
277  {
278  rd = rd1;
279  rd1 = rd2;
280  rd2 = rd;
281  }
282  momentummax = 0.0;
283  momentumsum = 0.0;
284  // daughter 0
285 
286  energy = rd2*(parentmass - sumofdaughtermass);
287  daughtermomentum[0] = std::sqrt(energy*energy + 2.0*energy* daughtermass[0]);
288  if ( daughtermomentum[0] >momentummax )momentummax = daughtermomentum[0];
289  momentumsum += daughtermomentum[0];
290 
291  // daughter 1
292  energy = (1.-rd1)*(parentmass - sumofdaughtermass);
293  daughtermomentum[1] = std::sqrt(energy*energy + 2.0*energy* daughtermass[1]);
294  if ( daughtermomentum[1] >momentummax )momentummax = daughtermomentum[1];
295  momentumsum += daughtermomentum[1];
296 
297  // daughter 2
298  energy = (rd1-rd2)*(parentmass - sumofdaughtermass);
299  daughtermomentum[2] = std::sqrt(energy*energy + 2.0*energy* daughtermass[2]);
300  if ( daughtermomentum[2] >momentummax )momentummax = daughtermomentum[2];
301  momentumsum += daughtermomentum[2];
302  } while (momentummax > momentumsum - momentummax );
303 
304  // output message
305  if (GetVerboseLevel()>1) {
306  G4cout << " daughter 0:" << daughtermomentum[0]/GeV << "[GeV/c]" <<G4endl;
307  G4cout << " daughter 1:" << daughtermomentum[1]/GeV << "[GeV/c]" <<G4endl;
308  G4cout << " daughter 2:" << daughtermomentum[2]/GeV << "[GeV/c]" <<G4endl;
309  G4cout << " momentum sum:" << momentumsum/GeV << "[GeV/c]" <<G4endl;
310  }
311 
312  //create daughter G4DynamicParticle
313  G4double costheta, sintheta, phi, sinphi, cosphi;
314  G4double costhetan, sinthetan, phin, sinphin, cosphin;
315  costheta = 2.*G4UniformRand()-1.0;
316  sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
317  phi = twopi*G4UniformRand()*rad;
318  sinphi = std::sin(phi);
319  cosphi = std::cos(phi);
320  G4ParticleMomentum direction0(sintheta*cosphi,sintheta*sinphi,costheta);
321  G4double Etotal=std::sqrt( daughtermass[0]*daughtermass[0] + daughtermomentum[0]*daughtermomentum[0]);
322  G4DynamicParticle * daughterparticle
323  = new G4DynamicParticle( G4MT_daughters[0], Etotal, direction0*daughtermomentum[0]);
324  products->PushProducts(daughterparticle);
325 
326  costhetan = (daughtermomentum[1]*daughtermomentum[1]-daughtermomentum[2]*daughtermomentum[2]-daughtermomentum[0]*daughtermomentum[0])/(2.0*daughtermomentum[2]*daughtermomentum[0]);
327  sinthetan = std::sqrt((1.0-costhetan)*(1.0+costhetan));
328  phin = twopi*G4UniformRand()*rad;
329  sinphin = std::sin(phin);
330  cosphin = std::cos(phin);
331  G4ParticleMomentum direction2;
332  direction2.setX( sinthetan*cosphin*costheta*cosphi - sinthetan*sinphin*sinphi + costhetan*sintheta*cosphi);
333  direction2.setY( sinthetan*cosphin*costheta*sinphi + sinthetan*sinphin*cosphi + costhetan*sintheta*sinphi);
334  direction2.setZ( -sinthetan*cosphin*sintheta + costhetan*costheta);
335  Etotal=std::sqrt( daughtermass[2]*daughtermass[2] + daughtermomentum[2]*daughtermomentum[2]/direction2.mag2());
336  daughterparticle = new G4DynamicParticle( G4MT_daughters[2],Etotal, direction2*(daughtermomentum[2]/direction2.mag()));
337  products->PushProducts(daughterparticle);
338  G4ThreeVector mom=(direction0*daughtermomentum[0] + direction2*(daughtermomentum[2]/direction2.mag()))*(-1.0);
339  Etotal= std::sqrt( daughtermass[1]*daughtermass[1] + mom.mag2() );
340  daughterparticle =
341  new G4DynamicParticle(G4MT_daughters[1], Etotal, mom);
342  products->PushProducts(daughterparticle);
343 
344  if (GetVerboseLevel()>1) {
345  G4cout << "G4GeneralPhaseSpaceDecay::ThreeBodyDecayIt ";
346  G4cout << " create decay products in rest frame " <<G4endl;
347  products->DumpInfo();
348  }
349  return products;
350 }
351 
353 // algorism of this code is originally written in FORTRAN by M.Asai
354 //*****************************************************************
355 // NBODY
356 // N-body phase space Monte-Carlo generator
357 // Makoto Asai
358 // Hiroshima Institute of Technology
359 // (asai@kekvax.kek.jp)
360 // Revised release : 19/Apr/1995
361 //
362 {
363  //return value
364  G4DecayProducts *products;
365 
366  if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay::ManyBodyDecayIt()"<<G4endl;
367 
368  //daughters'mass
369  G4double *daughtermass = new G4double[numberOfDaughters];
370  G4double sumofdaughtermass = 0.0;
371  for (G4int index=0; index<numberOfDaughters; index++){
372  daughtermass[index] = G4MT_daughters[index]->GetPDGMass();
373  sumofdaughtermass += daughtermass[index];
374  }
375 
376  //Calculate daughter momentum
377  G4double *daughtermomentum = new G4double[numberOfDaughters];
378  G4ParticleMomentum direction;
379  G4DynamicParticle **daughterparticle;
381  G4double tmas;
382  G4double weight = 1.0;
383  G4int numberOfTry = 0;
384  G4int index1;
385 
386  do {
387  //Generate rundom number in descending order
388  G4double temp;
390  rd[0] = 1.0;
391  for(index1 =1; index1 < numberOfDaughters -1; index1++)
392  rd[index1] = G4UniformRand();
393  rd[ numberOfDaughters -1] = 0.0;
394  for(index1 =1; index1 < numberOfDaughters -1; index1++) {
395  for(G4int index2 = index1+1; index2 < numberOfDaughters; index2++) {
396  if (rd[index1] < rd[index2]){
397  temp = rd[index1];
398  rd[index1] = rd[index2];
399  rd[index2] = temp;
400  }
401  }
402  }
403 
404  //calcurate virtual mass
405  tmas = parentmass - sumofdaughtermass;
406  temp = sumofdaughtermass;
407  for(index1 =0; index1 < numberOfDaughters; index1++) {
408  sm[index1] = rd[index1]*tmas + temp;
409  temp -= daughtermass[index1];
410  if (GetVerboseLevel()>1) {
411  G4cout << index1 << " rundom number:" << rd[index1];
412  G4cout << " virtual mass:" << sm[index1]/GeV << "[GeV/c/c]" <<G4endl;
413  }
414  }
415  delete [] rd;
416 
417  //Calculate daughter momentum
418  weight = 1.0;
419  index1 =numberOfDaughters-1;
420  daughtermomentum[index1]= Pmx( sm[index1-1],daughtermass[index1-1],sm[index1]);
421  if (GetVerboseLevel()>1) {
422  G4cout << " daughter " << index1 << ":" << *daughters_name[index1];
423  G4cout << " momentum:" << daughtermomentum[index1]/GeV << "[GeV/c]" <<G4endl;
424  }
425  for(index1 =numberOfDaughters-2; index1>=0; index1--) {
426  // calculate
427  daughtermomentum[index1]= Pmx( sm[index1],daughtermass[index1], sm[index1 +1]);
428  if(daughtermomentum[index1] < 0.0) {
429  // !!! illegal momentum !!!
430  if (GetVerboseLevel()>0) {
431  G4cout << "G4GeneralPhaseSpaceDecay::ManyBodyDecayIt ";
432  G4cout << " can not calculate daughter momentum " <<G4endl;
433  G4cout << " parent:" << *parent_name;
434  G4cout << " mass:" << parentmass/GeV << "[GeV/c/c]" <<G4endl;
435  G4cout << " daughter " << index1 << ":" << *daughters_name[index1];
436  G4cout << " mass:" << daughtermass[index1]/GeV << "[GeV/c/c]" ;
437  G4cout << " mass:" << daughtermomentum[index1]/GeV << "[GeV/c]" <<G4endl;
438  }
439  delete [] sm;
440  delete [] daughtermass;
441  delete [] daughtermomentum;
442  return NULL; // Error detection
443 
444  } else {
445  // calculate weight of this events
446  weight *= daughtermomentum[index1]/sm[index1];
447  if (GetVerboseLevel()>1) {
448  G4cout << " daughter " << index1 << ":" << *daughters_name[index1];
449  G4cout << " momentum:" << daughtermomentum[index1]/GeV << "[GeV/c]" <<G4endl;
450  }
451  }
452  }
453  if (GetVerboseLevel()>1) {
454  G4cout << " weight: " << weight <<G4endl;
455  }
456 
457  // exit if number of Try exceeds 100
458  if (numberOfTry++ >100) {
459  if (GetVerboseLevel()>0) {
460  G4cout << "G4GeneralPhaseSpaceDecay::ManyBodyDecayIt: ";
461  G4cout << " can not determine Decay Kinematics " << G4endl;
462  }
463  delete [] sm;
464  delete [] daughtermass;
465  delete [] daughtermomentum;
466  return NULL; // Error detection
467  }
468  } while ( weight > G4UniformRand());
469  if (GetVerboseLevel()>1) {
470  G4cout << "Start calulation of daughters momentum vector "<<G4endl;
471  }
472 
473  G4double costheta, sintheta, phi;
474  G4double beta;
475  daughterparticle = new G4DynamicParticle*[numberOfDaughters];
476 
477  index1 = numberOfDaughters -2;
478  costheta = 2.*G4UniformRand()-1.0;
479  sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
480  phi = twopi*G4UniformRand()*rad;
481  direction.setZ(costheta);
482  direction.setY(sintheta*std::sin(phi));
483  direction.setX(sintheta*std::cos(phi));
484  daughterparticle[index1] = new G4DynamicParticle( G4MT_daughters[index1], direction*daughtermomentum[index1] );
485  daughterparticle[index1+1] = new G4DynamicParticle( G4MT_daughters[index1+1], direction*(-1.0*daughtermomentum[index1]) );
486 
487  for (index1 = numberOfDaughters -3; index1 >= 0; index1--) {
488  //calculate momentum direction
489  costheta = 2.*G4UniformRand()-1.0;
490  sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
491  phi = twopi*G4UniformRand()*rad;
492  direction.setZ(costheta);
493  direction.setY(sintheta*std::sin(phi));
494  direction.setX(sintheta*std::cos(phi));
495 
496  // boost already created particles
497  beta = daughtermomentum[index1];
498  beta /= std::sqrt( daughtermomentum[index1]*daughtermomentum[index1] + sm[index1+1]*sm[index1+1] );
499  for (G4int index2 = index1+1; index2<numberOfDaughters; index2++) {
500  G4LorentzVector p4;
501  // make G4LorentzVector for secondaries
502  p4 = daughterparticle[index2]->Get4Momentum();
503 
504  // boost secondaries to new frame
505  p4.boost( direction.x()*beta, direction.y()*beta, direction.z()*beta);
506 
507  // change energy/momentum
508  daughterparticle[index2]->Set4Momentum(p4);
509  }
510  //create daughter G4DynamicParticle
511  daughterparticle[index1]= new G4DynamicParticle( G4MT_daughters[index1], direction*(-1.0*daughtermomentum[index1]));
512  }
513 
514  //create G4Decayproducts
515  G4DynamicParticle *parentparticle;
516  direction.setX(1.0); direction.setY(0.0); direction.setZ(0.0);
517  parentparticle = new G4DynamicParticle( G4MT_parent, direction, 0.0);
518  products = new G4DecayProducts(*parentparticle);
519  delete parentparticle;
520  for (index1 = 0; index1<numberOfDaughters; index1++) {
521  products->PushProducts(daughterparticle[index1]);
522  }
523  if (GetVerboseLevel()>1) {
524  G4cout << "G4GeneralPhaseSpaceDecay::ManyBodyDecayIt ";
525  G4cout << " create decay products in rest frame " << G4endl;
526  products->DumpInfo();
527  }
528 
529  delete [] daughterparticle;
530  delete [] daughtermomentum;
531  delete [] daughtermass;
532  delete [] sm;
533 
534  return products;
535 }
536 
537 
538 
539 
540 
CLHEP::Hep3Vector G4ThreeVector
G4int PushProducts(G4DynamicParticle *aParticle)
G4ParticleDefinition * G4MT_parent
G4ParticleDefinition ** G4MT_daughters
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
void DumpInfo() const
G4String * parent_name
static const double GeV
Definition: G4SIunits.hh:196
G4LorentzVector Get4Momentum() const
static const double rad
Definition: G4SIunits.hh:130
virtual G4DecayProducts * DecayIt(G4double mass=0.0)
void Set4Momentum(const G4LorentzVector &momentum)
G4int GetVerboseLevel() const
static G4double Pmx(G4double e, G4double p1, G4double p2)
G4double GetPDGMass() const
G4double energy(const ThreeVector &p, const G4double m)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4ThreeVector G4ParticleMomentum
G4String ** daughters_name
CLHEP::HepLorentzVector G4LorentzVector