Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VEnergyLoss Class Referenceabstract

#include <G4VEnergyLoss.hh>

+ Inheritance diagram for G4VEnergyLoss:

Public Member Functions

 G4VEnergyLoss (const G4String &, G4ProcessType aType=fNotDefined)
 
virtual ~G4VEnergyLoss ()
 
virtual G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &Step)=0
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &Step)=0
 
- Public Member Functions inherited from G4VContinuousDiscreteProcess
 G4VContinuousDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VContinuousDiscreteProcess (G4VContinuousDiscreteProcess &)
 
virtual ~G4VContinuousDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Static Public Member Functions

static void SetRndmStep (G4bool value)
 
static void SetEnlossFluc (G4bool value)
 
static void SetSubSec (G4bool value)
 
static void SetMinDeltaCutInRange (G4double value)
 
static void SetStepFunction (G4double c1, G4double c2)
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Protected Member Functions

G4double GetLossWithFluct (const G4DynamicParticle *aParticle, const G4MaterialCutsCouple *couple, G4double ChargeSquare, G4double MeanLoss, G4double step)
 
G4PhysicsTableBuildRangeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *theRangeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
G4PhysicsTableBuildLabTimeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *theLabTimeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
G4PhysicsTableBuildProperTimeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *ProperTimeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
G4PhysicsTableBuildRangeCoeffATable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffATable, G4double Tmin, G4double Tmax, G4int nbin)
 
G4PhysicsTableBuildRangeCoeffBTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffBTable, G4double Tmin, G4double Tmax, G4int nbin)
 
G4PhysicsTableBuildRangeCoeffCTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffCTable, G4double Tmin, G4double Tmax, G4int nbin)
 
G4PhysicsTableBuildInverseRangeTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theRangeCoeffATable, G4PhysicsTable *theRangeCoeffBTable, G4PhysicsTable *theRangeCoeffCTable, G4PhysicsTable *theInverseRangeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
G4bool CutsWhereModified ()
 
