G4ePolarizedIonisation Class Reference

#include <G4ePolarizedIonisation.hh>

Inheritance diagram for G4ePolarizedIonisation:

G4VEnergyLossProcess G4VContinuousDiscreteProcess G4VProcess

Public Member Functions

 G4ePolarizedIonisation (const G4String &name="pol-eIoni")
virtual ~G4ePolarizedIonisation ()
G4bool IsApplicable (const G4ParticleDefinition &p)
virtual void PrintInfo ()

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)

Detailed Description

Definition at line 66 of file G4ePolarizedIonisation.hh.


Constructor & Destructor Documentation

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

Definition at line 67 of file G4ePolarizedIonisation.cc.

References fIonisation, G4VProcess::SetProcessSubType(), and G4VProcess::verboseLevel.

00068   : G4VEnergyLossProcess(name),
00069     theElectron(G4Electron::Electron()),
00070     isElectron(true),
00071     isInitialised(false),
00072     theAsymmetryTable(NULL),
00073     theTransverseAsymmetryTable(NULL)
00074 {
00075   verboseLevel=0;
00076   //  SetDEDXBinning(120);
00077   //  SetLambdaBinning(120);
00078   //  numBinAsymmetryTable=78;
00079 
00080   //  SetMinKinEnergy(0.1*keV);
00081   //  SetMaxKinEnergy(100.0*TeV);
00082   //  PrintInfoDefinition();
00083   SetProcessSubType(fIonisation);
00084   flucModel = 0;
00085   emModel = 0; 
00086 }

G4ePolarizedIonisation::~G4ePolarizedIonisation (  )  [virtual]

Definition at line 90 of file G4ePolarizedIonisation.cc.

References G4PhysicsTable::clearAndDestroy().

00091 {
00092   if (theAsymmetryTable) {
00093     theAsymmetryTable->clearAndDestroy();
00094     delete theAsymmetryTable;
00095   }
00096   if (theTransverseAsymmetryTable) {
00097     theTransverseAsymmetryTable->clearAndDestroy();
00098     delete theTransverseAsymmetryTable;
00099   }
00100 }


Member Function Documentation

void G4ePolarizedIonisation::BuildPhysicsTable ( const G4ParticleDefinition  )  [protected, virtual]

Reimplemented from G4VEnergyLossProcess.

Definition at line 245 of file G4ePolarizedIonisation.cc.

References G4VEnergyLossProcess::BuildPhysicsTable(), G4PhysicsTable::clearAndDestroy(), ComputeAsymmetry(), G4PhysicsVector::Energy(), G4ProductionCutsTable::GetEnergyCutsVector(), G4ProductionCutsTable::GetMaterialCutsCouple(), G4ProductionCutsTable::GetProductionCutsTable(), G4ProductionCutsTable::GetTableSize(), G4PhysicsVector::GetVectorLength(), G4PhysicsTable::insertAt(), G4VEnergyLossProcess::LambdaPhysicsVector(), and G4PhysicsVector::PutValue().

00246 {
00247   // *** build DEDX and (unpolarized) cross section tables
00248   G4VEnergyLossProcess::BuildPhysicsTable(part);
00249   //  G4PhysicsTable* pt =
00250   //  BuildDEDXTable();
00251 
00252 
00253   // *** build asymmetry-table
00254   if (theAsymmetryTable) {
00255     theAsymmetryTable->clearAndDestroy(); delete theAsymmetryTable;}
00256   if (theTransverseAsymmetryTable) {
00257     theTransverseAsymmetryTable->clearAndDestroy(); delete theTransverseAsymmetryTable;}
00258 
00259   const G4ProductionCutsTable* theCoupleTable=
00260         G4ProductionCutsTable::GetProductionCutsTable();
00261   size_t numOfCouples = theCoupleTable->GetTableSize();
00262 
00263   theAsymmetryTable = new G4PhysicsTable(numOfCouples);
00264   theTransverseAsymmetryTable = new G4PhysicsTable(numOfCouples);
00265 
00266   for (size_t j=0 ; j < numOfCouples; j++ ) {
00267     // get cut value
00268     const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(j);
00269 
00270     G4double cut = (*theCoupleTable->GetEnergyCutsVector(1))[j];
00271 
00272     //create physics vectors then fill it (same parameters as lambda vector)
00273     G4PhysicsVector * ptrVectorA = LambdaPhysicsVector(couple,cut);
00274     G4PhysicsVector * ptrVectorB = LambdaPhysicsVector(couple,cut);
00275     size_t bins = ptrVectorA->GetVectorLength();
00276 
00277     for (size_t i = 0 ; i < bins ; i++ ) {
00278       G4double lowEdgeEnergy = ptrVectorA->Energy(i);
00279       G4double tasm=0.;
00280       G4double asym = ComputeAsymmetry(lowEdgeEnergy, couple, part, cut, tasm);
00281       ptrVectorA->PutValue(i,asym);
00282       ptrVectorB->PutValue(i,tasm);
00283     }
00284     theAsymmetryTable->insertAt( j , ptrVectorA ) ;
00285     theTransverseAsymmetryTable->insertAt( j , ptrVectorB ) ;
00286   }
00287 
00288 }

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

