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

#include <G4VeLowEnergyLoss.hh>

+ Inheritance diagram for G4VeLowEnergyLoss:

Public Member Functions

 G4VeLowEnergyLoss (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VeLowEnergyLoss (G4VeLowEnergyLoss &)
 
virtual ~G4VeLowEnergyLoss ()
 
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 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 MeanLoss, G4double step)
 
- 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 G4PhysicsTableBuildRangeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *theRangeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildLabTimeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *theLabTimeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildProperTimeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *ProperTimeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildRangeCoeffATable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffATable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildRangeCoeffBTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffBTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildRangeCoeffCTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffCTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildInverseRangeTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theRangeCoeffATable, G4PhysicsTable *theRangeCoeffBTable, G4PhysicsTable *theRangeCoeffCTable, G4PhysicsTable *theInverseRangeTable, G4double Tmin, G4double Tmax, G4int nbin)
 

Protected Attributes

const G4MateriallastMaterial
 
G4int imat
 
G4double f1Fluct
 
G4double f2Fluct
 
G4double e1Fluct
 
G4double e2Fluct
 
G4double rateFluct
 
G4double ipotFluct
 
G4double e1LogFluct
 
G4double e2LogFluct
 
G4double ipotLogFluct
 
const G4int nmaxCont1
 
const G4int nmaxCont2
 
- 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 ParticleMass
 
static G4double taulow
 
static G4double tauhigh
 
static G4double ltaulow
 
static G4double ltauhigh
 
static G4double dRoverRange = 20*perCent
 
static G4double finalRange = 200*micrometer
 
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
 

Detailed Description

Definition at line 77 of file G4VeLowEnergyLoss.hh.

Constructor & Destructor Documentation

◆ G4VeLowEnergyLoss() [1/2]

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

◆ G4VeLowEnergyLoss() [2/2]

G4VeLowEnergyLoss::G4VeLowEnergyLoss ( G4VeLowEnergyLoss right)

Definition at line 105 of file G4VeLowEnergyLoss.cc.

107 lastMaterial(0),
108 imat(-1),
109 f1Fluct(0),f2Fluct(0),e1Fluct(0),e2Fluct(0),rateFluct(0),
111 nmaxCont1(4),
112 nmaxCont2(16)
113{;}

◆ ~G4VeLowEnergyLoss()

G4VeLowEnergyLoss::~G4VeLowEnergyLoss ( )
virtual

Definition at line 99 of file G4VeLowEnergyLoss.cc.

100{
101}

Member Function Documentation

◆ AlongStepDoIt()

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

Reimplemented from G4VContinuousDiscreteProcess.

◆ BuildInverseRangeTable()

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

Definition at line 524 of file G4VeLowEnergyLoss.cc.

533{
534 G4bool b;
535
537
538 if(theInverseRangeTable)
539 { theInverseRangeTable->clearAndDestroy();
540 delete theInverseRangeTable; }
541 theInverseRangeTable = new G4PhysicsTable(numOfCouples);
542
543 // loop for materials
544 for (G4int i=0; i<numOfCouples; i++)
545 {
546
547 G4PhysicsVector* pv = (*theRangeTable)[i];
548 size_t nbins = pv->GetVectorLength();
549 G4double elow = pv->GetLowEdgeEnergy(0);
550 G4double ehigh = pv->GetLowEdgeEnergy(nbins-1);
551 G4double rlow = pv->GetValue(elow, b);
552 G4double rhigh = pv->GetValue(ehigh, b);
553
554 //rhigh *= std::exp(std::log(rhigh/rlow)/((G4double)(nbins-1)));
555
556 G4PhysicsLogVector* v = new G4PhysicsLogVector(rlow, rhigh, nbins-1);
557
558 v->PutValue(0,elow);
559 G4double energy1 = elow;
560 G4double range1 = rlow;
561 G4double energy2 = elow;
562 G4double range2 = rlow;
563 size_t ilow = 0;
564 size_t ihigh;
565
566 for (size_t j=1; j<nbins; j++) {
567
568 G4double range = v->GetLowEdgeEnergy(j);
569
570 for (ihigh=ilow+1; ihigh<nbins; ihigh++) {
571 energy2 = pv->GetLowEdgeEnergy(ihigh);
572 range2 = pv->GetValue(energy2, b);
573 if(range2 >= range || ihigh == nbins-1) {
574 ilow = ihigh - 1;
575 energy1 = pv->GetLowEdgeEnergy(ilow);
576 range1 = pv->GetValue(energy1, b);
577 break;
578 }
579 }
580
581 G4double e = std::log(energy1) + std::log(energy2/energy1)*std::log(range/range1)/std::log(range2/range1);
582
583 v->PutValue(j,std::exp(e));
584 }
585 theInverseRangeTable->insert(v);
586
587 }
588 return theInverseRangeTable ;
589}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
void clearAndDestroy()
void insert(G4PhysicsVector *)
size_t GetVectorLength() const
G4double GetValue(G4double theEnergy, G4bool &isOutRange)
virtual G4double GetLowEdgeEnergy(size_t binNumber) const
void PutValue(size_t index, G4double theValue)
static G4ProductionCutsTable * GetProductionCutsTable()

