80 G4cout <<
"Calling G4DNACPA100ExcitationModel::Initialise()" <<
G4endl;
84 if (p != fpParticle) {
85 std::ostringstream oss;
86 oss <<
" Model is not applied for this particle " << p->
GetParticleName();
87 G4Exception(
"G4DNACPA100ExcitationModel::G4DNACPA100ExcitationModel",
"CPA001",
93 if (path ==
nullptr) {
95 "G4LEDATA environment variable not set.");
100 if (fpG4_WATER !=
nullptr) {
101 index = fpG4_WATER->GetIndex();
106 if (fpGuanine !=
nullptr) {
107 index = fpGuanine->GetIndex();
112 if (fpDeoxyribose !=
nullptr) {
113 index = fpDeoxyribose->GetIndex();
118 if (fpCytosine !=
nullptr) {
119 index = fpCytosine->GetIndex();
124 if (fpThymine !=
nullptr) {
125 index = fpThymine->GetIndex();
130 if (fpAdenine !=
nullptr) {
131 index = fpAdenine->GetIndex();
136 if (fpPhosphate !=
nullptr) {
137 index = fpPhosphate->GetIndex();
138 AddCrossSectionData(index, p,
"dna/sigma_excitation_e_cpa100_phosphoric_acid", 1. * cm * cm);
150 if (dataModel ==
nullptr) {
151 G4cout <<
"G4DNACPA100ExcitationModel::CrossSectionPerVolume:: not good modelData" <<
G4endl;
154 fpModelData = dataModel;
157 isInitialised =
true;
175 lowLim = fpModelData->GetLowELimit(MatID, p);
178 highLim = fpModelData->GetHighELimit(MatID, p);
181 if (ekin >= lowLim && ekin < highLim) {
183 auto Data = fpModelData->GetData();
185 if ((*Data)[MatID][p] ==
nullptr) {
187 "No model is registered");
190 sigma = (*Data)[MatID][p]->FindValue(ekin);
195 G4cout <<
"__________________________________" <<
G4endl;
196 G4cout <<
"°°° G4DNACPA100ExcitationModel - XS INFO START" <<
G4endl;
197 G4cout <<
"°°° Kinetic energy(eV)=" << ekin / eV <<
" particle : " << particleName <<
G4endl;
198 G4cout <<
"°°° lowLim (eV) = " << lowLim / eV <<
" highLim (eV) : " << highLim / eV <<
G4endl;
200 G4cout <<
"°°° Cross section per " << MatID <<
" ID molecule (cm^2)=" << sigma / cm / cm
202 G4cout <<
"°°° Cross section per Phosphate molecule (cm^-1)="
203 << sigma * MolDensity / (1. / cm) <<
G4endl;
204 G4cout <<
"°°° G4DNACPA100ExcitationModel - XS INFO END" <<
G4endl;
210 return sigma * MolDensity;
223 G4double lowLim = fpModelData->GetLowELimit(materialID, particle);
224 G4double highLim = fpModelData->GetHighELimit(materialID, particle);
227 if (k >= lowLim && k < highLim) {
231 if (materialID == fpG4_WATER->GetIndex()) {
232 level = fpModelData->RandomSelectShell(k, particle, materialID);
233 excitationEnergy = eStructure.ExcitationEnergy(level, materialID);
237 level = eStructure.NumberOfLevels(materialID) *
G4UniformRand();
238 excitationEnergy = eStructure.ExcitationEnergy(level, materialID);
239 }
while ((k - eStructure.ExcitationEnergy(level, materialID)) < 0);
241 newEnergy = k - excitationEnergy;
243 if (k - newEnergy <= 0) {
244 G4cout <<
"k : " << k <<
" newEnergy : " << newEnergy <<
G4endl;
245 G4cout <<
"newEnergy : " << newEnergy <<
" k : " << k
246 <<
" excitationEnergy: " << excitationEnergy <<
G4endl;
247 G4cout <<
"G4DNACPA100ExcitationModel::level : " << eStructure.NumberOfLevels(materialID)
248 <<
" excitationEnergy : " << excitationEnergy <<
G4endl;
255 if (newEnergy >= 0) {
259 (excitationEnergy / k) / (1. + (k / (2 * electron_mass_c2)) * (1. - excitationEnergy / k));
261 cosTheta = std::sqrt(1. - cosTheta);
266 G4double CT1, ST1, CF1, SF1, CT2, ST2, CF2, SF2;
267 G4double sinTheta = std::sqrt(1 - cosTheta * cosTheta);
269 ST1 = std::sqrt(1. - CT1 * CT1);
271 ST1 != 0 ? CF1 = zVers.
x() / ST1 : CF1 = std::cos(2. * pi *
G4UniformRand());
272 ST1 != 0 ? SF1 = zVers.
y() / ST1 : SF1 = std::sqrt(1. - CF1 * CF1);
274 A3 = sinTheta * std::cos(phi);
275 A4 = A3 * CT1 + ST1 * cosTheta;
276 A5 = sinTheta * std::sin(phi);
277 A2 = A4 * SF1 + A5 * CF1;
278 A1 = A4 * CF1 - A5 * SF1;
280 CT2 = CT1 * cosTheta - ST1 * A3;
281 ST2 = std::sqrt(1. - CT2 * CT2);
290 fParticleChangeForGamma->ProposeMomentumDirection(zPrimeVers.
unit());
292 fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
295 fParticleChangeForGamma->SetProposedKineticEnergy(k);
298 fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy);
301 if (materialID == fpG4_WATER->GetIndex()) {
302 const G4Track* theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
308 G4cerr <<
"newEnergy : " << newEnergy <<
" k : " << k
309 <<
" excitationEnergy: " << excitationEnergy <<
G4endl;
310 G4cerr <<
"G4DNACPA100ExcitationModel::level : " << eStructure.NumberOfLevels(materialID)
311 <<
" excitationEnergy : " << excitationEnergy <<
G4endl;
316 "model is not registered for this energy");
320 G4cerr <<
"k : " << k <<
" lowLim : " << lowLim <<
" highLim : " << highLim <<
G4endl;
322 "model is not registered for this energy");
const char * G4FindDataDir(const char *)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
CLHEP::Hep3Vector G4ThreeVector
G4GLOB_DLL std::ostream G4cerr
G4GLOB_DLL std::ostream G4cout
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
Initialise Each model must implement an Initialize method.
G4DNACPA100ExcitationModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="DNACPA100ExcitationModel")
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
SampleSecondaries Each model must implement SampleSecondaries to decide if a particle will be created...
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
CrossSectionPerVolume Every model must implement its own CrossSectionPerVolume method....
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
static G4DNAMaterialManager * Instance()
void SetMasterDataModel(const DNAModelType &t, G4VEmModel *m)
G4VEmModel * GetModel(const DNAModelType &t)
const std::vector< G4double > * GetNumMolPerVolTableFor(const G4Material *) const
Retrieve a table of molecular densities (number of molecules per unit volume) in the G4 unit system f...
static G4DNAMolecularMaterial * Instance()
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * ElectronDefinition()
const G4Material * GetMaterial() const
std::size_t GetIndex() const
static G4MaterialTable * GetMaterialTable()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
const G4String & GetParticleName() const
G4String GetName()
GetName.
void LoadCrossSectionData(const G4ParticleDefinition *particleName)
LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresp...
G4VDNAModel(const G4String &nam, const G4String &applyToMaterial="")
G4VDNAModel Constructeur of the G4VDNAModel class.
void SetLowELimit(const size_t &materialID, const G4ParticleDefinition *particle, G4double lim)
SetLowEnergyLimit.
void SetHighELimit(const size_t &materialID, const G4ParticleDefinition *particle, G4double lim)
SetHighEnergyLimit.
void AddCrossSectionData(const size_t &materialName, const G4ParticleDefinition *particleName, const G4String &fileCS, const G4String &fileDiffCS, const G4double &scaleFactor)
AddCrossSectionData Method used during the initialization of the model class to add a new material....
G4ParticleChangeForGamma * GetParticleChangeForGamma()