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

#include <G4NeutrinoElectronProcess.hh>

+ Inheritance diagram for G4NeutrinoElectronProcess:

Public Member Functions

 G4NeutrinoElectronProcess (const G4String &anEnvelopeName, const G4String &procName="nuElectron")
 
 ~G4NeutrinoElectronProcess () override=default
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *) override
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
void ProcessDescription (std::ostream &outFile) const override
 
void SetLowestEnergy (G4double)
 
void SetBiasingFactors (G4double bfCc, G4double bfNc)
 
void SetBiasingFactor (G4double bf)
 
G4NeutrinoElectronProcessoperator= (const G4NeutrinoElectronProcess &right)=delete
 
 G4NeutrinoElectronProcess (const G4NeutrinoElectronProcess &)=delete
 
- Public Member Functions inherited from G4HadronicProcess
 G4HadronicProcess (const G4String &processName="Hadronic", G4ProcessType procType=fHadronic)
 
 G4HadronicProcess (const G4String &processName, G4HadronicProcessType subType)
 
 ~G4HadronicProcess () override
 
void RegisterMe (G4HadronicInteraction *a)
 
G4double GetElementCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
 
G4double GetMicroscopicCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
 
void StartTracking (G4Track *track) override
 
void PreparePhysicsTable (const G4ParticleDefinition &) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
void DumpPhysicsTable (const G4ParticleDefinition &p)
 
void AddDataSet (G4VCrossSectionDataSet *aDataSet)
 
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList ()
 
G4HadronicInteractionGetHadronicModel (const G4String &)
 
G4HadronicInteractionGetHadronicInteraction () const
 
const G4NucleusGetTargetNucleus () const
 
G4NucleusGetTargetNucleusPointer ()
 
const G4IsotopeGetTargetIsotope ()
 
G4double ComputeCrossSection (const G4ParticleDefinition *, const G4Material *, const G4double kinEnergy)
 
G4HadXSType CrossSectionType () const
 
void SetCrossSectionType (G4HadXSType val)
 
void BiasCrossSectionByFactor (G4double aScale)
 
void MultiplyCrossSectionBy (G4double factor)
 
G4double CrossSectionFactor () const
 
void SetIntegral (G4bool val)
 
void SetEpReportLevel (G4int level)
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
G4CrossSectionDataStoreGetCrossSectionDataStore ()
 
std::vector< G4TwoPeaksHadXS * > * TwoPeaksXS () const
 
std::vector< G4double > * EnergyOfCrossSectionMax () const
 
G4HadronicProcessoperator= (const G4HadronicProcess &right)=delete
 
 G4HadronicProcess (const G4HadronicProcess &)=delete
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
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 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 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 G4HadronicProcess
G4HadronicInteractionChooseHadronicInteraction (const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, const G4Material *aMaterial, const G4Element *anElement)
 
G4double GetLastCrossSection ()
 
void FillResult (G4HadFinalState *aR, const G4Track &aT)
 
void DumpState (const G4Track &, const G4String &, G4ExceptionDescription &)
 
G4HadFinalStateCheckResult (const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result)
 
void CheckEnergyMomentumConservation (const G4Track &, const G4Nucleus &)
 
- Protected Member Functions inherited from G4VDiscreteProcess
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4HadronicProcess
G4HadProjectile thePro
 
G4ParticleChangetheTotalResult
 
G4CrossSectionDataStoretheCrossSectionDataStore
 
G4double fWeight = 1.0
 
G4double aScaleFactor = 1.0
 
G4double theLastCrossSection = 0.0
 
G4double mfpKinEnergy = DBL_MAX
 
G4int epReportLevel = 0
 
G4HadXSType fXSType = fHadNoIntegral
 
- 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 50 of file G4NeutrinoElectronProcess.hh.

Constructor & Destructor Documentation

◆ G4NeutrinoElectronProcess() [1/2]

G4NeutrinoElectronProcess::G4NeutrinoElectronProcess ( const G4String & anEnvelopeName,
const G4String & procName = "nuElectron" )

Definition at line 64 of file G4NeutrinoElectronProcess.cc.

66{
67 lowestEnergy = 1.*keV;
68 fEnvelopeName = anEnvelopeName;
69 fTotXsc = new G4NeutrinoElectronTotXsc();
71 safetyHelper->InitialiseHelper();
72}
G4HadronicProcess(const G4String &processName="Hadronic", G4ProcessType procType=fHadronic)
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const

◆ ~G4NeutrinoElectronProcess()

G4NeutrinoElectronProcess::~G4NeutrinoElectronProcess ( )
overridedefault

◆ G4NeutrinoElectronProcess() [2/2]