◆ BuildLabTimeTable()

G4PhysicsTable * G4VeLowEnergyLoss::BuildLabTimeTable ( G4PhysicsTable theDEDXTable,
G4PhysicsTable theLabTimeTable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
)
staticprotected

Definition at line 268 of file G4VeLowEnergyLoss.cc.

273{
274
276
277 if(theLabTimeTable)
278 { theLabTimeTable->clearAndDestroy();
279 delete theLabTimeTable; }
280 theLabTimeTable = new G4PhysicsTable(numOfCouples);
281
282
283 for (G4int J=0; J<numOfCouples; J++)
284 {
285 G4PhysicsLogVector* aVector;
286
287 aVector = new G4PhysicsLogVector(lowestKineticEnergy,
288 highestKineticEnergy,TotBin);
289
290 BuildLabTimeVector(theDEDXTable,
291 lowestKineticEnergy,highestKineticEnergy,TotBin,J,aVector);
292 theLabTimeTable->insert(aVector);
293
294
295 }
296 return theLabTimeTable ;
297}

◆ BuildProperTimeTable()

G4PhysicsTable * G4VeLowEnergyLoss::BuildProperTimeTable ( G4PhysicsTable theDEDXTable,
G4PhysicsTable ProperTimeTable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
)
staticprotected

Definition at line 301 of file G4VeLowEnergyLoss.cc.

306{
307
309
310 if(theProperTimeTable)
311 { theProperTimeTable->clearAndDestroy();
312 delete theProperTimeTable; }
313 theProperTimeTable = new G4PhysicsTable(numOfCouples);
314
315
316 for (G4int J=0; J<numOfCouples; J++)
317 {
318 G4PhysicsLogVector* aVector;
319
320 aVector = new G4PhysicsLogVector(lowestKineticEnergy,
321 highestKineticEnergy,TotBin);
322
323 BuildProperTimeVector(theDEDXTable,
324 lowestKineticEnergy,highestKineticEnergy,TotBin,J,aVector);
325 theProperTimeTable->insert(aVector);
326
327
328 }
329 return theProperTimeTable ;
330}

◆ BuildRangeCoeffATable()

G4PhysicsTable * G4VeLowEnergyLoss::BuildRangeCoeffATable ( G4PhysicsTable theRangeTable,
G4PhysicsTable theCoeffATable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
)
staticprotected

Definition at line 649 of file G4VeLowEnergyLoss.cc.

