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

#include <G4ePolarizedIonisation.hh>

+ Inheritance diagram for G4ePolarizedIonisation:

Public Member Functions

 G4ePolarizedIonisation (const G4String &name="pol-eIoni")
 
virtual ~G4ePolarizedIonisation ()
 
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 G4bool IsApplicable (const G4ParticleDefinition &p)=0
 
virtual void PrintInfo ()=0
 
void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4PhysicsTableBuildDEDXTable (G4EmTableType tType=fRestricted)
 
G4PhysicsTableBuildLambdaTable (G4EmTableType tType=fRestricted)
 
void PrintInfoDefinition ()
 
void StartTracking (G4Track *)
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
G4double SampleSubCutSecondaries (std::vector< G4Track * > &, const G4Step &, G4VEmModel *model, G4int matIdx)
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
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)
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false)
 
G4int NumberOfModels ()
 
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 SetRandomStep (G4bool val)
 
void SetIntegral (G4bool val)
 
G4bool IsIntegral () const
 
void SetIonisation (G4bool val)
 
G4bool IsIonisationProcess () const
 
void SetLinearLossLimit (G4double val)
 
void SetMinSubRange (G4double val)
 
void SetLambdaFactor (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 SetLambdaBinning (G4int nbins)
 
void SetDEDXBinningForCSDARange (G4int nbins)
 
void SetMinKinEnergy (G4double e)
 
G4double MinKinEnergy () const
 
void SetMaxKinEnergy (G4double e)
 
G4double MaxKinEnergy () const
 
void SetMaxKinEnergyForCSDARange (G4double e)
 
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
 
G4PhysicsTableRangeTableForLoss () const
 
G4PhysicsTableInverseRangeTable () const
 
G4PhysicsTableLambdaTable ()
 
G4PhysicsTableSubLambdaTable ()
 
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 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
 

Protected Member Functions

virtual void InitialiseEnergyLossProcess (const G4ParticleDefinition *, const G4ParticleDefinition *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *, G4double cut)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4double ComputeAsymmetry (G4double energy, const G4MaterialCutsCouple *couple, const G4ParticleDefinition &particle, G4double cut, G4double &tasm)
 
const G4ParticleDefinitionDefineBaseParticle (const G4ParticleDefinition *p)
 
- Protected Member Functions inherited from G4VEnergyLossProcess
virtual void InitialiseEnergyLossProcess (const G4ParticleDefinition *, const G4ParticleDefinition *)=0
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *, G4double cut)
 
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *, G4double cut)
 
size_t CurrentMaterialCutsCoupleIndex () const
 
G4double GetCurrentRange () const
 
void SelectModel (G4double kinEnergy)
 
void SetParticle (const G4ParticleDefinition *p)
 
void SetSecondaryParticle (const G4ParticleDefinition *p)
 
- 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 ()
 

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
 
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
 

Detailed Description

Definition at line 66 of file G4ePolarizedIonisation.hh.

Constructor & Destructor Documentation

◆ G4ePolarizedIonisation()

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

Definition at line 67 of file G4ePolarizedIonisation.cc.

69 theElectron(G4Electron::Electron()),
70 isElectron(true),
71 isInitialised(false),
72 theAsymmetryTable(NULL),
73 theTransverseAsymmetryTable(NULL)
74{
76 // SetDEDXBinning(120);
77 // SetLambdaBinning(120);
78 // numBinAsymmetryTable=78;
79
80 // SetMinKinEnergy(0.1*keV);
81 // SetMaxKinEnergy(100.0*TeV);
82 // PrintInfoDefinition();
84 flucModel = 0;
85 emModel = 0;
86}
@ fIonisation
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403

◆ ~G4ePolarizedIonisation()

G4ePolarizedIonisation::~G4ePolarizedIonisation ( )
virtual

Definition at line 90 of file G4ePolarizedIonisation.cc.

91{
92 if (theAsymmetryTable) {
93 theAsymmetryTable->clearAndDestroy();
94 delete theAsymmetryTable;
95 }
96 if (theTransverseAsymmetryTable) {
97 theTransverseAsymmetryTable->clearAndDestroy();
98 delete theTransverseAsymmetryTable;
99 }
100}
void clearAndDestroy()

Member Function Documentation

◆ BuildPhysicsTable()