G4NeutrinoElectronProcess::G4NeutrinoElectronProcess ( const G4NeutrinoElectronProcess & )
delete

Member Function Documentation

◆ GetMeanFreePath()

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

Reimplemented from G4HadronicProcess.

Definition at line 102 of file G4NeutrinoElectronProcess.cc.

104{
105 //G4cout << "GetMeanFreePath " << aTrack.GetDefinition()->GetParticleName()
106 // << " Ekin= " << aTrack.GetKineticEnergy() << G4endl;
108 G4double totxsc =
110 aTrack.GetMaterial());
111
112 if ( rName == fEnvelopeName )
113 {
114 totxsc *= fNuEleTotXscBias;
115 }
116 G4double res = (totxsc>0.0) ? 1.0/totxsc : DBL_MAX;
117 //G4cout << " xsection= " << totxsc << G4endl;
118 return res;
119}
double G4double
Definition G4Types.hh:83
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
G4CrossSectionDataStore * GetCrossSectionDataStore()
G4Region * GetRegion() const
const G4String & GetName() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4StepPoint * GetPreStepPoint() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
const G4Step * GetStep() const
G4LogicalVolume * GetLogicalVolume() const
#define DBL_MAX
Definition templates.hh:62

◆ operator=()

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

◆ PostStepDoIt()

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

! is not needed for models inheriting G4NeutrinoElectron

Reimplemented from G4HadronicProcess.

Definition at line 133 of file G4NeutrinoElectronProcess.cc.