655{
656
658
659 if(theRangeCoeffATable)
660 { theRangeCoeffATable->clearAndDestroy();
661 delete theRangeCoeffATable; }
662 theRangeCoeffATable = new G4PhysicsTable(numOfCouples);
663
664 G4double RTable = std::exp(std::log(highestKineticEnergy/lowestKineticEnergy)/TotBin) ;
665 G4double R2 = RTable*RTable ;
666 G4double R1 = RTable+1.;
667 G4double w = R1*(RTable-1.)*(RTable-1.);
668 G4double w1 = RTable/w , w2 = -RTable*R1/w , w3 = R2/w ;
669 G4double Ti , Tim , Tip , Ri , Rim , Rip , Value ;
670 G4bool isOut;
671
672 // loop for materials
673 for (G4int J=0; J<numOfCouples; J++)
674 {
675 G4int binmax=TotBin ;
676 G4PhysicsLinearVector* aVector =
677 new G4PhysicsLinearVector(0.,binmax, TotBin);
678 Ti = lowestKineticEnergy ;
679 G4PhysicsVector* rangeVector= (*theRangeTable)[J];
680
681 for ( G4int i=0; i<=TotBin; i++)
682 {
683 Ri = rangeVector->GetValue(Ti,isOut) ;
684 if ( i==0 )
685 Rim = 0. ;
686 else
687 {
688 Tim = Ti/RTable ;
689 Rim = rangeVector->GetValue(Tim,isOut);
690 }
691 if ( i==TotBin)
692 Rip = Ri ;
693 else
694 {
695 Tip = Ti*RTable ;
696 Rip = rangeVector->GetValue(Tip,isOut);
697 }
698 Value = (w1*Rip + w2*Ri + w3*Rim)/(Ti*Ti) ;
699
700 aVector->PutValue(i,Value);
701 Ti = RTable*Ti ;
702 }
703
704 theRangeCoeffATable->insert(aVector);
705 }
706 return theRangeCoeffATable ;
707}

◆ BuildRangeCoeffBTable()

G4PhysicsTable * G4VeLowEnergyLoss::BuildRangeCoeffBTable ( G4PhysicsTable theRangeTable,
G4PhysicsTable theCoeffBTable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
)
staticprotected

Definition at line 711 of file G4VeLowEnergyLoss.cc.

717{
718
720
721 if(theRangeCoeffBTable)
722 { theRangeCoeffBTable->clearAndDestroy();
723 delete theRangeCoeffBTable; }
724 theRangeCoeffBTable = new G4PhysicsTable(numOfCouples);
725
726 G4double RTable = std::exp(std::log(highestKineticEnergy/lowestKineticEnergy)/TotBin) ;
727 G4double R2 = RTable*RTable ;
728 G4double R1 = RTable+1.;
729 G4double w = R1*(RTable-1.)*(RTable-1.);
730 G4double w1 = -R1/w , w2 = R1*(R2+1.)/w , w3 = -R2*R1/w ;
731 G4double Ti , Tim , Tip , Ri , Rim , Rip , Value ;
732 G4bool isOut;
733
734 // loop for materials
735 for (G4int J=0; J<numOfCouples; J++)
736 {
737 G4int binmax=TotBin ;
738 G4PhysicsLinearVector* aVector =
739 new G4PhysicsLinearVector(0.,binmax, TotBin);
740 Ti = lowestKineticEnergy ;
741 G4PhysicsVector* rangeVector= (*theRangeTable)[J];
742
743 for ( G4int i=0; i<=TotBin; i++)
744 {
745 Ri = rangeVector->GetValue(Ti,isOut) ;
746 if ( i==0 )
747 Rim = 0. ;
748 else
749 {
750 Tim = Ti/RTable ;
751 Rim = rangeVector->GetValue(Tim,isOut);
752 }
753 if ( i==TotBin)
754 Rip = Ri ;
755 else
756 {
757 Tip = Ti*RTable ;
758 Rip = rangeVector->GetValue(Tip,isOut);
759 }
760 Value = (w1*Rip + w2*Ri + w3*Rim)/Ti;
761
762 aVector->PutValue(i,Value);
763 Ti = RTable*Ti ;
764 }
765 theRangeCoeffBTable->insert(aVector);
766 }
767 return theRangeCoeffBTable ;
768}

◆ BuildRangeCoeffCTable()

G4PhysicsTable * G4VeLowEnergyLoss::BuildRangeCoeffCTable ( G4PhysicsTable theRangeTable,
G4PhysicsTable theCoeffCTable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
)
staticprotected