- Protected Member Functions inherited from G4VContinuousDiscreteProcess
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
virtual G4double GetContinuousStepLimit (const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
 
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Static Protected Member Functions

static G4bool EqualCutVectors (G4double *vec1, G4double *vec2)
 
static G4doubleCopyCutVectors (G4double *dest, G4double *source)
 

Protected Attributes

G4double ParticleMass
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Static Protected Attributes

static G4double dRoverRange = 20*perCent
 
static G4double finalRange = 1*mm
 
static G4double finalRangeRequested = -1*mm
 
static G4double c1lim = dRoverRange
 
static G4double c2lim = 2.*(1.-dRoverRange)*finalRange
 
static G4double c3lim = -(1.-dRoverRange)*finalRange*finalRange
 
static G4bool rndmStepFlag = false
 
static G4bool EnlossFlucFlag = true
 
static G4bool subSecFlag = false
 
static G4double MinDeltaCutInRange = 0.100*mm
 
static G4doubleMinDeltaEnergy = 0
 
static G4boolLowerLimitForced = 0
 
static G4bool setMinDeltaCutInRange = false
 
static G4EnergyLossMessengerELossMessenger = 0
 

Detailed Description

Definition at line 68 of file G4VEnergyLoss.hh.

Constructor & Destructor Documentation

◆ G4VEnergyLoss()

G4VEnergyLoss::G4VEnergyLoss ( const G4String aName,
G4ProcessType  aType = fNotDefined 
)

Definition at line 81 of file G4VEnergyLoss.cc.

82 : G4VContinuousDiscreteProcess(aName, aType),
83 lastMaterial(NULL),
84 nmaxCont1(4),
85 nmaxCont2(16)
86{
87 if(!ELossMessenger) {
88 G4cout << "### G4VEnergyLoss class is obsolete "
89 << "and will be removed for the next release." << G4endl;
91 }
92
93 imat = 0;
94 f1Fluct = f2Fluct = e1Fluct = e2Fluct = rateFluct = ipotFluct = e1LogFluct
95 = e2LogFluct = ipotLogFluct = taulow = tauhigh = ltaulow = ltauhigh
96 = ParticleMass = 0.0;
97}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4double ParticleMass
static G4EnergyLossMessenger * ELossMessenger

◆ ~G4VEnergyLoss()

G4VEnergyLoss::~G4VEnergyLoss ( )
virtual

Definition at line 101 of file G4VEnergyLoss.cc.

102{}

Member Function Documentation

◆ AlongStepDoIt()

virtual G4VParticleChange * G4VEnergyLoss::AlongStepDoIt ( const G4Track track,
const G4Step Step 
)
pure virtual

Reimplemented from G4VContinuousDiscreteProcess.

◆ BuildInverseRangeTable()

G4PhysicsTable * G4VEnergyLoss::BuildInverseRangeTable ( G4PhysicsTable theRangeTable,
G4PhysicsTable theRangeCoeffATable,
G4PhysicsTable theRangeCoeffBTable,
G4PhysicsTable theRangeCoeffCTable,
G4PhysicsTable theInverseRangeTable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
)
protected

Definition at line 618 of file G4VEnergyLoss.cc.

626{
627 G4double SmallestRange,BiggestRange ;
628 G4bool isOut ;
629 size_t numOfMaterials = theRangeTable->length();
630
631 if(theInverseRangeTable)
632 { theInverseRangeTable->clearAndDestroy();
633 delete theInverseRangeTable; }
634 theInverseRangeTable = new G4PhysicsTable(numOfMaterials);
635
636 // loop for materials
637 for (size_t J=0; J<numOfMaterials; J++)
638 {
639 SmallestRange = (*theRangeTable)(J)->
640 GetValue(LowestKineticEnergy,isOut) ;
641 BiggestRange = (*theRangeTable)(J)->
642 GetValue(HighestKineticEnergy,isOut) ;
643 G4PhysicsLogVector* aVector;
644 aVector = new G4PhysicsLogVector(SmallestRange,
645 BiggestRange,TotBin);
646
647 InvertRangeVector(theRangeTable,
648 theRangeCoeffATable,
649 theRangeCoeffBTable,
650 theRangeCoeffCTable,
651 LowestKineticEnergy,HighestKineticEnergy,TotBin,J, aVector);
652
653 theInverseRangeTable->insert(aVector);
654 }
655 return theInverseRangeTable ;
656}
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
size_t length() const
void clearAndDestroy()
void insert(G4PhysicsVector *)

◆ BuildLabTimeTable()

G4PhysicsTable * G4VEnergyLoss::BuildLabTimeTable ( G4PhysicsTable theDEDXTable,
G4PhysicsTable theLabTimeTable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
)
protected

Definition at line 367 of file G4VEnergyLoss.cc.

372{
373 size_t numOfMaterials = theDEDXTable->length();
374
375 if(theLabTimeTable)
376 { theLabTimeTable->clearAndDestroy();
377 delete theLabTimeTable; }
378 theLabTimeTable = new G4PhysicsTable(numOfMaterials);
379
380
381 for (size_t J=0; J<numOfMaterials; J++)
382 {
383 G4PhysicsLogVector* aVector;
384
385 aVector = new G4PhysicsLogVector(LowestKineticEnergy,
386 HighestKineticEnergy,TotBin);
387
388 BuildLabTimeVector(theDEDXTable,
389 LowestKineticEnergy,HighestKineticEnergy,TotBin,J,aVector);
390 theLabTimeTable->insert(aVector);
391
392
393 }
394 return theLabTimeTable ;
395}

◆ BuildProperTimeTable()

G4PhysicsTable * G4VEnergyLoss::BuildProperTimeTable ( G4PhysicsTable theDEDXTable,
G4PhysicsTable ProperTimeTable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
)
protected

Definition at line 399 of file G4VEnergyLoss.cc.

404{
405 size_t numOfMaterials = theDEDXTable->length();
406
407 if(theProperTimeTable)
408 { theProperTimeTable->clearAndDestroy();
409 delete theProperTimeTable; }
410 theProperTimeTable = new G4PhysicsTable(numOfMaterials);
411
412
413 for (size_t J=0; J<numOfMaterials; J++)
414 {
415 G4PhysicsLogVector* aVector;
416
417 aVector = new G4PhysicsLogVector(LowestKineticEnergy,
418 HighestKineticEnergy,TotBin);
419
420 BuildProperTimeVector(theDEDXTable,
421 LowestKineticEnergy,HighestKineticEnergy,TotBin,J,aVector);
422 theProperTimeTable->insert(aVector);
423
424
425 }
426 return theProperTimeTable ;
427}

◆ BuildRangeCoeffATable()

G4PhysicsTable * G4VEnergyLoss::BuildRangeCoeffATable ( G4PhysicsTable theRangeTable,
G4PhysicsTable theCoeffATable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
)
protected

Definition at line 716 of file G4VEnergyLoss.cc.

722{
723 G4int numOfMaterials = theRangeTable->length();
724
725 if(theRangeCoeffATable)
726 { theRangeCoeffATable->clearAndDestroy();
727 delete theRangeCoeffATable; }
728 theRangeCoeffATable = new G4PhysicsTable(numOfMaterials);
729
730 G4double RTable = std::exp(std::log(HighestKineticEnergy/LowestKineticEnergy)/TotBin) ;
731 G4double R2 = RTable*RTable ;
732 G4double R1 = RTable+1.;
733 G4double w = R1*(RTable-1.)*(RTable-1.);
734 G4double w1 = RTable/w , w2 = -RTable*R1/w , w3 = R2/w ;
735 G4double Ti , Tim , Tip , Ri , Rim , Rip , Value ;
736 G4bool isOut;
737
738 // loop for materials
739 for (G4int J=0; J<numOfMaterials; J++)
740 {
741 G4int binmax=TotBin ;
742 G4PhysicsLinearVector* aVector =
743 new G4PhysicsLinearVector(0.,binmax, TotBin);
744 Ti = LowestKineticEnergy ;
745 G4PhysicsVector* rangeVector= (*theRangeTable)[J];
746
747 for ( G4int i=0; i<TotBin; i++)
748 {
749 Ri = rangeVector->GetValue(Ti,isOut) ;
750 if ( i==0 )
751 Rim = 0. ;
752 else
753 {
754 Tim = Ti/RTable ;
755 Rim = rangeVector->GetValue(Tim,isOut);
756 }
757 if ( i==(TotBin-1))
758 Rip = Ri ;
759 else
760 {
761 Tip = Ti*RTable ;
762 Rip = rangeVector->GetValue(Tip,isOut);
763 }
764 Value = (w1*Rip + w2*Ri + w3*Rim)/(Ti*Ti) ;
765
766 aVector->PutValue(i,Value);
767 Ti = RTable*Ti ;
768 }
769
770 theRangeCoeffATable->insert(aVector);
771 }
772 return theRangeCoeffATable ;
773}
int G4int
Definition: G4Types.hh:66
G4double GetValue(G4double theEnergy, G4bool &isOutRange)
void PutValue(size_t index, G4double theValue)

◆ BuildRangeCoeffBTable()

G4PhysicsTable * G4VEnergyLoss::BuildRangeCoeffBTable ( G4PhysicsTable theRangeTable,
G4PhysicsTable theCoeffBTable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
)
protected

Definition at line 777 of file G4VEnergyLoss.cc.

783{
784 G4int numOfMaterials = theRangeTable->length();
785
786 if(theRangeCoeffBTable)
787 { theRangeCoeffBTable->clearAndDestroy();
788 delete theRangeCoeffBTable; }
789 theRangeCoeffBTable = new G4PhysicsTable(numOfMaterials);
790
791 G4double RTable = std::exp(std::log(HighestKineticEnergy/LowestKineticEnergy)/TotBin) ;
792 G4double R2 = RTable*RTable ;
793 G4double R1 = RTable+1.;
794 G4double w = R1*(RTable-1.)*(RTable-1.);
795 G4double w1 = -R1/w , w2 = R1*(R2+1.)/w , w3 = -R2*R1/w ;
796 G4double Ti , Tim , Tip , Ri , Rim , Rip , Value ;
797 G4bool isOut;
798
799 // loop for materials
800 for (G4int J=0; J<numOfMaterials; J++)
801 {
802 G4int binmax=TotBin ;
803 G4PhysicsLinearVector* aVector =
804 new G4PhysicsLinearVector(0.,binmax, TotBin);
805 Ti = LowestKineticEnergy ;
806 G4PhysicsVector* rangeVector= (*theRangeTable)[J];
807
808 for ( G4int i=0; i<TotBin; i++)
809 {
810 Ri = rangeVector->GetValue(Ti,isOut) ;
811 if ( i==0 )
812 Rim = 0. ;
813 else
814 {
815 Tim = Ti/RTable ;
816 Rim = rangeVector->GetValue(Tim,isOut);
817 }
818 if ( i==(TotBin-1))
819 Rip = Ri ;
820 else
821 {
822 Tip = Ti*RTable ;
823 Rip = rangeVector->GetValue(Tip,isOut);
824 }
825 Value = (w1*Rip + w2*Ri + w3*Rim)/Ti;
826
827 aVector->PutValue(i,Value);
828 Ti = RTable*Ti ;
829 }
830 theRangeCoeffBTable->insert(aVector);
831 }
832 return theRangeCoeffBTable ;
833}

◆ BuildRangeCoeffCTable()

G4PhysicsTable * G4VEnergyLoss::BuildRangeCoeffCTable ( G4PhysicsTable theRangeTable,
G4PhysicsTable theCoeffCTable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
)
protected

Definition at line 837 of file G4VEnergyLoss.cc.

843{
844 G4int numOfMaterials = theRangeTable->length();
845
846 if(theRangeCoeffCTable)
847 { theRangeCoeffCTable->clearAndDestroy();
848 delete theRangeCoeffCTable; }
849 theRangeCoeffCTable = new G4PhysicsTable(numOfMaterials);
850
851 G4double RTable = std::exp(std::log(HighestKineticEnergy/LowestKineticEnergy)/TotBin) ;
852 G4double R2 = RTable*RTable ;
853 G4double R1 = RTable+1.;
854 G4double w = R1*(RTable-1.)*(RTable-1.);
855 G4double w1 = 1./w , w2 = -RTable*R1/w , w3 = RTable*R2/w ;
856 G4double Ti , Tim , Tip , Ri , Rim , Rip , Value ;
857 G4bool isOut;
858
859 // loop for materials
860 for (G4int J=0; J<numOfMaterials; J++)
861 {
862 G4int binmax=TotBin ;
863 G4PhysicsLinearVector* aVector =
864 new G4PhysicsLinearVector(0.,binmax, TotBin);
865 Ti = LowestKineticEnergy ;
866 G4PhysicsVector* rangeVector= (*theRangeTable)[J];
867
868 for ( G4int i=0; i<TotBin; i++)
869 {
870 Ri = rangeVector->GetValue(Ti,isOut) ;
871 if ( i==0 )
872 Rim = 0. ;
873 else
874 {
875 Tim = Ti/RTable ;
876 Rim = rangeVector->GetValue(Tim,isOut);
877 }
878 if ( i==(TotBin-1))
879 Rip = Ri ;
880 else
881 {
882 Tip = Ti*RTable ;
883 Rip = rangeVector->GetValue(Tip,isOut);
884 }
885 Value = w1*Rip + w2*Ri + w3*Rim ;
886
887 aVector->PutValue(i,Value);
888 Ti = RTable*Ti ;
889 }
890 theRangeCoeffCTable->insert(aVector);
891 }
892 return theRangeCoeffCTable ;
893}

◆ BuildRangeTable()

G4PhysicsTable * G4VEnergyLoss::BuildRangeTable ( G4PhysicsTable theDEDXTable,
G4PhysicsTable theRangeTable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
)
protected

Definition at line 133 of file G4VEnergyLoss.cc.

137{
138 size_t numOfMaterials = theDEDXTable->length();
139
140 if(theRangeTable)
141 { theRangeTable->clearAndDestroy();
142 delete theRangeTable; }
143 theRangeTable = new G4PhysicsTable(numOfMaterials);
144
145 // loop for materials
146
147 for (size_t J=0; J<numOfMaterials; J++)
148 {
149 G4PhysicsLogVector* aVector;
150 aVector = new G4PhysicsLogVector(LowestKineticEnergy,
151 HighestKineticEnergy,TotBin);
152
153 BuildRangeVector(theDEDXTable,LowestKineticEnergy,HighestKineticEnergy,
154 TotBin,J,aVector);
155 theRangeTable->insert(aVector);
156 }
157 return theRangeTable ;
158}

◆ CopyCutVectors()

G4double * G4VEnergyLoss::CopyCutVectors ( G4double dest,
G4double source 
)
staticprotected

Definition at line 1131 of file G4VEnergyLoss.cc.

1132{
1133 if ( dest != 0) delete [] dest;
1135 for (size_t j=0; j<G4Material::GetNumberOfMaterials(); j++){
1136 dest[j] = source[j];
1137 }
1138 return dest;
1139}
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:569

◆ CutsWhereModified()

G4bool G4VEnergyLoss::CutsWhereModified ( )
protected

Definition at line 1143 of file G4VEnergyLoss.cc.

1144{
1145 G4bool wasModified = false;
1146 const G4ProductionCutsTable* theCoupleTable=
1148 size_t numOfCouples = theCoupleTable->GetTableSize();
1149
1150 for (size_t j=0; j<numOfCouples; j++){
1151 if (theCoupleTable->GetMaterialCutsCouple(j)->IsRecalcNeeded()) {
1152 wasModified = true;
1153 break;
1154 }
1155 }
1156 return wasModified;
1157}
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ProductionCutsTable * GetProductionCutsTable()

◆ EqualCutVectors()

G4bool G4VEnergyLoss::EqualCutVectors ( G4double vec1,
G4double vec2 
)
staticprotected

Definition at line 1116 of file G4VEnergyLoss.cc.

1117{
1118 if ( (vec1==0 ) || (vec2==0) ) return false;
1119
1120 G4bool flag = true;
1121
1122 for (size_t j=0; flag && j<G4Material::GetNumberOfMaterials(); j++){
1123 flag = (vec1[j] == vec2[j]);
1124 }
1125
1126 return flag;
1127}

◆ GetContinuousStepLimit()

virtual G4double G4VEnergyLoss::GetContinuousStepLimit ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety 
)
pure virtual

