#include <G4LivermorePolarizedRayleighModel.hh>
Inheritance diagram for G4LivermorePolarizedRayleighModel:
Public Member Functions | |
G4LivermorePolarizedRayleighModel (const G4ParticleDefinition *p=0, const G4String &nam="LivermorePolarizedRayleigh") | |
virtual | ~G4LivermorePolarizedRayleighModel () |
virtual void | Initialise (const G4ParticleDefinition *, const G4DataVector &) |
virtual G4double | ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) |
virtual void | SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) |
Protected Attributes | |
G4ParticleChangeForGamma * | fParticleChange |
Definition at line 43 of file G4LivermorePolarizedRayleighModel.hh.
G4LivermorePolarizedRayleighModel::G4LivermorePolarizedRayleighModel | ( | const G4ParticleDefinition * | p = 0 , |
|
const G4String & | nam = "LivermorePolarizedRayleigh" | |||
) |
Definition at line 53 of file G4LivermorePolarizedRayleighModel.cc.
References G4cout, G4endl, and G4VEmModel::SetHighEnergyLimit().
00055 :G4VEmModel(nam),fParticleChange(0),isInitialised(false), 00056 crossSectionHandler(0),formFactorData(0) 00057 { 00058 lowEnergyLimit = 250 * eV; 00059 highEnergyLimit = 100 * GeV; 00060 00061 //SetLowEnergyLimit(lowEnergyLimit); 00062 SetHighEnergyLimit(highEnergyLimit); 00063 // 00064 verboseLevel= 0; 00065 // Verbosity scale: 00066 // 0 = nothing 00067 // 1 = warning for energy non-conservation 00068 // 2 = details of energy budget 00069 // 3 = calculation of cross sections, file openings, sampling of atoms 00070 // 4 = entering in methods 00071 00072 if(verboseLevel > 0) { 00073 G4cout << "Livermore Polarized Rayleigh is constructed " << G4endl 00074 << "Energy range: " 00075 << lowEnergyLimit / eV << " eV - " 00076 << highEnergyLimit / GeV << " GeV" 00077 << G4endl; 00078 } 00079 }
G4LivermorePolarizedRayleighModel::~G4LivermorePolarizedRayleighModel | ( | ) | [virtual] |
Definition at line 83 of file G4LivermorePolarizedRayleighModel.cc.
00084 { 00085 if (crossSectionHandler) delete crossSectionHandler; 00086 if (formFactorData) delete formFactorData; 00087 }
G4double G4LivermorePolarizedRayleighModel::ComputeCrossSectionPerAtom | ( | const G4ParticleDefinition * | , | |
G4double | kinEnergy, | |||
G4double | Z, | |||
G4double | A = 0 , |
|||
G4double | cut = 0 , |
|||
G4double | emax = DBL_MAX | |||
) | [virtual] |
Reimplemented from G4VEmModel.
Definition at line 145 of file G4LivermorePolarizedRayleighModel.cc.
References G4VCrossSectionHandler::FindValue(), G4cout, and G4endl.
00150 { 00151 if (verboseLevel > 3) 00152 G4cout << "Calling CrossSectionPerAtom() of G4LivermorePolarizedRayleighModel" << G4endl; 00153 00154 if (GammaEnergy < lowEnergyLimit || GammaEnergy > highEnergyLimit) return 0.0; 00155 00156 G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy); 00157 return cs; 00158 }
void G4LivermorePolarizedRayleighModel::Initialise | ( | const G4ParticleDefinition * | , | |
const G4DataVector & | ||||
) | [virtual] |
Implements G4VEmModel.
Definition at line 91 of file G4LivermorePolarizedRayleighModel.cc.
References G4VCrossSectionHandler::Clear(), fParticleChange, G4cout, G4endl, G4VEmModel::GetParticleChangeForGamma(), G4VEmModel::HighEnergyLimit(), G4VEmModel::InitialiseElementSelectors(), G4VEMDataSet::LoadData(), G4VCrossSectionHandler::LoadData(), and G4VEmModel::LowEnergyLimit().
00093 { 00094 // Rayleigh process: The Quantum Theory of Radiation 00095 // W. Heitler, Oxford at the Clarendon Press, Oxford (1954) 00096 // Scattering function: A simple model of photon transport 00097 // D.E. Cullen, Nucl. Instr. Meth. in Phys. Res. B 101 (1995) 499-510 00098 // Polarization of the outcoming photon: Beam test of a prototype detector array for the PoGO astronomical hard X-ray/soft gamma-ray polarimeter 00099 // T. Mizuno et al., Nucl. Instr. Meth. in Phys. Res. A 540 (2005) 158-168 00100 00101 if (verboseLevel > 3) 00102 G4cout << "Calling G4LivermorePolarizedRayleighModel::Initialise()" << G4endl; 00103 00104 if (crossSectionHandler) 00105 { 00106 crossSectionHandler->Clear(); 00107 delete crossSectionHandler; 00108 } 00109 00110 // Read data files for all materials 00111 00112 crossSectionHandler = new G4CrossSectionHandler; 00113 crossSectionHandler->Clear(); 00114 G4String crossSectionFile = "rayl/re-cs-"; 00115 crossSectionHandler->LoadData(crossSectionFile); 00116 00117 G4VDataSetAlgorithm* ffInterpolation = new G4LogLogInterpolation; 00118 G4String formFactorFile = "rayl/re-ff-"; 00119 formFactorData = new G4CompositeEMDataSet(ffInterpolation,1.,1.); 00120 formFactorData->LoadData(formFactorFile); 00121 00122 InitialiseElementSelectors(particle,cuts); 00123 00124 // 00125 if (verboseLevel > 2) 00126 G4cout << "Loaded cross section files for Livermore Polarized Rayleigh model" << G4endl; 00127 00128 InitialiseElementSelectors(particle,cuts); 00129 00130 if (verboseLevel > 0) { 00131 G4cout << "Livermore Polarized Rayleigh model is initialized " << G4endl 00132 << "Energy range: " 00133 << LowEnergyLimit() / eV << " eV - " 00134 << HighEnergyLimit() / GeV << " GeV" 00135 << G4endl; 00136 } 00137 00138 if(isInitialised) return; 00139 fParticleChange = GetParticleChangeForGamma(); 00140 isInitialised = true; 00141 }
void G4LivermorePolarizedRayleighModel::SampleSecondaries | ( | std::vector< G4DynamicParticle * > * | , | |
const G4MaterialCutsCouple * | , | |||
const G4DynamicParticle * | , | |||
G4double | tmin, | |||
G4double | maxEnergy | |||
) | [virtual] |
Implements G4VEmModel.
Definition at line 162 of file G4LivermorePolarizedRayleighModel.cc.
References fParticleChange, fStopAndKill, G4cout, G4endl, G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4Element::GetZ(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4ParticleChangeForGamma::ProposePolarization(), G4VParticleChange::ProposeTrackStatus(), G4VEmModel::SelectRandomAtom(), and G4ParticleChangeForGamma::SetProposedKineticEnergy().
00167 { 00168 if (verboseLevel > 3) 00169 G4cout << "Calling SampleSecondaries() of G4LivermorePolarizedRayleighModel" << G4endl; 00170 00171 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy(); 00172 00173 if (photonEnergy0 <= lowEnergyLimit) 00174 { 00175 fParticleChange->ProposeTrackStatus(fStopAndKill); 00176 fParticleChange->SetProposedKineticEnergy(0.); 00177 fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0); 00178 return ; 00179 } 00180 00181 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection(); 00182 00183 // Select randomly one element in the current material 00184 // G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy0); 00185 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition(); 00186 const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0); 00187 G4int Z = (G4int)elm->GetZ(); 00188 00189 G4double outcomingPhotonCosTheta = GenerateCosTheta(photonEnergy0, Z); 00190 G4double outcomingPhotonPhi = GeneratePhi(outcomingPhotonCosTheta); 00191 G4double beta=GeneratePolarizationAngle(); 00192 00193 // incomingPhoton reference frame: 00194 // z = versor parallel to the incomingPhotonDirection 00195 // x = versor parallel to the incomingPhotonPolarization 00196 // y = defined as z^x 00197 00198 // outgoingPhoton reference frame: 00199 // z' = versor parallel to the outgoingPhotonDirection 00200 // x' = defined as x-x*z'z' normalized 00201 // y' = defined as z'^x' 00202 00203 G4ThreeVector z(aDynamicGamma->GetMomentumDirection().unit()); 00204 G4ThreeVector x(GetPhotonPolarization(*aDynamicGamma)); 00205 G4ThreeVector y(z.cross(x)); 00206 00207 // z' = std::cos(phi)*std::sin(theta) x + std::sin(phi)*std::sin(theta) y + std::cos(theta) z 00208 G4double xDir; 00209 G4double yDir; 00210 G4double zDir; 00211 zDir=outcomingPhotonCosTheta; 00212 xDir=std::sqrt(1-outcomingPhotonCosTheta*outcomingPhotonCosTheta); 00213 yDir=xDir; 00214 xDir*=std::cos(outcomingPhotonPhi); 00215 yDir*=std::sin(outcomingPhotonPhi); 00216 00217 G4ThreeVector zPrime((xDir*x + yDir*y + zDir*z).unit()); 00218 G4ThreeVector xPrime(x.perpPart(zPrime).unit()); 00219 G4ThreeVector yPrime(zPrime.cross(xPrime)); 00220 00221 // outgoingPhotonPolarization is directed as x' std::cos(beta) + y' std::sin(beta) 00222 G4ThreeVector outcomingPhotonPolarization(xPrime*std::cos(beta) + yPrime*std::sin(beta)); 00223 00224 fParticleChange->ProposeMomentumDirection(zPrime); 00225 fParticleChange->ProposePolarization(outcomingPhotonPolarization); 00226 fParticleChange->SetProposedKineticEnergy(photonEnergy0); 00227 00228 }
Definition at line 71 of file G4LivermorePolarizedRayleighModel.hh.
Referenced by Initialise(), and SampleSecondaries().