Geant4 10.7.0
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)
 
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 &)
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *) override
 
const G4NucleusGetTargetNucleus () const
 
const G4IsotopeGetTargetIsotope ()
 
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 ()
 
- 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 &)
 
- 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 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)
 
G4NucleusGetTargetNucleusPointer ()
 
void DumpState (const G4Track &, const G4String &, G4ExceptionDescription &)
 
G4HadronicInteractionGetHadronicInteraction () const
 
G4double GetLastCrossSection ()
 
void FillResult (G4HadFinalState *aR, const G4Track &aT)
 
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
 
G4double fWeight
 
G4int epReportLevel
 
- 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 69 of file G4ElNeutrinoNucleusProcess.cc.

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

◆ ~G4ElNeutrinoNucleusProcess()

G4ElNeutrinoNucleusProcess::~G4ElNeutrinoNucleusProcess ( )
virtual

Definition at line 84 of file G4ElNeutrinoNucleusProcess.cc.

85{
86 if( fTotXsc ) delete fTotXsc;
87}

Member Function Documentation

◆ GetMeanFreePath()

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

Reimplemented from G4HadronicProcess.

Definition at line 112 of file G4ElNeutrinoNucleusProcess.cc.

114{
115 //G4cout << "GetMeanFreePath " << aTrack.GetDefinition()->GetParticleName()
116 // << " Ekin= " << aTrack.GetKineticEnergy() << G4endl;
118 G4double totxsc(0.);
119
120 if( rName == fEnvelopeName && fNuNuclTotXscBias > 1.)
121 {
122 totxsc = fNuNuclTotXscBias*
124 aTrack.GetMaterial());
125 }
126 else
127 {
129 aTrack.GetMaterial());
130 }
131 G4double res = (totxsc>0.0) ? 1.0/totxsc : DBL_MAX;
132 //G4cout << " xsection= " << totxsc << G4endl;
133 return res;
134}
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 150 of file G4ElNeutrinoNucleusProcess.cc.

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

◆ PreparePhysicsTable()

void G4ElNeutrinoNucleusProcess::PreparePhysicsTable ( const G4ParticleDefinition part)
overridevirtual

Reimplemented from G4HadronicProcess.

Definition at line 421 of file G4ElNeutrinoNucleusProcess.cc.

422{
423 if(!isInitialised) {
424 isInitialised = true;
425 // if(G4Neutron::Neutron() == &part) { lowestEnergy = 1.e-6*eV; }
426 }
428}
void PreparePhysicsTable(const G4ParticleDefinition &) override

◆ ProcessDescription()

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

Reimplemented from G4HadronicProcess.

Definition at line 138 of file G4ElNeutrinoNucleusProcess.cc.

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

◆ SetBiasingFactor()

void G4ElNeutrinoNucleusProcess::SetBiasingFactor ( G4double  bf)

Definition at line 91 of file G4ElNeutrinoNucleusProcess.cc.

92{
93 fNuNuclTotXscBias = bf;
94
95 fTotXsc = new G4ElNeutrinoNucleusTotXsc();
96 fTotXsc->SetBiasingFactor(bf);
97}

Referenced by G4EmExtraPhysics::ConstructProcess().

◆ SetBiasingFactors()

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

Definition at line 101 of file G4ElNeutrinoNucleusProcess.cc.

102{
103 fNuNuclCcBias=bfCc;
104 fNuNuclNcBias=bfNc;
105
106 fTotXsc = new G4ElNeutrinoNucleusTotXsc();
107 // fTotXsc->SetBiasingFactors(bfCc, bfNc);
108}

◆ SetLowestEnergy()

void G4ElNeutrinoNucleusProcess::SetLowestEnergy ( G4double  val)
virtual

Definition at line 431 of file G4ElNeutrinoNucleusProcess.cc.

432{
433 lowestEnergy = val;
434}

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