Definition at line 772 of file G4VeLowEnergyLoss.cc.

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

◆ BuildRangeTable()

G4PhysicsTable * G4VeLowEnergyLoss::BuildRangeTable ( G4PhysicsTable theDEDXTable,
G4PhysicsTable theRangeTable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
)
staticprotected

Definition at line 134 of file G4VeLowEnergyLoss.cc.

140{
141
142 G4int numOfCouples = theDEDXTable->length();
143
144 if(theRangeTable)
145 { theRangeTable->clearAndDestroy();
146 delete theRangeTable; }
147 theRangeTable = new G4PhysicsTable(numOfCouples);
148
149 // loop for materials
150
151 for (G4int J=0; J<numOfCouples; J++)
152 {
153 G4PhysicsLogVector* aVector;
154 aVector = new G4PhysicsLogVector(lowestKineticEnergy,
155 highestKineticEnergy,TotBin);
156 BuildRangeVector(theDEDXTable,lowestKineticEnergy,highestKineticEnergy,
157 TotBin,J,aVector);
158 theRangeTable->insert(aVector);
159 }
160 return theRangeTable ;
161}
size_t length() const

◆ GetContinuousStepLimit()

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

◆ GetLossWithFluct()

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

Definition at line 833 of file G4VeLowEnergyLoss.cc.