Definition at line 291 of file G4ePolarizedIonisation.cc.

References G4VEmModel::CrossSection(), G4cout, G4PolarizedMollerBhabhaModel::SetBeamPolarization(), and G4PolarizedMollerBhabhaModel::SetTargetPolarization().

Referenced by BuildPhysicsTable().

00296 {
00297   G4double lAsymmetry = 0.0;
00298            tAsymmetry = 0.0;
00299   if (isElectron) {lAsymmetry = tAsymmetry = -1.0;}
00300 
00301   // calculate polarized cross section
00302   theTargetPolarization=G4ThreeVector(0.,0.,1.);
00303   emModel->SetTargetPolarization(theTargetPolarization);
00304   emModel->SetBeamPolarization(theTargetPolarization);
00305   G4double sigma2=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
00306 
00307   // calculate transversely polarized cross section
00308   theTargetPolarization=G4ThreeVector(1.,0.,0.);
00309   emModel->SetTargetPolarization(theTargetPolarization);
00310   emModel->SetBeamPolarization(theTargetPolarization);
00311   G4double sigma3=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
00312 
00313   // calculate unpolarized cross section
00314   theTargetPolarization=G4ThreeVector();
00315   emModel->SetTargetPolarization(theTargetPolarization);
00316   emModel->SetBeamPolarization(theTargetPolarization);
00317   G4double sigma0=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
00318   // determine assymmetries
00319   if (sigma0>0.) {
00320     lAsymmetry=sigma2/sigma0-1.;
00321     tAsymmetry=sigma3/sigma0-1.;
00322   }
00323   if (std::fabs(lAsymmetry)>1.) {
00324     G4cout<<" energy="<<energy<<"\n";
00325     G4cout<<"WARNING lAsymmetry= "<<lAsymmetry<<" ("<<std::fabs(lAsymmetry)-1.<<")\n";
00326   }
00327   if (std::fabs(tAsymmetry)>1.) {
00328     G4cout<<" energy="<<energy<<"\n";
00329     G4cout<<"WARNING tAsymmetry= "<<tAsymmetry<<" ("<<std::fabs(tAsymmetry)-1.<<")\n";
00330   }
00331 //   else {
00332 //     G4cout<<"        tAsymmetry= "<<tAsymmetry<<" ("<<std::fabs(tAsymmetry)-1.<<")\n";
00333 //   }
00334   return lAsymmetry;
00335 }

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

G4double G4ePolarizedIonisation::GetMeanFreePath ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
) [protected, virtual]

Reimplemented from G4VEnergyLossProcess.

Definition at line 139 of file G4ePolarizedIonisation.cc.

References G4VEnergyLossProcess::CurrentMaterialCutsCoupleIndex(), DBL_MAX, G4cout, G4endl, G4Track::GetDynamicParticle(), G4PolarizationManager::GetInstance(), G4VPhysicalVolume::GetLogicalVolume(), G4VEnergyLossProcess::GetMeanFreePath(), G4PolarizationHelper::GetParticleFrameX(), G4PolarizationHelper::GetParticleFrameY(), G4Track::GetPolarization(), G4Track::GetVolume(), G4PolarizationManager::GetVolumePolarization(), G4PolarizationManager::IsPolarized(), and G4StokesVector::IsZero().

