G4eSingleCoulombScatteringModel Class Reference

#include <G4eSingleCoulombScatteringModel.hh>

Inheritance diagram for G4eSingleCoulombScatteringModel:

G4VEmModel

Public Member Functions

 G4eSingleCoulombScatteringModel (const G4String &nam="eSingleCoulombScattering")
virtual ~G4eSingleCoulombScatteringModel ()
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
void SetRecoilThreshold (G4double eth)

Protected Member Functions

void DefineMaterial (const G4MaterialCutsCouple *)
void SetupParticle (const G4ParticleDefinition *)

Protected Attributes

G4ParticleTabletheParticleTable
G4ParticleChangeForGammafParticleChange
G4NistManagerfNistManager
G4ScreeningMottCrossSectionMottcross
const std::vector< G4double > * pCuts
const G4MaterialCutsCouplecurrentCouple
const G4MaterialcurrentMaterial
const G4ElementcurrentElement
G4int currentMaterialIndex
G4double cosThetaMin
G4double recoilThreshold
const G4ParticleDefinitionparticle
G4double mass
G4double lowEnergyLimit

Detailed Description

Definition at line 66 of file G4eSingleCoulombScatteringModel.hh.


Constructor & Destructor Documentation

G4eSingleCoulombScatteringModel::G4eSingleCoulombScatteringModel ( const G4String nam = "eSingleCoulombScattering"  ) 

Definition at line 70 of file G4eSingleCoulombScatteringModel.cc.

References currentCouple, currentElement, currentMaterial, currentMaterialIndex, fNistManager, fParticleChange, G4ParticleTable::GetParticleTable(), G4NistManager::Instance(), lowEnergyLimit, mass, Mottcross, particle, pCuts, recoilThreshold, and theParticleTable.

00071   : G4VEmModel(nam),
00072 
00073     cosThetaMin(1.0),
00074     isInitialised(false)
00075 {
00076         fNistManager = G4NistManager::Instance();
00077         theParticleTable = G4ParticleTable::GetParticleTable();
00078         fParticleChange = 0;
00079 
00080         pCuts=0;
00081         currentMaterial = 0;
00082         currentElement  = 0;
00083         currentCouple = 0;
00084 
00085         lowEnergyLimit  = 0*eV;
00086         recoilThreshold = 0.*eV;
00087         particle = 0;
00088         mass=0;
00089         currentMaterialIndex = -1;
00090 
00091         Mottcross = new G4ScreeningMottCrossSection(); 
00092 
00093 }

G4eSingleCoulombScatteringModel::~G4eSingleCoulombScatteringModel (  )  [virtual]

Definition at line 98 of file G4eSingleCoulombScatteringModel.cc.

References Mottcross.

00099 { delete  Mottcross;}


Member Function Documentation

G4double G4eSingleCoulombScatteringModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A,
G4double  cut,
G4double  emax 
) [virtual]

Reimplemented from G4VEmModel.

Definition at line 123 of file G4eSingleCoulombScatteringModel.cc.

References G4VEmModel::CurrentCouple(), DefineMaterial(), lowEnergyLimit, Mottcross, G4ScreeningMottCrossSection::NuclearCrossSection(), G4ScreeningMottCrossSection::SetupKinematic(), and SetupParticle().

00130 {
00131 
00132         SetupParticle(p);
00133  
00134         G4double cross =0.0;
00135         if(kinEnergy < lowEnergyLimit) return cross;
00136 
00137         DefineMaterial(CurrentCouple());
00138 
00139         //Total Cross section
00140         Mottcross->SetupKinematic(kinEnergy, Z);
00141         cross = Mottcross->NuclearCrossSection();
00142 
00143         //cout<< "....cross "<<G4BestUnit(cross,"Surface") << " cm2 "<< cross/cm2 <<endl;
00144   return cross;
00145 }

void G4eSingleCoulombScatteringModel::DefineMaterial ( const G4MaterialCutsCouple  )  [inline, protected]

Definition at line 148 of file G4eSingleCoulombScatteringModel.hh.

References currentCouple, currentMaterial, currentMaterialIndex, G4MaterialCutsCouple::GetIndex(), and G4MaterialCutsCouple::GetMaterial().

Referenced by ComputeCrossSectionPerAtom(), and SampleSecondaries().

00149 { 
00150         if(cup != currentCouple) {
00151                 currentCouple = cup;
00152                 currentMaterial = cup->GetMaterial();
00153                 currentMaterialIndex = currentCouple->GetIndex();
00154 
00155                 }
00156 }

void G4eSingleCoulombScatteringModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
) [virtual]

Implements G4VEmModel.

Definition at line 103 of file G4eSingleCoulombScatteringModel.cc.