◆ GetLossWithFluct()

G4double G4VEnergyLoss::GetLossWithFluct ( const G4DynamicParticle aParticle,
const G4MaterialCutsCouple couple,
G4double  ChargeSquare,
G4double  MeanLoss,
G4double  step 
)
protected

Definition at line 897 of file G4VEnergyLoss.cc.

904{
905 const G4double minLoss = 1.*eV ;
906 const G4double probLim = 0.01 ;
907 const G4double sumaLim = -std::log(probLim) ;
908 const G4double alim=10.;
909 const G4double kappa = 10. ;
910 const G4double factor = twopi_mc2_rcl2 ;
911 const G4Material* aMaterial = couple->GetMaterial();
912
913 // check if the material has changed ( cache mechanism)
914
915 if (aMaterial != lastMaterial)
916 {
917 lastMaterial = aMaterial;
918 imat = couple->GetIndex();
919 f1Fluct = aMaterial->GetIonisation()->GetF1fluct();
920 f2Fluct = aMaterial->GetIonisation()->GetF2fluct();
921 e1Fluct = aMaterial->GetIonisation()->GetEnergy1fluct();
922 e2Fluct = aMaterial->GetIonisation()->GetEnergy2fluct();
923 e1LogFluct = aMaterial->GetIonisation()->GetLogEnergy1fluct();
924 e2LogFluct = aMaterial->GetIonisation()->GetLogEnergy2fluct();
925 rateFluct = aMaterial->GetIonisation()->GetRateionexcfluct();
926 ipotFluct = aMaterial->GetIonisation()->GetMeanExcitationEnergy();
927 ipotLogFluct = aMaterial->GetIonisation()->GetLogMeanExcEnergy();
928 }
929 G4double threshold,w1,w2,C,
930 beta2,suma,e0,loss,lossc ,w,electronDensity;
931 G4double a1,a2,a3;
932 G4double p1,p2,p3 ;
933 G4int nb;
934 G4double Corrfac, na,alfa,rfac,namean,sa,alfa1,ea,sea;
935 G4double dp3;
936 G4double siga ;
937
938 // shortcut for very very small loss
939 if(MeanLoss < minLoss) return MeanLoss ;
940
941 // get particle data
942 G4double Tkin = aParticle->GetKineticEnergy();
943 ParticleMass = aParticle->GetMass() ;
944
946 ->GetEnergyCutsVector(1)))[imat];
947 G4double rmass = electron_mass_c2/ParticleMass;
948 G4double tau = Tkin/ParticleMass, tau1 = tau+1., tau2 = tau*(tau+2.);
949 G4double Tm = 2.*electron_mass_c2*tau2/(1.+2.*tau1*rmass+rmass*rmass);
950
951 if(Tm > threshold) Tm = threshold;
952
953 beta2 = tau2/(tau1*tau1);
954
955 // Gaussian fluctuation ?
956 if(MeanLoss >= kappa*Tm || Tm < kappa*ipotFluct)
957 {
958 electronDensity = aMaterial->GetElectronDensity() ;
959 siga = std::sqrt(Tm*(1.0-0.5*beta2)*step*
960 factor*electronDensity*ChargeSquare/beta2) ;
961 do {
962 loss = G4RandGauss::shoot(MeanLoss,siga) ;
963 } while (loss < 0. || loss > 2.*MeanLoss);
964 return loss;
965 }
966
967 w1 = Tm/ipotFluct;
968 w2 = std::log(2.*electron_mass_c2*tau2);
969
970 C = MeanLoss*(1.-rateFluct)/(w2-ipotLogFluct-beta2);
971
972 a1 = C*f1Fluct*(w2-e1LogFluct-beta2)/e1Fluct;
973 a2 = C*f2Fluct*(w2-e2LogFluct-beta2)/e2Fluct;
974 a3 = rateFluct*MeanLoss*(Tm-ipotFluct)/(ipotFluct*Tm*std::log(w1));
975 if(a1 < 0.) a1 = 0.;
976 if(a2 < 0.) a2 = 0.;
977 if(a3 < 0.) a3 = 0.;
978
979 suma = a1+a2+a3;
980
981 loss = 0. ;
982
983 if(suma < sumaLim) // very small Step
984 {
985 e0 = aMaterial->GetIonisation()->GetEnergy0fluct();
986
987 if(Tm == ipotFluct)
988 {
989 a3 = MeanLoss/e0;
990
991 if(a3>alim)
992 {
993 siga=std::sqrt(a3) ;
994 p3 = std::max(0.,G4RandGauss::shoot(a3,siga)+0.5);
995 }
996 else
997 p3 = G4float(G4Poisson(a3));
998
999 loss = p3*e0 ;
1000
1001 if(p3 > 0)
1002 loss += (1.-2.*G4UniformRand())*e0 ;
1003
1004 }
1005 else
1006 {
1007 Tm = Tm-ipotFluct+e0 ;
1008 a3 = MeanLoss*(Tm-e0)/(Tm*e0*std::log(Tm/e0));
1009
1010 if(a3>alim)
1011 {
1012 siga=std::sqrt(a3) ;
1013 p3 = std::max(0.,G4RandGauss::shoot(a3,siga)+0.5);
1014 }
1015 else
1016 p3 = G4float(G4Poisson(a3));
1017
1018 if(p3 > 0)
1019 {
1020 w = (Tm-e0)/Tm ;
1021 if(p3 > G4float(nmaxCont2))
1022 {
1023 dp3 = p3 ;
1024 Corrfac = dp3/G4float(nmaxCont2) ;
1025 p3 = G4float(nmaxCont2) ;
1026 }
1027 else
1028 Corrfac = 1. ;
1029
1030 for(G4int i=0; i<G4int(p3); i++) loss += 1./(1.-w*G4UniformRand()) ;
1031 loss *= e0*Corrfac ;
1032 }
1033 }
1034 }
1035
1036 else // not so small Step
1037 {
1038 // excitation type 1
1039 if(a1>alim)
1040 {
1041 siga=std::sqrt(a1) ;
1042 p1 = std::max(0.,G4RandGauss::shoot(a1,siga)+0.5);
1043 }
1044 else
1045 p1 = G4float(G4Poisson(a1));
1046
1047 // excitation type 2
1048 if(a2>alim)
1049 {
1050 siga=std::sqrt(a2) ;
1051 p2 = std::max(0.,G4RandGauss::shoot(a2,siga)+0.5);
1052 }
1053 else
1054 p2 = G4float(G4Poisson(a2));
1055
1056 loss = p1*e1Fluct+p2*e2Fluct;
1057 // smearing to avoid unphysical peaks
1058 if(p2 > 0)
1059 loss += (1.-2.*G4UniformRand())*e2Fluct;
1060 else if (loss>0.)
1061 loss += (1.-2.*G4UniformRand())*e1Fluct;
1062
1063 // ionisation .......................................
1064 if(a3 > 0.)
1065 {
1066 if(a3>alim)
1067 {
1068 siga=std::sqrt(a3) ;
1069 p3 = std::max(0.,G4RandGauss::shoot(a3,siga)+0.5);
1070 }
1071 else
1072 p3 = G4float(G4Poisson(a3));
1073
1074 lossc = 0.;
1075 if(p3 > 0)
1076 {
1077 na = 0.;
1078 alfa = 1.;
1079 if (p3 > G4float(nmaxCont2))
1080 {
1081 dp3 = p3;
1082 rfac = dp3/(G4float(nmaxCont2)+dp3);
1083 namean = p3*rfac;
1084 sa = G4float(nmaxCont1)*rfac;
1085 na = G4RandGauss::shoot(namean,sa);
1086 if (na > 0.)
1087 {
1088 alfa = w1*(G4float(nmaxCont2)+p3)/(w1*G4float(nmaxCont2)+p3);
1089 alfa1 = alfa*std::log(alfa)/(alfa-1.);
1090 ea = na*ipotFluct*alfa1;
1091 sea = ipotFluct*std::sqrt(na*(alfa-alfa1*alfa1));
1092 lossc += G4RandGauss::shoot(ea,sea);
1093 }
1094 }
1095
1096 nb = G4int(p3-na);
1097 if (nb > 0)
1098 {
1099 w2 = alfa*ipotFluct;
1100 w = (Tm-w2)/Tm;
1101 for (G4int k=0; k<nb; k++) lossc += w2/(1.-w*G4UniformRand());
1102 }
1103 }
1104 loss += lossc;
1105 }
1106 }
1107
1108 if( loss < 0.)
1109 loss = 0.;
1110
1111 return loss ;
1112}
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:50
float G4float
Definition: G4Types.hh:65
#define G4UniformRand()
Definition: Randomize.hh:53
G4double GetMass() const
G4double GetKineticEnergy() const
G4double GetF2fluct() const
G4double GetLogEnergy1fluct() const
G4double GetLogEnergy2fluct() const
G4double GetMeanExcitationEnergy() const
G4double GetF1fluct() const
G4double GetEnergy2fluct() const
G4double GetEnergy1fluct() const
G4double GetRateionexcfluct() const
G4double GetEnergy0fluct() const
G4double GetLogMeanExcEnergy() const
const G4Material * GetMaterial() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:225
G4double GetElectronDensity() const
Definition: G4Material.hh:216