00142 {
00143   // *** get unploarised mean free path from lambda table ***
00144   G4double mfp = G4VEnergyLossProcess::GetMeanFreePath(track, step, cond);
00145 
00146 
00147   // *** get asymmetry, if target is polarized ***
00148   G4VPhysicalVolume*  aPVolume  = track.GetVolume();
00149   G4LogicalVolume*    aLVolume  = aPVolume->GetLogicalVolume();
00150 
00151   G4PolarizationManager * polarizationManger = G4PolarizationManager::GetInstance();
00152   G4bool volumeIsPolarized = polarizationManger->IsPolarized(aLVolume);
00153   const G4StokesVector ePolarization = track.GetPolarization();
00154 
00155   if (mfp != DBL_MAX &&  volumeIsPolarized && !ePolarization.IsZero()) {
00156     const G4DynamicParticle* aDynamicElectron = track.GetDynamicParticle();
00157     G4double eEnergy = aDynamicElectron->GetKineticEnergy();
00158     const G4ParticleMomentum eDirection0 = aDynamicElectron->GetMomentumDirection();
00159 
00160     G4StokesVector volumePolarization = polarizationManger->GetVolumePolarization(aLVolume);
00161 
00162     G4bool isOutRange;
00163     size_t idx = CurrentMaterialCutsCoupleIndex();
00164     G4double lAsymmetry = (*theAsymmetryTable)(idx)->
00165                                   GetValue(eEnergy, isOutRange);
00166     G4double tAsymmetry = (*theTransverseAsymmetryTable)(idx)->
00167                                   GetValue(eEnergy, isOutRange);
00168 
00169     // calculate longitudinal spin component
00170     G4double polZZ = ePolarization.z()*
00171                         volumePolarization*eDirection0;
00172     // calculate transvers spin components
00173     G4double polXX = ePolarization.x()*
00174                         volumePolarization*G4PolarizationHelper::GetParticleFrameX(eDirection0);
00175     G4double polYY = ePolarization.y()*
00176                         volumePolarization*G4PolarizationHelper::GetParticleFrameY(eDirection0);
00177 
00178 
00179     G4double impact = 1. + polZZ*lAsymmetry + (polXX + polYY)*tAsymmetry;
00180     // determine polarization dependent mean free path
00181     mfp /= impact;
00182     if (mfp <=0.) {
00183      G4cout <<"PV impact ( "<<polXX<<" , "<<polYY<<" , "<<polZZ<<" )"<<G4endl;
00184      G4cout << " impact on MFP is "<< impact <<G4endl;
00185      G4cout<<" lAsymmetry= "<<lAsymmetry<<" ("<<std::fabs(lAsymmetry)-1.<<")\n";
00186      G4cout<<" tAsymmetry= "<<tAsymmetry<<" ("<<std::fabs(tAsymmetry)-1.<<")\n";
00187     }
00188   }
00189 
00190   return mfp;
00191 }

void G4ePolarizedIonisation::InitialiseEnergyLossProcess ( const G4ParticleDefinition ,
const G4ParticleDefinition  
) [protected, virtual]

Implements G4VEnergyLossProcess.

Definition at line 104 of file G4ePolarizedIonisation.cc.

References G4VEnergyLossProcess::AddEmModel(), G4VEnergyLossProcess::MaxKinEnergy(), G4VEnergyLossProcess::MinKinEnergy(), G4Positron::Positron(), G4VEmModel::SetHighEnergyLimit(), G4VEmModel::SetLowEnergyLimit(), and G4VEnergyLossProcess::SetSecondaryParticle().

00107 {
00108   if(!isInitialised) {
00109 
00110     if(part == G4Positron::Positron()) isElectron = false;
00111     SetSecondaryParticle(theElectron);
00112 
00113 
00114 
00115     flucModel = new G4UniversalFluctuation();
00116     //flucModel = new G4BohrFluctuations();
00117 
00118     //    G4VEmModel* em = new G4MollerBhabhaModel();
00119     emModel = new G4PolarizedMollerBhabhaModel;
00120     emModel->SetLowEnergyLimit(MinKinEnergy());
00121     emModel->SetHighEnergyLimit(MaxKinEnergy());
00122     AddEmModel(1, emModel, flucModel);
00123 
00124     isInitialised = true;
00125   }
00126 }

G4bool G4ePolarizedIonisation::IsApplicable ( const G4ParticleDefinition p  )  [inline, virtual]

Implements G4VEnergyLossProcess.

Definition at line 146 of file G4ePolarizedIonisation.hh.

References G4Electron::Electron(), and G4Positron::Positron().

00147 {
00148   return (&p == G4Electron::Electron() || &p == G4Positron::Positron());
00149 }