References cosThetaMin, currentCouple, currentMaterialIndex, fParticleChange, G4ProductionCutsTable::GetEnergyCutsVector(), G4VEmModel::GetParticleChangeForGamma(), G4ProductionCutsTable::GetProductionCutsTable(), G4ScreeningMottCrossSection::Initialise(), Mottcross, pCuts, G4VEmModel::PolarAngleLimit(), and SetupParticle().

00105 {
00106         SetupParticle(p);
00107         currentCouple = 0;
00108         currentMaterialIndex = -1;
00109         cosThetaMin = cos(PolarAngleLimit());
00110         Mottcross->Initialise(p,cosThetaMin);
00111  
00112         pCuts = G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
00113 
00114 
00115         if(!isInitialised) {
00116           isInitialised = true;
00117           fParticleChange = GetParticleChangeForGamma();
00118         }
00119 }

void G4eSingleCoulombScatteringModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle ,
G4double  tmin,
G4double  maxEnergy 
) [virtual]

Implements G4VEmModel.

Definition at line 149 of file G4eSingleCoulombScatteringModel.cc.

References currentElement, currentMaterialIndex, DefineMaterial(), fParticleChange, G4DynamicParticle::GetDefinition(), G4ParticleTable::GetIon(), G4DynamicParticle::GetKineticEnergy(), G4ScreeningMottCrossSection::GetMom2Lab(), G4DynamicParticle::GetMomentumDirection(), G4ScreeningMottCrossSection::GetNewDirection(), G4ScreeningMottCrossSection::GetTotalCross(), G4ScreeningMottCrossSection::GetTrec(), lowEnergyLimit, mass, Mottcross, particle, pCuts, G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4VParticleChange::ProposeNonIonizingEnergyDeposit(), recoilThreshold, G4VEmModel::SelectIsotopeNumber(), G4VEmModel::SelectRandomAtom(), G4ParticleChangeForGamma::SetProposedKineticEnergy(), SetupParticle(), and theParticleTable.

00155 {
00156         G4double kinEnergy = dp->GetKineticEnergy();
00157         //cout<<"--- kinEnergy "<<kinEnergy<<endl;
00158 
00159 
00160         if(kinEnergy < lowEnergyLimit) return;
00161         
00162         DefineMaterial(couple);
00163         SetupParticle(dp->GetDefinition());
00164 
00165         // Choose nucleus
00166         currentElement = SelectRandomAtom(couple,particle,
00167                                     kinEnergy,cutEnergy,kinEnergy);//last two :cutEnergy= min e kinEnergy=max
00168 
00169         G4double Z  = currentElement->GetZ();
00170         G4int iz    = G4int(Z);
00171         G4int ia = SelectIsotopeNumber(currentElement);
00172 
00173         //cout<<"Element "<<currentElement->GetName()<<endl;;   
00174 
00175         G4double cross= Mottcross->GetTotalCross();
00176 
00177 
00178         if(cross == 0.0)return;
00179                 
00180         G4ThreeVector dir = dp->GetMomentumDirection(); //old direction
00181         G4ThreeVector newDirection=Mottcross->GetNewDirection();//new direction
00182         newDirection.rotateUz(dir);   
00183   
00184         fParticleChange->ProposeMomentumDirection(newDirection);   
00185   
00186         //Recoil energy
00187         G4double trec= Mottcross->GetTrec();
00188         //Energy after scattering       
00189         G4double finalT = kinEnergy - trec;
00190 
00191 
00192         if(finalT <= lowEnergyLimit) { 
00193                 trec = kinEnergy;  
00194                 finalT = 0.0;
00195                 } 
00196     
00197         fParticleChange->SetProposedKineticEnergy(finalT);
00198 
00199         G4double tcut = recoilThreshold;
00200         if(pCuts) { tcut= std::min(tcut,(*pCuts)[currentMaterialIndex]); 
00201                 }
00202  
00203 
00204         if(trec > tcut) {
00205 
00206                 //cout<<"Trec "<<trec/eV<<endl;
00207                 G4ParticleDefinition* ion = theParticleTable->GetIon(iz, ia, 0.0);
00208 
00209                 //incident before scattering
00210                 G4double ptot=sqrt(Mottcross->GetMom2Lab());
00211                 //incident after scattering
00212                 G4double plab = sqrt(finalT*(finalT + 2.0*mass));
00213                 G4ThreeVector p2 = (ptot*dir - plab*newDirection).unit();
00214                 //secondary particle
00215                 G4DynamicParticle* newdp  = new G4DynamicParticle(ion, p2, trec);
00216                         fvect->push_back(newdp);
00217                         }
00218 
00219         else if(trec > 0.0) {
00220                 fParticleChange->ProposeNonIonizingEnergyDeposit(trec);
00221                 if(trec< tcut) fParticleChange->ProposeLocalEnergyDeposit(trec);
00222                 }
00223 
00224 
00225 }