839{
840 static const G4double minLoss = 1.*eV ;
841 static const G4double probLim = 0.01 ;
842 static const G4double sumaLim = -std::log(probLim) ;
843 static const G4double alim=10.;
844 static const G4double kappa = 10. ;
845 static const G4double factor = twopi_mc2_rcl2 ;
846 const G4Material* aMaterial = couple->GetMaterial();
847
848 // check if the material has changed ( cache mechanism)
849
850 if (aMaterial != lastMaterial)
851 {
852 lastMaterial = aMaterial;
853 imat = couple->GetIndex();
854 f1Fluct = aMaterial->GetIonisation()->GetF1fluct();
855 f2Fluct = aMaterial->GetIonisation()->GetF2fluct();
856 e1Fluct = aMaterial->GetIonisation()->GetEnergy1fluct();
857 e2Fluct = aMaterial->GetIonisation()->GetEnergy2fluct();
863 }
864 G4double threshold,w1,w2,C,
865 beta2,suma,e0,loss,lossc,w;
866 G4double a1,a2,a3;
867 G4int p1,p2,p3;
868 G4int nb;
869 G4double Corrfac, na,alfa,rfac,namean,sa,alfa1,ea,sea;
870 // G4double dp1;
871 G4double dp3;
872 G4double siga ;
873
874 // shortcut for very very small loss
875 if(MeanLoss < minLoss) return MeanLoss ;
876
877 // get particle data
878 G4double Tkin = aParticle->GetKineticEnergy();
879
880 // G4cout << "MGP -- Fluc Tkin " << Tkin/keV << " keV " << " MeanLoss = " << MeanLoss/keV << G4endl;
881
883 ->GetEnergyCutsVector(1)))[imat];
884 G4double rmass = electron_mass_c2/ParticleMass;
885 G4double tau = Tkin/ParticleMass, tau1 = tau+1., tau2 = tau*(tau+2.);
886 G4double Tm = 2.*electron_mass_c2*tau2/(1.+2.*tau1*rmass+rmass*rmass);
887
888 // G4cout << "MGP Particle mass " << ParticleMass/MeV << " Tm " << Tm << G4endl;
889
890 if(Tm > threshold) Tm = threshold;
891 beta2 = tau2/(tau1*tau1);
892
893 // Gaussian fluctuation ?
894 if(MeanLoss >= kappa*Tm || MeanLoss <= kappa*ipotFluct)
895 {
896 G4double electronDensity = aMaterial->GetElectronDensity() ;
897 siga = std::sqrt(Tm*(1.0-0.5*beta2)*step*
898 factor*electronDensity/beta2) ;
899 do {
900 loss = G4RandGauss::shoot(MeanLoss,siga) ;
901 } while (loss < 0. || loss > 2.0*MeanLoss);
902 return loss ;
903 }
904
905 w1 = Tm/ipotFluct;
906 w2 = std::log(2.*electron_mass_c2*tau2);
907
908 C = MeanLoss*(1.-rateFluct)/(w2-ipotLogFluct-beta2);
909
910 a1 = C*f1Fluct*(w2-e1LogFluct-beta2)/e1Fluct;
911 a2 = C*f2Fluct*(w2-e2LogFluct-beta2)/e2Fluct;
912 a3 = rateFluct*MeanLoss*(Tm-ipotFluct)/(ipotFluct*Tm*std::log(w1));
913
914 suma = a1+a2+a3;
915
916 loss = 0. ;
917
918 if(suma < sumaLim) // very small Step
919 {
920 e0 = aMaterial->GetIonisation()->GetEnergy0fluct();
921 // G4cout << "MGP e0 = " << e0/keV << G4endl;
922
923 if(Tm == ipotFluct)
924 {
925 a3 = MeanLoss/e0;
926
927 if(a3>alim)
928 {
929 siga=std::sqrt(a3) ;
930 p3 = std::max(0,G4int(G4RandGauss::shoot(a3,siga)+0.5));
931 }
932 else p3 = G4Poisson(a3);
933
934 loss = p3*e0 ;
935
936 if(p3 > 0) loss += (1.-2.*G4UniformRand())*e0 ;
937 // G4cout << "MGP very small step " << loss/keV << G4endl;
938 }
939 else
940 {
941 // G4cout << "MGP old Tm = " << Tm << " " << ipotFluct << " " << e0 << G4endl;
942 Tm = Tm-ipotFluct+e0 ;
943
944 // MGP ---- workaround to avoid log argument<0, TO BE CHECKED
945 if (Tm <= 0.)
946 {
947 loss = MeanLoss;
948 p3 = 0;
949 // G4cout << "MGP correction loss = MeanLoss " << loss/keV << G4endl;
950 }
951 else
952 {
953 a3 = MeanLoss*(Tm-e0)/(Tm*e0*std::log(Tm/e0));
954
955 // G4cout << "MGP new Tm = " << Tm << " " << ipotFluct << " " << e0 << " a3= " << a3 << G4endl;
956
957 if(a3>alim)
958 {
959 siga=std::sqrt(a3) ;
960 p3 = std::max(0,G4int(G4RandGauss::shoot(a3,siga)+0.5));
961 }
962 else
963 p3 = G4Poisson(a3);
964 //G4cout << "MGP p3 " << p3 << G4endl;
965
966 }
967
968 if(p3 > 0)
969 {
970 w = (Tm-e0)/Tm ;
971 if(p3 > nmaxCont2)
972 {
973 // G4cout << "MGP dp3 " << dp3 << " p3 " << p3 << " " << nmaxCont2 << G4endl;
974 dp3 = G4double(p3) ;
975 Corrfac = dp3/G4double(nmaxCont2) ;
976 p3 = nmaxCont2 ;
977 }
978 else
979 Corrfac = 1. ;
980
981 for(G4int i=0; i<p3; i++) loss += 1./(1.-w*G4UniformRand()) ;
982 loss *= e0*Corrfac ;
983 // G4cout << "MGP Corrfac = " << Corrfac << " e0 = " << e0/keV << " loss = " << loss/keV << G4endl;
984 }
985 }
986 }
987
988 else // not so small Step
989 {
990 // excitation type 1
991 if(a1>alim)
992 {
993 siga=std::sqrt(a1) ;
994 p1 = std::max(0,int(G4RandGauss::shoot(a1,siga)+0.5));
995 }
996 else
997 p1 = G4Poisson(a1);
998
999 // excitation type 2
1000 if(a2>alim)
1001 {
1002 siga=std::sqrt(a2) ;
1003 p2 = std::max(0,int(G4RandGauss::shoot(a2,siga)+0.5));
1004 }
1005 else
1006 p2 = G4Poisson(a2);
1007
1008 loss = p1*e1Fluct+p2*e2Fluct;
1009
1010 // smearing to avoid unphysical peaks
1011 if(p2 > 0)
1012 loss += (1.-2.*G4UniformRand())*e2Fluct;
1013 else if (loss>0.)
1014 loss += (1.-2.*G4UniformRand())*e1Fluct;
1015
1016 // ionisation .......................................
1017 if(a3 > 0.)
1018 {
1019 if(a3>alim)
1020 {
1021 siga=std::sqrt(a3) ;
1022 p3 = std::max(0,int(G4RandGauss::shoot(a3,siga)+0.5));
1023 }
1024 else
1025 p3 = G4Poisson(a3);
1026
1027 lossc = 0.;
1028 if(p3 > 0)
1029 {
1030 na = 0.;
1031 alfa = 1.;
1032 if (p3 > nmaxCont2)
1033 {
1034 dp3 = G4double(p3);
1035 rfac = dp3/(G4double(nmaxCont2)+dp3);
1036 namean = G4double(p3)*rfac;
1037 sa = G4double(nmaxCont1)*rfac;
1038 na = G4RandGauss::shoot(namean,sa);
1039 if (na > 0.)
1040 {
1041 alfa = w1*G4double(nmaxCont2+p3)/(w1*G4double(nmaxCont2)+G4double(p3));
1042 alfa1 = alfa*std::log(alfa)/(alfa-1.);
1043 ea = na*ipotFluct*alfa1;
1044 sea = ipotFluct*std::sqrt(na*(alfa-alfa1*alfa1));
1045 lossc += G4RandGauss::shoot(ea,sea);
1046 }
1047 }
1048
1049 nb = G4int(G4double(p3)-na);
1050 if (nb > 0)
1051 {
1052 w2 = alfa*ipotFluct;
1053 w = (Tm-w2)/Tm;
1054 for (G4int k=0; k<nb; k++) lossc += w2/(1.-w*G4UniformRand());
1055 }
1056 }
1057
1058 loss += lossc;
1059 }
1060 }
1061
1062 return loss ;
1063}
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:50
#define G4UniformRand()
Definition: Randomize.hh:53
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
static G4double ParticleMass

