Geant4 11.1.1
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-inelastic")
 
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-inelastic")

Definition at line 53 of file G4NeutrinoElectronCcModel.cc.

55{
56 SetMinEnergy( 0.0*GeV );
58 SetMinEnergy(1.e-6*eV);
59
60 theNeutrinoE = G4NeutrinoE::NeutrinoE();
61 theAntiNeutrinoE = G4AntiNeutrinoE::AntiNeutrinoE();
62
63 theNeutrinoMu = G4NeutrinoMu::NeutrinoMu();
64 theAntiNeutrinoMu = G4AntiNeutrinoMu::AntiNeutrinoMu();
65
66 theNeutrinoTau = G4NeutrinoTau::NeutrinoTau();
67 theAntiNeutrinoTau = G4AntiNeutrinoTau::AntiNeutrinoTau();
68
69 theMuonMinus = G4MuonMinus::MuonMinus();
70 theTauMinus = G4TauMinus::TauMinus();
71
72 // PDG2016: sin^2 theta Weinberg
73
74 fSin2tW = 0.23129; // 0.2312;
75
76 fCutEnergy = 0.; // default value
77
78 // Creator model ID
79 secID = G4PhysicsModelCatalog::GetModelID( "model_" + GetModelName() );
80}
static G4AntiNeutrinoE * AntiNeutrinoE()
static G4AntiNeutrinoMu * AntiNeutrinoMu()
static G4AntiNeutrinoTau * AntiNeutrinoTau()
void SetMinEnergy(G4double anEnergy)
const G4String & GetModelName() const
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 G4int GetModelID(const G4int modelIndex)
static G4TauMinus * TauMinus()
Definition: G4TauMinus.cc:134

◆ ~G4NeutrinoElectronCcModel()

G4NeutrinoElectronCcModel::~G4NeutrinoElectronCcModel ( )
virtual

Definition at line 83 of file G4NeutrinoElectronCcModel.cc.

84{}

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 126 of file G4NeutrinoElectronCcModel.cc.

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

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

◆ ModelDescription()

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

Reimplemented from G4HadronicInteraction.

Definition at line 87 of file G4NeutrinoElectronCcModel.cc.

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

◆ SampleCosCMS()

G4double G4NeutrinoElectronCcModel::SampleCosCMS ( const G4HadProjectile aParticle)

Definition at line 206 of file G4NeutrinoElectronCcModel.cc.

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

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: