60 :
G4VEmModel(
"G4JAEAElasticScatteringModel"),isInitialised(false)
74 G4cout <<
"G4JAEAElasticScatteringModel is constructed " <<
G4endl;
103 G4cout <<
"Calling Initialise() of G4JAEAElasticScatteringModel." <<
G4endl
115 char* path = std::getenv(
"G4LEDATA");
120 for(
G4int i=0; i<numOfCouples; ++i)
128 for (
G4int j=0; j<nelm; ++j)
158 G4cout <<
"Calling ReadData() of G4JAEAElasticScatteringModel"
164 const char* datadir = path;
168 datadir = std::getenv(
"G4LEDATA");
171 G4Exception(
"G4JAEAElasticScatteringModel::ReadData()",
"em0006",
173 "Environment variable G4LEDATA not defined");
183 std::ostringstream ostCS;
184 ostCS << datadir <<
"/JAEAESData/amp_Z_" <<
Z ;
185 std::ifstream ES_Data_Buffer(ostCS.str().c_str(),ios::binary);
186 if( !ES_Data_Buffer.is_open() )
189 ed <<
"G4JAEAElasticScattertingModel data file <" << ostCS.str().c_str()
190 <<
"> is not opened!" <<
G4endl;
193 "G4LEDATA version should be G4EMLOW7.11 or later. Elastic Scattering Data are not loaded");
199 G4cout <<
"File " << ostCS.str()
200 <<
" is opened by G4JAEAElasticScatteringModel" <<
G4endl;
207 while (ES_Data_Buffer.read(
reinterpret_cast<char*
>(&buffer_var),
sizeof(
float)))
219 for (
G4int i=0;i<300;++i)
237 G4cout <<
"G4JAEAElasticScatteringModel::ComputeCrossSectionPerAtom()"
247 if(intZ < 1 || intZ >
maxZ) {
return xs; }
256 if(!pv) {
return xs; }
264 }
else if(e >= pv->
Energy(0)) {
270 G4cout <<
"****** DEBUG: tcs value for Z=" <<
Z <<
" at energy (MeV)="
272 G4cout <<
" cs (Geant4 internal unit)=" << xs <<
G4endl;
273 G4cout <<
" -> first E*E*cs value in CS data file (iu) =" << (*pv)[0]
275 G4cout <<
" -> last E*E*cs value in CS data file (iu) =" << (*pv)[
n]
277 G4cout <<
"*********************************************************"
287 std::vector<G4DynamicParticle*>*,
293 G4cout <<
"Calling SampleSecondaries() of G4JAEAElasticScatteringModel."
311 Xi1=gammaPolarization0.
x();
312 Xi2=gammaPolarization0.
y();
313 Xi3=gammaPolarization0.
z();
315 G4double polarization_magnitude=Xi1*Xi1+Xi2*Xi2+Xi3*Xi3;
316 if ((polarization_magnitude)>0 || (Xi1*Xi1>0) || (Xi2*Xi2>0) || (Xi3*Xi3>0))
318 G4cout<<
"WARNING: G4JAEAElasticScatteringModel is only compatible with non-polarized photons."
329 G4int energyindex=round(100*photonEnergy0)-1;
336 for (
G4int i=0;i<=180;++i)
338 a1=
ES_Data[
Z]->at(4*i+300+181*4*(energyindex));
339 a2=
ES_Data[
Z]->at(4*i+1+300+181*4*(energyindex));
340 a3=
ES_Data[
Z]->at(4*i+2+300+181*4*(energyindex));
341 a4=
ES_Data[
Z]->at(4*i+3+300+181*4*(energyindex));
347 for (
G4int i =0;i<=180;++i)
351 for (
G4int i=0; i<=180;++i)
353 cdfsum=cdfsum+
pdf[i];
361 G4double theta = (cdfindex+cdfinv)/180.;
375 G4ThreeVector finaldirection(sin(theta)*cos(phi), sin(theta)*sin(phi), cos(theta));
std::vector< const G4Element * > G4ElementVector
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
static constexpr double keV
static constexpr double eV
static constexpr double GeV
#define G4MUTEX_INITIALIZER
G4GLOB_DLL std::ostream G4cout
Hep3Vector & rotateUz(const Hep3Vector &)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void InitialiseForElement(const G4ParticleDefinition *, G4int Z) override
virtual ~G4JAEAElasticScatteringModel()
static G4DataVector * ES_Data[maxZ+1]
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
static G4PhysicsFreeVector * dataCS[maxZ+1]
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
G4JAEAElasticScatteringModel()
void ReadData(size_t Z, const char *path=0)
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double distribution[181]
G4ParticleChangeForGamma * fParticleChange
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double Energy(const std::size_t index) const
G4double Value(const G4double energy, std::size_t &lastidx) const
std::size_t GetVectorLength() const
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
std::vector< G4EmElementSelector * > * GetElementSelectors()
G4double HighEnergyLimit() const
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static constexpr double twopi
static constexpr double pi
G4Mutex G4JAEAElasticScatteringModelMutex