◆ GetMeanFreePath()

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

◆ PostStepDoIt()

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

Reimplemented from G4VContinuousDiscreteProcess.

◆ SetEnlossFluc()

void G4VeLowEnergyLoss::SetEnlossFluc ( G4bool  value)
static

Definition at line 120 of file G4VeLowEnergyLoss.cc.

121{
122 EnlossFlucFlag = value;
123}
static G4bool EnlossFlucFlag

◆ SetRndmStep()

void G4VeLowEnergyLoss::SetRndmStep ( G4bool  value)
static

Definition at line 115 of file G4VeLowEnergyLoss.cc.

116{
117 rndmStepFlag = value;
118}
static G4bool rndmStepFlag

◆ SetStepFunction()

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

Definition at line 125 of file G4VeLowEnergyLoss.cc.

126{
127 dRoverRange = c1;
128 finalRange = c2;
132}
static G4double c2lim
static G4double finalRange
static G4double c1lim
static G4double dRoverRange
static G4double c3lim

Member Data Documentation

◆ c1lim

G4double G4VeLowEnergyLoss::c1lim = dRoverRange
staticprotected

Definition at line 237 of file G4VeLowEnergyLoss.hh.

Referenced by SetStepFunction().

◆ c2lim

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

Definition at line 237 of file G4VeLowEnergyLoss.hh.

Referenced by SetStepFunction().

◆ c3lim

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

Definition at line 237 of file G4VeLowEnergyLoss.hh.

Referenced by SetStepFunction().

◆ dRoverRange

G4double G4VeLowEnergyLoss::dRoverRange = 20*perCent
staticprotected

Definition at line 234 of file G4VeLowEnergyLoss.hh.

Referenced by SetStepFunction().

◆ e1Fluct

G4double G4VeLowEnergyLoss::e1Fluct
protected

Definition at line 125 of file G4VeLowEnergyLoss.hh.