void G4ePolarizedIonisation::BuildPhysicsTable ( const G4ParticleDefinition part)
protectedvirtual

Reimplemented from G4VProcess.

Definition at line 245 of file G4ePolarizedIonisation.cc.

246{
247 // *** build DEDX and (unpolarized) cross section tables
249 // G4PhysicsTable* pt =
250 // BuildDEDXTable();
251
252
253 // *** build asymmetry-table
254 if (theAsymmetryTable) {
255 theAsymmetryTable->clearAndDestroy(); delete theAsymmetryTable;}
256 if (theTransverseAsymmetryTable) {
257 theTransverseAsymmetryTable->clearAndDestroy(); delete theTransverseAsymmetryTable;}
258
259 const G4ProductionCutsTable* theCoupleTable=
261 size_t numOfCouples = theCoupleTable->GetTableSize();
262
263 theAsymmetryTable = new G4PhysicsTable(numOfCouples);
264 theTransverseAsymmetryTable = new G4PhysicsTable(numOfCouples);
265
266 for (size_t j=0 ; j < numOfCouples; j++ ) {
267 // get cut value
268 const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(j);
269
270 G4double cut = (*theCoupleTable->GetEnergyCutsVector(1))[j];
271
272 //create physics vectors then fill it (same parameters as lambda vector)
273 G4PhysicsVector * ptrVectorA = LambdaPhysicsVector(couple,cut);
274 G4PhysicsVector * ptrVectorB = LambdaPhysicsVector(couple,cut);
275 size_t bins = ptrVectorA->GetVectorLength();
276
277 for (size_t i = 0 ; i < bins ; i++ ) {
278 G4double lowEdgeEnergy = ptrVectorA->Energy(i);
279 G4double tasm=0.;
280 G4double asym = ComputeAsymmetry(lowEdgeEnergy, couple, part, cut, tasm);
281 ptrVectorA->PutValue(i,asym);
282 ptrVectorB->PutValue(i,tasm);
283 }
284 theAsymmetryTable->insertAt( j , ptrVectorA ) ;
285 theTransverseAsymmetryTable->insertAt( j , ptrVectorB ) ;
286 }
287
288}
double G4double
Definition: G4Types.hh:64
void insertAt(size_t, G4PhysicsVector *)
size_t GetVectorLength() const
G4double Energy(size_t index) const
void PutValue(size_t index, G4double theValue)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *, G4double cut)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double ComputeAsymmetry(G4double energy, const G4MaterialCutsCouple *couple, const G4ParticleDefinition &particle, G4double cut, G4double &tasm)

◆ ComputeAsymmetry()

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

Definition at line 291 of file G4ePolarizedIonisation.cc.

296{
297 G4double lAsymmetry = 0.0;
298 tAsymmetry = 0.0;
299 if (isElectron) {lAsymmetry = tAsymmetry = -1.0;}
300
301 // calculate polarized cross section
302 theTargetPolarization=G4ThreeVector(0.,0.,1.);
303 emModel->SetTargetPolarization(theTargetPolarization);
304 emModel->SetBeamPolarization(theTargetPolarization);
305 G4double sigma2=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
306
307 // calculate transversely polarized cross section
308 theTargetPolarization=G4ThreeVector(1.,0.,0.);
309 emModel->SetTargetPolarization(theTargetPolarization);
310 emModel->SetBeamPolarization(theTargetPolarization);
311 G4double sigma3=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
312
313 // calculate unpolarized cross section
314 theTargetPolarization=G4ThreeVector();
315 emModel->SetTargetPolarization(theTargetPolarization);
316 emModel->SetBeamPolarization(theTargetPolarization);
317 G4double sigma0=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
318 // determine assymmetries
319 if (sigma0>0.) {
320 lAsymmetry=sigma2/sigma0-1.;
321 tAsymmetry=sigma3/sigma0-1.;
322 }
323 if (std::fabs(lAsymmetry)>1.) {
324 G4cout<<" energy="<<energy<<"\n";
325 G4cout<<"WARNING lAsymmetry= "<<lAsymmetry<<" ("<<std::fabs(lAsymmetry)-1.<<")\n";
326 }
327 if (std::fabs(tAsymmetry)>1.) {
328 G4cout<<" energy="<<energy<<"\n";
329 G4cout<<"WARNING tAsymmetry= "<<tAsymmetry<<" ("<<std::fabs(tAsymmetry)-1.<<")\n";
330 }
331// else {
332// G4cout<<" tAsymmetry= "<<tAsymmetry<<" ("<<std::fabs(tAsymmetry)-1.<<")\n";
333// }
334 return lAsymmetry;
335}
CLHEP::Hep3Vector G4ThreeVector
G4DLLIMPORT std::ostream G4cout
void SetBeamPolarization(const G4ThreeVector &pBeam)
void SetTargetPolarization(const G4ThreeVector &pTarget)
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:418