◆ GetMeanFreePath()

virtual G4double G4VEnergyLoss::GetMeanFreePath ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
pure virtual

◆ PostStepDoIt()

virtual G4VParticleChange * G4VEnergyLoss::PostStepDoIt ( const G4Track track,
const G4Step Step 
)
pure virtual

Reimplemented from G4VContinuousDiscreteProcess.

◆ SetEnlossFluc()

void G4VEnergyLoss::SetEnlossFluc ( G4bool  value)
static

Definition at line 110 of file G4VEnergyLoss.cc.

110{EnlossFlucFlag = value;}
static G4bool EnlossFlucFlag

◆ SetMinDeltaCutInRange()

void G4VEnergyLoss::SetMinDeltaCutInRange ( G4double  value)
static

Definition at line 118 of file G4VEnergyLoss.cc.

static G4bool setMinDeltaCutInRange
static G4double MinDeltaCutInRange

◆ SetRndmStep()

void G4VEnergyLoss::SetRndmStep ( G4bool  value)
static

Definition at line 106 of file G4VEnergyLoss.cc.

106{rndmStepFlag = value;}
static G4bool rndmStepFlag

◆ SetStepFunction()

void G4VEnergyLoss::SetStepFunction ( G4double  c1,
G4double  c2 
)
static

