Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4NeutrinoElectronCcModel.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 // $Id: G4NeutrinoElectronCcModel.cc 91806 2015-08-06 12:20:45Z gcosmo $
27 //
28 // Geant4 Header : G4NeutrinoElectronCcModel
29 //
30 // Author : V.Grichine 26.4.17
31 //
32 
34 #include "G4SystemOfUnits.hh"
35 #include "G4ParticleTable.hh"
36 #include "G4ParticleDefinition.hh"
37 #include "G4IonTable.hh"
38 #include "Randomize.hh"
39 #include "G4NeutrinoE.hh"
40 #include "G4AntiNeutrinoE.hh"
41 #include "G4MuonMinus.hh"
42 #include "G4TauMinus.hh"
43 
44 using namespace std;
45 using namespace CLHEP;
46 
48  : G4HadronicInteraction(name)
49 {
50  SetMinEnergy( 0.0*GeV );
51  SetMaxEnergy( 100.*TeV );
52  SetMinEnergy(1.e-6*eV);
53 
58 
59  // PDG2016: sin^2 theta Weinberg
60 
61  fSin2tW = 0.23129; // 0.2312;
62 
63  fCutEnergy = 0.; // default value
64 
65 }
66 
67 
69 {}
70 
71 
72 void G4NeutrinoElectronCcModel::ModelDescription(std::ostream& outFile) const
73 {
74 
75  outFile << "G4NeutrinoElectronCcModel is a neutrino-electron (neutral current) elastic scattering\n"
76  << "model which uses the standard model \n"
77  << "transfer parameterization. The model is fully relativistic\n";
78 
79 }
80 
82 
84  G4Nucleus & targetNucleus)
85 {
86  G4bool result = false;
87  G4String pName = aPart.GetDefinition()->GetParticleName();
88  G4double minEnergy = 0., energy = aPart.GetTotalEnergy();
89  G4double fmass, emass = electron_mass_c2;
90 
91  if( pName == "nu_mu" || pName == "anti_nu_mu" ) fmass = theMuonMinus->GetPDGMass();
92  else if( pName == "nu_tau" || pName == "anti_nu_tau" ) fmass = theTauMinus->GetPDGMass();
93  else fmass = emass;
94 
95  minEnergy = (fmass-emass)*(fmass+emass)/emass;
96  SetMinEnergy( minEnergy );
97 
98  if( ( pName == "nu_mu" || pName == "anti_nu_mu" ||
99  pName == "nu_tau" || pName == "anti_nu_tau" ) &&
100  energy > minEnergy )
101  {
102  result = true;
103  }
104  G4int Z = targetNucleus.GetZ_asInt();
105  Z *= 1;
106 
107  return result;
108 }
109 
111 //
112 //
113 
115  const G4HadProjectile& aTrack, G4Nucleus& targetNucleus)
116 {
118 
119  const G4HadProjectile* aParticle = &aTrack;
120  G4double energy = aParticle->GetTotalEnergy();
121 
122  G4String pName = aParticle->GetDefinition()->GetParticleName();
123  G4double minEnergy(0.), fmass(0.), emass = electron_mass_c2;
124 
125  if( pName == "nu_mu" || pName == "anti_nu_mu" ) fmass = theMuonMinus->GetPDGMass();
126  else if( pName == "nu_tau" || pName == "anti_nu_tau" ) fmass = theTauMinus->GetPDGMass();
127  else fmass = emass;
128 
129  minEnergy = (fmass-emass)*(fmass+emass)/emass;
130 
131  if( energy <= minEnergy )
132  {
135  return &theParticleChange;
136  }
137  G4double massf(0.), massf2(0.); // , emass = electron_mass_c2;
138  G4double sTot = 2.*energy*emass + emass*emass;
139 
140  G4LorentzVector lvp1 = aParticle->Get4Momentum();
141  G4LorentzVector lvt1(0.,0.,0.,electron_mass_c2);
142  G4LorentzVector lvsum = lvp1+lvt1;
143  G4ThreeVector bst = lvsum.boostVector();
144 
145  // sample and make final state in CMS frame
146 
147  G4double cost = SampleCosCMS( aParticle );
148  G4double sint = std::sqrt( (1.0 - cost)*(1.0 + cost) );
150 
151  G4ThreeVector eP( sint*std::cos(phi), sint*std::sin(phi), cost );
152 
153  if( pName == "nu_mu" || pName == "anti_nu_mu" ) massf = theMuonMinus->GetPDGMass();
154  else if( pName == "nu_tau" || pName == "anti_nu_tau") massf = theTauMinus->GetPDGMass();
155 
156  massf2 = massf*massf;
157 
158  G4double epf = 0.5*(sTot - massf2)/sqrt(sTot);
159  // G4double etf = epf*(sTot + massf2)/(sTot - massf2);
160 
161  eP *= epf;
162  G4LorentzVector lvp2( eP, epf );
163  lvp2.boost(bst); // back to lab frame
164 
165  G4LorentzVector lvt2 = lvsum - lvp2; // ?
166 
167  G4DynamicParticle* aNu = nullptr;
168  G4DynamicParticle* aLept = nullptr;
169 
170  if( pName == "nu_mu" || pName == "nu_tau")
171  {
172  aNu = new G4DynamicParticle( theNeutrinoE, lvp2 );
173  }
174  else if( pName == "anti_nu_mu" || pName == "anti_nu_tau")
175  {
176  aNu = new G4DynamicParticle( theAntiNeutrinoE, lvp2 );
177  }
178  if( pName == "nu_mu" || pName == "anti_nu_mu")
179  {
180  aLept = new G4DynamicParticle( theMuonMinus, lvt2 );
181  }
182  else if( pName == "nu_tau" || pName == "anti_nu_tau")
183  {
184  aLept = new G4DynamicParticle( theTauMinus, lvt2 );
185  }
186  if(aNu) { theParticleChange.AddSecondary( aNu ); }
187  if(aLept) { theParticleChange.AddSecondary( aLept ); }
188 
189  G4int Z = targetNucleus.GetZ_asInt();
190  Z *= 1;
191 
192  return &theParticleChange;
193 }
194 
196 //
197 // sample recoil electron energy in lab frame
198 
200 {
201  G4double result = 0., cofL, cofR, cofLR, massf2, sTot, emass = electron_mass_c2, emass2;
202 
203  G4double energy = aParticle->GetTotalEnergy();
204 
205  if( energy == 0.) return result; // vmg: < th?? as in xsc
206 
207  G4String pName = aParticle->GetDefinition()->GetParticleName();
208 
209  if( pName == "nu_mu" || pName == "nu_tau")
210  {
211  return 2.*G4UniformRand()-1.; // uniform scattering cos in CMS
212  }
213  else if( pName == "anti_nu_mu" || pName == "anti_nu_tau")
214  {
215  emass2 = emass*emass;
216  sTot = 2.*energy*emass + emass2;
217 
218  cofL = (sTot-emass2)/(sTot+emass2);
219 
220  if(pName == "anti_nu_mu") massf2 = theMuonMinus->GetPDGMass()*theMuonMinus->GetPDGMass();
221  else massf2 = theTauMinus->GetPDGMass()*theTauMinus->GetPDGMass();
222 
223  cofR = (sTot-massf2)/(sTot+massf2);
224 
225  cofLR = cofL*cofR/3.;
226 
227  // cofs of cos 3rd equation
228 
229  G4double a = cofLR;
230  G4double b = 0.5*(cofR+cofL);
231  G4double c = 1.;
232 
233  G4double d = -G4UniformRand()*2.*(1.+ cofLR);
234  d += c - b + a;
235 
236  // G4cout<<a<<" "<<b<<" "<<c<<" "<<d<<G4endl<<G4endl;
237 
238  // cofs of the incomplete 3rd equation
239 
240  G4double p = c/a;
241  p -= b*b/a/a/3.;
242 
243  G4double q = d/a;
244  q -= b*c/a/a/3.;
245  q += 2*b*b*b/a/a/a/27.;
246 
247 
248  // cofs for the incomplete colutions
249 
250  G4double D = p*p*p/3./3./3.;
251  D += q*q/2./2.;
252 
253  // G4cout<<"D = "<<D<<G4endl;
254  if(D < 0.) D = -D;
255  // G4complex A1 = G4complex(- q/2., std::sqrt(-D) );
256  // G4complex A = std::pow(A1,1./3.);
257 
258  // G4complex B1 = G4complex(- q/2., -std::sqrt(-D) );
259  // G4complex B = std::pow(B1,1./3.);
260 
261  G4double A, B;
262 
263  G4double A1 = - q/2. + std::sqrt(D);
264  if (A1 < 0.) A1 = -A1;
265  A = std::pow(A1,1./3.);
266  if (A1 < 0.) A = -A;
267 
268  G4double B1 = - q/2. - std::sqrt(D);
269  // G4double B = std::pow(-B1,1./3.);
270  if(B1 < 0.) B1 = -B1;
271  B = std::pow(B1,1./3.);
272  if(B1 < 0.) B = -B;
273  // G4cout<<"A1 = "<<A1<<"; A = "<<A<<"; B1 = "<<B1<<"; B = "<<B<<G4endl;
274  // roots of the incomplete 3rd equation
275 
276  G4complex y1 = A + B;
277  // G4complex y2 = -0.5*(A + B) + 0.5*std::sqrt(3.)*(A - B)*G4complex(0.,1.);
278  // G4complex y3 = -0.5*(A + B) - 0.5*std::sqrt(3.)*(A - B)*G4complex(0.,1.);
279 
280  G4complex x1 = y1 - b/a/3.;
281  // G4complex x2 = y2 - b/a/3.;
282  // G4complex x3 = y3 - b/a/3.;
283  // G4cout<<"re_x1 = "<<real(x1)<<" + i*"<<imag(x1)<<G4endl;
284  // G4cout<<"re_x1 = "<<real(x1)<<"; re_x2 = "<<real(x2)<<"; re_x3 = "<<real(x3)<<G4endl;
285  // G4cout<<"im_x1 = "<<imag(x1)<<"; im_x2 = "<<imag(x2)<<"; im_x3 = "<<imag(x3)<<G4endl<<G4endl;
286 
287  result = real(x1);
288  }
289  else
290  {
291  return result;
292  }
293  return result;
294 }
295 
296 //
297 //
const XML_Char * name
Definition: expat.h:151
void SetMinEnergy(G4double anEnergy)
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
void SetMaxEnergy(const G4double anEnergy)
Float_t y1[n_points_granero]
Definition: compare.C:5
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
void SetMomentumChange(const G4ThreeVector &aV)
Float_t x1[n_points_granero]
Definition: compare.C:5
double D(double temp)
const char * p
Definition: xmltok.h:285
const G4String & GetParticleName() const
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
G4NeutrinoElectronCcModel(const G4String &name="nu-e-elastic")
G4double GetTotalEnergy() const
virtual void ModelDescription(std::ostream &) const
G4double GetPDGMass() const
G4ParticleDefinition * theMuonMinus
G4double SampleCosCMS(const G4HadProjectile *aParticle)
Float_t Z
void SetEnergyChange(G4double anEnergy)
static constexpr double TeV
Definition: G4SIunits.hh:218
std::complex< G4double > G4complex
Definition: G4Types.hh:81
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4ParticleDefinition * theAntiNeutrinoE
static G4AntiNeutrinoE * AntiNeutrinoE()
G4ParticleDefinition * theNeutrinoE
static constexpr double electron_mass_c2
virtual G4bool IsApplicable(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
#define G4UniformRand()
Definition: Randomize.hh:53
double energy
Definition: plottest35.C:25
double A(double temperature)
Float_t d
static constexpr double eV
Definition: G4SIunits.hh:215
const G4LorentzVector & Get4Momentum() const
Hep3Vector unit() const
G4double G4ParticleHPJENDLHEData::G4double result
static constexpr double twopi
Definition: SystemOfUnits.h:55
const G4ParticleDefinition * GetDefinition() const
static G4TauMinus * TauMinus()
Definition: G4TauMinus.cc:135
int G4int
Definition: G4Types.hh:78
static G4NeutrinoE * NeutrinoE()
Definition: G4NeutrinoE.cc:85
G4ParticleDefinition * theTauMinus
Hep3Vector boostVector() const
Hep3Vector vect() const
double B(double temperature)
static constexpr double GeV
Definition: G4SIunits.hh:217
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
HepLorentzVector & boost(double, double, double)