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

#include <G4NuVacOscProcess.hh>

+ Inheritance diagram for G4NuVacOscProcess:

Public Member Functions

 G4NuVacOscProcess (const G4String &anEnvelopeName, const G4String &procName="nuVacOscillation")
 
 ~G4NuVacOscProcess () override=default
 
void InitParameters ()
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *) override
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
G4int NuVacProbability (G4int aa, G4double Enu, G4double Lnu)
 
void ProcessDescription (std::ostream &outFile) const override
 
void SetLowestEnergy (G4double e)
 
void SetBiasingFactor (G4double bf)
 
G4NuVacOscProcessoperator= (const G4NuVacOscProcess &right)=delete
 
 G4NuVacOscProcess (const G4NuVacOscProcess &)=delete
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VDiscreteProcess
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 49 of file G4NuVacOscProcess.hh.

Constructor & Destructor Documentation

◆ G4NuVacOscProcess() [1/2]

G4NuVacOscProcess::G4NuVacOscProcess ( const G4String & anEnvelopeName,
const G4String & procName = "nuVacOscillation" )

Definition at line 62 of file G4NuVacOscProcess.cc.

64{
66 fLowestEnergy = 1.*eV;
67 fEnvelopeName = eName;
68 theNuE = G4NeutrinoE::NeutrinoE();
69 theAntiNuE = G4AntiNeutrinoE::AntiNeutrinoE();
70 theNuMu = G4NeutrinoMu::NeutrinoMu();
72 theNuTau = G4NeutrinoTau::NeutrinoTau();
74
76}
@ fHadronic
static G4AntiNeutrinoE * AntiNeutrinoE()
static G4AntiNeutrinoMu * AntiNeutrinoMu()
static G4AntiNeutrinoTau * AntiNeutrinoTau()
static G4NeutrinoE * NeutrinoE()
static G4NeutrinoMu * NeutrinoMu()
static G4NeutrinoTau * NeutrinoTau()
void SetProcessSubType(G4int)

◆ ~G4NuVacOscProcess()

G4NuVacOscProcess::~G4NuVacOscProcess ( )
overridedefault

◆ G4NuVacOscProcess() [2/2]

G4NuVacOscProcess::G4NuVacOscProcess ( const G4NuVacOscProcess & )
delete

Member Function Documentation

◆ GetMeanFreePath()

G4double G4NuVacOscProcess::GetMeanFreePath ( const G4Track & aTrack,
G4double ,
G4ForceCondition *  )
overridevirtual

Implements G4VDiscreteProcess.

Definition at line 147 of file G4NuVacOscProcess.cc.

149{
150 const G4String rName =
152 G4double lambda(0.);
154
155 lambda = 0.4*CLHEP::hbarc*energy/( fDsm32 + fDsm21 );
156
157 if( rName == fEnvelopeName && fNuNuclTotXscBias > 1.)
158 {
159 lambda /= fNuNuclTotXscBias;
160 }
161 return lambda;
162}
double G4double
Definition G4Types.hh:83
G4Region * GetRegion() const
const G4String & GetName() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4StepPoint * GetPreStepPoint() const
G4double GetKineticEnergy() const
const G4Step * GetStep() const
G4LogicalVolume * GetLogicalVolume() const
G4double energy(const ThreeVector &p, const G4double m)

◆ InitParameters()

void G4NuVacOscProcess::InitParameters ( )

Definition at line 82 of file G4NuVacOscProcess.cc.