Definition at line 123 of file G4VEnergyLoss.cc.

124{
129}
static G4double dRoverRange
static G4double c3lim
static G4double c1lim
static G4double finalRangeRequested
static G4double c2lim

◆ SetSubSec()

void G4VEnergyLoss::SetSubSec ( G4bool  value)
static

Definition at line 114 of file G4VEnergyLoss.cc.

114{subSecFlag = value;}
static G4bool subSecFlag

Member Data Documentation

◆ c1lim

G4double G4VEnergyLoss::c1lim = dRoverRange
staticprotected

Definition at line 231 of file G4VEnergyLoss.hh.

Referenced by SetStepFunction().

◆ c2lim

G4double G4VEnergyLoss::c2lim = 2.*(1.-dRoverRange)*finalRange
staticprotected

Definition at line 231 of file G4VEnergyLoss.hh.

Referenced by SetStepFunction().

◆ c3lim

G4double G4VEnergyLoss::c3lim = -(1.-dRoverRange)*finalRange*finalRange
staticprotected

Definition at line 231 of file G4VEnergyLoss.hh.

Referenced by SetStepFunction().

◆ dRoverRange

G4double G4VEnergyLoss::dRoverRange = 20*perCent
staticprotected

Definition at line 227 of file G4VEnergyLoss.hh.

