Geant4-11
G4KineticTrack.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// -----------------------------------------------------------------------------
28// GEANT 4 class implementation file
29//
30// History: first implementation, A. Feliciello, 20th May 1998
31// -----------------------------------------------------------------------------
32
33#include "globals.hh"
34#include "G4ios.hh"
35//#include <cmath>
36
37#include "Randomize.hh"
39#include "G4ThreeVector.hh"
40#include "G4LorentzVector.hh"
41#include "G4KineticTrack.hh"
44#include "G4DecayTable.hh"
46#include "G4DecayProducts.hh"
47#include "G4LorentzRotation.hh"
48#include "G4SampleResonance.hh"
49#include "G4Integrator.hh"
50#include "G4KaonZero.hh"
51#include "G4KaonZeroShort.hh"
52#include "G4KaonZeroLong.hh"
53#include "G4AntiKaonZero.hh"
54
55#include "G4HadTmpUtil.hh"
56
57//
58// Some static clobal for integration
59//
60
62
63//
64// Default constructor
65//
66
68 theDefinition(0),
69 theFormationTime(0),
70 thePosition(0),
71 the4Momentum(0),
72 theFermi3Momentum(0),
73 theTotal4Momentum(0),
74 theNucleon(0),
75 nChannels(0),
76 theActualMass(0),
77 theActualWidth(0),
78 theDaughterMass(0),
79 theDaughterWidth(0),
80 theStateToNucleus(undefined),
81 theProjectilePotential(0),
82 theCreatorModel(-1)
83{
85// DEBUG //
87
88/*
89 G4cerr << G4endl << G4endl << G4endl;
90 G4cerr << " G4KineticTrack default constructor invoked! \n";
91 G4cerr << " =========================================== \n" << G4endl;
92*/
93}
94
95
96
97//
98// Copy constructor
99//
100
102{
105 thePosition = right.GetPosition();
110 nChannels = right.GetnChannels();
113 for (G4int i = 0; i < nChannels; i++)
114 {
115 theActualWidth[i] = right.theActualWidth[i];
116 }
117 theDaughterMass = 0;
123// DEBUG //
125
126/*
127 G4cerr << G4endl << G4endl << G4endl;
128 G4cerr << " G4KineticTrack copy constructor invoked! \n";
129 G4cerr << " ======================================== \n" <<G4endl;
130*/
131}
132
133
134//
135// By argument constructor
136//
137
139 G4double aFormationTime,
140 const G4ThreeVector& aPosition,
141 const G4LorentzVector& a4Momentum) :
142 theDefinition(aDefinition),
143 theFormationTime(aFormationTime),
144 thePosition(aPosition),
145 the4Momentum(a4Momentum),
146 theFermi3Momentum(0),
147 theTotal4Momentum(a4Momentum),
148 theNucleon(0),
149 theStateToNucleus(undefined),
150 theProjectilePotential(0),
151 theCreatorModel(-1)
152{
155 {
156 if(G4UniformRand()<0.5)
157 {
159 }
160 else
161 {
163 }
164 }
165
166//
167// Get the number of decay channels
168//
169
170 G4DecayTable* theDecayTable = theDefinition->GetDecayTable();
171 if (theDecayTable != 0)
172 {
173 nChannels = theDecayTable->entries();
174
175 }
176 else
177 {
178 nChannels = 0;
179 }
180
181//
182// Get the actual mass value
183//
184
186
187//
188// Create an array to Store the actual partial widths
189// of the decay channels
190//
191
192 theDaughterMass = 0;
194 theActualWidth = 0;
195 G4bool * theDaughterIsShortLived = 0;
196
198
199 // cout << " ****CONSTR*** ActualMass ******* " << theActualMass << G4endl;
200 G4int index;
201 for (index = nChannels - 1; index >= 0; --index)
202 {
203 G4VDecayChannel* theChannel = theDecayTable->GetDecayChannel(index);
204 G4int nDaughters = theChannel->GetNumberOfDaughters();
205 G4double theMotherWidth;
206 if (nDaughters == 2 || nDaughters == 3)
207 {
208 G4double thePoleMass = theDefinition->GetPDGMass();
209 theMotherWidth = theDefinition->GetPDGWidth();
210 G4double thePoleWidth = theChannel->GetBR()*theMotherWidth;
211 const G4ParticleDefinition* aDaughter;
212 theDaughterMass = new G4double[nDaughters];
213 theDaughterWidth = new G4double[nDaughters];
214 theDaughterIsShortLived = new G4bool[nDaughters];
215 for (G4int n = 0; n < nDaughters; ++n)
216 {
217 aDaughter = theChannel->GetDaughter(n);
218 theDaughterMass[n] = aDaughter->GetPDGMass();
219 theDaughterWidth[n] = aDaughter->GetPDGWidth();
220 theDaughterIsShortLived[n] = aDaughter->IsShortLived();
221 }
222
223//
224// Check whether both the decay products are stable
225//
226
227 G4double theActualMom = 0.0;
228 G4double thePoleMom = 0.0;
229 G4SampleResonance aSampler;
230 if (nDaughters==2)
231 {
232 if ( !theDaughterIsShortLived[0] && !theDaughterIsShortLived[1] )
233 {
234
235 // G4cout << G4endl << "Both the " << nDaughters <<
236 // " decay products are stable!";
237 // cout << " LB: Both decay products STABLE !" << G4endl;
238 // cout << " parent: " << theChannel->GetParentName() << G4endl;
239 // cout << " particle1: " << theChannel->GetDaughterName(0) << G4endl;
240 // cout << " particle2: " << theChannel->GetDaughterName(1) << G4endl;
241
242 theActualMom = EvaluateCMMomentum(theActualMass,
244 thePoleMom = EvaluateCMMomentum(thePoleMass,
246 // cout << G4endl;
247 // cout << " LB: ActualMass/DaughterMass " << theActualMass << " " << theDaughterMass << G4endl;
248 // cout << " LB: ActualMom " << theActualMom << G4endl;
249 // cout << " LB: PoleMom " << thePoleMom << G4endl;
250 // cout << G4endl;
251 }
252 else if ( !theDaughterIsShortLived[0] && theDaughterIsShortLived[1] )
253 {
254
255 // G4cout << G4endl << "Only the first of the " << nDaughters <<" decay products is stable!";
256 // cout << " LB: only the first decay product is STABLE !" << G4endl;
257 // cout << " parent: " << theChannel->GetParentName() << G4endl;
258 // cout << " particle1: " << theChannel->GetDaughterName(0) << G4endl;
259 // cout << " particle2: " << theChannel->GetDaughterName(1) << G4endl;
260
261// global variable definition
262 G4double lowerLimit = aSampler.GetMinimumMass(theChannel->GetDaughter(1));
263 theActualMom = IntegrateCMMomentum(lowerLimit);
264 thePoleMom = IntegrateCMMomentum(lowerLimit, thePoleMass);
265 // cout << " LB Parent Mass = " << G4KineticTrack_Gmass << G4endl;
266 // cout << " LB Actual Mass = " << theActualMass << G4endl;
267 // cout << " LB Daughter1 Mass = " << G4KineticTrack_Gmass1 << G4endl;
268 // cout << " LB Daughter2 Mass = " << G4KineticTrack_Gmass2 << G4endl;
269 // cout << " The Actual Momentum = " << theActualMom << G4endl;
270 // cout << " The Pole Momentum = " << thePoleMom << G4endl;
271 // cout << G4endl;
272
273 }
274 else if ( theDaughterIsShortLived[0] && !theDaughterIsShortLived[1] )
275 {
276
277 // G4cout << G4endl << "Only the second of the " << nDaughters <<
278 // " decay products is stable!";
279 // cout << " LB: only the second decay product is STABLE !" << G4endl;
280 // cout << " parent: " << theChannel->GetParentName() << G4endl;
281 // cout << " particle1: " << theChannel->GetDaughterName(0) << G4endl;
282 // cout << " particle2: " << theChannel->GetDaughterName(1) << G4endl;
283
284//
285// Swap the content of the theDaughterMass and theDaughterWidth arrays!!!
286//
287
290
291// global variable definition
292 G4double lowerLimit = aSampler.GetMinimumMass(theChannel->GetDaughter(0));
293 theActualMom = IntegrateCMMomentum(lowerLimit);
294 thePoleMom = IntegrateCMMomentum(lowerLimit, thePoleMass);
295 // cout << " LB Parent Mass = " << G4KineticTrack_Gmass << G4endl;
296 // cout << " LB Actual Mass = " << theActualMass << G4endl;
297 // cout << " LB Daughter1 Mass = " << G4KineticTrack_Gmass1 << G4endl;
298 // cout << " LB Daughter2 Mass = " << G4KineticTrack_Gmass2 << G4endl;
299 // cout << " The Actual Momentum = " << theActualMom << G4endl;
300 // cout << " The Pole Momentum = " << thePoleMom << G4endl;
301 // cout << G4endl;
302
303 }
304 else if ( theDaughterIsShortLived[0] && theDaughterIsShortLived[1] )
305 {
306
307// G4cout << G4endl << "Both the " << nDaughters <<
308// " decay products are resonances!";
309 // cout << " LB: both decay products are RESONANCES !" << G4endl;
310 // cout << " parent: " << theChannel->GetParentName() << G4endl;
311 // cout << " particle1: " << theChannel->GetDaughterName(0) << G4endl;
312 // cout << " particle2: " << theChannel->GetDaughterName(1) << G4endl;
313
314// global variable definition
316 theActualMom = IntegrateCMMomentum2();
317 G4KineticTrack_Gmass = thePoleMass;
318 thePoleMom = IntegrateCMMomentum2();
319 // cout << " LB Parent Mass = " << G4KineticTrack_Gmass << G4endl;
320 // cout << " LB Daughter1 Mass = " << G4KineticTrack_Gmass1 << G4endl;
321 // cout << " LB Daughter2 Mass = " << G4KineticTrack_Gmass2 << G4endl;
322 // cout << " The Actual Momentum = " << theActualMom << G4endl;
323 // cout << " The Pole Momentum = " << thePoleMom << G4endl;
324 // cout << G4endl;
325
326 }
327 }
328 else // (nDaughter==3)
329 {
330
331 G4int nShortLived = 0;
332 if ( theDaughterIsShortLived[0] )
333 {
334 ++nShortLived;
335 }
336 if ( theDaughterIsShortLived[1] )
337 {
338 ++nShortLived;
341 }
342 if ( theDaughterIsShortLived[2] )
343 {
344 ++nShortLived;
347 }
348 if ( nShortLived == 0 )
349 {
351 theActualMom = EvaluateCMMomentum(theActualMass,
353 thePoleMom = EvaluateCMMomentum(thePoleMass,
355 }
356// else if ( nShortLived == 1 )
357 else if ( nShortLived >= 1 )
358 {
359 // need the shortlived particle in slot 1! (very bad style...)
363 theActualMom = IntegrateCMMomentum(0.0);
364 thePoleMom = IntegrateCMMomentum(0.0, thePoleMass);
365 }
366// else
367// {
368// throw G4HadronicException(__FILE__, __LINE__, ("can't handle more than one shortlived in 3 particle output channel");
369// }
370
371 }
372
373 //if(nDaughters<3) theChannel->GetAngularMomentum();
374 G4double theMassRatio = thePoleMass / theActualMass;
375 G4double theMomRatio = theActualMom / thePoleMom;
376 // VI 11.06.2015: for l=0 one not need use pow
377 //G4double l=0;
378 //theActualWidth[index] = thePoleWidth * theMassRatio *
379 // std::pow(theMomRatio, (2 * l + 1)) *
380 // (1.2 / (1+ 0.2*std::pow(theMomRatio, (2 * l))));
381 theActualWidth[index] = thePoleWidth * theMassRatio *
382 theMomRatio;
383 delete [] theDaughterMass;
384 theDaughterMass = 0;
385 delete [] theDaughterWidth;
387 delete [] theDaughterIsShortLived;
388 theDaughterIsShortLived = 0;
389 }
390
391 else // nDaughter = 1 ( e.g. K0 decays 50% to Kshort, 50% Klong
392 {
393 theMotherWidth = theDefinition->GetPDGWidth();
394 theActualWidth[index] = theChannel->GetBR()*theMotherWidth;
395 }
396 }
397
399// DEBUG //
401
402// for (G4int y = nChannels - 1; y >= 0; --y)
403// {
404// G4cout << G4endl << theActualWidth[y];
405// }
406// G4cout << G4endl << G4endl << G4endl;
407
408 /*
409 G4cerr << G4endl << G4endl << G4endl;
410 G4cerr << " G4KineticTrack by argument constructor invoked! \n";
411 G4cerr << " =============================================== \n" << G4endl;
412 */
413
414}
415
417 const G4ThreeVector& aPosition,
418 const G4LorentzVector& a4Momentum)
419 : theDefinition(nucleon->GetDefinition()),
420 theFormationTime(0),
421 thePosition(aPosition),
422 the4Momentum(a4Momentum),
423 theFermi3Momentum(nucleon->GetMomentum()),
424 theNucleon(nucleon),
425 nChannels(0),
426 theActualMass(nucleon->GetDefinition()->GetPDGMass()),
427 theActualWidth(0),
428 theDaughterMass(0),
429 theDaughterWidth(0),
430 theStateToNucleus(undefined),
431 theProjectilePotential(0),
432 theCreatorModel(-1)
433{
435 Set4Momentum(a4Momentum);
436}
437
438
440{
441 if (theActualWidth != 0) delete [] theActualWidth;
442 if (theDaughterMass != 0) delete [] theDaughterMass;
443 if (theDaughterWidth != 0) delete [] theDaughterWidth;
444}
445
446
447
449{
450 if (this != &right)
451 {
454 the4Momentum = right.the4Momentum;
460 if (theActualWidth != 0) delete [] theActualWidth;
461 nChannels = right.GetnChannels();
463 for (G4int i = 0; i < nChannels; ++i) theActualWidth[i] = right.theActualWidth[i];
465 }
466 return *this;
467}
468
469
470
472{
473 return (this == & right);
474}
475
476
477
479{
480 return (this != & right);
481}
482
483
484
486{
487//
488// Select a possible decay channel
489//
490/*
491 G4int index1;
492 for (index1 = nChannels - 1; index1 >= 0; --index1)
493 G4cout << "DECAY Actual Width IND/ActualW " << index1 << " " << theActualWidth[index1] << G4endl;
494 G4cout << "DECAY Actual Mass " << theActualMass << G4endl;
495*/
496 const G4ParticleDefinition* thisDefinition = this->GetDefinition();
497 if(!thisDefinition)
498 {
499 G4cerr << "Error condition encountered in G4KineticTrack::Decay()"<<G4endl;
500 G4cerr << " track has no particle definition associated."<<G4endl;
501 return 0;
502 }
503 G4DecayTable* theDecayTable = thisDefinition->GetDecayTable();
504 if(!theDecayTable)
505 {
506 G4cerr << "Error condition encountered in G4KineticTrack::Decay()"<<G4endl;
507 G4cerr << " particle definition has no decay table associated."<<G4endl;
508 G4cerr << " particle was "<<thisDefinition->GetParticleName()<<G4endl;
509 return 0;
510 }
511
512 G4int chargeBalance = G4lrint(theDefinition->GetPDGCharge() );
513 G4int baryonBalance = G4lrint(theDefinition->GetBaryonNumber() );
514 G4LorentzVector energyMomentumBalance(Get4Momentum());
515 G4double theTotalActualWidth = this->EvaluateTotalActualWidth();
516 if (theTotalActualWidth !=0)
517 {
518
519 //AR-16Aug2016 : Repeat the sampling of the decay channel until is
520 // kinematically above threshold or a max number of attempts is reached
521 G4bool isChannelBelowThreshold = true;
522 const G4int maxNumberOfLoops = 10000;
523 G4int loopCounter = 0;
524
525 G4int chosench;
526 G4String theParentName;
527 G4double theParentMass;
528 G4double theBR;
529 G4int theNumberOfDaughters;
530 G4String theDaughtersName1;
531 G4String theDaughtersName2;
532 G4String theDaughtersName3;
533 G4String theDaughtersName4;
534 G4double masses[4]={0.,0.,0.,0.};
535
536 do {
537
538 G4double theSumActualWidth = 0.0;
539 G4double* theCumActualWidth = new G4double[nChannels]{};
540 for (G4int index = nChannels - 1; index >= 0; --index)
541 {
542 theSumActualWidth += theActualWidth[index];
543 theCumActualWidth[index] = theSumActualWidth;
544 // cout << "DECAY Cum. Width " << index << " " << theCumActualWidth[index] << G4endl;
545 }
546 // cout << "DECAY Total Width " << theSumActualWidth << G4endl;
547 // cout << "DECAY Total Width " << theTotalActualWidth << G4endl;
548 G4double r = theTotalActualWidth * G4UniformRand();
549 G4VDecayChannel* theDecayChannel(0);
550 chosench=-1;
551 for (G4int index = nChannels - 1; index >= 0; --index)
552 {
553 if (r < theCumActualWidth[index])
554 {
555 theDecayChannel = theDecayTable->GetDecayChannel(index);
556 // cout << "DECAY SELECTED CHANNEL" << index << G4endl;
557 chosench=index;
558 break;
559 }
560 }
561
562 delete [] theCumActualWidth;
563
564 if(!theDecayChannel)
565 {
566 G4cerr << "Error condition encountered in G4KineticTrack::Decay()"<<G4endl;
567 G4cerr << " decay channel has 0x0 channel associated."<<G4endl;
568 G4cerr << " particle was "<<thisDefinition->GetParticleName()<<G4endl;
569 G4cerr << " channel index "<< chosench << "of "<<nChannels<<"channels"<<G4endl;
570 return 0;
571 }
572 theParentName = theDecayChannel->GetParentName();
573 theParentMass = this->GetActualMass();
574 theBR = theActualWidth[chosench];
575 // cout << "**BR*** DECAYNEW " << theBR << G4endl;
576 theNumberOfDaughters = theDecayChannel->GetNumberOfDaughters();
577 theDaughtersName1 = "";
578 theDaughtersName2 = "";
579 theDaughtersName3 = "";
580 theDaughtersName4 = "";
581
582 for (G4int i=0; i < 4; ++i) masses[i]=0.;
583 G4int shortlivedDaughters[4];
584 G4int numberOfShortliveds(0);
585 G4double SumLongLivedMass(0);
586 for (G4int aD=0; aD < theNumberOfDaughters ; ++aD)
587 {
588 const G4ParticleDefinition* aDaughter = theDecayChannel->GetDaughter(aD);
589 masses[aD] = aDaughter->GetPDGMass();
590 if ( aDaughter->IsShortLived() )
591 {
592 shortlivedDaughters[numberOfShortliveds]=aD;
593 ++numberOfShortliveds;
594 } else {
595 SumLongLivedMass += aDaughter->GetPDGMass();
596 }
597
598 }
599 switch (theNumberOfDaughters)
600 {
601 case 0:
602 break;
603 case 1:
604 theDaughtersName1 = theDecayChannel->GetDaughterName(0);
605 theDaughtersName2 = "";
606 theDaughtersName3 = "";
607 theDaughtersName4 = "";
608 break;
609 case 2:
610 theDaughtersName1 = theDecayChannel->GetDaughterName(0);
611 theDaughtersName2 = theDecayChannel->GetDaughterName(1);
612 theDaughtersName3 = "";
613 theDaughtersName4 = "";
614 if ( numberOfShortliveds == 1)
615 { G4SampleResonance aSampler;
616 G4double massmax=theParentMass - SumLongLivedMass;
617 const G4ParticleDefinition * aDaughter=theDecayChannel->GetDaughter(shortlivedDaughters[0]);
618 masses[shortlivedDaughters[0]]= aSampler.SampleMass(aDaughter,massmax);
619 } else if ( numberOfShortliveds == 2) {
620 // choose masses one after the other, start with randomly choosen
621 G4int zero= (G4UniformRand() > 0.5) ? 0 : 1;
622 G4int one = 1-zero;
623 G4SampleResonance aSampler;
624 G4double massmax=theParentMass - aSampler.GetMinimumMass(theDecayChannel->GetDaughter(shortlivedDaughters[one]));
625 const G4ParticleDefinition * aDaughter=theDecayChannel->GetDaughter(shortlivedDaughters[zero]);
626 masses[shortlivedDaughters[zero]]=aSampler.SampleMass(aDaughter,massmax);
627 massmax=theParentMass - masses[shortlivedDaughters[zero]];
628 aDaughter=theDecayChannel->GetDaughter(shortlivedDaughters[one]);
629 masses[shortlivedDaughters[one]]=aSampler.SampleMass(aDaughter,massmax);
630 }
631 break;
632 case 3:
633 theDaughtersName1 = theDecayChannel->GetDaughterName(0);
634 theDaughtersName2 = theDecayChannel->GetDaughterName(1);
635 theDaughtersName3 = theDecayChannel->GetDaughterName(2);
636 theDaughtersName4 = "";
637 if ( numberOfShortliveds == 1)
638 { G4SampleResonance aSampler;
639 G4double massmax=theParentMass - SumLongLivedMass;
640 const G4ParticleDefinition * aDaughter=theDecayChannel->GetDaughter(shortlivedDaughters[0]);
641 masses[shortlivedDaughters[0]]= aSampler.SampleMass(aDaughter,massmax);
642 }
643 break;
644 default:
645 theDaughtersName1 = theDecayChannel->GetDaughterName(0);
646 theDaughtersName2 = theDecayChannel->GetDaughterName(1);
647 theDaughtersName3 = theDecayChannel->GetDaughterName(2);
648 theDaughtersName4 = theDecayChannel->GetDaughterName(3);
649 if ( numberOfShortliveds == 1)
650 { G4SampleResonance aSampler;
651 G4double massmax=theParentMass - SumLongLivedMass;
652 const G4ParticleDefinition * aDaughter=theDecayChannel->GetDaughter(shortlivedDaughters[0]);
653 masses[shortlivedDaughters[0]]= aSampler.SampleMass(aDaughter,massmax);
654 }
655 if ( theNumberOfDaughters > 4 ) {
657 ed << "More than 4 decay daughters: kept only the first 4" << G4endl;
658 G4Exception( "G4KineticTrack::Decay()", "KINTRK5", JustWarning, ed );
659 }
660 break;
661 }
662
663 //AR-16Aug2016 : Check whether the sum of the masses of the daughters is smaller than the parent mass.
664 // If this is still not the case, but the max number of attempts has been reached,
665 // then the subsequent call thePhaseSpaceDecayChannel.DecayIt() will throw an exception.
666 G4double sumDaughterMasses = 0.0;
667 for (G4int i=0; i < 4; ++i) sumDaughterMasses += masses[i];
668 if ( theParentMass - sumDaughterMasses > 0.0 ) isChannelBelowThreshold = false;
669
670 } while ( isChannelBelowThreshold && ++loopCounter < maxNumberOfLoops ); /* Loop checking, 16.08.2016, A.Ribon */
671
672//
673// Get the decay products List
674//
675
676 G4GeneralPhaseSpaceDecay thePhaseSpaceDecayChannel(theParentName,
677 theParentMass,
678 theBR,
679 theNumberOfDaughters,
680 theDaughtersName1,
681 theDaughtersName2,
682 theDaughtersName3,
683 theDaughtersName4,
684 masses);
685 G4DecayProducts* theDecayProducts = thePhaseSpaceDecayChannel.DecayIt();
686 if(!theDecayProducts)
687 {
689 ed << "Error condition encountered: phase-space decay failed." << G4endl
690 << "\t the decaying particle is: " << thisDefinition->GetParticleName() << G4endl
691 << "\t the channel index is: "<< chosench << " of "<< nChannels << "channels" << G4endl
692 << "\t " << theNumberOfDaughters << " daughter particles: "
693 << theDaughtersName1 << " " << theDaughtersName2 << " " << theDaughtersName3 << " "
694 << theDaughtersName4 << G4endl;
695 G4Exception( "G4KineticTrack::Decay ", "HAD_KINTRACK_001", JustWarning, ed );
696 return 0;
697 }
698
699//
700// Create the kinetic track List associated to the decay products
701//
702// For the decay products of hadronic resonances, we assign as creator model ID
703// the same as their parent
704 G4LorentzRotation toMoving(Get4Momentum().boostVector());
705 G4DynamicParticle* theDynamicParticle;
706 G4double formationTime = 0.0;
708 G4LorentzVector momentum;
709 G4LorentzVector momentumBalanceCMS(0);
710 G4KineticTrackVector* theDecayProductList = new G4KineticTrackVector;
711 G4int dEntries = theDecayProducts->entries();
712 const G4ParticleDefinition * aProduct = 0;
713 for (G4int i=dEntries; i > 0; --i)
714 {
715 theDynamicParticle = theDecayProducts->PopProducts();
716 aProduct = theDynamicParticle->GetDefinition();
717 chargeBalance -= G4lrint(aProduct->GetPDGCharge() );
718 baryonBalance -= G4lrint(aProduct->GetBaryonNumber() );
719 momentumBalanceCMS += theDynamicParticle->Get4Momentum();
720 momentum = toMoving*theDynamicParticle->Get4Momentum();
721 energyMomentumBalance -= momentum;
722 G4KineticTrack* aDaughter = new G4KineticTrack (aProduct,
723 formationTime,
724 position,
725 momentum);
726 if (aDaughter != nullptr) aDaughter->SetCreatorModelID(GetCreatorModelID());
727 theDecayProductList->push_back(aDaughter);
728 delete theDynamicParticle;
729 }
730 delete theDecayProducts;
731 if(std::getenv("DecayEnergyBalanceCheck"))
732 std::cout << "DEBUGGING energy balance in cms and lab, charge baryon balance : "
733 << momentumBalanceCMS << " "
734 <<energyMomentumBalance << " "
735 <<chargeBalance<<" "
736 <<baryonBalance<<" "
737 <<G4endl;
738 return theDecayProductList;
739 }
740 else
741 {
742 return 0;
743 }
744}
745
747{
748 G4double mass = theActualMass; /* the actual mass value */
749 G4double mass1 = theDaughterMass[0];
750 G4double mass2 = theDaughterMass[1];
751 G4double gamma2 = theDaughterWidth[1];
752
753 G4double result = (1. / (2 * mass)) *
754 std::sqrt(std::max((((mass * mass) - (mass1 + xmass) * (mass1 + xmass)) *
755 ((mass * mass) - (mass1 - xmass) * (mass1 - xmass))),0.0)) *
756 BrWig(gamma2, mass2, xmass);
757 return result;
758}
759
761{
762 G4double mass = theDefinition->GetPDGMass(); /* the pole mass value */
763 G4double mass1 = theDaughterMass[0];
764 G4double mass2 = theDaughterMass[1];
765 G4double gamma2 = theDaughterWidth[1];
766 G4double result = (1. / (2 * mass)) *
767 std::sqrt(std::max((((mass * mass) - (mass1 + xmass) * (mass1 + xmass)) *
768 ((mass * mass) - (mass1 - xmass) * (mass1 - xmass))),0.0)) *
769 BrWig(gamma2, mass2, xmass);
770 return result;
771}
772
774{
775 const G4double mass = G4KineticTrack_Gmass; /* the actual mass value */
776// const G4double mass1 = theDaughterMass[0];
777 const G4double mass2 = theDaughterMass[1];
778 const G4double gamma2 = theDaughterWidth[1];
779
780 const G4double result = (1. / (2 * mass)) *
781 std::sqrt(((mass * mass) - (G4KineticTrack_xmass1 + xmass) * (G4KineticTrack_xmass1 + xmass)) *
782 ((mass * mass) - (G4KineticTrack_xmass1 - xmass) * (G4KineticTrack_xmass1 - xmass))) *
783 BrWig(gamma2, mass2, xmass);
784 return result;
785}
786
788{
789 const G4double mass = G4KineticTrack_Gmass;
790 const G4double mass1 = theDaughterMass[0];
791 const G4double gamma1 = theDaughterWidth[0];
792// const G4double mass2 = theDaughterMass[1];
793
794 G4KineticTrack_xmass1 = xmass;
795
796 const G4double theLowerLimit = 0.0;
797 const G4double theUpperLimit = mass - xmass;
798 const G4int nIterations = 100;
799
800 G4Integrator<const G4KineticTrack, G4double(G4KineticTrack::*)(G4double) const> integral;
801 G4double result = BrWig(gamma1, mass1, xmass)*
802 integral.Simpson(this, &G4KineticTrack::IntegrandFunction3, theLowerLimit, theUpperLimit, nIterations);
803
804 return result;
805}
806
808{
809 const G4double theUpperLimit = theActualMass - theDaughterMass[0];
810 const G4int nIterations = 100;
811
812 if (theLowerLimit>=theUpperLimit) return 0.0;
813
814 G4Integrator<const G4KineticTrack, G4double(G4KineticTrack::*)(G4double) const> integral;
815 G4double theIntegralOverMass2 = integral.Simpson(this, &G4KineticTrack::IntegrandFunction1,
816 theLowerLimit, theUpperLimit, nIterations);
817 return theIntegralOverMass2;
818}
819
820G4double G4KineticTrack::IntegrateCMMomentum(const G4double theLowerLimit, const G4double poleMass) const
821{
822 const G4double theUpperLimit = poleMass - theDaughterMass[0];
823 const G4int nIterations = 100;
824
825 if (theLowerLimit>=theUpperLimit) return 0.0;
826
827 G4Integrator<const G4KineticTrack, G4double(G4KineticTrack::*)(G4double) const> integral;
828 const G4double theIntegralOverMass2 = integral.Simpson(this, &G4KineticTrack::IntegrandFunction2,
829 theLowerLimit, theUpperLimit, nIterations);
830 return theIntegralOverMass2;
831}
832
833
835{
836 const G4double theLowerLimit = 0.0;
837 const G4double theUpperLimit = theActualMass;
838 const G4int nIterations = 100;
839
840 if (theLowerLimit>=theUpperLimit) return 0.0;
841
842 G4Integrator<const G4KineticTrack, G4double(G4KineticTrack::*)(G4double) const> integral;
843 G4double theIntegralOverMass2 = integral.Simpson(this, &G4KineticTrack::IntegrandFunction4,
844 theLowerLimit, theUpperLimit, nIterations);
845 return theIntegralOverMass2;
846}
847
848
849
850
851
852
853
854
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
static G4ThreadLocal G4double G4KineticTrack_xmass1
static G4ThreadLocal G4double G4KineticTrack_Gmass
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
#define G4UniformRand()
Definition: Randomize.hh:52
static G4AntiKaonZero * AntiKaonZero()
G4int entries() const
G4DynamicParticle * PopProducts()
G4VDecayChannel * GetDecayChannel(G4int index) const
G4int entries() const
G4ParticleDefinition * GetDefinition() const
G4LorentzVector Get4Momentum() const
virtual G4DecayProducts * DecayIt(G4double mass=0.0)
static G4KaonZeroLong * KaonZeroLong()
static G4KaonZeroShort * KaonZeroShort()
static G4KaonZero * KaonZero()
Definition: G4KaonZero.cc:103
G4double GetFormationTime() const
G4double * theDaughterMass
G4double IntegrandFunction1(G4double xmass) const
G4double IntegrateCMMomentum2() const
G4KineticTrackVector * Decay()
G4KineticTrack & operator=(const G4KineticTrack &right)
void SetCreatorModelID(G4int id)
G4int GetCreatorModelID() const
G4LorentzVector the4Momentum
G4bool operator!=(const G4KineticTrack &right) const
G4ThreeVector thePosition
void Set4Momentum(const G4LorentzVector &a4Momentum)
G4LorentzVector theTotal4Momentum
G4double * theActualWidth
G4int GetnChannels() const
const G4ThreeVector & GetPosition() const
G4Nucleon * theNucleon
const G4ParticleDefinition * GetDefinition() const
G4bool operator==(const G4KineticTrack &right) const
G4double theFormationTime
const G4LorentzVector & GetTrackingMomentum() const
G4double * theDaughterWidth
G4double theProjectilePotential
G4double IntegrateCMMomentum(const G4double lowerLimit) const
G4double IntegrandFunction3(G4double xmass) const
G4double IntegrandFunction4(G4double xmass) const
G4double IntegrandFunction2(G4double xmass) const
G4double EvaluateCMMomentum(const G4double mass, const G4double *m_ij) const
G4double BrWig(const G4double Gamma, const G4double rmass, const G4double mass) const
G4double theActualMass
CascadeState theStateToNucleus
G4LorentzVector theFermi3Momentum
const G4LorentzVector & Get4Momentum() const
const G4ParticleDefinition * theDefinition
G4double GetActualMass() const
G4double EvaluateTotalActualWidth()
G4double GetPDGWidth() const
G4double GetPDGCharge() const
G4DecayTable * GetDecayTable() const
const G4String & GetParticleName() const
G4double SampleMass(const G4double poleMass, const G4double gamma, const G4double minMass, const G4double maxMass) const
G4double GetMinimumMass(const G4ParticleDefinition *p) const
G4double GetBR() const
const G4String & GetParentName() const
G4int GetNumberOfDaughters() const
G4ParticleDefinition * GetDaughter(G4int anIndex)
const G4String & GetDaughterName(G4int anIndex) const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4bool nucleon(G4int ityp)
static const G4LorentzVector zero(0., 0., 0., 0.)
void G4SwapObj(T *a, T *b)
Definition: templates.hh:112
int G4lrint(double ad)
Definition: templates.hh:134
#define G4ThreadLocal
Definition: tls.hh:77