Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutrinoElectronCcModel Class Reference

#include <G4NeutrinoElectronCcModel.hh>

+ Inheritance diagram for G4NeutrinoElectronCcModel:

Public Member Functions

 G4NeutrinoElectronCcModel (const G4String &name="nu-e-elastic")
 
virtual ~G4NeutrinoElectronCcModel ()
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4double SampleCosCMS (const G4HadProjectile *aParticle)
 
void SetCutEnergy (G4double ec)
 
G4double GetCutEnergy ()
 
virtual void ModelDescription (std::ostream &) const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void InitialiseModel ()
 
 G4HadronicInteraction (const G4HadronicInteraction &right)=delete
 
const G4HadronicInteractionoperator= (const G4HadronicInteraction &right)=delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 49 of file G4NeutrinoElectronCcModel.hh.

Constructor & Destructor Documentation

◆ G4NeutrinoElectronCcModel()

G4NeutrinoElectronCcModel::G4NeutrinoElectronCcModel ( const G4String name = "nu-e-elastic")

Definition at line 52 of file G4NeutrinoElectronCcModel.cc.

54{
55 SetMinEnergy( 0.0*GeV );
57 SetMinEnergy(1.e-6*eV);
58
59 theNeutrinoE = G4NeutrinoE::NeutrinoE();
60 theAntiNeutrinoE = G4AntiNeutrinoE::AntiNeutrinoE();
61
62 theNeutrinoMu = G4NeutrinoMu::NeutrinoMu();
63 theAntiNeutrinoMu = G4AntiNeutrinoMu::AntiNeutrinoMu();
64
65 theNeutrinoTau = G4NeutrinoTau::NeutrinoTau();
66 theAntiNeutrinoTau = G4AntiNeutrinoTau::AntiNeutrinoTau();
67
68 theMuonMinus = G4MuonMinus::MuonMinus();
69 theTauMinus = G4TauMinus::TauMinus();
70
71 // PDG2016: sin^2 theta Weinberg
72
73 fSin2tW = 0.23129; // 0.2312;
74
75 fCutEnergy = 0.; // default value
76
77}
static G4AntiNeutrinoE * AntiNeutrinoE()
static G4AntiNeutrinoMu * AntiNeutrinoMu()
static G4AntiNeutrinoTau * AntiNeutrinoTau()
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
static G4HadronicParameters * Instance()
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:99
static G4NeutrinoE * NeutrinoE()
Definition: G4NeutrinoE.cc:84
static G4NeutrinoMu * NeutrinoMu()
Definition: G4NeutrinoMu.cc:84
static G4NeutrinoTau * NeutrinoTau()
static G4TauMinus * TauMinus()
Definition: G4TauMinus.cc:134

◆ ~G4NeutrinoElectronCcModel()

G4NeutrinoElectronCcModel::~G4NeutrinoElectronCcModel ( )
virtual

Definition at line 80 of file G4NeutrinoElectronCcModel.cc.

81{}

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4NeutrinoElectronCcModel::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
)
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 125 of file G4NeutrinoElectronCcModel.cc.