Referenced by SetStepFunction().

◆ ELossMessenger

G4EnergyLossMessenger * G4VEnergyLoss::ELossMessenger = 0
staticprotected

Definition at line 243 of file G4VEnergyLoss.hh.

Referenced by G4VEnergyLoss().

◆ EnlossFlucFlag

G4bool G4VEnergyLoss::EnlossFlucFlag = true
staticprotected

Definition at line 234 of file G4VEnergyLoss.hh.

Referenced by SetEnlossFluc().

◆ finalRange

G4double G4VEnergyLoss::finalRange = 1*mm
staticprotected

Definition at line 229 of file G4VEnergyLoss.hh.

◆ finalRangeRequested

G4double G4VEnergyLoss::finalRangeRequested = -1*mm
staticprotected

Definition at line 230 of file G4VEnergyLoss.hh.

Referenced by SetStepFunction().

◆ LowerLimitForced

G4bool * G4VEnergyLoss::LowerLimitForced = 0
staticprotected

Definition at line 239 of file G4VEnergyLoss.hh.

◆ MinDeltaCutInRange

G4double G4VEnergyLoss::MinDeltaCutInRange = 0.100*mm
staticprotected

Definition at line 237 of file G4VEnergyLoss.hh.

Referenced by SetMinDeltaCutInRange().

◆ MinDeltaEnergy

G4double * G4VEnergyLoss::MinDeltaEnergy = 0
staticprotected

Definition at line 238 of file G4VEnergyLoss.hh.

◆ ParticleMass

G4double G4VEnergyLoss::ParticleMass
protected

Definition at line 175 of file G4VEnergyLoss.hh.

Referenced by G4VEnergyLoss(), and GetLossWithFluct().

◆ rndmStepFlag

G4bool G4VEnergyLoss::rndmStepFlag = false
staticprotected

Definition at line 233 of file G4VEnergyLoss.hh.

Referenced by SetRndmStep().

◆ setMinDeltaCutInRange

G4bool G4VEnergyLoss::setMinDeltaCutInRange = false
staticprotected

Definition at line 241 of file G4VEnergyLoss.hh.

Referenced by SetMinDeltaCutInRange().

◆ subSecFlag

G4bool G4VEnergyLoss::subSecFlag = false
staticprotected

Definition at line 235 of file G4VEnergyLoss.hh.

Referenced by SetSubSec().


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