Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4GeneralPhaseSpaceDecay.cc
이 파일의 문서화 페이지로 가기
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 95909 2016-03-02 11:01:47Z 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  G4double theParentMass,
123  G4double theBR,
124  G4int theNumberOfDaughters,
125  const G4String& theDaughterName1,
126  const G4String& theDaughterName2,
127  const G4String& theDaughterName3,
128  const G4String& theDaughterName4,
129  const G4double *masses) :
130  G4VDecayChannel("Phase Space",
131  theParentName,theBR,
132  theNumberOfDaughters,
133  theDaughterName1,
134  theDaughterName2,
135  theDaughterName3,
136  theDaughterName4),
137  parentmass(theParentMass),
138  theDaughterMasses(masses)
139 {
140  if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay:: constructor " << G4endl;
141 }
142 
144 {
145 }
146 
148 {
149  if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay::DecayIt ";
150  G4DecayProducts * products = NULL;
151 
154 
155  switch (numberOfDaughters){
156  case 0:
157  if (GetVerboseLevel()>0) {
158  G4cout << "G4GeneralPhaseSpaceDecay::DecayIt ";
159  G4cout << " daughters not defined " <<G4endl;
160  }
161  break;
162  case 1:
163  products = OneBodyDecayIt();
164  break;
165  case 2:
166  products = TwoBodyDecayIt();
167  break;
168  case 3:
169  products = ThreeBodyDecayIt();
170  break;
171  default:
172  products = ManyBodyDecayIt();
173  break;
174  }
175  if ((products == NULL) && (GetVerboseLevel()>0)) {
176  G4cout << "G4GeneralPhaseSpaceDecay::DecayIt ";
177  G4cout << *parent_name << " can not decay " << G4endl;
178  DumpInfo();
179  }
180  return products;
181 }
182 
184 {
185  if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay::OneBodyDecayIt()"<<G4endl;
186 
187 // G4double daughtermass = daughters[0]->GetPDGMass();
188 
189  //create parent G4DynamicParticle at rest
190  G4ParticleMomentum dummy;
191  G4DynamicParticle * parentparticle = new G4DynamicParticle(G4MT_parent, dummy, 0.0);
192 
193  //create G4Decayproducts
194  G4DecayProducts *products = new G4DecayProducts(*parentparticle);
195  delete parentparticle;
196 
197  //create daughter G4DynamicParticle at rest
198  G4DynamicParticle * daughterparticle = new G4DynamicParticle(G4MT_daughters[0], dummy, 0.0);
199  products->PushProducts(daughterparticle);
200 
201  if (GetVerboseLevel()>1)
202  {
203  G4cout << "G4GeneralPhaseSpaceDecay::OneBodyDecayIt ";
204  G4cout << " create decay products in rest frame " <<G4endl;
205  products->DumpInfo();
206  }
207  return products;
208 }
209 
211 {
212  if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay::TwoBodyDecayIt()"<<G4endl;
213 
214  //daughters'mass
215  G4double daughtermass[2];
216  G4double daughtermomentum;
217  if ( theDaughterMasses )
218  {
219  daughtermass[0]= *(theDaughterMasses);
220  daughtermass[1] = *(theDaughterMasses+1);
221  } else {
222  daughtermass[0] = G4MT_daughters[0]->GetPDGMass();
223  daughtermass[1] = G4MT_daughters[1]->GetPDGMass();
224  }
225 
226 // G4double sumofdaughtermass = daughtermass[0] + daughtermass[1];
227 
228  //create parent G4DynamicParticle at rest
229  G4ParticleMomentum dummy;
230  G4DynamicParticle * parentparticle = new G4DynamicParticle( G4MT_parent, dummy, 0.0);
231 
232  //create G4Decayproducts @@GF why dummy parentparticle?
233  G4DecayProducts *products = new G4DecayProducts(*parentparticle);
234  delete parentparticle;
235 
236  //calculate daughter momentum
237  daughtermomentum = Pmx(parentmass,daughtermass[0],daughtermass[1]);
238  G4double costheta = 2.*G4UniformRand()-1.0;
239  G4double sintheta = std::sqrt((1.0 - costheta)*(1.0 + costheta));
241  G4ParticleMomentum direction(sintheta*std::cos(phi),sintheta*std::sin(phi),costheta);
242 
243  //create daughter G4DynamicParticle
244  G4double Etotal= std::sqrt(daughtermass[0]*daughtermass[0] + daughtermomentum*daughtermomentum);
245  G4DynamicParticle * daughterparticle = new G4DynamicParticle( G4MT_daughters[0],Etotal, direction*daughtermomentum);
246  products->PushProducts(daughterparticle);
247  Etotal= std::sqrt(daughtermass[1]*daughtermass[1] + daughtermomentum*daughtermomentum);
248  daughterparticle = new G4DynamicParticle( G4MT_daughters[1],Etotal, direction*(-1.0*daughtermomentum));
249  products->PushProducts(daughterparticle);
250 
251  if (GetVerboseLevel()>1)
252  {
253  G4cout << "G4GeneralPhaseSpaceDecay::TwoBodyDecayIt ";
254  G4cout << " create decay products in rest frame " <<G4endl;
255  products->DumpInfo();
256  }
257  return products;
258 }
259 
261 // algorism of this code is originally written in GDECA3 of GEANT3
262 {
263  if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay::ThreeBodyDecayIt()"<<G4endl;
264 
265  //daughters'mass
266  G4double daughtermass[3];
267  G4double sumofdaughtermass = 0.0;
268  for (G4int index=0; index<3; index++)
269  {
270  if ( theDaughterMasses )
271  {
272  daughtermass[index]= *(theDaughterMasses+index);
273  } else {
274  daughtermass[index] = G4MT_daughters[index]->GetPDGMass();
275  }
276  sumofdaughtermass += daughtermass[index];
277  }
278 
279  //create parent G4DynamicParticle at rest
280  G4ParticleMomentum dummy;
281  G4DynamicParticle * parentparticle = new G4DynamicParticle( G4MT_parent, dummy, 0.0);
282 
283  //create G4Decayproducts
284  G4DecayProducts *products = new G4DecayProducts(*parentparticle);
285  delete parentparticle;
286 
287  //calculate daughter momentum
288  // Generate two
289  G4double rd1, rd2, rd;
290  G4double daughtermomentum[3];
291  G4double momentummax=0.0, momentumsum = 0.0;
293  const G4int maxNumberOfLoops = 10000;
294  G4int loopCounter = 0;
295 
296  do
297  {
298  rd1 = G4UniformRand();
299  rd2 = G4UniformRand();
300  if (rd2 > rd1)
301  {
302  rd = rd1;
303  rd1 = rd2;
304  rd2 = rd;
305  }
306  momentummax = 0.0;
307  momentumsum = 0.0;
308  // daughter 0
309 
310  energy = rd2*(parentmass - sumofdaughtermass);
311  daughtermomentum[0] = std::sqrt(energy*energy + 2.0*energy* daughtermass[0]);
312  if ( daughtermomentum[0] >momentummax )momentummax = daughtermomentum[0];
313  momentumsum += daughtermomentum[0];
314 
315  // daughter 1
316  energy = (1.-rd1)*(parentmass - sumofdaughtermass);
317  daughtermomentum[1] = std::sqrt(energy*energy + 2.0*energy* daughtermass[1]);
318  if ( daughtermomentum[1] >momentummax )momentummax = daughtermomentum[1];
319  momentumsum += daughtermomentum[1];
320 
321  // daughter 2
322  energy = (rd1-rd2)*(parentmass - sumofdaughtermass);
323  daughtermomentum[2] = std::sqrt(energy*energy + 2.0*energy* daughtermass[2]);
324  if ( daughtermomentum[2] >momentummax )momentummax = daughtermomentum[2];
325  momentumsum += daughtermomentum[2];
326  } while ( ( momentummax > momentumsum - momentummax ) && /* Loop checking, 02.11.2015, A.Ribon */
327  ++loopCounter < maxNumberOfLoops );
328  if ( loopCounter >= maxNumberOfLoops ) {
330  ed << " Failed sampling after maxNumberOfLoops attempts : forced exit" << G4endl;
331  G4Exception( " G4GeneralPhaseSpaceDecay::ThreeBodyDecayIt ", "HAD_PHASESPACE_001", FatalException, ed );
332  }
333 
334  // output message
335  if (GetVerboseLevel()>1) {
336  G4cout << " daughter 0:" << daughtermomentum[0]/GeV << "[GeV/c]" <<G4endl;
337  G4cout << " daughter 1:" << daughtermomentum[1]/GeV << "[GeV/c]" <<G4endl;
338  G4cout << " daughter 2:" << daughtermomentum[2]/GeV << "[GeV/c]" <<G4endl;
339  G4cout << " momentum sum:" << momentumsum/GeV << "[GeV/c]" <<G4endl;
340  }
341 
342  //create daughter G4DynamicParticle
343  G4double costheta, sintheta, phi, sinphi, cosphi;
344  G4double costhetan, sinthetan, phin, sinphin, cosphin;
345  costheta = 2.*G4UniformRand()-1.0;
346  sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
347  phi = twopi*G4UniformRand()*rad;
348  sinphi = std::sin(phi);
349  cosphi = std::cos(phi);
350  G4ParticleMomentum direction0(sintheta*cosphi,sintheta*sinphi,costheta);
351  G4double Etotal=std::sqrt( daughtermass[0]*daughtermass[0] + daughtermomentum[0]*daughtermomentum[0]);
352  G4DynamicParticle * daughterparticle
353  = new G4DynamicParticle( G4MT_daughters[0], Etotal, direction0*daughtermomentum[0]);
354  products->PushProducts(daughterparticle);
355 
356  costhetan = (daughtermomentum[1]*daughtermomentum[1]-daughtermomentum[2]*daughtermomentum[2]-daughtermomentum[0]*daughtermomentum[0])/(2.0*daughtermomentum[2]*daughtermomentum[0]);
357  sinthetan = std::sqrt((1.0-costhetan)*(1.0+costhetan));
358  phin = twopi*G4UniformRand()*rad;
359  sinphin = std::sin(phin);
360  cosphin = std::cos(phin);
361  G4ParticleMomentum direction2;
362  direction2.setX( sinthetan*cosphin*costheta*cosphi - sinthetan*sinphin*sinphi + costhetan*sintheta*cosphi);
363  direction2.setY( sinthetan*cosphin*costheta*sinphi + sinthetan*sinphin*cosphi + costhetan*sintheta*sinphi);
364  direction2.setZ( -sinthetan*cosphin*sintheta + costhetan*costheta);
365  Etotal=std::sqrt( daughtermass[2]*daughtermass[2] + daughtermomentum[2]*daughtermomentum[2]/direction2.mag2());
366  daughterparticle = new G4DynamicParticle( G4MT_daughters[2],Etotal, direction2*(daughtermomentum[2]/direction2.mag()));
367  products->PushProducts(daughterparticle);
368  G4ThreeVector mom=(direction0*daughtermomentum[0] + direction2*(daughtermomentum[2]/direction2.mag()))*(-1.0);
369  Etotal= std::sqrt( daughtermass[1]*daughtermass[1] + mom.mag2() );
370  daughterparticle =
371  new G4DynamicParticle(G4MT_daughters[1], Etotal, mom);
372  products->PushProducts(daughterparticle);
373 
374  if (GetVerboseLevel()>1) {
375  G4cout << "G4GeneralPhaseSpaceDecay::ThreeBodyDecayIt ";
376  G4cout << " create decay products in rest frame " <<G4endl;
377  products->DumpInfo();
378  }
379  return products;
380 }
381 
383 // algorism of this code is originally written in FORTRAN by M.Asai
384 //*****************************************************************
385 // NBODY
386 // N-body phase space Monte-Carlo generator
387 // Makoto Asai
388 // Hiroshima Institute of Technology
389 // (asai@kekvax.kek.jp)
390 // Revised release : 19/Apr/1995
391 //
392 {
393  //return value
394  G4DecayProducts *products;
395 
396  if (GetVerboseLevel()>1) G4cout << "G4GeneralPhaseSpaceDecay::ManyBodyDecayIt()"<<G4endl;
397 
398  //daughters'mass
399  G4double *daughtermass = new G4double[numberOfDaughters];
400  G4double sumofdaughtermass = 0.0;
401  for (G4int index=0; index<numberOfDaughters; index++){
402  daughtermass[index] = G4MT_daughters[index]->GetPDGMass();
403  sumofdaughtermass += daughtermass[index];
404  }
405 
406  //Calculate daughter momentum
407  G4double *daughtermomentum = new G4double[numberOfDaughters];
408  G4ParticleMomentum direction;
409  G4DynamicParticle **daughterparticle;
411  G4double tmas;
412  G4double weight = 1.0;
413  G4int numberOfTry = 0;
414  G4int index1;
415 
416  do {
417  //Generate rundom number in descending order
418  G4double temp;
420  rd[0] = 1.0;
421  for(index1 =1; index1 < numberOfDaughters -1; index1++)
422  rd[index1] = G4UniformRand();
423  rd[ numberOfDaughters -1] = 0.0;
424  for(index1 =1; index1 < numberOfDaughters -1; index1++) {
425  for(G4int index2 = index1+1; index2 < numberOfDaughters; index2++) {
426  if (rd[index1] < rd[index2]){
427  temp = rd[index1];
428  rd[index1] = rd[index2];
429  rd[index2] = temp;
430  }
431  }
432  }
433 
434  //calcurate virtual mass
435  tmas = parentmass - sumofdaughtermass;
436  temp = sumofdaughtermass;
437  for(index1 =0; index1 < numberOfDaughters; index1++) {
438  sm[index1] = rd[index1]*tmas + temp;
439  temp -= daughtermass[index1];
440  if (GetVerboseLevel()>1) {
441  G4cout << index1 << " rundom number:" << rd[index1];
442  G4cout << " virtual mass:" << sm[index1]/GeV << "[GeV/c/c]" <<G4endl;
443  }
444  }
445  delete [] rd;
446 
447  //Calculate daughter momentum
448  weight = 1.0;
449  index1 =numberOfDaughters-1;
450  daughtermomentum[index1]= Pmx( sm[index1-1],daughtermass[index1-1],sm[index1]);
451  if (GetVerboseLevel()>1) {
452  G4cout << " daughter " << index1 << ":" << *daughters_name[index1];
453  G4cout << " momentum:" << daughtermomentum[index1]/GeV << "[GeV/c]" <<G4endl;
454  }
455  for(index1 =numberOfDaughters-2; index1>=0; index1--) {
456  // calculate
457  daughtermomentum[index1]= Pmx( sm[index1],daughtermass[index1], sm[index1 +1]);
458  if(daughtermomentum[index1] < 0.0) {
459  // !!! illegal momentum !!!
460  if (GetVerboseLevel()>0) {
461  G4cout << "G4GeneralPhaseSpaceDecay::ManyBodyDecayIt ";
462  G4cout << " can not calculate daughter momentum " <<G4endl;
463  G4cout << " parent:" << *parent_name;
464  G4cout << " mass:" << parentmass/GeV << "[GeV/c/c]" <<G4endl;
465  G4cout << " daughter " << index1 << ":" << *daughters_name[index1];
466  G4cout << " mass:" << daughtermass[index1]/GeV << "[GeV/c/c]" ;
467  G4cout << " mass:" << daughtermomentum[index1]/GeV << "[GeV/c]" <<G4endl;
468  }
469  delete [] sm;
470  delete [] daughtermass;
471  delete [] daughtermomentum;
472  return NULL; // Error detection
473 
474  } else {
475  // calculate weight of this events
476  weight *= daughtermomentum[index1]/sm[index1];
477  if (GetVerboseLevel()>1) {
478  G4cout << " daughter " << index1 << ":" << *daughters_name[index1];
479  G4cout << " momentum:" << daughtermomentum[index1]/GeV << "[GeV/c]" <<G4endl;
480  }
481  }
482  }
483  if (GetVerboseLevel()>1) {
484  G4cout << " weight: " << weight <<G4endl;
485  }
486 
487  // exit if number of Try exceeds 100
488  if (numberOfTry++ >100) {
489  if (GetVerboseLevel()>0) {
490  G4cout << "G4GeneralPhaseSpaceDecay::ManyBodyDecayIt: ";
491  G4cout << " can not determine Decay Kinematics " << G4endl;
492  }
493  delete [] sm;
494  delete [] daughtermass;
495  delete [] daughtermomentum;
496  return NULL; // Error detection
497  }
498  } while ( weight > G4UniformRand()); /* Loop checking, 02.11.2015, A.Ribon */
499  if (GetVerboseLevel()>1) {
500  G4cout << "Start calulation of daughters momentum vector "<<G4endl;
501  }
502 
503  G4double costheta, sintheta, phi;
504  G4double beta;
505  daughterparticle = new G4DynamicParticle*[numberOfDaughters];
506 
507  index1 = numberOfDaughters -2;
508  costheta = 2.*G4UniformRand()-1.0;
509  sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
510  phi = twopi*G4UniformRand()*rad;
511  direction.setZ(costheta);
512  direction.setY(sintheta*std::sin(phi));
513  direction.setX(sintheta*std::cos(phi));
514  daughterparticle[index1] = new G4DynamicParticle( G4MT_daughters[index1], direction*daughtermomentum[index1] );
515  daughterparticle[index1+1] = new G4DynamicParticle( G4MT_daughters[index1+1], direction*(-1.0*daughtermomentum[index1]) );
516 
517  for (index1 = numberOfDaughters -3; index1 >= 0; index1--) {
518  //calculate momentum direction
519  costheta = 2.*G4UniformRand()-1.0;
520  sintheta = std::sqrt((1.0-costheta)*(1.0+costheta));
521  phi = twopi*G4UniformRand()*rad;
522  direction.setZ(costheta);
523  direction.setY(sintheta*std::sin(phi));
524  direction.setX(sintheta*std::cos(phi));
525 
526  // boost already created particles
527  beta = daughtermomentum[index1];
528  beta /= std::sqrt( daughtermomentum[index1]*daughtermomentum[index1] + sm[index1+1]*sm[index1+1] );
529  for (G4int index2 = index1+1; index2<numberOfDaughters; index2++) {
530  G4LorentzVector p4;
531  // make G4LorentzVector for secondaries
532  p4 = daughterparticle[index2]->Get4Momentum();
533 
534  // boost secondaries to new frame
535  p4.boost( direction.x()*beta, direction.y()*beta, direction.z()*beta);
536 
537  // change energy/momentum
538  daughterparticle[index2]->Set4Momentum(p4);
539  }
540  //create daughter G4DynamicParticle
541  daughterparticle[index1]= new G4DynamicParticle( G4MT_daughters[index1], direction*(-1.0*daughtermomentum[index1]));
542  }
543 
544  //create G4Decayproducts
545  G4DynamicParticle *parentparticle;
546  direction.setX(1.0); direction.setY(0.0); direction.setZ(0.0);
547  parentparticle = new G4DynamicParticle( G4MT_parent, direction, 0.0);
548  products = new G4DecayProducts(*parentparticle);
549  delete parentparticle;
550  for (index1 = 0; index1<numberOfDaughters; index1++) {
551  products->PushProducts(daughterparticle[index1]);
552  }
553  if (GetVerboseLevel()>1) {
554  G4cout << "G4GeneralPhaseSpaceDecay::ManyBodyDecayIt ";
555  G4cout << " create decay products in rest frame " << G4endl;
556  products->DumpInfo();
557  }
558 
559  delete [] daughterparticle;
560  delete [] daughtermomentum;
561  delete [] daughtermass;
562  delete [] sm;
563 
564  return products;
565 }
566 
567 
568 
569 
570 
G4int PushProducts(G4DynamicParticle *aParticle)
G4String ** daughters_name
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:45
static G4double Pmx(G4double e, G4double p1, G4double p2)
#define G4endl
Definition: G4ios.hh:61
G4String * parent_name
double z() const
Double_t beta
G4double GetPDGMass() const
void setX(double)
double G4double
Definition: G4Types.hh:76
void setZ(double)
G4ParticleDefinition * G4MT_parent
void DumpInfo() const
#define G4UniformRand()
Definition: Randomize.hh:53
static constexpr double twopi
Definition: G4SIunits.hh:76
double energy
Definition: plottest35.C:25
void Set4Momentum(const G4LorentzVector &momentum)
G4LorentzVector Get4Momentum() const
G4ParticleDefinition ** G4MT_daughters
double mag2() const
static constexpr double rad
Definition: G4SIunits.hh:149
G4int GetVerboseLevel() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
double mag() const
double weight
Definition: plottest35.C:25
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
double x() const
double y() const
virtual G4DecayProducts * DecayIt(G4double mass=0.0)
static constexpr double GeV
Definition: G4SIunits.hh:217
void setY(double)
HepLorentzVector & boost(double, double, double)
void CheckAndFillDaughters()