127{
129
130 const G4HadProjectile* aParticle = &aTrack;
131 G4double energy = aParticle->GetTotalEnergy();
132
133 G4String pName = aParticle->GetDefinition()->GetParticleName();
134 G4double minEnergy(0.), fmass(0.), emass = electron_mass_c2;
135
136 if( pName == "nu_mu" ) fmass = theMuonMinus->GetPDGMass();
137 else if( pName == "nu_tau" ) fmass = theTauMinus->GetPDGMass();
138 else fmass = emass;
139
140 minEnergy = (fmass-emass)*(fmass+emass)/emass;
141
142 if( energy <= minEnergy )
143 {
146 return &theParticleChange;
147 }
148 G4double massf(0.), massf2(0.); // , emass = electron_mass_c2;
149 G4double sTot = 2.*energy*emass + emass*emass;
150
151 G4LorentzVector lvp1 = aParticle->Get4Momentum();
152 G4LorentzVector lvt1(0.,0.,0.,electron_mass_c2);
153 G4LorentzVector lvsum = lvp1+lvt1;
154 G4ThreeVector bst = lvsum.boostVector();
155
156 // sample and make final state in CMS frame
157
158 G4double cost = SampleCosCMS( aParticle );
159 G4double sint = std::sqrt( (1.0 - cost)*(1.0 + cost) );
160 G4double phi = G4UniformRand()*CLHEP::twopi;
161
162 G4ThreeVector eP( sint*std::cos(phi), sint*std::sin(phi), cost );
163
164 if( pName == "nu_mu" ) massf = theMuonMinus->GetPDGMass();
165 else if( pName == "nu_tau" ) massf = theTauMinus->GetPDGMass();
166
167 massf2 = massf*massf;
168
169 G4double epf = 0.5*(sTot - massf2)/sqrt(sTot);
170 // G4double etf = epf*(sTot + massf2)/(sTot - massf2);
171
172 eP *= epf;
173 G4LorentzVector lvp2( eP, epf );
174 lvp2.boost(bst); // back to lab frame
175
176 G4LorentzVector lvt2 = lvsum - lvp2; // ?
177
178 G4DynamicParticle* aNu = nullptr;
179 G4DynamicParticle* aLept = nullptr;
180
181 if( pName == "nu_mu" || pName == "nu_tau")
182 {
183 aNu = new G4DynamicParticle( theNeutrinoE, lvp2 );
184 }
185 else if( pName == "anti_nu_e" ) aNu = new G4DynamicParticle( theAntiNeutrinoMu, lvp2 ); // s-channel for mu (tau later)
186
187 if( pName == "nu_mu" || pName == "anti_nu_e")
188 {
189 aLept = new G4DynamicParticle( theMuonMinus, lvt2 );
190 }
191 else if( pName == "nu_tau" ) // || pName == "anti_nu_tau")
192 {
193 aLept = new G4DynamicParticle( theTauMinus, lvt2 );
194 }
195 if(aNu) { theParticleChange.AddSecondary( aNu ); }
196 if(aLept) { theParticleChange.AddSecondary( aLept ); }
197
198 G4int Z = targetNucleus.GetZ_asInt();
199 Z *= 1;
200
201 return &theParticleChange;
202}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector boostVector() const
Hep3Vector vect() const
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & Get4Momentum() const
G4double GetTotalEnergy() const
G4double SampleCosCMS(const G4HadProjectile *aParticle)
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
const G4String & GetParticleName() const
G4double energy(const ThreeVector &p, const G4double m)

◆ GetCutEnergy()

G4double G4NeutrinoElectronCcModel::GetCutEnergy ( )
inline

Definition at line 68 of file G4NeutrinoElectronCcModel.hh.

68{return fCutEnergy;};

◆ IsApplicable()

G4bool G4NeutrinoElectronCcModel::IsApplicable ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
)
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 95 of file G4NeutrinoElectronCcModel.cc.

97{
98 G4bool result = false;
99 G4String pName = aPart.GetDefinition()->GetParticleName();
100 if(pName == "anti_nu_mu" || pName == "anti_nu_tau") return result; // no cc for anti_nu_(mu,tau)
101 G4double minEnergy = 0., energy = aPart.GetTotalEnergy();
102 G4double fmass, emass = electron_mass_c2;
103
104 if( pName == "nu_mu" ) fmass = theMuonMinus->GetPDGMass();
105 else if( pName == "nu_tau" ) fmass = theTauMinus->GetPDGMass();
106 else fmass = emass;
107
108 minEnergy = (fmass-emass)*(fmass+emass)/emass;
109 SetMinEnergy( minEnergy );
110
111 if( ( pName == "nu_mu" || pName == "nu_tau" || pName == "anti_nu_e" ) && energy > minEnergy )
112 {
113 result = true;
114 }
115 G4int Z = targetNucleus.GetZ_asInt();
116 Z *= 1;
117
118 return result;
119}
bool G4bool
Definition: G4Types.hh:86

◆ ModelDescription()

void G4NeutrinoElectronCcModel::ModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 84 of file G4NeutrinoElectronCcModel.cc.

85{
86
87 outFile << "G4NeutrinoElectronCcModel is a neutrino-electron (neutral current) elastic scattering\n"
88 << "model which uses the standard model \n"
89 << "transfer parameterization. The model is fully relativistic\n";
90
91}

◆ SampleCosCMS()

G4double G4NeutrinoElectronCcModel::SampleCosCMS ( const G4HadProjectile aParticle)

Definition at line 208 of file G4NeutrinoElectronCcModel.cc.

