Geant4  10.02.p03
G4ePolarizedIonisation Class Reference

#include <G4ePolarizedIonisation.hh>

Inheritance diagram for G4ePolarizedIonisation:
Collaboration diagram for G4ePolarizedIonisation:

Public Member Functions

 G4ePolarizedIonisation (const G4String &name="pol-eIoni")
 
virtual ~G4ePolarizedIonisation ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p)
 
virtual void PrintInfo ()
 
- Public Member Functions inherited from G4VEnergyLossProcess
 G4VEnergyLossProcess (const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
 
virtual ~G4VEnergyLossProcess ()
 
virtual void ProcessDescription (std::ostream &outFile) const
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
G4PhysicsTableBuildDEDXTable (G4EmTableType tType=fRestricted)
 
G4PhysicsTableBuildLambdaTable (G4EmTableType tType=fRestricted)
 
void PrintInfoDefinition (const G4ParticleDefinition &part)
 
virtual void StartTracking (G4Track *)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
virtual G4VParticleChange * AlongStepDoIt (const G4Track &, const G4Step &)
 
G4double SampleSubCutSecondaries (std::vector< G4Track *> &, const G4Step &, G4VEmModel *model, G4int matIdx)
 
virtual G4VParticleChange * PostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii)
 
G4double GetDEDXDispersion (const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)
 