83{
84 if( fNormOrd ) // normal mass ordering
85 {
86 fSin2t12 = 0.31;
87 fSin2t23 = 0.558;
88 fSin2t13 = 0.02241;
89 fDsm21 = 7.390e-5*CLHEP::eV*CLHEP::eV;
90 fDsm32 = 2.449e-3*CLHEP::eV*CLHEP::eV;
91 fdcp = CLHEP::degree * 222.; // 270.; // 90.; // 120.; //
92 }
93 else
94 {
95 fSin2t12 = 0.31;
96 fSin2t23 = 0.563;
97 fSin2t13 = 0.02261;
98 fDsm21 = 7.3900e-5*CLHEP::eV*CLHEP::eV;
99 fDsm32 = -2.509e-3*CLHEP::eV*CLHEP::eV;
100 fdcp = CLHEP::degree * 285.; // 120. //
101 }
102 G4double c12(1.), s12(0.), c13(1.), s13(0.), c23(1.), s23(0.);
103
104 s12 = std::sqrt( fSin2t12 );
105 s23 = std::sqrt( fSin2t23 );
106 s13 = std::sqrt( fSin2t13 );
107
108 c12 = std::sqrt( 1. - fSin2t12 );
109 c23 = std::sqrt( 1. - fSin2t23 );
110 c13 = std::sqrt( 1. - fSin2t13 );
111
112 G4complex expdcp = G4complex( std::cos(fdcp), std::sin(fdcp) ); // exp(i*deltaCP)
113
114 G4complex u11, u12, u13, u21, u22, u23, u31, u32, u33;
115
116 u11 = c12*c13; u12 = c13*s12; u13 = s13*conj(expdcp);
117
118 u21 = -s12*c23 - s13*s23*c12*expdcp; u22 = c12*c23 - s12*s23*s13*expdcp; u23 = c13*s23;
119
120 u31 = s12*s23 - s13*c12*c23*expdcp; u32 = -c12*s23 - s12*s13*c23*expdcp; u33 = c13*c23;
121
122 // fUdcp[3][3] = { { u11, u12, u13 }, { u21, u22, u23 }, { u31, u32, u33 } };
123
124 // fUdcp[3][3] = { u11, u12, u13, u21, u22, u23, u31, u32, u33 };
125
126 fUdcp[0][0] = u11; fUdcp[0][1] = u12; fUdcp[0][2] = u13;
127 fUdcp[1][0] = u21; fUdcp[1][1] = u22; fUdcp[1][2] = u23;
128 fUdcp[2][0] = u31; fUdcp[2][1] = u32; fUdcp[2][2] = u33;
129
130 G4double m12, m13, m21, m23, m31, m32; //m11(0.), m22(0.),, m33(0.)
131
132 m12 = -fDsm21; m13 = -fDsm21-fDsm32;
133 m21 = -m12; m23 = -fDsm32;
134 m31 = -m13; m32 = -m23;
135
136 fDms[0][0] = fDms[1][1] = fDms[2][2] = 0.; // asymmetric
137 fDms[0][1] = m12; fDms[0][2] = m13;
138 fDms[1][0] = m21; fDms[1][2] = m23;
139 fDms[2][0] = m31; fDms[2][1] = m32;
140}
std::complex< G4double > G4complex
Definition G4Types.hh:88

Referenced by G4NuVacOscProcess().

◆ NuVacProbability()

G4int G4NuVacOscProcess::NuVacProbability ( G4int aa,
G4double Enu,
G4double Lnu )

Definition at line 238 of file G4NuVacOscProcess.cc.

239{
240 G4double probab(0.), probac(0.), probaa(0.), rr(0.), elCof(0.), delta[3][3];
241
242 G4int bb(0), cc(0);
243
244 if ( aa == 0 ) { bb = 1; cc = 2; }
245 else if( aa == 1 ) { bb = 0; cc = 2; }
246 else if( aa == 2 ) { bb = 0; cc = 1; }
247
248 elCof = 0.5*Lnu/Enu/CLHEP::hbarc;
249
250 G4complex tmp(0.,0.), sum1(0.,0.), sum2(0.,0.), expdel;
251
252 for( G4int i = 0; i < 3; ++i )
253 {
254 for( G4int j = 0; j < 3; ++j ) delta[i][j] = fDms[i][j]*elCof;
255 }
256 if( !fAnti )
257 {
258 for( G4int j = 0; j < 3; ++j )
259 {
260 for( G4int k = j+1; k < 3; ++k )
261 {
262 expdel = G4complex( std::cos( delta[k][j] ), -std::sin( delta[k][j] ) );
263
264 tmp = conj( fUdcp[bb][k] ) * fUdcp[aa][k] * fUdcp[bb][j] * conj( fUdcp[aa][j] );
265
266 sum1 += tmp * std::sin( delta[k][j]*0.5 ) * std::sin( delta[k][j]*0.5 );
267 sum2 += tmp * std::sin( delta[k][j] );
268 }
269 }
270 probab = 2.*imag(sum2) - 4.*real(sum1);
271
272 sum1 = sum2 = G4complex( 0., 0. );
273
274 for( G4int j = 0; j < 3; ++j )
275 {
276 for( G4int k = j+1; k < 3; ++k )
277 {
278 expdel = G4complex( std::cos( delta[k][j] ), -std::sin( delta[k][j] ) );
279
280 tmp = conj( fUdcp[cc][k] ) * fUdcp[aa][k] * fUdcp[cc][j] * conj( fUdcp[aa][j] );
281
282 sum1 += tmp * std::sin( delta[k][j]*0.5 ) * std::sin( delta[k][j]*0.5 );
283 sum2 += tmp * std::sin( delta[k][j] );
284 }
285 }
286 probac = 2.*imag(sum2) - 4.*real(sum1);
287 }
288 else // anti CP: exp(-i*delta)
289 {
290 for( G4int j = 0; j < 3; ++j )
291 {
292 for( G4int k = j+1; k < 3; ++k )
293 {
294 expdel = G4complex( std::cos( delta[k][j] ), -std::sin( delta[k][j] ) );
295
296 tmp = fUdcp[bb][k] * conj( fUdcp[aa][k] ) *conj( fUdcp[bb][j] ) * fUdcp[aa][j];
297
298 sum1 += tmp * std::sin( delta[k][j]*0.5 ) * std::sin( delta[k][j]*0.5 );
299 sum2 += tmp * std::sin( delta[k][j] );
300 }
301 }
302 probab = 2.*imag(sum2) - 4.*real(sum1);
303 sum1 = sum2 = G4complex( 0., 0. );
304
305 for( G4int j = 0; j < 3; ++j )
306 {
307 for( G4int k = j+1; k < 3; ++k )
308 {
309 expdel = G4complex( std::cos( delta[k][j] ), -std::sin( delta[k][j] ) );
310
311 tmp = fUdcp[cc][k] * conj( fUdcp[aa][k] ) * conj( fUdcp[cc][j] ) * fUdcp[aa][j];
312
313 sum1 += tmp * std::sin( delta[k][j]*0.5 ) * std::sin( delta[k][j]*0.5 );
314 sum2 += tmp * std::sin( delta[k][j] );
315 }
316 }
317 probac = 2.*imag(sum2) - 4.*real(sum1);
318 }
319 probaa = 1. - probab - probac;
320
321 if ( probaa < 0.)
322 {
323 G4cout<<" sum neutrino disappearance > 1. "<<G4endl;
324
325 rr = G4UniformRand()*( probab + probac );
326
327 if( rr <= probab ) return bb;
328 else return cc;
329 }
330 else
331 {
332 rr = G4UniformRand();
333
334 if ( rr <= probab ) return bb;
335 else if( rr > probab && rr <= probab + probac ) return cc;
336 else return aa;
337 }
338}
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52