Referenced by BuildPhysicsTable().

◆ DefineBaseParticle()

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

◆ GetMeanFreePath()

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

Implements G4VContinuousDiscreteProcess.

Definition at line 139 of file G4ePolarizedIonisation.cc.

142{
143 // *** get unploarised mean free path from lambda table ***
144 G4double mfp = G4VEnergyLossProcess::GetMeanFreePath(track, step, cond);
145
146
147 // *** get asymmetry, if target is polarized ***
148 G4VPhysicalVolume* aPVolume = track.GetVolume();
149 G4LogicalVolume* aLVolume = aPVolume->GetLogicalVolume();
150
152 G4bool volumeIsPolarized = polarizationManger->IsPolarized(aLVolume);
153 const G4StokesVector ePolarization = track.GetPolarization();
154
155 if (mfp != DBL_MAX && volumeIsPolarized && !ePolarization.IsZero()) {
156 const G4DynamicParticle* aDynamicElectron = track.GetDynamicParticle();
157 G4double eEnergy = aDynamicElectron->GetKineticEnergy();
158 const G4ParticleMomentum eDirection0 = aDynamicElectron->GetMomentumDirection();
159
160 G4StokesVector volumePolarization = polarizationManger->GetVolumePolarization(aLVolume);
161
162 G4bool isOutRange;
163 size_t idx = CurrentMaterialCutsCoupleIndex();
164 G4double lAsymmetry = (*theAsymmetryTable)(idx)->
165 GetValue(eEnergy, isOutRange);
166 G4double tAsymmetry = (*theTransverseAsymmetryTable)(idx)->
167 GetValue(eEnergy, isOutRange);
168
169 // calculate longitudinal spin component
170 G4double polZZ = ePolarization.z()*
171 volumePolarization*eDirection0;
172 // calculate transvers spin components
173 G4double polXX = ePolarization.x()*
174 volumePolarization*G4PolarizationHelper::GetParticleFrameX(eDirection0);
175 G4double polYY = ePolarization.y()*
176 volumePolarization*G4PolarizationHelper::GetParticleFrameY(eDirection0);
177
178
179 G4double impact = 1. + polZZ*lAsymmetry + (polXX + polYY)*tAsymmetry;
180 // determine polarization dependent mean free path
181 mfp /= impact;
182 if (mfp <=0.) {
183 G4cout <<"PV impact ( "<<polXX<<" , "<<polYY<<" , "<<polZZ<<" )"<<G4endl;
184 G4cout << " impact on MFP is "<< impact <<G4endl;
185 G4cout<<" lAsymmetry= "<<lAsymmetry<<" ("<<std::fabs(lAsymmetry)-1.<<")\n";
186 G4cout<<" tAsymmetry= "<<tAsymmetry<<" ("<<std::fabs(tAsymmetry)-1.<<")\n";
187 }
188 }
189
190 return mfp;
191}
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
double z() const
double x() const
double y() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4ThreeVector GetParticleFrameY(const G4ThreeVector &)
static G4ThreeVector GetParticleFrameX(const G4ThreeVector &)
bool IsPolarized(G4LogicalVolume *lVol) const
static G4PolarizationManager * GetInstance()
const G4ThreeVector & GetVolumePolarization(G4LogicalVolume *lVol) const
G4bool IsZero() const
G4VPhysicalVolume * GetVolume() const
const G4DynamicParticle * GetDynamicParticle() const
const G4ThreeVector & GetPolarization() const
size_t CurrentMaterialCutsCoupleIndex() const
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4LogicalVolume * GetLogicalVolume() const
#define DBL_MAX
Definition: templates.hh:83

◆ InitialiseEnergyLossProcess()

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

Implements G4VEnergyLossProcess.