G4double CrossSectionPerVolume (G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double MeanFreePath (const G4Track &track)
 
G4double ContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
 
G4VEmModelSelectModelForMaterial (G4double kinEnergy, size_t &idx) const
 
void AddEmModel (G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
void SetEmModel (G4VEmModel *, G4int index=1)
 
G4VEmModelEmModel (G4int index=1) const
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) const
 
G4int NumberOfModels () const
 
void SetFluctModel (G4VEmFluctuationModel *)
 
G4VEmFluctuationModelFluctModel ()
 
void SetBaseParticle (const G4ParticleDefinition *p)
 
const G4ParticleDefinitionParticle () const
 
const G4ParticleDefinitionBaseParticle () const
 
const G4ParticleDefinitionSecondaryParticle () const
 
void ActivateSubCutoff (G4bool val, const G4Region *region=0)
 
void SetCrossSectionBiasingFactor (G4double f, G4bool flag=true)
 
void ActivateForcedInteraction (G4double length=0.0, const G4String &region="", G4bool flag=true)
 
void ActivateSecondaryBiasing (const G4String &region, G4double factor, G4double energyLimit)
 
void AddCollaborativeProcess (G4VEnergyLossProcess *)
 
void SetLossFluctuations (G4bool val)
 
void SetIntegral (G4bool val)
 
G4bool IsIntegral () const
 
void SetIonisation (G4bool val)
 
G4bool IsIonisationProcess () const
 
void SetLinearLossLimit (G4double val)
 
void SetStepFunction (G4double v1, G4double v2)
 
void SetLowestEnergyLimit (G4double)
 
G4int NumberOfSubCutoffRegions () const
 
void SetDEDXTable (G4PhysicsTable *p, G4EmTableType tType)
 
void SetCSDARangeTable (G4PhysicsTable *pRange)
 
void SetRangeTableForLoss (G4PhysicsTable *p)
 
void SetSecondaryRangeTable (G4PhysicsTable *p)
 
void SetInverseRangeTable (G4PhysicsTable *p)
 
void SetLambdaTable (G4PhysicsTable *p)
 
void SetSubLambdaTable (G4PhysicsTable *p)
 
void SetDEDXBinning (G4int nbins)
 
void SetMinKinEnergy (G4double e)
 
G4double MinKinEnergy () const
 
void SetMaxKinEnergy (G4double e)
 
G4double MaxKinEnergy () const
 
G4double CrossSectionBiasingFactor () const
 
G4double GetDEDX (G4double &kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetDEDXForSubsec (G4double &kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetRange (G4double &kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetCSDARange (G4double &kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetRangeForLoss (G4double &kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetKineticEnergy (G4double &range, const G4MaterialCutsCouple *)
 
G4double GetLambda (G4double &kineticEnergy, const G4MaterialCutsCouple *)
 
G4bool TablesAreBuilt () const
 
G4PhysicsTableDEDXTable () const
 
G4PhysicsTableDEDXTableForSubsec () const
 
G4PhysicsTableDEDXunRestrictedTable () const
 
G4PhysicsTableIonisationTable () const
 
G4PhysicsTableIonisationTableForSubsec () const
 
G4PhysicsTableCSDARangeTable () const
 
G4PhysicsTableSecondaryRangeTable () const
 
G4PhysicsTableRangeTableForLoss () const
 
G4PhysicsTableInverseRangeTable () const
 
G4PhysicsTableLambdaTable () const
 
G4PhysicsTableSubLambdaTable () const
 
const G4ElementGetCurrentElement () const
 
void SetDynamicMassCharge (G4double massratio, G4double charge2ratio)
 
- Public Member Functions inherited from G4VContinuousDiscreteProcess
 G4VContinuousDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VContinuousDiscreteProcess (G4VContinuousDiscreteProcess &)
 
virtual ~G4VContinuousDiscreteProcess ()
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChange * AtRestDoIt (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
 
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)
 
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 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 &)
 

Protected Member Functions

virtual void InitialiseEnergyLossProcess (const G4ParticleDefinition *, const G4ParticleDefinition *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *, G4double cut)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
const G4ParticleDefinitionDefineBaseParticle (const G4ParticleDefinition *p)
 
- Protected Member Functions inherited from G4VEnergyLossProcess
virtual G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *, G4double cut)
 
size_t CurrentMaterialCutsCoupleIndex () const
 
void SelectModel (G4double kinEnergy)
 
void SetParticle (const G4ParticleDefinition *p)
 
void SetSecondaryParticle (const G4ParticleDefinition *p)
 
- Protected Member Functions inherited from G4VContinuousDiscreteProcess
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Private Member Functions

void CleanTables ()
 
void BuildAsymmetryTables (const G4ParticleDefinition &part)
 
G4double ComputeAsymmetry (G4double energy, const G4MaterialCutsCouple *couple, const G4ParticleDefinition &particle, G4double cut, G4double &tasm)
 
G4double ComputeSaturationFactor (const G4Track &aTrack)
 
G4ePolarizedIonisationoperator= (const G4ePolarizedIonisation &right)
 
 G4ePolarizedIonisation (const G4ePolarizedIonisation &)
 

Private Attributes

G4ParticleDefinitiontheElectron
 
G4VEmFluctuationModelflucModel
 
G4PolarizedMollerBhabhaModelemModel
 
G4bool isElectron
 
G4bool isInitialised
 
G4ThreeVector theTargetPolarization
 
G4PhysicsTabletheAsymmetryTable
 
G4PhysicsTabletheTransverseAsymmetryTable
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VEnergyLossProcess
G4ParticleChangeForLoss fParticleChange
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChange * pParticleChange
 
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
 

Detailed Description

Definition at line 66 of file G4ePolarizedIonisation.hh.

Constructor & Destructor Documentation

◆ G4ePolarizedIonisation() [1/2]

G4ePolarizedIonisation::G4ePolarizedIonisation ( const G4String name = "pol-eIoni")

Definition at line 70 of file G4ePolarizedIonisation.cc.

71  : G4VEnergyLossProcess(name),
73  isElectron(true),
74  isInitialised(false),
75  theTargetPolarization(0.,0.,0.),
76  theAsymmetryTable(nullptr),
78 {
79  verboseLevel=0;
82  flucModel = nullptr;
83  emModel = nullptr;
84 }
G4VEnergyLossProcess(const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
G4int verboseLevel
Definition: G4VProcess.hh:368
G4PolarizedMollerBhabhaModel * emModel
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
void SetSecondaryParticle(const G4ParticleDefinition *p)
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4ParticleDefinition * theElectron
G4VEmFluctuationModel * flucModel
G4PhysicsTable * theTransverseAsymmetryTable
Here is the call graph for this function:

◆ ~G4ePolarizedIonisation()

G4ePolarizedIonisation::~G4ePolarizedIonisation ( )
virtual

Definition at line 88 of file G4ePolarizedIonisation.cc.

89 {
90  CleanTables();
91 }
Here is the call graph for this function:

◆ G4ePolarizedIonisation() [2/2]

G4ePolarizedIonisation::G4ePolarizedIonisation ( const G4ePolarizedIonisation )
private

Member Function Documentation

◆ BuildAsymmetryTables()

void G4ePolarizedIonisation::BuildAsymmetryTables ( const G4ParticleDefinition part)
private

Definition at line 287 of file G4ePolarizedIonisation.cc.

289 {
290  // cleanup old, initialise new table
291  CleanTables();
296 
297  const G4ProductionCutsTable* theCoupleTable=
299  size_t numOfCouples = theCoupleTable->GetTableSize();
300 
301  for (size_t j=0 ; j < numOfCouples; j++ ) {
302  // get cut value
303  const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(j);
304 
305  G4double cut = (*theCoupleTable->GetEnergyCutsVector(1))[j];
306 
307  //create physics vectors then fill it (same parameters as lambda vector)
308  G4PhysicsVector * ptrVectorA = LambdaPhysicsVector(couple,cut);
309  G4PhysicsVector * ptrVectorB = LambdaPhysicsVector(couple,cut);
310  size_t bins = ptrVectorA->GetVectorLength();
311 
312  for (size_t i = 0 ; i < bins ; i++ ) {
313  G4double lowEdgeEnergy = ptrVectorA->Energy(i);
314  G4double tasm=0.;
315  G4double asym = ComputeAsymmetry(lowEdgeEnergy, couple, part, cut, tasm);
316  ptrVectorA->PutValue(i,asym);
317  ptrVectorB->PutValue(i,tasm);
318  }
319  theAsymmetryTable->insertAt( j , ptrVectorA ) ;
320  theTransverseAsymmetryTable->insertAt( j , ptrVectorB ) ;
321  }
322 }
G4double ComputeAsymmetry(G4double energy, const G4MaterialCutsCouple *couple, const G4ParticleDefinition &particle, G4double cut, G4double &tasm)
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *, G4double cut)
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
void PutValue(size_t index, G4double theValue)
size_t GetVectorLength() const
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
void insertAt(size_t, G4PhysicsVector *)
G4double Energy(size_t index) const
double G4double
Definition: G4Types.hh:76
G4PhysicsTable * theTransverseAsymmetryTable
Here is the call graph for this function:
Here is the caller graph for this function:

◆ BuildPhysicsTable()

void G4ePolarizedIonisation::BuildPhysicsTable ( const G4ParticleDefinition part)
protectedvirtual

Reimplemented from G4VEnergyLossProcess.

Definition at line 273 of file G4ePolarizedIonisation.cc.

275 {
276  // *** build DEDX and (unpolarized) cross section tables
278  G4bool master = true;
279  const G4ePolarizedIonisation* masterProcess =
280  static_cast<const G4ePolarizedIonisation*>(GetMasterProcess());
281  if(masterProcess && masterProcess != this) { master = false; }
282  if(master) { BuildAsymmetryTables(part); }
283 }
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:538
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
void BuildAsymmetryTables(const G4ParticleDefinition &part)
bool G4bool
Definition: G4Types.hh:79
Here is the call graph for this function:

◆ CleanTables()

void G4ePolarizedIonisation::CleanTables ( )
private

Definition at line 95 of file G4ePolarizedIonisation.cc.

96 {
97  if(theAsymmetryTable) {
99  delete theAsymmetryTable;
100  theAsymmetryTable = nullptr;
101  }
105  theTransverseAsymmetryTable = nullptr;
106  }
107 }
void clearAndDestroy()
G4PhysicsTable * theTransverseAsymmetryTable
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeAsymmetry()

G4double G4ePolarizedIonisation::ComputeAsymmetry ( G4double  energy,
const G4MaterialCutsCouple couple,
const G4ParticleDefinition particle,
G4double  cut,
G4double tasm 
)
private

Definition at line 327 of file G4ePolarizedIonisation.cc.

332 {
333  G4double lAsymmetry = 0.0;
334  tAsymmetry = 0.0;
335  if (isElectron) { lAsymmetry = tAsymmetry = -1.0; }
336 
337  // calculate polarized cross section
341  G4double sigma2=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
342 
343  // calculate transversely polarized cross section
347  G4double sigma3=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
348 
349  // calculate unpolarized cross section
353  G4double sigma0 = emModel->CrossSection(couple,&aParticle,energy,cut,energy);
354  // determine assymmetries
355  if (sigma0 > 0.) {
356  lAsymmetry=sigma2/sigma0 - 1.;
357  tAsymmetry=sigma3/sigma0 - 1.;
358  }
359  if (std::fabs(lAsymmetry)>1.) {
360  G4cout<<"G4ePolarizedIonisation::ComputeAsymmetry WARNING: E(MeV)= "
361  << energy << " lAsymmetry= "<<lAsymmetry
362  <<" ("<<std::fabs(lAsymmetry)-1.<<")\n";
363  }
364  if (std::fabs(tAsymmetry)>1.) {
365  G4cout<<" energy="<<energy<<"\n";
366  G4cout<<"G4ePolarizedIonisation::ComputeAsymmetry WARNING: E(MeV)= "
367  << energy << " tAsymmetry= "<<tAsymmetry
368  <<" ("<<std::fabs(tAsymmetry)-1.<<")\n";
369  }
370  return lAsymmetry;
371 }
CLHEP::Hep3Vector G4ThreeVector
G4PolarizedMollerBhabhaModel * emModel
G4GLOB_DLL std::ostream G4cout
double energy
Definition: plottest35.C:25
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:501
void SetTargetPolarization(const G4ThreeVector &pTarget)
double G4double
Definition: G4Types.hh:76
void SetBeamPolarization(const G4ThreeVector &pBeam)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeSaturationFactor()

G4double G4ePolarizedIonisation::ComputeSaturationFactor ( const G4Track &  aTrack)
private

Definition at line 203 of file G4ePolarizedIonisation.cc.

204 {
205  G4Material* aMaterial = track.GetMaterial();
206  G4VPhysicalVolume* aPVolume = track.GetVolume();
207  G4LogicalVolume* aLVolume = aPVolume->GetLogicalVolume();
208 
210 
211  const G4bool volumeIsPolarized = polarizationManger->IsPolarized(aLVolume);
212  G4StokesVector volPolarization = polarizationManger->GetVolumePolarization(aLVolume);
213 
214  G4double factor = 1.0;
215 
216  if (volumeIsPolarized && !volPolarization.IsZero()) {
217 
218  // *** get asymmetry, if target is polarized ***
219  const G4DynamicParticle* aDynamicPart = track.GetDynamicParticle();
220  const G4double energy = aDynamicPart->GetKineticEnergy();
221  const G4StokesVector polarization = track.GetPolarization();
222  const G4ParticleMomentum direction0 = aDynamicPart->GetMomentumDirection();
223 
224  if (verboseLevel>=2) {
225  G4cout << "G4ePolarizedIonisation::ComputeSaturationFactor: " << G4endl;
226  G4cout << " Energy(MeV) " << energy/MeV << G4endl;
227  G4cout << " Direction " << direction0 << G4endl;
228  G4cout << " Polarization " << polarization << G4endl;
229  G4cout << " MaterialPol. " << volPolarization << G4endl;
230  G4cout << " Phys. Volume " << aPVolume->GetName() << G4endl;
231  G4cout << " Log. Volume " << aLVolume->GetName() << G4endl;
232  G4cout << " Material " << aMaterial << G4endl;
233  }
234 
235  size_t midx = CurrentMaterialCutsCoupleIndex();
236  const G4PhysicsVector* aVector = nullptr;
237  const G4PhysicsVector* bVector = nullptr;
238  if(midx < theAsymmetryTable->size()) {
239  aVector = (*theAsymmetryTable)(midx);
240  }
241  if(midx < theTransverseAsymmetryTable->size()) {
242  bVector = (*theTransverseAsymmetryTable)(midx);
243  }
244  if(aVector && bVector) {
245  G4double lAsymmetry = aVector->Value(energy);
246  G4double tAsymmetry = bVector->Value(energy);
247  G4double polZZ = polarization.z()*(volPolarization*direction0);
248  G4double polXX = polarization.x()*
249  (volPolarization*G4PolarizationHelper::GetParticleFrameX(direction0));
250  G4double polYY = polarization.y()*
251  (volPolarization*G4PolarizationHelper::GetParticleFrameY(direction0));
252 
253  factor /= (1. + polZZ*lAsymmetry + (polXX + polYY)*tAsymmetry);
254 
255  if (verboseLevel>=2) {
256  G4cout << " Asymmetry: " << lAsymmetry << ", " << tAsymmetry << G4endl;
257  G4cout << " PolProduct: " << polXX << ", " << polYY << ", " << polZZ << G4endl;
258  G4cout << " Factor: " << factor << G4endl;
259  }
260  } else {
262  ed << "Problem with asymmetry tables: material index " << midx
263  << " is out of range or tables are not filled";
264  G4Exception("G4ePolarizedIonisation::ComputeSaturationFactor","em0048",
265  JustWarning, ed, "");
266  }
267  }
268  return factor;
269 }
static const double MeV
Definition: G4SIunits.hh:211
G4int verboseLevel
Definition: G4VProcess.hh:368
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const G4ThreeVector & GetVolumePolarization(G4LogicalVolume *lVol) const
static G4PolarizationManager * GetInstance()
G4double GetKineticEnergy() const
bool IsPolarized(G4LogicalVolume *lVol) const
G4GLOB_DLL std::ostream G4cout
double energy
Definition: plottest35.C:25
const G4String & GetName() const
bool G4bool
Definition: G4Types.hh:79
G4double Value(G4double theEnergy, size_t &lastidx) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
double x() const
size_t CurrentMaterialCutsCoupleIndex() const
const G4String & GetName() const
static const G4double factor
G4bool IsZero() const
double y() const
const G4ThreeVector & GetMomentumDirection() const
double z() const
static G4ThreeVector GetParticleFrameY(const G4ThreeVector &)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4LogicalVolume * GetLogicalVolume() const
static G4ThreeVector GetParticleFrameX(const G4ThreeVector &)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DefineBaseParticle()

const G4ParticleDefinition* G4ePolarizedIonisation::DefineBaseParticle ( const G4ParticleDefinition p)
protected

◆ GetMeanFreePath()

G4double G4ePolarizedIonisation::GetMeanFreePath ( const G4Track &  track,
G4double  previousStepSize,
G4ForceCondition *  condition 
)
protectedvirtual

Reimplemented from G4VEnergyLossProcess.

Definition at line 157 of file G4ePolarizedIonisation.cc.

160 {
161  // *** get unploarised mean free path from lambda table ***
165  }
166  if (verboseLevel>=2) {
167  G4cout << "G4ePolarizedIonisation::MeanFreePath: "
168  << mfp / mm << " mm " << G4endl;
169  }
170  return mfp;
171 }
G4int verboseLevel
Definition: G4VProcess.hh:368
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4GLOB_DLL std::ostream G4cout
G4double ComputeSaturationFactor(const G4Track &aTrack)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
static const double mm
Definition: G4SIunits.hh:114
G4PhysicsTable * theTransverseAsymmetryTable
Here is the call graph for this function:

◆ InitialiseEnergyLossProcess()

void G4ePolarizedIonisation::InitialiseEnergyLossProcess ( const G4ParticleDefinition part,
const G4ParticleDefinition  
)
protectedvirtual

Implements G4VEnergyLossProcess.

Definition at line 128 of file G4ePolarizedIonisation.cc.

131 {
132  if(!isInitialised) {
133 
134  if(part == G4Positron::Positron()) { isElectron = false; }
135 
137  flucModel = FluctModel();
138 
140  SetEmModel(emModel, 1);
145 
146  isInitialised = true;
147  }
148 }
void SetFluctModel(G4VEmFluctuationModel *)
G4PolarizedMollerBhabhaModel * emModel
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
G4VEmFluctuationModel * FluctModel()
G4double MinKinEnergy() const
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4double MaxKinEnergy() const
static G4EmParameters * Instance()
void SetEmModel(G4VEmModel *, G4int index=1)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4VEmFluctuationModel * flucModel
Here is the call graph for this function:

◆ IsApplicable()

G4bool G4ePolarizedIonisation::IsApplicable ( const G4ParticleDefinition p)
virtual

Implements G4VEnergyLossProcess.

Definition at line 122 of file G4ePolarizedIonisation.cc.

123 {
124  return (&p == G4Electron::Electron() || &p == G4Positron::Positron());
125 }
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4Electron * Electron()
Definition: G4Electron.cc:94
Here is the call graph for this function:

◆ MinPrimaryEnergy()

G4double G4ePolarizedIonisation::MinPrimaryEnergy ( const G4ParticleDefinition ,
const G4Material ,
G4double  cut 
)
protectedvirtual

Reimplemented from G4VEnergyLossProcess.

Definition at line 112 of file G4ePolarizedIonisation.cc.

114 {
115  G4double x = cut;
116  if(isElectron) { x += cut; }
117  return x;
118 }
double G4double
Definition: G4Types.hh:76

◆ operator=()

G4ePolarizedIonisation& G4ePolarizedIonisation::operator= ( const G4ePolarizedIonisation right)
private

◆ PostStepGetPhysicalInteractionLength()

G4double G4ePolarizedIonisation::PostStepGetPhysicalInteractionLength ( const G4Track &  track,
G4double  previousStepSize,
G4ForceCondition *  condition 
)
protectedvirtual

Reimplemented from G4VEnergyLossProcess.

Definition at line 175 of file G4ePolarizedIonisation.cc.

178 {
179  // save previous value
181 
182  // *** get unploarised mean free path from lambda table ***
184 
187  if(nLength > 0.0) {
189  std::max(nLength - step/curLength, 0.0);
190  }
191  x = theNumberOfInteractionLengthLeft * curLength;
192  }
193  if (verboseLevel>=2) {
194  G4cout << "G4ePolarizedIonisation::PostStepGetPhysicalInteractionLength: "
195  << x/mm << " mm " << G4endl;
196  }
197  return x;
198 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4GLOB_DLL std::ostream G4cout
G4double currentInteractionLength
Definition: G4VProcess.hh:297
G4double ComputeSaturationFactor(const G4Track &aTrack)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
static const double mm
Definition: G4SIunits.hh:114
G4PhysicsTable * theTransverseAsymmetryTable
Here is the call graph for this function:

◆ PrintInfo()

void G4ePolarizedIonisation::PrintInfo ( )
virtual

Implements G4VEnergyLossProcess.

Definition at line 152 of file G4ePolarizedIonisation.cc.

153 {}

Member Data Documentation

◆ emModel

G4PolarizedMollerBhabhaModel* G4ePolarizedIonisation::emModel
private

Definition at line 123 of file G4ePolarizedIonisation.hh.

◆ flucModel

G4VEmFluctuationModel* G4ePolarizedIonisation::flucModel
private

Definition at line 122 of file G4ePolarizedIonisation.hh.

◆ isElectron

G4bool G4ePolarizedIonisation::isElectron
private

Definition at line 125 of file G4ePolarizedIonisation.hh.

◆ isInitialised

G4bool G4ePolarizedIonisation::isInitialised
private

Definition at line 126 of file G4ePolarizedIonisation.hh.

◆ theAsymmetryTable

G4PhysicsTable* G4ePolarizedIonisation::theAsymmetryTable
private

Definition at line 131 of file G4ePolarizedIonisation.hh.

◆ theElectron

G4ParticleDefinition* G4ePolarizedIonisation::theElectron
private

Definition at line 121 of file G4ePolarizedIonisation.hh.

◆ theTargetPolarization

G4ThreeVector G4ePolarizedIonisation::theTargetPolarization
private

Definition at line 129 of file G4ePolarizedIonisation.hh.

◆ theTransverseAsymmetryTable

G4PhysicsTable* G4ePolarizedIonisation::theTransverseAsymmetryTable
private

Definition at line 132 of file G4ePolarizedIonisation.hh.


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