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

#include <G4TauNeutrinoNucleusProcess.hh>

+ Inheritance diagram for G4TauNeutrinoNucleusProcess:

Public Member Functions

 G4TauNeutrinoNucleusProcess (const G4String &anEnvelopeName, const G4String &procName="tau-neutrino-nucleus")
 
 ~G4TauNeutrinoNucleusProcess () 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 SetLowestEnergy (G4double)
 
void ProcessDescription (std::ostream &outFile) const override
 
void SetBiasingFactors (G4double bfCc, G4double bfNc)
 
void SetBiasingFactor (G4double bf)
 
G4TauNeutrinoNucleusProcessoperator= (const G4TauNeutrinoNucleusProcess &right)=delete
 
 G4TauNeutrinoNucleusProcess (const G4TauNeutrinoNucleusProcess &)=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 48 of file G4TauNeutrinoNucleusProcess.hh.

Constructor & Destructor Documentation

◆ G4TauNeutrinoNucleusProcess() [1/2]

G4TauNeutrinoNucleusProcess::G4TauNeutrinoNucleusProcess ( const G4String & anEnvelopeName,
const G4String & procName = "tau-neutrino-nucleus" )

Definition at line 62 of file G4TauNeutrinoNucleusProcess.cc.

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

◆ ~G4TauNeutrinoNucleusProcess()

G4TauNeutrinoNucleusProcess::~G4TauNeutrinoNucleusProcess ( )
overridedefault

◆ G4TauNeutrinoNucleusProcess() [2/2]

G4TauNeutrinoNucleusProcess::G4TauNeutrinoNucleusProcess ( const G4TauNeutrinoNucleusProcess & )
delete

Member Function Documentation

◆ GetMeanFreePath()

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

Reimplemented from G4HadronicProcess.

Definition at line 100 of file G4TauNeutrinoNucleusProcess.cc.

102{
103 //G4cout << "GetMeanFreePath " << aTrack.GetDefinition()->GetParticleName()
104 // << " Ekin= " << aTrack.GetKineticEnergy() << G4endl;
106 G4double totxsc =
108 aTrack.GetMaterial());
109
110 if( rName == fEnvelopeName )
111 {
112 totxsc *= fNuNuclTotXscBias;
113 }
114 G4double res = (totxsc>0.0) ? 1.0/totxsc : DBL_MAX;
115 //G4cout << " xsection= " << totxsc << G4endl;
116 return res;
117}
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=()

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

◆ PostStepDoIt()

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

! is not needed for models inheriting G4TauNeutrinoNucleus

Reimplemented from G4HadronicProcess.

Definition at line 133 of file G4TauNeutrinoNucleusProcess.cc.

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

Reimplemented from G4HadronicProcess.

Definition at line 90 of file G4TauNeutrinoNucleusProcess.cc.

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

◆ ProcessDescription()

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

Reimplemented from G4HadronicProcess.

Definition at line 121 of file G4TauNeutrinoNucleusProcess.cc.

122{
123
124 outFile << "G4TauNeutrinoNucleusProcess handles the inelastic scattering of \n"
125 << "tau-neutrino on nucleus by invoking the following model(s) and \n"
126 << "cross section(s).\n";
127
128}

◆ SetBiasingFactor()

void G4TauNeutrinoNucleusProcess::SetBiasingFactor ( G4double bf)

Definition at line 74 of file G4TauNeutrinoNucleusProcess.cc.

75{
76 fNuNuclTotXscBias = bf;
77}

◆ SetBiasingFactors()

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

Definition at line 81 of file G4TauNeutrinoNucleusProcess.cc.

82{
83 fNuNuclCcBias = bfCc;
84 fNuNuclNcBias = bfNc;
85 fNuNuclTotXscBias = std::max(fNuNuclCcBias, fNuNuclNcBias);
86}

◆ SetLowestEnergy()

void G4TauNeutrinoNucleusProcess::SetLowestEnergy ( G4double val)

Definition at line 392 of file G4TauNeutrinoNucleusProcess.cc.

393{
394 lowestEnergy = val;
395}

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