Referenced by PostStepDoIt().

◆ operator=()

G4NuVacOscProcess & G4NuVacOscProcess::operator= ( const G4NuVacOscProcess & right)
delete

◆ PostStepDoIt()

G4VParticleChange * G4NuVacOscProcess::PostStepDoIt ( const G4Track & aTrack,
const G4Step & aStep )
overridevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 176 of file G4NuVacOscProcess.cc.

177{
180 if ( track.GetTrackStatus() != fAlive )
181 {
182 return &aParticleChange;
183 }
184 G4double weight = track.GetWeight();
186 G4double kineticEnergy = track.GetKineticEnergy();
187 if ( kineticEnergy <= fLowestEnergy )
188 {
189 return &aParticleChange;
190 }
191 const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
192 const G4ParticleDefinition* part = dynParticle->GetDefinition();
193 G4LorentzVector lv1 = dynParticle->Get4Momentum();
194
195 G4double ll = track.GetTrackLength(); // total track length
196 const G4String rName =
197 step.GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetRegion()->GetName();
198 if(rName == fEnvelopeName && fNuNuclTotXscBias > 1.) ll *= fNuNuclTotXscBias;
199 G4DynamicParticle* aLept = nullptr;
200
201 fAnti = (part == theAntiNuE || part == theAntiNuMu || part == theAntiNuTau);
202
203 // neutrino flavors aa and bb
204 G4int aa = 2;
205 if (part == theNuE || part == theAntiNuE) { aa = 0; }
206 else if(part == theNuMu || part == theAntiNuMu ) { aa = 1; }
207 G4int bb = NuVacProbability( aa, kineticEnergy, ll); // oscillation engine
208
209 if( bb == aa ) // no change
210 {
211 return &aParticleChange;
212 }
213 else if( bb == 0 ) // new flavor (anti)neutrino - kill initial & add new
214 {
215 if( !fAnti ) aLept = new G4DynamicParticle( theNuE, lv1 );
216 else aLept = new G4DynamicParticle( theAntiNuE, lv1 );
217 }
218 else if( bb == 1 )
219 {
220 if( !fAnti ) aLept = new G4DynamicParticle( theNuMu, lv1 );
221 else aLept = new G4DynamicParticle( theAntiNuMu, lv1 );
222 }
223 else if( bb == 2 )
224 {
225 if( !fAnti ) aLept = new G4DynamicParticle( theNuTau, lv1 );
226 else aLept = new G4DynamicParticle( theAntiNuTau, lv1 );
227 }
230
231 return &aParticleChange;
232}
@ fAlive
@ fStopAndKill
G4ParticleDefinition * GetDefinition() const
G4LorentzVector Get4Momentum() const
G4int NuVacProbability(G4int aa, G4double Enu, G4double Lnu)
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &) override
void ProposeTrackStatus(G4TrackStatus status)
void ProposeWeight(G4double finalWeight)
G4ParticleChange aParticleChange

◆ ProcessDescription()

void G4NuVacOscProcess::ProcessDescription ( std::ostream & outFile) const
overridevirtual

Reimplemented from G4VProcess.

Definition at line 166 of file G4NuVacOscProcess.cc.

167{
168 outFile << "G4NuVacOscProcess handles the oscillation of \n"
169 << "three flavor neutrinos on electrons by invoking the following model(s) and \n"
170 << "mean pathe much smaller than the oscillation period.\n";
171}

◆ SetBiasingFactor()

void G4NuVacOscProcess::SetBiasingFactor ( G4double bf)
inline

Definition at line 71 of file G4NuVacOscProcess.hh.

71 {
72 if(bf > 1.0) { fNuNuclTotXscBias = bf; fBiased = true; }
73 }

◆ SetLowestEnergy()

void G4NuVacOscProcess::SetLowestEnergy ( G4double e)
inline

Definition at line 69 of file G4NuVacOscProcess.hh.

69{ fLowestEnergy = e; }

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