Definition at line 104 of file G4ePolarizedIonisation.cc.

107{
108 if(!isInitialised) {
109
110 if(part == G4Positron::Positron()) isElectron = false;
111 SetSecondaryParticle(theElectron);
112
113
114
115 flucModel = new G4UniversalFluctuation();
116 //flucModel = new G4BohrFluctuations();
117
118 // G4VEmModel* em = new G4MollerBhabhaModel();
119 emModel = new G4PolarizedMollerBhabhaModel;
122 AddEmModel(1, emModel, flucModel);
123
124 isInitialised = true;
125 }
126}
static G4Positron * Positron()
Definition: G4Positron.cc:94
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
void SetSecondaryParticle(const G4ParticleDefinition *p)

◆ IsApplicable()

G4bool G4ePolarizedIonisation::IsApplicable ( const G4ParticleDefinition p)
inlinevirtual

Implements G4VEnergyLossProcess.

Definition at line 146 of file G4ePolarizedIonisation.hh.

147{
148 return (&p == G4Electron::Electron() || &p == G4Positron::Positron());
149}

◆ MinPrimaryEnergy()

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

Reimplemented from G4VEnergyLossProcess.

Definition at line 135 of file G4ePolarizedIonisation.hh.

138{
139 G4double x = cut;
140 if(isElectron) x += cut;
141 return x;
142}

◆ PostStepGetPhysicalInteractionLength()

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

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 193 of file G4ePolarizedIonisation.cc.

196{
197 // *** get unploarised mean free path from lambda table ***
199
200
201 // *** get asymmetry, if target is polarized ***
202 G4VPhysicalVolume* aPVolume = track.GetVolume();
203 G4LogicalVolume* aLVolume = aPVolume->GetLogicalVolume();
204
206 G4bool volumeIsPolarized = polarizationManger->IsPolarized(aLVolume);
207 const G4StokesVector ePolarization = track.GetPolarization();
208
209 if (mfp != DBL_MAX && volumeIsPolarized && !ePolarization.IsZero()) {
210 const G4DynamicParticle* aDynamicElectron = track.GetDynamicParticle();
211 G4double eEnergy = aDynamicElectron->GetKineticEnergy();
212 const G4ParticleMomentum eDirection0 = aDynamicElectron->GetMomentumDirection();
213
214 G4StokesVector volumePolarization = polarizationManger->GetVolumePolarization(aLVolume);
215
216 size_t idx = CurrentMaterialCutsCoupleIndex();
217 G4double lAsymmetry = (*theAsymmetryTable)(idx)->Value(eEnergy);
218 G4double tAsymmetry = (*theTransverseAsymmetryTable)(idx)->Value(eEnergy);
219
220 // calculate longitudinal spin component
221 G4double polZZ = ePolarization.z()*
222 volumePolarization*eDirection0;
223 // calculate transvers spin components
224 G4double polXX = ePolarization.x()*
225 volumePolarization*G4PolarizationHelper::GetParticleFrameX(eDirection0);
226 G4double polYY = ePolarization.y()*
227 volumePolarization*G4PolarizationHelper::GetParticleFrameY(eDirection0);
228
229
230 G4double impact = 1. + polZZ*lAsymmetry + (polXX + polYY)*tAsymmetry;
231 // determine polarization dependent mean free path
232 mfp /= impact;
233 if (mfp <=0.) {
234 G4cout <<"PV impact ( "<<polXX<<" , "<<polYY<<" , "<<polZZ<<" )"<<G4endl;
235 G4cout << " impact on MFP is "<< impact <<G4endl;
236 G4cout<<" lAsymmetry= "<<lAsymmetry<<" ("<<std::fabs(lAsymmetry)-1.<<")\n";
237 G4cout<<" tAsymmetry= "<<tAsymmetry<<" ("<<std::fabs(tAsymmetry)-1.<<")\n";
238 }
239 }
240
241 return mfp;
242}
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)

◆ PrintInfo()

void G4ePolarizedIonisation::PrintInfo ( )
virtual

Implements G4VEnergyLossProcess.

Definition at line 130 of file G4ePolarizedIonisation.cc.

131{
132 G4cout << " Delta cross sections from Moller+Bhabha, "
133 << "good description from 1 KeV to 100 GeV."
134 << G4endl;
135}

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