G4double G4ePolarizedIonisation::MinPrimaryEnergy ( const G4ParticleDefinition ,
const G4Material ,
G4double  cut 
) [inline, protected, virtual]

Reimplemented from G4VEnergyLossProcess.

Definition at line 135 of file G4ePolarizedIonisation.hh.

00138 {
00139   G4double x = cut;
00140   if(isElectron) x += cut;
00141   return x;
00142 }

G4double G4ePolarizedIonisation::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
) [protected, virtual]

Reimplemented from G4VEnergyLossProcess.

Definition at line 193 of file G4ePolarizedIonisation.cc.

References G4VEnergyLossProcess::CurrentMaterialCutsCoupleIndex(), DBL_MAX, G4cout, G4endl, G4Track::GetDynamicParticle(), G4PolarizationManager::GetInstance(), G4VPhysicalVolume::GetLogicalVolume(), G4PolarizationHelper::GetParticleFrameX(), G4PolarizationHelper::GetParticleFrameY(), G4Track::GetPolarization(), G4Track::GetVolume(), G4PolarizationManager::GetVolumePolarization(), G4PolarizationManager::IsPolarized(), G4StokesVector::IsZero(), and G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength().

00196 {
00197   // *** get unploarised mean free path from lambda table ***
00198   G4double mfp = G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(track, step, cond);
00199 
00200 
00201   // *** get asymmetry, if target is polarized ***
00202   G4VPhysicalVolume*  aPVolume  = track.GetVolume();
00203   G4LogicalVolume*    aLVolume  = aPVolume->GetLogicalVolume();
00204 
00205   G4PolarizationManager * polarizationManger = G4PolarizationManager::GetInstance();
00206   G4bool volumeIsPolarized = polarizationManger->IsPolarized(aLVolume);
00207   const G4StokesVector ePolarization = track.GetPolarization();
00208 
00209   if (mfp != DBL_MAX &&  volumeIsPolarized && !ePolarization.IsZero()) {
00210     const G4DynamicParticle* aDynamicElectron = track.GetDynamicParticle();
00211     G4double eEnergy = aDynamicElectron->GetKineticEnergy();
00212     const G4ParticleMomentum eDirection0 = aDynamicElectron->GetMomentumDirection();
00213 
00214     G4StokesVector volumePolarization = polarizationManger->GetVolumePolarization(aLVolume);
00215 
00216     size_t idx = CurrentMaterialCutsCoupleIndex();
00217     G4double lAsymmetry = (*theAsymmetryTable)(idx)->Value(eEnergy);
00218     G4double tAsymmetry = (*theTransverseAsymmetryTable)(idx)->Value(eEnergy);
00219 
00220     // calculate longitudinal spin component
00221     G4double polZZ = ePolarization.z()*
00222                         volumePolarization*eDirection0;
00223     // calculate transvers spin components
00224     G4double polXX = ePolarization.x()*
00225                         volumePolarization*G4PolarizationHelper::GetParticleFrameX(eDirection0);
00226     G4double polYY = ePolarization.y()*
00227                         volumePolarization*G4PolarizationHelper::GetParticleFrameY(eDirection0);
00228 
00229 
00230     G4double impact = 1. + polZZ*lAsymmetry + (polXX + polYY)*tAsymmetry;
00231     // determine polarization dependent mean free path
00232     mfp /= impact;
00233     if (mfp <=0.) {
00234      G4cout <<"PV impact ( "<<polXX<<" , "<<polYY<<" , "<<polZZ<<" )"<<G4endl;
00235      G4cout << " impact on MFP is "<< impact <<G4endl;
00236      G4cout<<" lAsymmetry= "<<lAsymmetry<<" ("<<std::fabs(lAsymmetry)-1.<<")\n";
00237      G4cout<<" tAsymmetry= "<<tAsymmetry<<" ("<<std::fabs(tAsymmetry)-1.<<")\n";
00238     }
00239   }
00240 
00241   return mfp;
00242 }

void G4ePolarizedIonisation::PrintInfo (  )  [virtual]

Implements G4VEnergyLossProcess.

Definition at line 130 of file G4ePolarizedIonisation.cc.

References G4cout, and G4endl.

00131 {
00132   G4cout << "      Delta cross sections from Moller+Bhabha, "
00133          << "good description from 1 KeV to 100 GeV."
00134          << G4endl;
00135 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:51:54 2013 for Geant4 by  doxygen 1.4.7