Geant4 11.1.1
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 (G4String anEnvelopeName, const G4String &procName="neutrino-electron")
 
virtual ~G4NeutrinoElectronProcess ()
 
virtual G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual void SetLowestEnergy (G4double)
 
virtual void ProcessDescription (std::ostream &outFile) const
 
void SetBiasingFactors (G4double bfCc, G4double bfNc)
 
void SetBiasingFactor (G4double bf)
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *)
 
- 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
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double, G4ForceCondition *) override
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) 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
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *) override
 
const G4NucleusGetTargetNucleus () const
 
G4NucleusGetTargetNucleusPointer ()
 
const G4IsotopeGetTargetIsotope ()
 
G4double ComputeCrossSection (const G4ParticleDefinition *, const G4Material *, const G4double kinEnergy)
 
G4HadXSType CrossSectionType () const
 
void SetCrossSectionType (G4HadXSType val)
 
void ProcessDescription (std::ostream &outFile) const override
 
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 PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
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
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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
 
virtual void ProcessDescription (std::ostream &outfile) 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 &)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- 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
 
G4long 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 52 of file G4NeutrinoElectronProcess.hh.

Constructor & Destructor Documentation

◆ G4NeutrinoElectronProcess()

G4NeutrinoElectronProcess::G4NeutrinoElectronProcess ( G4String  anEnvelopeName,
const G4String procName = "neutrino-electron" 
)

Definition at line 66 of file G4NeutrinoElectronProcess.cc.

67 : G4HadronicProcess( pName, fHadronElastic ), isInitialised(false), fBiased(true) // fHadronElastic???
68{
69 lowestEnergy = 1.*keV;
70 fEnvelope = nullptr;
71 fEnvelopeName = anEnvelopeName;
72 fTotXsc = nullptr; // new G4NeutrinoElectronTotXsc();
73 fNuEleCcBias=1.;
74 fNuEleNcBias=1.;
75 fNuEleTotXscBias=1.;
77 safetyHelper->InitialiseHelper();
78}
@ fHadronElastic
void InitialiseHelper()
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const

◆ ~G4NeutrinoElectronProcess()

G4NeutrinoElectronProcess::~G4NeutrinoElectronProcess ( )
virtual

Definition at line 80 of file G4NeutrinoElectronProcess.cc.

81{
82 // if( fTotXsc ) delete fTotXsc;
83}

Member Function Documentation

◆ GetMeanFreePath()

G4double G4NeutrinoElectronProcess::GetMeanFreePath ( const G4Track aTrack,
G4double  ,
G4ForceCondition  
)
virtual

Reimplemented from G4HadronicProcess.

Definition at line 108 of file G4NeutrinoElectronProcess.cc.

110{
111 //G4cout << "GetMeanFreePath " << aTrack.GetDefinition()->GetParticleName()
112 // << " Ekin= " << aTrack.GetKineticEnergy() << G4endl;
114 G4double totxsc(0.);
115 try
116 {
117 if( rName == fEnvelopeName && fNuEleTotXscBias > 1.)
118 {
119 totxsc = fNuEleTotXscBias*
121 aTrack.GetMaterial());
122 }
123 else
124 {
126 aTrack.GetMaterial());
127 }
128 }
129 catch(G4HadronicException & aR)
130 {
132 aR.Report(ed);
133 DumpState(aTrack,"GetMeanFreePath",ed);
134 ed << " Cross section is not available" << G4endl;
135 G4Exception("G4NeutrinoElectronProcess::GetMeanFreePath", "had002", FatalException,
136 ed);
137 }
138 G4double res = (totxsc>0.0) ? 1.0/totxsc : DBL_MAX;
139 //G4cout << " xsection= " << totxsc << G4endl;
140 return res;
141}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
void Report(std::ostream &aS) const
G4CrossSectionDataStore * GetCrossSectionDataStore()
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
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

◆ PostStepDoIt()

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

! is not needed for models inheriting G4NeutrinoElectron

Reimplemented from G4HadronicProcess.

Definition at line 157 of file G4NeutrinoElectronProcess.cc.

158{
159 // track.GetVolume()->GetLogicalVolume()->GetName()
160 // if( track.GetVolume()->GetLogicalVolume() != fEnvelope )
161
162 G4String rName = track.GetStep()->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetRegion()->GetName();
163
164 if( rName != fEnvelopeName )
165 {
166 if( verboseLevel > 0 )
167 {
168 G4cout<<"Go out from G4NeutrinoElectronProcess::PostStepDoIt: wrong volume "<<G4endl;
169 }
170 return G4VDiscreteProcess::PostStepDoIt( track, step );
171 }
174 G4double weight = track.GetWeight();
176
177 if( track.GetTrackStatus() != fAlive )
178 {
179 return theTotalResult;
180 }
181 // Next check for illegal track status
182 //
183 if (track.GetTrackStatus() != fAlive &&
184 track.GetTrackStatus() != fSuspend)
185 {
186 if (track.GetTrackStatus() == fStopAndKill ||
187 track.GetTrackStatus() == fKillTrackAndSecondaries ||
188 track.GetTrackStatus() == fPostponeToNextEvent)
189 {
191 ed << "G4HadronicProcess: track in unusable state - "
192 << track.GetTrackStatus() << G4endl;
193 ed << "G4HadronicProcess: returning unchanged track " << G4endl;
194 DumpState(track,"PostStepDoIt",ed);
195 G4Exception("G4HadronicProcess::PostStepDoIt", "had004", JustWarning, ed);
196 }
197 // No warning for fStopButAlive which is a legal status here
198 return theTotalResult;
199 }
200
201 // For elastic scattering, _any_ result is considered an interaction
203
204 G4double kineticEnergy = track.GetKineticEnergy();
205 const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
206 const G4ParticleDefinition* part = dynParticle->GetDefinition();
207
208 // NOTE: Very low energy scatters were causing numerical (FPE) errors
209 // in earlier releases; these limits have not been changed since.
210
211 if ( kineticEnergy <= lowestEnergy ) return theTotalResult;
212
213 const G4Material* material = track.GetMaterial();
214 G4Nucleus* targNucleus = GetTargetNucleusPointer();
215
216 //////////////// uniform random spread of the neutrino interaction point ////////////
217
218 const G4StepPoint* pPostStepPoint = step.GetPostStepPoint();
219 const G4DynamicParticle* aParticle = track.GetDynamicParticle();
220 G4ThreeVector position = pPostStepPoint->GetPosition(), newPosition=position;
221 G4ParticleMomentum direction = aParticle->GetMomentumDirection();
222
223 if( fNuEleCcBias > 1.0 || fNuEleNcBias > 1.0) // = true, if fBiasingfactor != 1., i.e. xsc is biased
224 {
225 const G4RotationMatrix* rotM = pPostStepPoint->GetTouchable()->GetRotation();
226 G4ThreeVector transl = pPostStepPoint->GetTouchable()->GetTranslation();
227 G4AffineTransform transform = G4AffineTransform(rotM,transl);
228 transform.Invert();
229
230 G4ThreeVector localP = transform.TransformPoint(position);
231 G4ThreeVector localV = transform.TransformAxis(direction);
232
233 G4double forward = track.GetVolume()->GetLogicalVolume()->GetSolid()->DistanceToOut(localP, localV);
234 G4double backward = track.GetVolume()->GetLogicalVolume()->GetSolid()->DistanceToOut(localP, -localV);
235
236 G4double distance = forward+backward;
237
238 // G4cout<<distance/cm<<", ";
239
240 // uniform sampling of nu-e interaction point
241 // along neutrino direction in current volume
242
243 G4double range = -backward+G4UniformRand()*distance;
244
245 newPosition = position + range*direction;
246
247 safetyHelper->ReLocateWithinVolume(newPosition);
248
249 theTotalResult->ProposePosition(newPosition); // G4Exception : GeomNav1002
250 }
251 G4HadProjectile theProj( track );
252 G4HadronicInteraction* hadi = nullptr;
253 G4HadFinalState* result = nullptr;
254
255 // Select element
256 const G4Element* elm = nullptr;
257
258 try
259 {
260 elm = GetCrossSectionDataStore()->SampleZandA(dynParticle, material,
261 *targNucleus);
262 }
263 catch( G4HadronicException & aR )
264 {
266 aR.Report(ed);
267 DumpState(track,"SampleZandA",ed);
268 ed << " PostStepDoIt failed on element selection" << G4endl;
269 G4Exception("G4NeutrinoElectronProcess::PostStepDoIt", "had003",
270 FatalException, ed);
271 }
272
273 G4double ccTotRatio = fTotXsc->GetCcRatio();
274
275 if( G4UniformRand() < ccTotRatio ) // Cc-model
276 {
277 // Initialize the hadronic projectile from the track
278 thePro.Initialise(track);
279
280 hadi = (GetHadronicInteractionList())[0];
281
282 result = hadi->ApplyYourself( thePro, *targNucleus);
283
285
287
288 FillResult(result, track);
289 }
290 else // Nc-model, like 'elastic', 2->2 scattering
291 {
292
293 hadi = (GetHadronicInteractionList())[1];
294
295 size_t idx = track.GetMaterialCutsCouple()->GetIndex();
296
298
299 hadi->SetRecoilEnergyThreshold(tcut);
300
301 if( verboseLevel > 1 )
302 {
303 G4cout << "G4NeutrinoElectronProcess::PostStepDoIt for "
304 << part->GetParticleName()
305 << " in " << material->GetName()
306 << " Target Z= " << targNucleus->GetZ_asInt()
307 << " A= " << targNucleus->GetA_asInt() << G4endl;
308 }
309 try
310 {
311 result = hadi->ApplyYourself( theProj, *targNucleus);
312 }
313 catch(G4HadronicException & aR)
314 {
316 aR.Report(ed);
317 ed << "Call for " << hadi->GetModelName() << G4endl;
318 ed << "Target element "<< elm->GetName()<<" Z= "
319 << targNucleus->GetZ_asInt()
320 << " A= " << targNucleus->GetA_asInt() << G4endl;
321 DumpState(track,"ApplyYourself",ed);
322 ed << " ApplyYourself failed" << G4endl;
323 G4Exception("G4NeutrinoElectronProcess::PostStepDoIt", "had006",
324 FatalException, ed);
325 }
326 // directions
327
328 G4ThreeVector indir = track.GetMomentumDirection();
329 G4double phi = CLHEP::twopi*G4UniformRand();
330 G4ThreeVector it(0., 0., 1.);
331 G4ThreeVector outdir = result->GetMomentumChange();
332
333 if(verboseLevel>1)
334 {
335 G4cout << "Efin= " << result->GetEnergyChange()
336 << " de= " << result->GetLocalEnergyDeposit()
337 << " nsec= " << result->GetNumberOfSecondaries()
338 << " dir= " << outdir
339 << G4endl;
340 }
341 // energies
342
343 G4double edep = result->GetLocalEnergyDeposit();
344 G4double efinal = result->GetEnergyChange();
345
346 if(efinal < 0.0) { efinal = 0.0; }
347 if(edep < 0.0) { edep = 0.0; }
348
349 // NOTE: Very low energy scatters were causing numerical (FPE) errors
350 // in earlier releases; these limits have not been changed since.
351
352 if(efinal <= lowestEnergy)
353 {
354 edep += efinal;
355 efinal = 0.0;
356 }
357 // primary change
358
360
361 G4TrackStatus status = track.GetTrackStatus();
362
363 if(efinal > 0.0)
364 {
365 outdir.rotate(phi, it);
366 outdir.rotateUz(indir);
368 }
369 else
370 {
371 if( part->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
372 {
373 status = fStopButAlive;
374 }
375 else
376 {
377 status = fStopAndKill;
378 }
380 }
381 //G4cout << "Efinal= " << efinal << " TrackStatus= " << status << G4endl;
382
384
385 // recoil
386
387 if( result->GetNumberOfSecondaries() > 0 )
388 {
389 G4DynamicParticle* p = result->GetSecondary(0)->GetParticle();
390
391 if(p->GetKineticEnergy() > tcut)
392 {
395
396 // G4cout << "recoil " << pdir << G4endl;
397 //!! is not needed for models inheriting G4NeutrinoElectron
398
399 pdir.rotate(phi, it);
400 pdir.rotateUz(indir);
401
402 // G4cout << "recoil rotated " << pdir << G4endl;
403
404 p->SetMomentumDirection(pdir);
405
406 // in elastic scattering time and weight are not changed
407
408 G4Track* t = new G4Track(p, track.GetGlobalTime(),
409 track.GetPosition());
410 t->SetWeight(weight);
411 t->SetTouchableHandle(track.GetTouchableHandle());
413 }
414 else
415 {
416 edep += p->GetKineticEnergy();
417 delete p;
418 }
419 }
422 result->Clear();
423 }
424 return theTotalResult;
425}
@ JustWarning
G4TrackStatus
@ fKillTrackAndSecondaries
@ fSuspend
@ fAlive
@ fStopAndKill
@ fStopButAlive
@ fPostponeToNextEvent
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:24
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:127
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()
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()
const G4String & GetName() const
Definition: G4Material.hh:172
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:691
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
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)
G4TrackStatus GetTrackStatus() const
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:428
G4int verboseLevel
Definition: G4VProcess.hh:360
virtual const G4ThreeVector & GetTranslation(G4int depth=0) const =0
virtual const G4RotationMatrix * GetRotation(G4int depth=0) const =0

◆ PreparePhysicsTable()

void G4NeutrinoElectronProcess::PreparePhysicsTable ( const G4ParticleDefinition part)
virtual

Reimplemented from G4HadronicProcess.

Definition at line 428 of file G4NeutrinoElectronProcess.cc.

429{
430 if(!isInitialised) {
431 isInitialised = true;
432 if(G4Neutron::Neutron() == &part) { lowestEnergy = 1.e-6*eV; }
433 }
435}
void PreparePhysicsTable(const G4ParticleDefinition &) override
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103

◆ ProcessDescription()

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

Reimplemented from G4HadronicProcess.

Definition at line 145 of file G4NeutrinoElectronProcess.cc.

146{
147
148 outFile << "G4NeutrinoElectronProcess handles the scattering of \n"
149 << "neutrino on electrons by invoking the following model(s) and \n"
150 << "cross section(s).\n";
151
152}

◆ SetBiasingFactor()

void G4NeutrinoElectronProcess::SetBiasingFactor ( G4double  bf)

Definition at line 87 of file G4NeutrinoElectronProcess.cc.

88{
89 fNuEleTotXscBias = bf;
90
91 fTotXsc = new G4NeutrinoElectronTotXsc();
92 // fTotXsc->SetBiasingFactor(bf);
93}

Referenced by G4EmExtraPhysics::ConstructProcess().

◆ SetBiasingFactors()

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

Definition at line 97 of file G4NeutrinoElectronProcess.cc.

98{
99 fNuEleCcBias=bfCc;
100 fNuEleNcBias=bfNc;
101
102 fTotXsc = new G4NeutrinoElectronTotXsc();
103 fTotXsc->SetBiasingFactors(bfCc, bfNc);
104}
void SetBiasingFactors(G4double bfCc, G4double bfNc)

Referenced by G4EmExtraPhysics::ConstructProcess().

◆ SetLowestEnergy()

void G4NeutrinoElectronProcess::SetLowestEnergy ( G4double  val)
virtual

Definition at line 438 of file G4NeutrinoElectronProcess.cc.

439{
440 lowestEnergy = val;
441}

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