Referenced by GetLossWithFluct().

◆ e1LogFluct

G4double G4VeLowEnergyLoss::e1LogFluct
protected

Definition at line 126 of file G4VeLowEnergyLoss.hh.

Referenced by GetLossWithFluct().

◆ e2Fluct

G4double G4VeLowEnergyLoss::e2Fluct
protected

Definition at line 125 of file G4VeLowEnergyLoss.hh.

Referenced by GetLossWithFluct().

◆ e2LogFluct

G4double G4VeLowEnergyLoss::e2LogFluct
protected

Definition at line 126 of file G4VeLowEnergyLoss.hh.

Referenced by GetLossWithFluct().

◆ EnlossFlucFlag

G4bool G4VeLowEnergyLoss::EnlossFlucFlag = true
staticprotected

Definition at line 240 of file G4VeLowEnergyLoss.hh.

Referenced by SetEnlossFluc().

◆ f1Fluct

G4double G4VeLowEnergyLoss::f1Fluct
protected

Definition at line 125 of file G4VeLowEnergyLoss.hh.

Referenced by GetLossWithFluct().

◆ f2Fluct

G4double G4VeLowEnergyLoss::f2Fluct
protected

Definition at line 125 of file G4VeLowEnergyLoss.hh.

Referenced by GetLossWithFluct().

◆ finalRange

G4double G4VeLowEnergyLoss::finalRange = 200*micrometer
staticprotected

Definition at line 236 of file G4VeLowEnergyLoss.hh.

Referenced by SetStepFunction().

◆ imat

G4int G4VeLowEnergyLoss::imat
protected

Definition at line 124 of file G4VeLowEnergyLoss.hh.

Referenced by GetLossWithFluct().

◆ ipotFluct

G4double G4VeLowEnergyLoss::ipotFluct
protected

Definition at line 125 of file G4VeLowEnergyLoss.hh.

Referenced by GetLossWithFluct().

◆ ipotLogFluct

G4double G4VeLowEnergyLoss::ipotLogFluct
protected

Definition at line 126 of file G4VeLowEnergyLoss.hh.

Referenced by GetLossWithFluct().

◆ lastMaterial

const G4Material* G4VeLowEnergyLoss::lastMaterial
protected

Definition at line 123 of file G4VeLowEnergyLoss.hh.

Referenced by GetLossWithFluct().

◆ ltauhigh

G4double G4VeLowEnergyLoss::ltauhigh
staticprotected

Definition at line 231 of file G4VeLowEnergyLoss.hh.

◆ ltaulow

G4double G4VeLowEnergyLoss::ltaulow
staticprotected

Definition at line 231 of file G4VeLowEnergyLoss.hh.

◆ nmaxCont1

const G4int G4VeLowEnergyLoss::nmaxCont1
protected

Definition at line 128 of file G4VeLowEnergyLoss.hh.

Referenced by GetLossWithFluct().

◆ nmaxCont2

const G4int G4VeLowEnergyLoss::nmaxCont2
protected

Definition at line 128 of file G4VeLowEnergyLoss.hh.

Referenced by GetLossWithFluct().

◆ ParticleMass

G4double G4VeLowEnergyLoss::ParticleMass
staticprotected

Definition at line 231 of file G4VeLowEnergyLoss.hh.

Referenced by GetLossWithFluct().

◆ rateFluct

G4double G4VeLowEnergyLoss::rateFluct
protected

Definition at line 125 of file G4VeLowEnergyLoss.hh.

Referenced by GetLossWithFluct().

◆ rndmStepFlag

G4bool G4VeLowEnergyLoss::rndmStepFlag = false
staticprotected

Definition at line 239 of file G4VeLowEnergyLoss.hh.

Referenced by SetRndmStep().

◆ tauhigh

G4double G4VeLowEnergyLoss::tauhigh
staticprotected

Definition at line 231 of file G4VeLowEnergyLoss.hh.

◆ taulow

G4double G4VeLowEnergyLoss::taulow
staticprotected

Definition at line 231 of file G4VeLowEnergyLoss.hh.


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