134{
135 G4String rName = track.GetStep()->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetRegion()->GetName();
136
137 if( rName != fEnvelopeName )
138 {
139 if( verboseLevel > 0 )
140 {
141 G4cout<<"Go out from G4NeutrinoElectronProcess::PostStepDoIt: wrong volume "<<G4endl;
142 }
143 return G4VDiscreteProcess::PostStepDoIt( track, step );
144 }
147 G4double weight = track.GetWeight();
149
150 if( track.GetTrackStatus() != fAlive )
151 {
152 return theTotalResult;
153 }
154 // Next check for illegal track status
155 //
156 if (track.GetTrackStatus() != fAlive &&
157 track.GetTrackStatus() != fSuspend)
158 {
159 if (track.GetTrackStatus() == fStopAndKill ||
160 track.GetTrackStatus() == fKillTrackAndSecondaries ||
161 track.GetTrackStatus() == fPostponeToNextEvent)
162 {
164 ed << "G4HadronicProcess: track in unusable state - "
165 << track.GetTrackStatus() << G4endl;
166 ed << "G4HadronicProcess: returning unchanged track " << G4endl;
167 DumpState(track,"PostStepDoIt",ed);
168 G4Exception("G4HadronicProcess::PostStepDoIt", "had004", JustWarning, ed);
169 }
170 // No warning for fStopButAlive which is a legal status here
171 return theTotalResult;
172 }
173
174 // For elastic scattering, _any_ result is considered an interaction
176
177 G4double kineticEnergy = track.GetKineticEnergy();
178 const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
179 const G4ParticleDefinition* part = dynParticle->GetDefinition();
180
181 // NOTE: Very low energy scatters were causing numerical (FPE) errors
182 // in earlier releases; these limits have not been changed since.
183
184 if ( kineticEnergy <= lowestEnergy ) return theTotalResult;
185
186 const G4Material* material = track.GetMaterial();
187 G4Nucleus* targNucleus = GetTargetNucleusPointer();
188
189 //////////////// uniform random spread of the neutrino interaction point ////////////
190
191 const G4StepPoint* pPostStepPoint = step.GetPostStepPoint();
192 const G4DynamicParticle* aParticle = track.GetDynamicParticle();
193 G4ThreeVector position = pPostStepPoint->GetPosition(), newPosition=position;
194 G4ParticleMomentum direction = aParticle->GetMomentumDirection();
195
196 if( fNuEleCcBias > 1.0 || fNuEleNcBias > 1.0) // = true, if fBiasingfactor != 1., i.e. xsc is biased
197 {
198 const G4RotationMatrix* rotM = pPostStepPoint->GetTouchable()->GetRotation();
199 G4ThreeVector transl = pPostStepPoint->GetTouchable()->GetTranslation();
200 G4AffineTransform transform = G4AffineTransform(rotM,transl);
201 transform.Invert();
202
203 G4ThreeVector localP = transform.TransformPoint(position);
204 G4ThreeVector localV = transform.TransformAxis(direction);
205
206 G4double forward = track.GetVolume()->GetLogicalVolume()->GetSolid()->DistanceToOut(localP, localV);
207 G4double backward = track.GetVolume()->GetLogicalVolume()->GetSolid()->DistanceToOut(localP, -localV);
208
209 G4double distance = forward+backward;
210
211 // G4cout<<distance/cm<<", ";
212
213 // uniform sampling of nu-e interaction point
214 // along neutrino direction in current volume
215
216 G4double range = -backward+G4UniformRand()*distance;
217
218 newPosition = position + range*direction;
219
220 safetyHelper->ReLocateWithinVolume(newPosition);
221
222 theTotalResult->ProposePosition(newPosition); // G4Exception : GeomNav1002
223 }
224 G4HadProjectile theProj( track );
225 G4HadronicInteraction* hadi = nullptr;
226 G4HadFinalState* result = nullptr;
227
228 // Select element
229 const G4Element* elm = nullptr;
230
231 try
232 {
233 elm = GetCrossSectionDataStore()->SampleZandA(dynParticle, material,
234 *targNucleus);
235 }
236 catch( G4HadronicException & aR )
237 {
239 aR.Report(ed);
240 DumpState(track,"SampleZandA",ed);
241 ed << " PostStepDoIt failed on element selection" << G4endl;
242 G4Exception("G4NeutrinoElectronProcess::PostStepDoIt", "had003",
243 FatalException, ed);
244 }
245
246 G4double ccTotRatio = fTotXsc->GetCcRatio();
247
248 if( G4UniformRand() < ccTotRatio ) // Cc-model
249 {
250 // Initialize the hadronic projectile from the track
251 thePro.Initialise(track);
252
253 hadi = (GetHadronicInteractionList())[0];
254
255 result = hadi->ApplyYourself( thePro, *targNucleus);
256
258
260
261 FillResult(result, track);
262 }
263 else // Nc-model, like 'elastic', 2->2 scattering
264 {
265
266 hadi = (GetHadronicInteractionList())[1];
267
268 size_t idx = track.GetMaterialCutsCouple()->GetIndex();
269
271
272 hadi->SetRecoilEnergyThreshold(tcut);
273
274 if( verboseLevel > 1 )
275 {
276 G4cout << "G4NeutrinoElectronProcess::PostStepDoIt for "
277 << part->GetParticleName()
278 << " in " << material->GetName()
279 << " Target Z= " << targNucleus->GetZ_asInt()
280 << " A= " << targNucleus->GetA_asInt() << G4endl;
281 }
282 try
283 {
284 result = hadi->ApplyYourself( theProj, *targNucleus);
285 }
286 catch(G4HadronicException & aR)
287 {
289 aR.Report(ed);
290 ed << "Call for " << hadi->GetModelName() << G4endl;
291 ed << "Target element "<< elm->GetName()<<" Z= "
292 << targNucleus->GetZ_asInt()
293 << " A= " << targNucleus->GetA_asInt() << G4endl;
294 DumpState(track,"ApplyYourself",ed);
295 ed << " ApplyYourself failed" << G4endl;
296 G4Exception("G4NeutrinoElectronProcess::PostStepDoIt", "had006",
297 FatalException, ed);
298 }
299 // directions
300
301 G4ThreeVector indir = track.GetMomentumDirection();
302 G4double phi = CLHEP::twopi*G4UniformRand();
303 G4ThreeVector it(0., 0., 1.);
304 G4ThreeVector outdir = result->GetMomentumChange();
305
306 if(verboseLevel>1)
307 {
308 G4cout << "Efin= " << result->GetEnergyChange()
309 << " de= " << result->GetLocalEnergyDeposit()
310 << " nsec= " << result->GetNumberOfSecondaries()
311 << " dir= " << outdir
312 << G4endl;
313 }
314 // energies
315
316 G4double edep = result->GetLocalEnergyDeposit();
317 G4double efinal = result->GetEnergyChange();
318
319 if(efinal < 0.0) { efinal = 0.0; }
320 if(edep < 0.0) { edep = 0.0; }
321
322 // NOTE: Very low energy scatters were causing numerical (FPE) errors
323 // in earlier releases; these limits have not been changed since.
324
325 if(efinal <= lowestEnergy)
326 {
327 edep += efinal;
328 efinal = 0.0;
329 }
330 // primary change
331
333
334 G4TrackStatus status = track.GetTrackStatus();
335
336 if(efinal > 0.0)
337 {
338 outdir.rotate(phi, it);
339 outdir.rotateUz(indir);
341 }
342 else
343 {
344 if( part->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
345 {
346 status = fStopButAlive;
347 }
348 else
349 {
350 status = fStopAndKill;
351 }
353 }
354 //G4cout << "Efinal= " << efinal << " TrackStatus= " << status << G4endl;
355
357
358 // recoil
359
360 if( result->GetNumberOfSecondaries() > 0 )
361 {
362 G4DynamicParticle* p = result->GetSecondary(0)->GetParticle();
363
364 if(p->GetKineticEnergy() > tcut)
365 {
368
369 // G4cout << "recoil " << pdir << G4endl;
370 //!! is not needed for models inheriting G4NeutrinoElectron
371
372 pdir.rotate(phi, it);
373 pdir.rotateUz(indir);
374
375 // G4cout << "recoil rotated " << pdir << G4endl;
376
377 p->SetMomentumDirection(pdir);
378
379 // in elastic scattering time and weight are not changed
380
381 G4Track* t = new G4Track(p, track.GetGlobalTime(),
382 track.GetPosition());
383 t->SetWeight(weight);
384 t->SetTouchableHandle(track.GetTouchableHandle());
386 }
387 else
388 {
389 edep += p->GetKineticEnergy();
390 delete p;
391 }
392 }
395 result->Clear();
396 }
397 return theTotalResult;
398}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4TrackStatus
@ fKillTrackAndSecondaries
@ fSuspend
@ fAlive
@ fStopAndKill
@ fStopButAlive
@ fPostponeToNextEvent
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector & rotateUz(const Hep3Vector &)
Hep3Vector & rotate(double, const Hep3Vector &)
G4AffineTransform & Invert()
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
const G4Element * SampleZandA(const G4DynamicParticle *, const G4Material *, G4Nucleus &target)
void SetMomentumDirection(const G4ThreeVector &aDirection)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4String & GetName() const
Definition G4Element.hh:115
G4double GetEnergyChange() const
void SetTrafoToLab(const G4LorentzRotation &aT)
G4double GetLocalEnergyDeposit() const
const G4ThreeVector & GetMomentumChange() const
std::size_t GetNumberOfSecondaries() const
G4HadSecondary * GetSecondary(size_t i)
void Initialise(const G4Track &aT)
G4LorentzRotation & GetTrafoToLab()
G4DynamicParticle * GetParticle()
void Report(std::ostream &aS) const
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
const G4String & GetModelName() const
void SetRecoilEnergyThreshold(G4double val)
void FillResult(G4HadFinalState *aR, const G4Track &aT)
G4HadProjectile thePro
G4Nucleus * GetTargetNucleusPointer()
G4ParticleChange * theTotalResult
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList()
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
const G4String & GetName() const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
G4int GetA_asInt() const
Definition G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition G4Nucleus.hh:105
void AddSecondary(G4Track *aSecondary)
void ProposePosition(G4double x, G4double y, G4double z)
void Initialize(const G4Track &) override
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
void ReLocateWithinVolume(const G4ThreeVector &pGlobalPoint)
const G4VTouchable * GetTouchable() const
const G4ThreeVector & GetPosition() const
virtual const G4RotationMatrix * GetRotation(G4int depth=0) const
virtual const G4ThreeVector & GetTranslation(G4int depth=0) const
void SetWeight(G4double aValue)
void SetTouchableHandle(const G4TouchableHandle &apValue)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
void ProposeWeight(G4double finalWeight)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
void ClearNumberOfInteractionLengthLeft()
G4int verboseLevel

◆ PostStepGetPhysicalInteractionLength()

G4double G4NeutrinoElectronProcess::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Reimplemented from G4HadronicProcess.

Definition at line 92 of file G4NeutrinoElectronProcess.cc.

95{
97 previousStepSize, condition );
98}
G4double condition(const G4ErrorSymMatrix &m)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)

◆ ProcessDescription()

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

Reimplemented from G4HadronicProcess.

Definition at line 123 of file G4NeutrinoElectronProcess.cc.

124{
125 outFile << "G4NeutrinoElectronProcess handles the scattering of \n"
126 << "neutrino on electrons by invoking the following model(s) and \n"
127 << "cross section(s).\n";
128}

◆ SetBiasingFactor()

void G4NeutrinoElectronProcess::SetBiasingFactor ( G4double bf)

Definition at line 76 of file G4NeutrinoElectronProcess.cc.

77{
78 fNuEleTotXscBias = bf;
79}

◆ SetBiasingFactors()

void G4NeutrinoElectronProcess::SetBiasingFactors ( G4double bfCc,
G4double bfNc )

Definition at line 83 of file G4NeutrinoElectronProcess.cc.

84{
85 fNuEleCcBias = bfCc;
86 fNuEleNcBias = bfNc;
87 fNuEleTotXscBias = std::max( fNuEleCcBias, fNuEleNcBias );
88}

◆ SetLowestEnergy()

void G4NeutrinoElectronProcess::SetLowestEnergy ( G4double val)

Definition at line 401 of file G4NeutrinoElectronProcess.cc.

402{
403 lowestEnergy = val;
404}

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