Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4NeutrinoElectronCcXsc.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 
29 #include "G4PhysicalConstants.hh"
30 #include "G4SystemOfUnits.hh"
31 #include "G4DynamicParticle.hh"
32 #include "G4ParticleTable.hh"
33 #include "G4IonTable.hh"
34 #include "G4HadTmpUtil.hh"
35 #include "G4NistManager.hh"
36 
37 #include "G4MuonMinus.hh"
38 #include "G4TauMinus.hh"
39 
40 using namespace std;
41 using namespace CLHEP;
42 
44  : G4VCrossSectionDataSet("NuElectronCcXsc")
45 {
46  // PDG2016: Gf=1.1663787(6)e-5*(hc)^3/GeV^2
47  // fCofXsc = Gf*Gf*MeC2*2/pi
48 
49  fCofXsc = 1.36044e-22;
51  fCofXsc /= halfpi;
52 
53  // G4cout<<"fCofXsc = "<<fCofXsc*GeV/cm2<<" cm2/GeV"<<G4endl;
54 
55  // G4cout<<"hbarc = "<<hbarc/MeV/fermi<<" MeV*fermi"<<G4endl;
56 
57  // PDG2016: sin^2 theta Weinberg
58 
59  fSin2tW = 0.23129; // 0.2312;
60 
61  fCutEnergy = 0.; // default value
62 
63  fBiasingFactor = 1.; // default as physics
64 
67 }
68 
70 {}
71 
73 
74 G4bool
76 {
77  G4bool result = false;
78  G4String pName = aPart->GetDefinition()->GetParticleName();
79  G4double minEnergy = 0., energy = aPart->GetTotalEnergy();
80  G4double fmass, emass = electron_mass_c2;
81 
82  if( pName == "nu_mu" || pName == "anti_nu_mu" ) fmass = theMuonMinus->GetPDGMass();
83  else if( pName == "nu_tau" || pName == "anti_nu_tau" ) fmass = theTauMinus->GetPDGMass();
84  else fmass = emass;
85 
86  minEnergy = (fmass-emass)*(fmass+emass)/emass;
87 
88  if( ( pName == "nu_mu" || pName == "anti_nu_mu" ||
89  pName == "nu_tau" || pName == "anti_nu_tau" ) &&
90  energy > minEnergy )
91  {
92  result = true;
93  }
94  return result;
95 }
96 
98 
101  const G4Material*)
102 {
103  G4double result = 0., totS, fmass, fmass2, emass=electron_mass_c2, emass2;
104 
105  G4double energy = aPart->GetTotalEnergy();
106  G4String pName = aPart->GetDefinition()->GetParticleName();
107 
108  emass2 = emass*emass;
109  totS = 2.*energy*emass + emass2;
110 
111  if( pName == "nu_mu")
112  {
113  fmass = theMuonMinus->GetPDGMass();
114  fmass2 = fmass*fmass;
115  result = (1. - fmass2/totS)*(1. - fmass2/totS);
116  }
117  else if( pName == "anti_nu_mu")
118  {
119  fmass = theMuonMinus->GetPDGMass();
120  fmass2 = fmass*fmass;
121 
122  result = (1.+ emass2/totS)*(1.+ fmass2/totS);
123  result += (1.- emass2/totS)*(1.- fmass2/totS)/3.;
124  result *= 0.25*(1. - fmass2/totS)*(1. - fmass2/totS);
125  }
126  else if( pName == "nu_tau")
127  {
128  fmass = theTauMinus->GetPDGMass();
129  fmass2 = fmass*fmass;
130  result = (1. - fmass2/totS)*(1. - fmass2/totS);
131  }
132  else if( pName == "anti_nu_tau")
133  {
134  fmass = theTauMinus->GetPDGMass();
135  fmass2 = fmass*fmass;
136 
137  result = (1.+ emass2/totS)*(1.+ fmass2/totS);
138  result += (1.- emass2/totS)*(1.- fmass2/totS)/3.;
139  result *= 0.25*(1. - fmass2/totS)*(1. - fmass2/totS);
140  }
141  else
142  {
143  return result;
144  }
145  // if( energy <= electron_mass_c2 ) return result;
146 
147  result *= fCofXsc; //*energy;
148  result *= energy + 0.5*emass;
149  result *= ZZ; // incoherent sum over all element electrons
150 
151  result *= fBiasingFactor; // biasing up, if set >1
152 
153  return result;
154 }
155 
G4ParticleDefinition * theTauMinus
G4ParticleDefinition * theMuonMinus
static constexpr double hbarc
const G4String & GetParticleName() const
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *)
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
G4double GetPDGMass() const
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4ParticleDefinition * GetDefinition() const
static constexpr double electron_mass_c2
double energy
Definition: plottest35.C:25
static constexpr double halfpi
Definition: G4SIunits.hh:77
G4double G4ParticleHPJENDLHEData::G4double result
static G4TauMinus * TauMinus()
Definition: G4TauMinus.cc:135
int G4int
Definition: G4Types.hh:78
G4double GetTotalEnergy() const