Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DalitzDecayChannel.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// G4DalitzDecayChannel class implementation
27//
28// Author: H.Kurashige, 30 May 1997
29// --------------------------------------------------------------------
30
32#include "G4SystemOfUnits.hh"
34#include "G4DecayProducts.hh"
35#include "G4VDecayChannel.hh"
38#include "Randomize.hh"
39#include "G4LorentzVector.hh"
40#include "G4LorentzRotation.hh"
41
44{
45}
46
47G4DalitzDecayChannel::
48G4DalitzDecayChannel(const G4String& theParentName,
49 G4double theBR,
50 const G4String& theLeptonName,
51 const G4String& theAntiLeptonName)
52 : G4VDecayChannel("Dalitz Decay", 1)
53{
54 // set names for daughter particles
55 SetParent(theParentName);
56 SetBR(theBR);
58 G4String gammaName = "gamma";
59 SetDaughter(idGamma, gammaName);
60 SetDaughter(idLepton, theLeptonName);
61 SetDaughter(idAntiLepton, theAntiLeptonName);
62}
63
65{
66}
67
69 : G4VDecayChannel(right)
70{
71}
72
75{
76 if (this != &right)
77 {
80 rbranch = right.rbranch;
81
82 // copy parent name
83 parent_name = new G4String(*right.parent_name);
84
85 // clear daughters_name array
87
88 // recreate array
90 if ( numberOfDaughters >0 )
91 {
92 if (daughters_name != nullptr) ClearDaughtersName();
94 //copy daughters name
95 for (G4int index=0; index < numberOfDaughters; ++index)
96 {
97 daughters_name[index] = new G4String(*right.daughters_name[index]);
98 }
99 }
100 }
101 return *this;
102}
103
105{
106#ifdef G4VERBOSE
107 if (GetVerboseLevel()>1) G4cout << "G4DalitzDecayChannel::DecayIt ";
108#endif
111
112 // parent mass
113 G4double parentmass = G4MT_parent->GetPDGMass();
114
115 // create parent G4DynamicParticle at rest
116 G4ThreeVector dummy;
117 G4DynamicParticle* parentparticle
118 = new G4DynamicParticle( G4MT_parent, dummy, 0.0);
119
120 //daughters' mass
121 G4double leptonmass = G4MT_daughters[idLepton]->GetPDGMass();
122
123 // Generate t ( = std::exp(x):mass Square of (l+ l-) system)
124 G4double xmin = 2.0*std::log(2.0*leptonmass);
125 G4double xmax = 2.0*std::log(parentmass);
126 G4double wmax = 1.5;
127 G4double x, w, ww, w1, w2, w3, t;
128 const std::size_t MAX_LOOP = 10000;
129 for (std::size_t loop_counter=0; loop_counter<MAX_LOOP; ++loop_counter)
130 {
131 x = G4UniformRand()*(xmax-xmin) + xmin;
132 w = G4UniformRand()*wmax;
133 t = std::exp(x);
134 w1 = (1.0-4.0*leptonmass*leptonmass/t);
135 if ( w1 > 0.0)
136 {
137 w2 = ( 1.0 + 2.0*leptonmass*leptonmass/t );
138 w3 = ( 1.0 - t/parentmass/parentmass );
139 w3 = w3 * w3 * w3;
140 ww = w3 * w2 * std::sqrt(w1);
141 }
142 else
143 {
144 ww = 0.0;
145 }
146 if (w <= ww) break;
147 }
148
149 // calculate gamma momentum
150 G4double Pgamma =
151 G4PhaseSpaceDecayChannel::Pmx(parentmass, 0.0, std::sqrt(t));
152 G4double costheta = 2.*G4UniformRand()-1.0;
153 G4double sintheta = std::sqrt((1.0 - costheta)*(1.0 + costheta));
154 G4double phi = twopi*G4UniformRand()*rad;
155 G4ThreeVector gdirection(sintheta*std::cos(phi),
156 sintheta*std::sin(phi),
157 costheta);
158
159 // create G4DynamicParticle for gamma
160 G4DynamicParticle* gammaparticle
161 = new G4DynamicParticle(G4MT_daughters[idGamma] , gdirection, Pgamma);
162
163 // calculate beta of (l+ l-)system
164 G4double beta = Pgamma/(parentmass-Pgamma);
165
166 // calculate lepton momentum in the rest frame of (l+ l-)system
167 G4double Plepton =
168 G4PhaseSpaceDecayChannel::Pmx(std::sqrt(t),leptonmass, leptonmass);
169 G4double Elepton = std::sqrt(Plepton*Plepton + leptonmass*leptonmass );
170 costheta = 2.*G4UniformRand()-1.0;
171 sintheta = std::sqrt((1.0 - costheta)*(1.0 + costheta));
172 phi = twopi*G4UniformRand()*rad;
173 G4ThreeVector ldirection(sintheta*std::cos(phi),
174 sintheta*std::sin(phi),
175 costheta);
176 // create G4DynamicParticle for leptons in the rest frame of (l+ l-)system
177 G4DynamicParticle* leptonparticle
178 = new G4DynamicParticle(G4MT_daughters[idLepton] ,
179 ldirection, Elepton-leptonmass );
180 G4DynamicParticle* antileptonparticle
181 = new G4DynamicParticle(G4MT_daughters[idAntiLepton] ,
182 -1.0*ldirection, Elepton-leptonmass );
183 //boost leptons in the rest frame of the parent
184 G4LorentzVector p4 = leptonparticle->Get4Momentum();
185 p4.boost( -1.0*gdirection.x()*beta,
186 -1.0*gdirection.y()*beta,
187 -1.0*gdirection.z()*beta);
188 leptonparticle->Set4Momentum(p4);
189 p4 = antileptonparticle->Get4Momentum();
190 p4.boost( -1.0*gdirection.x()*beta,
191 -1.0*gdirection.y()*beta,
192 -1.0*gdirection.z()*beta);
193 antileptonparticle->Set4Momentum(p4);
194
195 //create G4Decayproducts
196 G4DecayProducts* products = new G4DecayProducts(*parentparticle);
197 delete parentparticle;
198 products->PushProducts(gammaparticle);
199 products->PushProducts(leptonparticle);
200 products->PushProducts(antileptonparticle);
201
202#ifdef G4VERBOSE
203 if (GetVerboseLevel()>1)
204 {
205 G4cout << "G4DalitzDecayChannel::DecayIt ";
206 G4cout << " create decay products in rest frame " <<G4endl;
207 products->DumpInfo();
208 }
209#endif
210 return products;
211}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
double x() const
double y() const
HepLorentzVector & boost(double, double, double)
G4DalitzDecayChannel(const G4String &theParentName, G4double theBR, const G4String &theLeptonName, const G4String &theAntiLeptonName)
G4DalitzDecayChannel & operator=(const G4DalitzDecayChannel &)
virtual G4DecayProducts * DecayIt(G4double)
void DumpInfo() const
G4int PushProducts(G4DynamicParticle *aParticle)
G4LorentzVector Get4Momentum() const
void Set4Momentum(const G4LorentzVector &momentum)
static G4double Pmx(G4double e, G4double p1, G4double p2)
G4ParticleDefinition ** G4MT_daughters
G4String * parent_name
void SetBR(G4double value)
G4String ** daughters_name
G4int GetVerboseLevel() const
void SetNumberOfDaughters(G4int value)
G4ParticleDefinition * G4MT_parent
void CheckAndFillDaughters()
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
G4String kinematics_name
void SetParent(const G4ParticleDefinition *particle_type)