209{
210 G4double result = 0., cofL, cofR, cofLR, massf2, sTot, emass = electron_mass_c2, emass2;
211
212 G4double energy = aParticle->GetTotalEnergy();
213
214 if( energy == 0.) return result; // vmg: < th?? as in xsc
215
216 G4String pName = aParticle->GetDefinition()->GetParticleName();
217
218 if( pName == "nu_mu" || pName == "nu_tau")
219 {
220 return 2.*G4UniformRand()-1.; // uniform scattering cos in CMS
221 }
222 else if( pName == "anti_nu_mu" || pName == "anti_nu_tau")
223 {
224 emass2 = emass*emass;
225 sTot = 2.*energy*emass + emass2;
226
227 cofL = (sTot-emass2)/(sTot+emass2);
228
229 if(pName == "anti_nu_mu") massf2 = theMuonMinus->GetPDGMass()*theMuonMinus->GetPDGMass();
230 else massf2 = theTauMinus->GetPDGMass()*theTauMinus->GetPDGMass();
231
232 cofR = (sTot-massf2)/(sTot+massf2);
233
234 cofLR = cofL*cofR/3.;
235
236 // cofs of cos 3rd equation
237
238 G4double a = cofLR;
239 G4double b = 0.5*(cofR+cofL);
240 G4double c = 1.;
241
242 G4double d = -G4UniformRand()*2.*(1.+ cofLR);
243 d += c - b + a;
244
245 // G4cout<<a<<" "<<b<<" "<<c<<" "<<d<<G4endl<<G4endl;
246
247 // cofs of the incomplete 3rd equation
248
249 G4double p = c/a;
250 p -= b*b/a/a/3.;
251
252 G4double q = d/a;
253 q -= b*c/a/a/3.;
254 q += 2*b*b*b/a/a/a/27.;
255
256
257 // cofs for the incomplete colutions
258
259 G4double D = p*p*p/3./3./3.;
260 D += q*q/2./2.;
261
262 // G4cout<<"D = "<<D<<G4endl;
263 if(D < 0.) D = -D;
264 // G4complex A1 = G4complex(- q/2., std::sqrt(-D) );
265 // G4complex A = std::pow(A1,1./3.);
266
267 // G4complex B1 = G4complex(- q/2., -std::sqrt(-D) );
268 // G4complex B = std::pow(B1,1./3.);
269
270 G4double A, B;
271
272 G4double A1 = - q/2. + std::sqrt(D);
273 if (A1 < 0.) A1 = -A1;
274 A = std::pow(A1,1./3.);
275 if (A1 < 0.) A = -A;
276
277 G4double B1 = - q/2. - std::sqrt(D);
278 // G4double B = std::pow(-B1,1./3.);
279 if(B1 < 0.) B1 = -B1;
280 B = std::pow(B1,1./3.);
281 if(B1 < 0.) B = -B;
282 // G4cout<<"A1 = "<<A1<<"; A = "<<A<<"; B1 = "<<B1<<"; B = "<<B<<G4endl;
283 // roots of the incomplete 3rd equation
284
285 G4complex y1 = A + B;
286 // G4complex y2 = -0.5*(A + B) + 0.5*std::sqrt(3.)*(A - B)*G4complex(0.,1.);
287 // G4complex y3 = -0.5*(A + B) - 0.5*std::sqrt(3.)*(A - B)*G4complex(0.,1.);
288
289 G4complex x1 = y1 - b/a/3.;
290 // G4complex x2 = y2 - b/a/3.;
291 // G4complex x3 = y3 - b/a/3.;
292 // G4cout<<"re_x1 = "<<real(x1)<<" + i*"<<imag(x1)<<G4endl;
293 // G4cout<<"re_x1 = "<<real(x1)<<"; re_x2 = "<<real(x2)<<"; re_x3 = "<<real(x3)<<G4endl;
294 // G4cout<<"im_x1 = "<<imag(x1)<<"; im_x2 = "<<imag(x2)<<"; im_x3 = "<<imag(x3)<<G4endl<<G4endl;
295
296 result = real(x1);
297 }
298 else
299 {
300 return result;
301 }
302 return result;
303}
double B(double temperature)
double A(double temperature)
double D(double temp)
std::complex< G4double > G4complex
Definition: G4Types.hh:88

Referenced by ApplyYourself().

◆ SetCutEnergy()

void G4NeutrinoElectronCcModel::SetCutEnergy ( G4double  ec)
inline

Definition at line 67 of file G4NeutrinoElectronCcModel.hh.

67{fCutEnergy=ec;};

The documentation for this class was generated from the following files: