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

#include <G4ElNeutrinoNucleusProcess.hh>

+ Inheritance diagram for G4ElNeutrinoNucleusProcess:

Public Member Functions

 G4ElNeutrinoNucleusProcess (G4String anEnvelopeName, const G4String &procName="mu-neutrino-nucleus")
 
virtual ~G4ElNeutrinoNucleusProcess ()
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
void PreparePhysicsTable (const G4ParticleDefinition &) override
 
virtual void SetLowestEnergy (G4double)
 
void ProcessDescription (std::ostream &outFile) const override
 
void SetBiasingFactors (G4double bfCc, G4double bfNc)
 
void SetBiasingFactor (G4double bf)
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *) override
 
- 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 G4ElNeutrinoNucleusProcess.hh.

Constructor & Destructor Documentation

◆ G4ElNeutrinoNucleusProcess()

G4ElNeutrinoNucleusProcess::G4ElNeutrinoNucleusProcess ( G4String  anEnvelopeName,
const G4String procName = "mu-neutrino-nucleus" 
)

Definition at line 67 of file G4ElNeutrinoNucleusProcess.cc.

68 : G4HadronicProcess( pName, fHadronInelastic ), isInitialised(false), fBiased(true) // fHadronElastic???
69{
70 lowestEnergy = 1.*keV;
71 fEnvelope = nullptr;
72 fEnvelopeName = anEnvelopeName;
73 fTotXsc = nullptr; // new G4ElNeutrinoNucleusTotXsc();
74 fNuNuclCcBias=1.;
75 fNuNuclNcBias=1.;
76 fNuNuclTotXscBias=1.;
78 safetyHelper->InitialiseHelper();
79}
@ fHadronInelastic
void InitialiseHelper()
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const

◆ ~G4ElNeutrinoNucleusProcess()

G4ElNeutrinoNucleusProcess::~G4ElNeutrinoNucleusProcess ( )
virtual

Definition at line 81 of file G4ElNeutrinoNucleusProcess.cc.

82{
83 if( fTotXsc ) delete fTotXsc;
84}

Member Function Documentation

◆ GetMeanFreePath()

G4double G4ElNeutrinoNucleusProcess::GetMeanFreePath ( const G4Track aTrack,
G4double  ,
G4ForceCondition  
)
overridevirtual

Reimplemented from G4HadronicProcess.

Definition at line 109 of file G4ElNeutrinoNucleusProcess.cc.

111{
112 //G4cout << "GetMeanFreePath " << aTrack.GetDefinition()->GetParticleName()
113 // << " Ekin= " << aTrack.GetKineticEnergy() << G4endl;
115 G4double totxsc(0.);
116
117 if( rName == fEnvelopeName && fNuNuclTotXscBias > 1.)
118 {
119 totxsc = fNuNuclTotXscBias*
121 aTrack.GetMaterial());
122 }
123 else
124 {
126 aTrack.GetMaterial());
127 }
128 G4double res = (totxsc>0.0) ? 1.0/totxsc : DBL_MAX;
129 //G4cout << " xsection= " << totxsc << G4endl;
130 return res;
131}
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

◆ PostStepDoIt()

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

! is not needed for models inheriting G4ElNeutrinoNucleus

Reimplemented from G4HadronicProcess.

Definition at line 147 of file G4ElNeutrinoNucleusProcess.cc.

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

Reimplemented from G4HadronicProcess.

Definition at line 403 of file G4ElNeutrinoNucleusProcess.cc.

404{
405 if(!isInitialised) {
406 isInitialised = true;
407 // if(G4Neutron::Neutron() == &part) { lowestEnergy = 1.e-6*eV; }
408 }
410}
void PreparePhysicsTable(const G4ParticleDefinition &) override

◆ ProcessDescription()

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

Reimplemented from G4HadronicProcess.

Definition at line 135 of file G4ElNeutrinoNucleusProcess.cc.

136{
137
138 outFile << "G4ElNeutrinoNucleusProcess handles the scattering of \n"
139 << "neutrino on electrons by invoking the following model(s) and \n"
140 << "cross section(s).\n";
141
142}

◆ SetBiasingFactor()

void G4ElNeutrinoNucleusProcess::SetBiasingFactor ( G4double  bf)

Definition at line 88 of file G4ElNeutrinoNucleusProcess.cc.

89{
90 fNuNuclTotXscBias = bf;
91
92 fTotXsc = new G4ElNeutrinoNucleusTotXsc();
93 fTotXsc->SetBiasingFactor(bf);
94}

Referenced by G4EmExtraPhysics::ConstructProcess().

◆ SetBiasingFactors()

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

Definition at line 98 of file G4ElNeutrinoNucleusProcess.cc.

99{
100 fNuNuclCcBias=bfCc;
101 fNuNuclNcBias=bfNc;
102
103 fTotXsc = new G4ElNeutrinoNucleusTotXsc();
104 // fTotXsc->SetBiasingFactors(bfCc, bfNc);
105}

◆ SetLowestEnergy()

void G4ElNeutrinoNucleusProcess::SetLowestEnergy ( G4double  val)
virtual

Definition at line 413 of file G4ElNeutrinoNucleusProcess.cc.

414{
415 lowestEnergy = val;
416}

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