void G4eSingleCoulombScatteringModel::SetRecoilThreshold ( G4double  eth  )  [inline]

Definition at line 171 of file G4eSingleCoulombScatteringModel.hh.

References recoilThreshold.

00172 {
00173         recoilThreshold = eth;
00174 }

void G4eSingleCoulombScatteringModel::SetupParticle ( const G4ParticleDefinition  )  [inline, protected]

Definition at line 161 of file G4eSingleCoulombScatteringModel.hh.

References G4ParticleDefinition::GetPDGMass(), mass, Mottcross, particle, and G4ScreeningMottCrossSection::SetupParticle().

Referenced by ComputeCrossSectionPerAtom(), Initialise(), and SampleSecondaries().

00162 {
00163         if(p != particle) {
00164                 particle = p;
00165                 mass = particle->GetPDGMass();
00166                 Mottcross->SetupParticle(p);
00167                 }
00168 }


Field Documentation

G4double G4eSingleCoulombScatteringModel::cosThetaMin [protected]

Definition at line 131 of file G4eSingleCoulombScatteringModel.hh.

Referenced by Initialise().

const G4MaterialCutsCouple* G4eSingleCoulombScatteringModel::currentCouple [protected]

Definition at line 125 of file G4eSingleCoulombScatteringModel.hh.

Referenced by DefineMaterial(), G4eSingleCoulombScatteringModel(), and Initialise().

const G4Element* G4eSingleCoulombScatteringModel::currentElement [protected]

Definition at line 127 of file G4eSingleCoulombScatteringModel.hh.

Referenced by G4eSingleCoulombScatteringModel(), and SampleSecondaries().

const G4Material* G4eSingleCoulombScatteringModel::currentMaterial [protected]

Definition at line 126 of file G4eSingleCoulombScatteringModel.hh.

Referenced by DefineMaterial(), and G4eSingleCoulombScatteringModel().

G4int G4eSingleCoulombScatteringModel::currentMaterialIndex [protected]

Definition at line 128 of file G4eSingleCoulombScatteringModel.hh.

Referenced by DefineMaterial(), G4eSingleCoulombScatteringModel(), Initialise(), and SampleSecondaries().

G4NistManager* G4eSingleCoulombScatteringModel::fNistManager [protected]

Definition at line 121 of file G4eSingleCoulombScatteringModel.hh.

Referenced by G4eSingleCoulombScatteringModel().

G4ParticleChangeForGamma* G4eSingleCoulombScatteringModel::fParticleChange [protected]

Definition at line 120 of file G4eSingleCoulombScatteringModel.hh.

Referenced by G4eSingleCoulombScatteringModel(), Initialise(), and SampleSecondaries().

G4double G4eSingleCoulombScatteringModel::lowEnergyLimit [protected]

Definition at line 138 of file G4eSingleCoulombScatteringModel.hh.

Referenced by ComputeCrossSectionPerAtom(), G4eSingleCoulombScatteringModel(), and SampleSecondaries().

G4double G4eSingleCoulombScatteringModel::mass [protected]

Definition at line 137 of file G4eSingleCoulombScatteringModel.hh.

Referenced by G4eSingleCoulombScatteringModel(), SampleSecondaries(), and SetupParticle().

G4ScreeningMottCrossSection* G4eSingleCoulombScatteringModel::Mottcross [protected]

Definition at line 122 of file G4eSingleCoulombScatteringModel.hh.

Referenced by ComputeCrossSectionPerAtom(), G4eSingleCoulombScatteringModel(), Initialise(), SampleSecondaries(), SetupParticle(), and ~G4eSingleCoulombScatteringModel().

const G4ParticleDefinition* G4eSingleCoulombScatteringModel::particle [protected]

Definition at line 136 of file G4eSingleCoulombScatteringModel.hh.

Referenced by G4eSingleCoulombScatteringModel(), SampleSecondaries(), and SetupParticle().

const std::vector<G4double>* G4eSingleCoulombScatteringModel::pCuts [protected]

Definition at line 124 of file G4eSingleCoulombScatteringModel.hh.

Referenced by G4eSingleCoulombScatteringModel(), Initialise(), and SampleSecondaries().

G4double G4eSingleCoulombScatteringModel::recoilThreshold [protected]

Definition at line 132 of file G4eSingleCoulombScatteringModel.hh.

Referenced by G4eSingleCoulombScatteringModel(), SampleSecondaries(), and SetRecoilThreshold().

G4ParticleTable* G4eSingleCoulombScatteringModel::theParticleTable [protected]

Definition at line 119 of file G4eSingleCoulombScatteringModel.hh.

Referenced by G4eSingleCoulombScatteringModel(), and SampleSecondaries().


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