| 1 | //
|
|---|
| 2 | // ********************************************************************
|
|---|
| 3 | // * License and Disclaimer *
|
|---|
| 4 | // * *
|
|---|
| 5 | // * The Geant4 software is copyright of the Copyright Holders of *
|
|---|
| 6 | // * the Geant4 Collaboration. It is provided under the terms and *
|
|---|
| 7 | // * conditions of the Geant4 Software License, included in the file *
|
|---|
| 8 | // * LICENSE and available at http://cern.ch/geant4/license . These *
|
|---|
| 9 | // * include a list of copyright holders. *
|
|---|
| 10 | // * *
|
|---|
| 11 | // * Neither the authors of this software system, nor their employing *
|
|---|
| 12 | // * institutes,nor the agencies providing financial support for this *
|
|---|
| 13 | // * work make any representation or warranty, express or implied, *
|
|---|
| 14 | // * regarding this software system or assume any liability for its *
|
|---|
| 15 | // * use. Please see the license in the file LICENSE and URL above *
|
|---|
| 16 | // * for the full disclaimer and the limitation of liability. *
|
|---|
| 17 | // * *
|
|---|
| 18 | // * This code implementation is the result of the scientific and *
|
|---|
| 19 | // * technical work of the GEANT4 collaboration. *
|
|---|
| 20 | // * By using, copying, modifying or distributing the software (or *
|
|---|
| 21 | // * any work based on the software) you agree to acknowledge its *
|
|---|
| 22 | // * use in resulting scientific publications, and indicate your *
|
|---|
| 23 | // * acceptance of all terms of the Geant4 Software license. *
|
|---|
| 24 | // ********************************************************************
|
|---|
| 25 | //
|
|---|
| 26 | // $Id: G4PenelopeIonisation.cc,v 1.19 2006/06/29 19:40:49 gunter Exp $
|
|---|
| 27 | // GEANT4 tag $Name: geant4-09-02-ref-02 $
|
|---|
| 28 | //
|
|---|
| 29 | // --------------------------------------------------------------
|
|---|
| 30 | //
|
|---|
| 31 | // File name: G4PenelopeIonisation
|
|---|
| 32 | //
|
|---|
| 33 | // Author: Luciano Pandola
|
|---|
| 34 | //
|
|---|
| 35 | // Creation date: March 2003
|
|---|
| 36 | //
|
|---|
| 37 | // Modifications:
|
|---|
| 38 | //
|
|---|
| 39 | // 25.03.03 L.Pandola First implementation
|
|---|
| 40 | // 03.06.03 L.Pandola Added continuous part
|
|---|
| 41 | // 30.06.03 L.Pandola Added positrons
|
|---|
| 42 | // 01.07.03 L.Pandola Changed cross section files for e- and e+
|
|---|
| 43 | // Interface with PenelopeCrossSectionHandler
|
|---|
| 44 | // 18.01.04 M.Mendenhall (Vanderbilt University) [bug report 568]
|
|---|
| 45 | // Changed returns in CalculateDiscreteForElectrons()
|
|---|
| 46 | // to eliminate leaks
|
|---|
| 47 | // 20.01.04 L.Pandola Changed returns in CalculateDiscreteForPositrons()
|
|---|
| 48 | // to eliminate the same bug
|
|---|
| 49 | // 10.03.04 L.Pandola Bug fixed with reference system of delta rays
|
|---|
| 50 | // 17.03.04 L.Pandola Removed unnecessary calls to std::pow(a,b)
|
|---|
| 51 | // 18.03.04 L.Pandola Bug fixed in the destructor
|
|---|
| 52 | // 01.06.04 L.Pandola StopButAlive for positrons on PostStepDoIt
|
|---|
| 53 | // 10.03.05 L.Pandola Fix of bug report 729. The solution works but it is
|
|---|
| 54 | // quite un-elegant. Something better to be found.
|
|---|
| 55 | // --------------------------------------------------------------
|
|---|
| 56 |
|
|---|
| 57 | #include "G4PenelopeIonisation.hh"
|
|---|
| 58 | #include "G4PenelopeCrossSectionHandler.hh"
|
|---|
| 59 | #include "G4AtomicTransitionManager.hh"
|
|---|
| 60 | #include "G4AtomicShell.hh"
|
|---|
| 61 | #include "G4eIonisationSpectrum.hh"
|
|---|
| 62 | #include "G4VDataSetAlgorithm.hh"
|
|---|
| 63 | #include "G4SemiLogInterpolation.hh"
|
|---|
| 64 | #include "G4LogLogInterpolation.hh"
|
|---|
| 65 | #include "G4EMDataSet.hh"
|
|---|
| 66 | #include "G4VEMDataSet.hh"
|
|---|
| 67 | #include "G4CompositeEMDataSet.hh"
|
|---|
| 68 | #include "G4EnergyLossTables.hh"
|
|---|
| 69 | #include "G4UnitsTable.hh"
|
|---|
| 70 | #include "G4Electron.hh"
|
|---|
| 71 | #include "G4Gamma.hh"
|
|---|
| 72 | #include "G4Positron.hh"
|
|---|
| 73 | #include "G4ProductionCutsTable.hh"
|
|---|
| 74 | #include "G4ProcessManager.hh"
|
|---|
| 75 |
|
|---|
| 76 | G4PenelopeIonisation::G4PenelopeIonisation(const G4String& nam)
|
|---|
| 77 | : G4eLowEnergyLoss(nam),
|
|---|
| 78 | crossSectionHandler(0),
|
|---|
| 79 | theMeanFreePath(0),
|
|---|
| 80 | kineticEnergy1(0.0),
|
|---|
| 81 | cosThetaPrimary(1.0),
|
|---|
| 82 | energySecondary(0.0),
|
|---|
| 83 | cosThetaSecondary(0.0),
|
|---|
| 84 | iOsc(-1)
|
|---|
| 85 | {
|
|---|
| 86 | cutForPhotons = 250.0*eV;
|
|---|
| 87 | cutForElectrons = 250.0*eV;
|
|---|
| 88 | verboseLevel = 0;
|
|---|
| 89 | ionizationEnergy = new std::map<G4int,G4DataVector*>;
|
|---|
| 90 | resonanceEnergy = new std::map<G4int,G4DataVector*>;
|
|---|
| 91 | occupationNumber = new std::map<G4int,G4DataVector*>;
|
|---|
| 92 | shellFlag = new std::map<G4int,G4DataVector*>;
|
|---|
| 93 | ReadData(); //Read data from file
|
|---|
| 94 |
|
|---|
| 95 | }
|
|---|
| 96 |
|
|---|
| 97 |
|
|---|
| 98 | G4PenelopeIonisation::~G4PenelopeIonisation()
|
|---|
| 99 | {
|
|---|
| 100 | delete crossSectionHandler;
|
|---|
| 101 | delete theMeanFreePath;
|
|---|
| 102 | for (G4int Z=1;Z<100;Z++)
|
|---|
| 103 | {
|
|---|
| 104 | if (ionizationEnergy->count(Z)) delete (ionizationEnergy->find(Z)->second);
|
|---|
| 105 | if (resonanceEnergy->count(Z)) delete (resonanceEnergy->find(Z)->second);
|
|---|
| 106 | if (occupationNumber->count(Z)) delete (occupationNumber->find(Z)->second);
|
|---|
| 107 | if (shellFlag->count(Z)) delete (shellFlag->find(Z)->second);
|
|---|
| 108 | }
|
|---|
| 109 | delete ionizationEnergy;
|
|---|
| 110 | delete resonanceEnergy;
|
|---|
| 111 | delete occupationNumber;
|
|---|
| 112 | delete shellFlag;
|
|---|
| 113 | }
|
|---|
| 114 |
|
|---|
| 115 |
|
|---|
| 116 | void G4PenelopeIonisation::BuildPhysicsTable(const G4ParticleDefinition& aParticleType)
|
|---|
| 117 | {
|
|---|
| 118 | if(verboseLevel > 0) {
|
|---|
| 119 | G4cout << "G4PenelopeIonisation::BuildPhysicsTable start" << G4endl;
|
|---|
| 120 | }
|
|---|
| 121 |
|
|---|
| 122 | cutForDelta.clear();
|
|---|
| 123 |
|
|---|
| 124 | // Create and fill G4CrossSectionHandler once
|
|---|
| 125 | if ( crossSectionHandler != 0 ) delete crossSectionHandler;
|
|---|
| 126 | G4VDataSetAlgorithm* interpolation = new G4LogLogInterpolation();
|
|---|
| 127 | G4double lowKineticEnergy = GetLowerBoundEloss();
|
|---|
| 128 | G4double highKineticEnergy = GetUpperBoundEloss();
|
|---|
| 129 | G4int totBin = GetNbinEloss();
|
|---|
| 130 | crossSectionHandler = new G4PenelopeCrossSectionHandler(this,aParticleType,
|
|---|
| 131 | interpolation,
|
|---|
| 132 | lowKineticEnergy,
|
|---|
| 133 | highKineticEnergy,
|
|---|
| 134 | totBin);
|
|---|
| 135 |
|
|---|
| 136 | if (&aParticleType == G4Electron::Electron())
|
|---|
| 137 | {
|
|---|
| 138 | crossSectionHandler->LoadData("penelope/ion-cs-el-");
|
|---|
| 139 | }
|
|---|
| 140 | else if (&aParticleType == G4Positron::Positron())
|
|---|
| 141 | {
|
|---|
| 142 | crossSectionHandler->LoadData("penelope/ion-cs-po-");
|
|---|
| 143 | }
|
|---|
| 144 |
|
|---|
| 145 | if (verboseLevel > 0) {
|
|---|
| 146 | G4cout << GetProcessName() << " is created." << G4endl;
|
|---|
| 147 | }
|
|---|
| 148 |
|
|---|
| 149 | // Build loss table for Ionisation
|
|---|
| 150 |
|
|---|
| 151 | BuildLossTable(aParticleType);
|
|---|
| 152 |
|
|---|
| 153 | if(verboseLevel > 0) {
|
|---|
| 154 | G4cout << "The loss table is built"
|
|---|
| 155 | << G4endl;
|
|---|
| 156 | }
|
|---|
| 157 |
|
|---|
| 158 | if (&aParticleType==G4Electron::Electron()) {
|
|---|
| 159 |
|
|---|
| 160 | RecorderOfElectronProcess[CounterOfElectronProcess] = (*this).theLossTable;
|
|---|
| 161 | CounterOfElectronProcess++;
|
|---|
| 162 | PrintInfoDefinition();
|
|---|
| 163 |
|
|---|
| 164 | } else {
|
|---|
| 165 |
|
|---|
| 166 | RecorderOfPositronProcess[CounterOfPositronProcess] = (*this).theLossTable;
|
|---|
| 167 | CounterOfPositronProcess++;
|
|---|
| 168 | PrintInfoDefinition();
|
|---|
| 169 | }
|
|---|
| 170 |
|
|---|
| 171 | // Build mean free path data using cut values
|
|---|
| 172 |
|
|---|
| 173 | if( theMeanFreePath ) delete theMeanFreePath;
|
|---|
| 174 | theMeanFreePath = crossSectionHandler->
|
|---|
| 175 | BuildMeanFreePathForMaterials(&cutForDelta);
|
|---|
| 176 |
|
|---|
| 177 | if(verboseLevel > 0) {
|
|---|
| 178 | G4cout << "The MeanFreePath table is built"
|
|---|
| 179 | << G4endl;
|
|---|
| 180 | if(verboseLevel > 1) theMeanFreePath->PrintData();
|
|---|
| 181 | }
|
|---|
| 182 |
|
|---|
| 183 | // Build common DEDX table for all ionisation processes
|
|---|
| 184 |
|
|---|
| 185 | BuildDEDXTable(aParticleType);
|
|---|
| 186 |
|
|---|
| 187 | if (verboseLevel > 0) {
|
|---|
| 188 | G4cout << "G4PenelopeIonisation::BuildPhysicsTable end"
|
|---|
| 189 | << G4endl;
|
|---|
| 190 | }
|
|---|
| 191 | }
|
|---|
| 192 |
|
|---|
| 193 |
|
|---|
| 194 | void G4PenelopeIonisation::BuildLossTable(
|
|---|
| 195 | const G4ParticleDefinition& aParticleType)
|
|---|
| 196 | {
|
|---|
| 197 | // Build table for energy loss due to soft brems
|
|---|
| 198 | // the tables are built for *MATERIALS* binning is taken from LowEnergyLoss
|
|---|
| 199 |
|
|---|
| 200 | G4double lowKineticEnergy = GetLowerBoundEloss();
|
|---|
| 201 | G4double highKineticEnergy = GetUpperBoundEloss();
|
|---|
| 202 | size_t totBin = GetNbinEloss();
|
|---|
| 203 |
|
|---|
| 204 | // create table
|
|---|
| 205 |
|
|---|
| 206 | if (theLossTable) {
|
|---|
| 207 | theLossTable->clearAndDestroy();
|
|---|
| 208 | delete theLossTable;
|
|---|
| 209 | }
|
|---|
| 210 | const G4ProductionCutsTable* theCoupleTable=
|
|---|
| 211 | G4ProductionCutsTable::GetProductionCutsTable();
|
|---|
| 212 | size_t numOfCouples = theCoupleTable->GetTableSize();
|
|---|
| 213 | theLossTable = new G4PhysicsTable(numOfCouples);
|
|---|
| 214 |
|
|---|
| 215 | // Clean up the vector of cuts
|
|---|
| 216 |
|
|---|
| 217 | cutForDelta.clear();
|
|---|
| 218 |
|
|---|
| 219 | // Loop for materials
|
|---|
| 220 |
|
|---|
| 221 | for (size_t m=0; m<numOfCouples; m++) {
|
|---|
| 222 |
|
|---|
| 223 | // create physics vector and fill it
|
|---|
| 224 | G4PhysicsLogVector* aVector = new G4PhysicsLogVector(lowKineticEnergy,
|
|---|
| 225 | highKineticEnergy,
|
|---|
| 226 | totBin);
|
|---|
| 227 |
|
|---|
| 228 | // get material parameters needed for the energy loss calculation
|
|---|
| 229 | const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(m);
|
|---|
| 230 | const G4Material* material= couple->GetMaterial();
|
|---|
| 231 |
|
|---|
| 232 | // the cut cannot be below lowest limit
|
|---|
| 233 | G4double tCut = 0.0;
|
|---|
| 234 | tCut = (*(theCoupleTable->GetEnergyCutsVector(1)))[m];
|
|---|
| 235 | tCut = std::min(tCut,highKineticEnergy);
|
|---|
| 236 |
|
|---|
| 237 | cutForDelta.push_back(tCut);
|
|---|
| 238 |
|
|---|
| 239 | const G4ElementVector* theElementVector = material->GetElementVector();
|
|---|
| 240 | size_t NumberOfElements = material->GetNumberOfElements() ;
|
|---|
| 241 | const G4double* theAtomicNumDensityVector =
|
|---|
| 242 | material->GetAtomicNumDensityVector();
|
|---|
| 243 | const G4double electronVolumeDensity =
|
|---|
| 244 | material->GetTotNbOfElectPerVolume(); //electron density
|
|---|
| 245 | if(verboseLevel > 0) {
|
|---|
| 246 | G4cout << "Energy loss for material # " << m
|
|---|
| 247 | << " tCut(keV)= " << tCut/keV
|
|---|
| 248 | << G4endl;
|
|---|
| 249 | }
|
|---|
| 250 |
|
|---|
| 251 | // now comes the loop for the kinetic energy values
|
|---|
| 252 | for (size_t i = 0; i<totBin; i++) {
|
|---|
| 253 |
|
|---|
| 254 | G4double lowEdgeEnergy = aVector->GetLowEdgeEnergy(i);
|
|---|
| 255 | G4double ionloss = 0.;
|
|---|
| 256 |
|
|---|
| 257 | // loop for elements in the material
|
|---|
| 258 | for (size_t iel=0; iel<NumberOfElements; iel++ ) {
|
|---|
| 259 |
|
|---|
| 260 | G4int Z = (G4int)((*theElementVector)[iel]->GetZ());
|
|---|
| 261 | ionloss +=
|
|---|
| 262 | CalculateContinuous(lowEdgeEnergy,tCut,Z,electronVolumeDensity,
|
|---|
| 263 | aParticleType) * theAtomicNumDensityVector[iel];
|
|---|
| 264 |
|
|---|
| 265 | if(verboseLevel > 1) {
|
|---|
| 266 | G4cout << "Z= " << Z
|
|---|
| 267 | << " E(keV)= " << lowEdgeEnergy/keV
|
|---|
| 268 | << " loss= " << ionloss
|
|---|
| 269 | << " rho= " << theAtomicNumDensityVector[iel]
|
|---|
| 270 | << G4endl;
|
|---|
| 271 | }
|
|---|
| 272 | }
|
|---|
| 273 | aVector->PutValue(i,ionloss);
|
|---|
| 274 | }
|
|---|
| 275 | theLossTable->insert(aVector);
|
|---|
| 276 | }
|
|---|
| 277 | }
|
|---|
| 278 |
|
|---|
| 279 |
|
|---|
| 280 | G4VParticleChange* G4PenelopeIonisation::PostStepDoIt(const G4Track& track,
|
|---|
| 281 | const G4Step& step)
|
|---|
| 282 | {
|
|---|
| 283 | aParticleChange.Initialize(track);
|
|---|
| 284 |
|
|---|
| 285 | const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
|
|---|
| 286 | const G4DynamicParticle* incidentElectron = track.GetDynamicParticle();
|
|---|
| 287 | const G4Material* material = couple->GetMaterial();
|
|---|
| 288 | const G4double electronVolumeDensity =
|
|---|
| 289 | material->GetTotNbOfElectPerVolume(); //electron density
|
|---|
| 290 | const G4ParticleDefinition* aParticleType = track.GetDefinition();
|
|---|
| 291 | G4double kineticEnergy0 = incidentElectron->GetKineticEnergy();
|
|---|
| 292 | G4ParticleMomentum electronDirection0 = incidentElectron->GetMomentumDirection();
|
|---|
| 293 |
|
|---|
| 294 | //Inizialisation of variables
|
|---|
| 295 | kineticEnergy1=kineticEnergy0;
|
|---|
| 296 | cosThetaPrimary=1.0;
|
|---|
| 297 | energySecondary=0.0;
|
|---|
| 298 | cosThetaSecondary=1.0;
|
|---|
| 299 |
|
|---|
| 300 | G4int Z = crossSectionHandler->SelectRandomAtom(couple, kineticEnergy0);
|
|---|
| 301 | G4int index = couple->GetIndex();
|
|---|
| 302 | G4double tCut = cutForDelta[index];
|
|---|
| 303 |
|
|---|
| 304 | if (aParticleType==G4Electron::Electron()){
|
|---|
| 305 | CalculateDiscreteForElectrons(kineticEnergy0,tCut,Z,electronVolumeDensity);
|
|---|
| 306 | }
|
|---|
| 307 | else if (aParticleType==G4Positron::Positron()){
|
|---|
| 308 | CalculateDiscreteForPositrons(kineticEnergy0,tCut,Z,electronVolumeDensity);
|
|---|
| 309 | }
|
|---|
| 310 | // the method CalculateDiscrete() sets the private variables:
|
|---|
| 311 | // kineticEnergy1 = energy of the primary electron after the interaction
|
|---|
| 312 | // cosThetaPrimary = std::cos(theta) of the primary after the interaction
|
|---|
| 313 | // energySecondary = energy of the secondary electron
|
|---|
| 314 | // cosThetaSecondary = std::cos(theta) of the secondary
|
|---|
| 315 |
|
|---|
| 316 | if(energySecondary == 0.0)
|
|---|
| 317 | {
|
|---|
| 318 | return G4VContinuousDiscreteProcess::PostStepDoIt(track, step);
|
|---|
| 319 | }
|
|---|
| 320 |
|
|---|
| 321 | //Update the primary particle
|
|---|
| 322 | G4double sint = std::sqrt(1. - cosThetaPrimary*cosThetaPrimary);
|
|---|
| 323 | G4double phi = twopi * G4UniformRand();
|
|---|
| 324 | G4double dirx = sint * std::cos(phi);
|
|---|
| 325 | G4double diry = sint * std::sin(phi);
|
|---|
| 326 | G4double dirz = cosThetaPrimary;
|
|---|
| 327 |
|
|---|
| 328 | G4ThreeVector electronDirection1(dirx,diry,dirz);
|
|---|
| 329 | electronDirection1.rotateUz(electronDirection0);
|
|---|
| 330 | aParticleChange.ProposeMomentumDirection(electronDirection1) ;
|
|---|
| 331 |
|
|---|
| 332 | if (kineticEnergy1 > 0.)
|
|---|
| 333 | {
|
|---|
| 334 | aParticleChange.ProposeEnergy(kineticEnergy1) ;
|
|---|
| 335 | }
|
|---|
| 336 | else
|
|---|
| 337 | {
|
|---|
| 338 | aParticleChange.ProposeEnergy(0.) ;
|
|---|
| 339 | if (aParticleType->GetProcessManager()->GetAtRestProcessVector()->size())
|
|---|
| 340 | //In this case there is at least one AtRest process
|
|---|
| 341 | {
|
|---|
| 342 | aParticleChange.ProposeTrackStatus(fStopButAlive);
|
|---|
| 343 | }
|
|---|
| 344 | else
|
|---|
| 345 | {
|
|---|
| 346 | aParticleChange.ProposeTrackStatus(fStopAndKill);
|
|---|
| 347 | }
|
|---|
| 348 | }
|
|---|
| 349 |
|
|---|
| 350 | //Generate the delta day
|
|---|
| 351 | G4int iosc2 = 0;
|
|---|
| 352 | G4double ioniEnergy = 0.0;
|
|---|
| 353 | if (iOsc > 0) {
|
|---|
| 354 | ioniEnergy=(*(ionizationEnergy->find(Z)->second))[iOsc];
|
|---|
| 355 | iosc2 = (ionizationEnergy->find(Z)->second->size()) - iOsc; //they are in reversed order
|
|---|
| 356 | }
|
|---|
| 357 |
|
|---|
| 358 | const G4AtomicTransitionManager* transitionManager = G4AtomicTransitionManager::Instance();
|
|---|
| 359 | G4double bindingEnergy = 0.0;
|
|---|
| 360 | G4int shellId = 0;
|
|---|
| 361 | if (iOsc > 0){
|
|---|
| 362 | const G4AtomicShell* shell = transitionManager->Shell(Z,iosc2-1); // Modified by Alf
|
|---|
| 363 | bindingEnergy = shell->BindingEnergy();
|
|---|
| 364 | shellId = shell->ShellId();
|
|---|
| 365 | }
|
|---|
| 366 |
|
|---|
| 367 | G4double ionEnergy = bindingEnergy; //energy spent to ionise the atom according to G4dabatase
|
|---|
| 368 | G4double eKineticEnergy = energySecondary;
|
|---|
| 369 |
|
|---|
| 370 | //This is an awful thing: Penelope generates the fluorescence only for L and K shells
|
|---|
| 371 | //(i.e. Osc = 1 --> 4). For high-Z, the other shells can be quite relevant. In this case
|
|---|
| 372 | //one MUST ensure ''by hand'' the energy conservation. Then there is the other problem that
|
|---|
| 373 | //the fluorescence database of Penelope doesn not match that of Geant4.
|
|---|
| 374 |
|
|---|
| 375 | G4double energyBalance = kineticEnergy0 - kineticEnergy1 - energySecondary; //Penelope Balance
|
|---|
| 376 |
|
|---|
| 377 | if (std::abs(energyBalance) < 1*eV)
|
|---|
| 378 | {
|
|---|
| 379 | //in this case Penelope didn't subtract the fluorescence energy: do here by hand
|
|---|
| 380 | eKineticEnergy = energySecondary - bindingEnergy;
|
|---|
| 381 | }
|
|---|
| 382 | else
|
|---|
| 383 | {
|
|---|
| 384 | //Penelope subtracted the fluorescence, but one has to match the databases
|
|---|
| 385 | eKineticEnergy = energySecondary+ioniEnergy-bindingEnergy;
|
|---|
| 386 | }
|
|---|
| 387 |
|
|---|
| 388 | //Now generates the various secondaries
|
|---|
| 389 | size_t nTotPhotons=0;
|
|---|
| 390 | G4int nPhotons=0;
|
|---|
| 391 |
|
|---|
| 392 | const G4ProductionCutsTable* theCoupleTable=
|
|---|
| 393 | G4ProductionCutsTable::GetProductionCutsTable();
|
|---|
| 394 | size_t indx = couple->GetIndex();
|
|---|
| 395 | G4double cutg = (*(theCoupleTable->GetEnergyCutsVector(0)))[indx];
|
|---|
| 396 | cutg = std::min(cutForPhotons,cutg);
|
|---|
| 397 |
|
|---|
| 398 | G4double cute = (*(theCoupleTable->GetEnergyCutsVector(1)))[indx];
|
|---|
| 399 | cute = std::min(cutForPhotons,cute);
|
|---|
| 400 |
|
|---|
| 401 | std::vector<G4DynamicParticle*>* photonVector=0;
|
|---|
| 402 | G4DynamicParticle* aPhoton;
|
|---|
| 403 | if (Z>5 && (ionEnergy > cutg || ionEnergy > cute))
|
|---|
| 404 | {
|
|---|
| 405 | photonVector = deexcitationManager.GenerateParticles(Z,shellId);
|
|---|
| 406 | nTotPhotons = photonVector->size();
|
|---|
| 407 | for (size_t k=0;k<nTotPhotons;k++){
|
|---|
| 408 | aPhoton = (*photonVector)[k];
|
|---|
| 409 | if (aPhoton)
|
|---|
| 410 | {
|
|---|
| 411 | G4double itsCut = cutg;
|
|---|
| 412 | if (aPhoton->GetDefinition() == G4Electron::Electron()) itsCut = cute;
|
|---|
| 413 | G4double itsEnergy = aPhoton->GetKineticEnergy();
|
|---|
| 414 | if (itsEnergy > itsCut && itsEnergy <= ionEnergy)
|
|---|
| 415 | {
|
|---|
| 416 | nPhotons++;
|
|---|
| 417 | ionEnergy -= itsEnergy;
|
|---|
| 418 | }
|
|---|
| 419 | else
|
|---|
| 420 | {
|
|---|
| 421 | delete aPhoton;
|
|---|
| 422 | (*photonVector)[k]=0;
|
|---|
| 423 | }
|
|---|
| 424 | }
|
|---|
| 425 | }
|
|---|
| 426 | }
|
|---|
| 427 | G4double energyDeposit=ionEnergy; //il deposito locale e' quello che rimane
|
|---|
| 428 | G4int nbOfSecondaries=nPhotons;
|
|---|
| 429 |
|
|---|
| 430 |
|
|---|
| 431 | // Generate the delta ray
|
|---|
| 432 | G4double sin2 = std::sqrt(1. - cosThetaSecondary*cosThetaSecondary);
|
|---|
| 433 | G4double phi2 = twopi * G4UniformRand();
|
|---|
| 434 | G4DynamicParticle* electron = 0;
|
|---|
| 435 |
|
|---|
| 436 | G4double xEl = sin2 * std::cos(phi2);
|
|---|
| 437 | G4double yEl = sin2 * std::sin(phi2);
|
|---|
| 438 | G4double zEl = cosThetaSecondary;
|
|---|
| 439 | G4ThreeVector eDirection(xEl,yEl,zEl); //electron direction
|
|---|
| 440 | eDirection.rotateUz(electronDirection0);
|
|---|
| 441 |
|
|---|
| 442 | electron = new G4DynamicParticle (G4Electron::Electron(),
|
|---|
| 443 | eDirection,eKineticEnergy) ;
|
|---|
| 444 | nbOfSecondaries++;
|
|---|
| 445 |
|
|---|
| 446 | aParticleChange.SetNumberOfSecondaries(nbOfSecondaries);
|
|---|
| 447 | if (electron) aParticleChange.AddSecondary(electron);
|
|---|
| 448 |
|
|---|
| 449 | G4double energySumTest = kineticEnergy1 + eKineticEnergy;
|
|---|
| 450 |
|
|---|
| 451 | for (size_t ll=0;ll<nTotPhotons;ll++)
|
|---|
| 452 | {
|
|---|
| 453 | aPhoton = (*photonVector)[ll];
|
|---|
| 454 | if (aPhoton) {
|
|---|
| 455 | aParticleChange.AddSecondary(aPhoton);
|
|---|
| 456 | energySumTest += aPhoton->GetKineticEnergy();
|
|---|
| 457 | }
|
|---|
| 458 | }
|
|---|
| 459 | delete photonVector;
|
|---|
| 460 | if (energyDeposit < 0)
|
|---|
| 461 | {
|
|---|
| 462 | G4cout << "WARNING-"
|
|---|
| 463 | << "G4PenelopeIonisaition::PostStepDoIt - Negative energy deposit"
|
|---|
| 464 | << G4endl;
|
|---|
| 465 | energyDeposit=0;
|
|---|
| 466 | }
|
|---|
| 467 | energySumTest += energyDeposit;
|
|---|
| 468 | if (std::abs(energySumTest-kineticEnergy0)>1*eV)
|
|---|
| 469 | {
|
|---|
| 470 | G4cout << "WARNING-"
|
|---|
| 471 | << "G4PenelopeIonisaition::PostStepDoIt - Energy non conservation"
|
|---|
| 472 | << G4endl;
|
|---|
| 473 | G4cout << "Final energy - initial energy = " <<
|
|---|
| 474 | (energySumTest-kineticEnergy0)/eV << " eV" << G4endl;
|
|---|
| 475 | }
|
|---|
| 476 | aParticleChange.ProposeLocalEnergyDeposit(energyDeposit);
|
|---|
| 477 | return G4VContinuousDiscreteProcess::PostStepDoIt(track, step);
|
|---|
| 478 | }
|
|---|
| 479 |
|
|---|
| 480 |
|
|---|
| 481 | void G4PenelopeIonisation::PrintInfoDefinition()
|
|---|
| 482 | {
|
|---|
| 483 | G4String comments = "Total cross sections from EEDL database.";
|
|---|
| 484 | comments += "\n Delta energy sampled from a parametrised formula.";
|
|---|
| 485 | comments += "\n Implementation of the continuous dE/dx part.";
|
|---|
| 486 | comments += "\n At present it can be used for electrons and positrons ";
|
|---|
| 487 | comments += "in the energy range [250eV,100GeV].";
|
|---|
| 488 | comments += "\n The process must work with G4PenelopeBremsstrahlung.";
|
|---|
| 489 |
|
|---|
| 490 | G4cout << G4endl << GetProcessName() << ": " << comments << G4endl;
|
|---|
| 491 | }
|
|---|
| 492 |
|
|---|
| 493 | G4bool G4PenelopeIonisation::IsApplicable(const G4ParticleDefinition& particle)
|
|---|
| 494 | {
|
|---|
| 495 | return ( (&particle == G4Electron::Electron()) || (
|
|---|
| 496 | &particle == G4Positron::Positron()) );
|
|---|
| 497 | }
|
|---|
| 498 |
|
|---|
| 499 | G4double G4PenelopeIonisation::GetMeanFreePath(const G4Track& track,
|
|---|
| 500 | G4double, // previousStepSize
|
|---|
| 501 | G4ForceCondition* cond)
|
|---|
| 502 | {
|
|---|
| 503 | *cond = NotForced;
|
|---|
| 504 | G4int index = (track.GetMaterialCutsCouple())->GetIndex();
|
|---|
| 505 | const G4VEMDataSet* data = theMeanFreePath->GetComponent(index);
|
|---|
| 506 | G4double meanFreePath = data->FindValue(track.GetKineticEnergy());
|
|---|
| 507 | return meanFreePath;
|
|---|
| 508 | }
|
|---|
| 509 |
|
|---|
| 510 | void G4PenelopeIonisation::SetCutForLowEnSecPhotons(G4double cut)
|
|---|
| 511 | {
|
|---|
| 512 | cutForPhotons = cut;
|
|---|
| 513 | deexcitationManager.SetCutForSecondaryPhotons(cut);
|
|---|
| 514 | }
|
|---|
| 515 |
|
|---|
| 516 | void G4PenelopeIonisation::SetCutForLowEnSecElectrons(G4double cut)
|
|---|
| 517 | {
|
|---|
| 518 | cutForElectrons = cut;
|
|---|
| 519 | deexcitationManager.SetCutForAugerElectrons(cut);
|
|---|
| 520 | }
|
|---|
| 521 |
|
|---|
| 522 | void G4PenelopeIonisation::ActivateAuger(G4bool val)
|
|---|
| 523 | {
|
|---|
| 524 | deexcitationManager.ActivateAugerElectronProduction(val);
|
|---|
| 525 | }
|
|---|
| 526 |
|
|---|
| 527 |
|
|---|
| 528 | void G4PenelopeIonisation::CalculateDiscreteForElectrons(G4double ene,G4double cutoff,
|
|---|
| 529 | G4int Z,G4double electronVolumeDensity)
|
|---|
| 530 | {
|
|---|
| 531 | kineticEnergy1=ene;
|
|---|
| 532 | cosThetaPrimary=1.0;
|
|---|
| 533 | energySecondary=0.0;
|
|---|
| 534 | cosThetaSecondary=1.0;
|
|---|
| 535 | iOsc=-1;
|
|---|
| 536 | //constants
|
|---|
| 537 | G4double rb=ene+2.0*electron_mass_c2;
|
|---|
| 538 | G4double gamma = 1.0+ene/electron_mass_c2;
|
|---|
| 539 | G4double gamma2 = gamma*gamma;
|
|---|
| 540 | G4double beta2 = (gamma2-1.0)/gamma2;
|
|---|
| 541 | G4double amol = (gamma-1.0)*(gamma-1.0)/gamma2;
|
|---|
| 542 | G4double cps = ene*rb;
|
|---|
| 543 | G4double cp = std::sqrt(cps);
|
|---|
| 544 |
|
|---|
| 545 | G4double delta = CalculateDeltaFermi(ene,Z,electronVolumeDensity);
|
|---|
| 546 | G4double distantTransvCS0 = std::max(std::log(gamma2)-beta2-delta,0.0);
|
|---|
| 547 |
|
|---|
| 548 | G4double rl,rl1;
|
|---|
| 549 |
|
|---|
| 550 | if (cutoff > ene) return; //delta rays are not generated
|
|---|
| 551 |
|
|---|
| 552 | G4DataVector* qm = new G4DataVector();
|
|---|
| 553 | G4DataVector* cumulHardCS = new G4DataVector();
|
|---|
| 554 | G4DataVector* typeOfInteraction = new G4DataVector();
|
|---|
| 555 | G4DataVector* nbOfLevel = new G4DataVector();
|
|---|
| 556 |
|
|---|
| 557 | //Hard close collisions with outer shells
|
|---|
| 558 | G4double wmaxc = 0.5*ene;
|
|---|
| 559 | G4double closeCS0 = 0.0;
|
|---|
| 560 | G4double closeCS = 0.0;
|
|---|
| 561 | if (cutoff>0.1*eV)
|
|---|
| 562 | {
|
|---|
| 563 | rl=cutoff/ene;
|
|---|
| 564 | rl1=1.0-rl;
|
|---|
| 565 | if (rl < 0.5)
|
|---|
| 566 | closeCS0 = (amol*(0.5-rl)+(1.0/rl)-(1.0/rl1)+(1.0-amol)*std::log(rl/rl1))/ene;
|
|---|
| 567 | }
|
|---|
| 568 |
|
|---|
| 569 | // Cross sections for the different oscillators
|
|---|
| 570 |
|
|---|
| 571 | // totalHardCS contains the cumulative hard interaction cross section for the different
|
|---|
| 572 | // excitable levels and the different interaction channels (close, distant, etc.),
|
|---|
| 573 | // i.e.
|
|---|
| 574 | // cumulHardCS[0] = 0.0
|
|---|
| 575 | // cumulHardCS[1] = 1st excitable level (distant longitudinal only)
|
|---|
| 576 | // cumulHardCS[2] = 1st excitable level (distant longitudinal + transverse)
|
|---|
| 577 | // cumulHardCS[3] = 1st excitable level (distant longitudinal + transverse + close)
|
|---|
| 578 | // cumulHardCS[4] = 1st excitable level (all channels) + 2nd excitable level (distant long only)
|
|---|
| 579 | // etc.
|
|---|
| 580 | // This is used for sampling the atomic level which is ionised and the channel of the
|
|---|
| 581 | // interaction.
|
|---|
| 582 | //
|
|---|
| 583 | // For each index iFill of the cumulHardCS vector,
|
|---|
| 584 | // nbOfLevel[iFill] contains the current excitable atomic level and
|
|---|
| 585 | // typeOfInteraction[iFill] contains the current interaction channel, with the legenda:
|
|---|
| 586 | // 1 = distant longitudinal interaction
|
|---|
| 587 | // 2 = distant transverse interaction
|
|---|
| 588 | // 3 = close collision
|
|---|
| 589 | // 4 = close collision with outer shells (in this case nbOfLevel < 0 --> no binding energy)
|
|---|
| 590 |
|
|---|
| 591 |
|
|---|
| 592 | G4int nOscil = ionizationEnergy->find(Z)->second->size();
|
|---|
| 593 | G4double totalHardCS = 0.0;
|
|---|
| 594 | G4double involvedElectrons = 0.0;
|
|---|
| 595 | for (G4int i=0;i<nOscil;i++){
|
|---|
| 596 | G4double wi = (*(resonanceEnergy->find(Z)->second))[i];
|
|---|
| 597 | G4int occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i];
|
|---|
| 598 | //Distant excitations
|
|---|
| 599 | if (wi>cutoff && wi<ene)
|
|---|
| 600 | {
|
|---|
| 601 | if (wi>(1e-6*ene)){
|
|---|
| 602 | G4double cpp=std::sqrt((ene-wi)*(ene-wi+2.0*electron_mass_c2));
|
|---|
| 603 | qm->push_back(std::sqrt((cp-cpp)*(cp-cpp)+electron_mass_c2*electron_mass_c2)-electron_mass_c2);
|
|---|
| 604 | }
|
|---|
| 605 | else
|
|---|
| 606 | {
|
|---|
| 607 | qm->push_back((wi*wi)/(beta2*2.0*electron_mass_c2));
|
|---|
| 608 | }
|
|---|
| 609 | //verificare che quando arriva qui il vettore ha SEMPRE l'i-esimo elemento
|
|---|
| 610 | if ((*qm)[i] < wi)
|
|---|
| 611 | {
|
|---|
| 612 |
|
|---|
| 613 | G4double distantLongitCS = occupNb*std::log(wi*((*qm)[i]+2.0*electron_mass_c2)/
|
|---|
| 614 | ((*qm)[i]*(wi+2.0*electron_mass_c2)))/wi;
|
|---|
| 615 | cumulHardCS->push_back(totalHardCS);
|
|---|
| 616 | typeOfInteraction->push_back(1.0); //distant longitudinal
|
|---|
| 617 | nbOfLevel->push_back((G4double) i); //only excitable level are counted
|
|---|
| 618 | totalHardCS += distantLongitCS;
|
|---|
| 619 |
|
|---|
| 620 | G4double distantTransvCS = occupNb*distantTransvCS0/wi;
|
|---|
| 621 |
|
|---|
| 622 | cumulHardCS->push_back(totalHardCS);
|
|---|
| 623 | typeOfInteraction->push_back(2.0); //distant tranverse
|
|---|
| 624 | nbOfLevel->push_back((G4double) i);
|
|---|
| 625 | totalHardCS += distantTransvCS;
|
|---|
| 626 | }
|
|---|
| 627 | }
|
|---|
| 628 | else
|
|---|
| 629 | {
|
|---|
| 630 | qm->push_back(wi);
|
|---|
| 631 | }
|
|---|
| 632 | //close collisions
|
|---|
| 633 | if(wi < wmaxc){
|
|---|
| 634 | if (wi < cutoff) {
|
|---|
| 635 | involvedElectrons += occupNb;
|
|---|
| 636 | }
|
|---|
| 637 | else
|
|---|
| 638 | {
|
|---|
| 639 | rl=wi/ene;
|
|---|
| 640 | rl1=1.0-rl;
|
|---|
| 641 | closeCS = occupNb*(amol*(0.5-rl)+(1.0/rl)-(1.0/rl1)+(1.0-amol)*std::log(rl/rl1))/ene;
|
|---|
| 642 | cumulHardCS->push_back(totalHardCS);
|
|---|
| 643 | typeOfInteraction->push_back(3.0); //close
|
|---|
| 644 | nbOfLevel->push_back((G4double) i);
|
|---|
| 645 | totalHardCS += closeCS;
|
|---|
| 646 | }
|
|---|
| 647 | }
|
|---|
| 648 | } // loop on the levels
|
|---|
| 649 |
|
|---|
| 650 | cumulHardCS->push_back(totalHardCS);
|
|---|
| 651 | typeOfInteraction->push_back(4.0); //close interaction with outer shells
|
|---|
| 652 | nbOfLevel->push_back(-1.0);
|
|---|
| 653 | totalHardCS += involvedElectrons*closeCS0;
|
|---|
| 654 | cumulHardCS->push_back(totalHardCS); //this is the final value of the totalHardCS
|
|---|
| 655 |
|
|---|
| 656 | if (totalHardCS < 1e-30) {
|
|---|
| 657 | kineticEnergy1=ene;
|
|---|
| 658 | cosThetaPrimary=1.0;
|
|---|
| 659 | energySecondary=0.0;
|
|---|
| 660 | cosThetaSecondary=0.0;
|
|---|
| 661 | iOsc=-1;
|
|---|
| 662 | delete qm;
|
|---|
| 663 | delete cumulHardCS;
|
|---|
| 664 | delete typeOfInteraction;
|
|---|
| 665 | delete nbOfLevel;
|
|---|
| 666 | return;
|
|---|
| 667 | }
|
|---|
| 668 |
|
|---|
| 669 |
|
|---|
| 670 | //Selection of the active oscillator on the basis of the cumulative cross sections
|
|---|
| 671 | G4double TST = totalHardCS*G4UniformRand();
|
|---|
| 672 | G4int is=0;
|
|---|
| 673 | G4int js= nbOfLevel->size();
|
|---|
| 674 | do{
|
|---|
| 675 | G4int it=(is+js)/2;
|
|---|
| 676 | if (TST > (*cumulHardCS)[it]) is=it;
|
|---|
| 677 | if (TST <= (*cumulHardCS)[it]) js=it;
|
|---|
| 678 | }while((js-is) > 1);
|
|---|
| 679 |
|
|---|
| 680 | G4double UII=0.0;
|
|---|
| 681 | G4double rkc=cutoff/ene;
|
|---|
| 682 | G4double dde;
|
|---|
| 683 | G4int kks;
|
|---|
| 684 |
|
|---|
| 685 | G4double sampledInteraction = (*typeOfInteraction)[is];
|
|---|
| 686 | iOsc = (G4int) (*nbOfLevel)[is];
|
|---|
| 687 |
|
|---|
| 688 | //Generates the final state according to the sampled level and
|
|---|
| 689 | //interaction channel
|
|---|
| 690 |
|
|---|
| 691 | if (sampledInteraction == 1.0) //Hard distant longitudinal collisions
|
|---|
| 692 | {
|
|---|
| 693 | dde= (*(resonanceEnergy->find(Z)->second))[iOsc];
|
|---|
| 694 | kineticEnergy1=ene-dde;
|
|---|
| 695 | G4double qs=(*qm)[iOsc]/(1.0+((*qm)[iOsc]/(2.0*electron_mass_c2)));
|
|---|
| 696 | G4double q=qs/(std::pow((qs/dde)*(1.0+(0.5*dde/electron_mass_c2)),G4UniformRand())-(0.5*qs/electron_mass_c2));
|
|---|
| 697 | G4double qtrev = q*(q+2.0*electron_mass_c2);
|
|---|
| 698 | G4double cpps = kineticEnergy1*(kineticEnergy1+2.0*electron_mass_c2);
|
|---|
| 699 | cosThetaPrimary = (cpps+cps-qtrev)/(2.0*cp*std::sqrt(cpps));
|
|---|
| 700 | if (cosThetaPrimary>1.0) cosThetaPrimary=1.0;
|
|---|
| 701 | //Energy and emission angle of the delta ray
|
|---|
| 702 | kks = (G4int) (*(shellFlag->find(Z)->second))[iOsc];
|
|---|
| 703 | if (kks>4)
|
|---|
| 704 | {
|
|---|
| 705 | energySecondary=dde;
|
|---|
| 706 | }
|
|---|
| 707 | else
|
|---|
| 708 | {
|
|---|
| 709 | energySecondary=dde-(*(ionizationEnergy->find(Z)->second))[iOsc];
|
|---|
| 710 | }
|
|---|
| 711 | cosThetaSecondary = 0.5*(dde*(ene+rb-dde)+qtrev)/std::sqrt(cps*qtrev);
|
|---|
| 712 | if (cosThetaSecondary>1.0) cosThetaSecondary=1.0;
|
|---|
| 713 | }
|
|---|
| 714 |
|
|---|
| 715 | else if (sampledInteraction == 2.0) //Hard distant transverse collisions
|
|---|
| 716 | {
|
|---|
| 717 | dde=(*(resonanceEnergy->find(Z)->second))[iOsc];
|
|---|
| 718 | kineticEnergy1=ene-dde;
|
|---|
| 719 | cosThetaPrimary=1.0;
|
|---|
| 720 | //Energy and emission angle of the delta ray
|
|---|
| 721 | kks = (G4int) (*(shellFlag->find(Z)->second))[iOsc];
|
|---|
| 722 | if (kks>4)
|
|---|
| 723 | {
|
|---|
| 724 | energySecondary=dde;
|
|---|
| 725 | }
|
|---|
| 726 | else
|
|---|
| 727 | {
|
|---|
| 728 | energySecondary=dde-(*(ionizationEnergy->find(Z)->second))[iOsc];
|
|---|
| 729 | }
|
|---|
| 730 | cosThetaSecondary = 1.0;
|
|---|
| 731 | }
|
|---|
| 732 |
|
|---|
| 733 | else if (sampledInteraction == 3.0 || sampledInteraction == 4.0) //Close interaction
|
|---|
| 734 | {
|
|---|
| 735 | if (sampledInteraction == 4.0) //interaction with inner shells
|
|---|
| 736 | {
|
|---|
| 737 | UII=0.0;
|
|---|
| 738 | rkc = cutoff/ene;
|
|---|
| 739 | iOsc = -1;
|
|---|
| 740 | }
|
|---|
| 741 | else
|
|---|
| 742 | {
|
|---|
| 743 | kks = (G4int) (*(shellFlag->find(Z)->second))[iOsc];
|
|---|
| 744 | if (kks > 4) {
|
|---|
| 745 | UII=0.0;
|
|---|
| 746 | }
|
|---|
| 747 | else
|
|---|
| 748 | {
|
|---|
| 749 | UII = (*(ionizationEnergy->find(Z)->second))[iOsc];
|
|---|
| 750 | }
|
|---|
| 751 | rkc = (*(resonanceEnergy->find(Z)->second))[iOsc]/ene;
|
|---|
| 752 | }
|
|---|
| 753 | G4double A = 0.5*amol;
|
|---|
| 754 | G4double arkc = A*0.5*rkc;
|
|---|
| 755 | G4double phi,rk2,rk,rkf;
|
|---|
| 756 | do{
|
|---|
| 757 | G4double fb = (1.0+arkc)*G4UniformRand();
|
|---|
| 758 | if (fb<1.0)
|
|---|
| 759 | {
|
|---|
| 760 | rk=rkc/(1.0-fb*(1.0-(rkc*2.0)));
|
|---|
| 761 | }
|
|---|
| 762 | else{
|
|---|
| 763 | rk = rkc+(fb-1.0)*(0.5-rkc)/arkc;
|
|---|
| 764 | }
|
|---|
| 765 | rk2 = rk*rk;
|
|---|
| 766 | rkf = rk/(1.0-rk);
|
|---|
| 767 | phi = 1.0+(rkf*rkf)-rkf+amol*(rk2+rkf);
|
|---|
| 768 | }while ((G4UniformRand()*(1.0+A*rk2)) > phi);
|
|---|
| 769 | //Energy and scattering angle (primary electron);
|
|---|
| 770 | kineticEnergy1 = ene*(1.0-rk);
|
|---|
| 771 | cosThetaPrimary = std::sqrt(kineticEnergy1*rb/(ene*(rb-(rk*ene))));
|
|---|
| 772 | //Energy and scattering angle of the delta ray
|
|---|
| 773 | energySecondary = ene-kineticEnergy1-UII;
|
|---|
| 774 | cosThetaSecondary = std::sqrt(rk*ene*rb/(ene*(rk*ene+2.0*electron_mass_c2)));
|
|---|
| 775 | }
|
|---|
| 776 |
|
|---|
| 777 | else
|
|---|
| 778 |
|
|---|
| 779 | {
|
|---|
| 780 | G4String excep = "G4PenelopeIonisation - Error in the calculation of the final state";
|
|---|
| 781 | G4Exception(excep);
|
|---|
| 782 | }
|
|---|
| 783 |
|
|---|
| 784 | delete qm;
|
|---|
| 785 | delete cumulHardCS;
|
|---|
| 786 | delete typeOfInteraction;
|
|---|
| 787 | delete nbOfLevel;
|
|---|
| 788 |
|
|---|
| 789 | return;
|
|---|
| 790 | }
|
|---|
| 791 |
|
|---|
| 792 | void G4PenelopeIonisation::ReadData()
|
|---|
| 793 | {
|
|---|
| 794 | char* path = getenv("G4LEDATA");
|
|---|
| 795 | if (!path)
|
|---|
| 796 | {
|
|---|
| 797 | G4String excep = "G4PenelopeIonisation - G4LEDATA environment variable not set!";
|
|---|
| 798 | G4Exception(excep);
|
|---|
| 799 | }
|
|---|
| 800 | G4String pathString(path);
|
|---|
| 801 | G4String pathFile = pathString + "/penelope/ion-pen.dat";
|
|---|
| 802 | std::ifstream file(pathFile);
|
|---|
| 803 | std::filebuf* lsdp = file.rdbuf();
|
|---|
| 804 |
|
|---|
| 805 | if (!(lsdp->is_open()))
|
|---|
| 806 | {
|
|---|
| 807 | G4String excep = "G4PenelopeIonisation - data file " + pathFile + " not found!";
|
|---|
| 808 | G4Exception(excep);
|
|---|
| 809 | }
|
|---|
| 810 |
|
|---|
| 811 | G4int k1,test,test1,k2,k3;
|
|---|
| 812 | G4double a1,a2,a3,a4;
|
|---|
| 813 | G4int Z=1,nLevels=0;
|
|---|
| 814 | G4DataVector* x1;
|
|---|
| 815 | G4DataVector* x2;
|
|---|
| 816 | G4DataVector* x3;
|
|---|
| 817 | G4DataVector* x4;
|
|---|
| 818 |
|
|---|
| 819 | do{
|
|---|
| 820 | x1 = new G4DataVector;
|
|---|
| 821 | x2 = new G4DataVector;
|
|---|
| 822 | x3 = new G4DataVector;
|
|---|
| 823 | x4 = new G4DataVector;
|
|---|
| 824 | file >> Z >> nLevels;
|
|---|
| 825 | for (G4int h=0;h<nLevels;h++){
|
|---|
| 826 | //index,occup number,ion energy,res energy,fj0,kz,shell flag
|
|---|
| 827 | file >> k1 >> a1 >> a2 >> a3 >> a4 >> k2 >> k3;
|
|---|
| 828 | x1->push_back(a1);
|
|---|
| 829 | x2->push_back(a2);
|
|---|
| 830 | x3->push_back(a3);
|
|---|
| 831 | x4->push_back((G4double) k3);
|
|---|
| 832 | }
|
|---|
| 833 | occupationNumber->insert(std::make_pair(Z,x1));
|
|---|
| 834 | ionizationEnergy->insert(std::make_pair(Z,x2));
|
|---|
| 835 | resonanceEnergy->insert(std::make_pair(Z,x3));
|
|---|
| 836 | shellFlag->insert(std::make_pair(Z,x4));
|
|---|
| 837 | file >> test >> test1; //-1 -1 close the data for each Z
|
|---|
| 838 | if (test > 0) {
|
|---|
| 839 | G4String excep = "G4PenelopeIonisation - data file corrupted!";
|
|---|
| 840 | G4Exception(excep);
|
|---|
| 841 | }
|
|---|
| 842 | }while (test != -2); //the very last Z is closed with -2 instead of -1
|
|---|
| 843 | }
|
|---|
| 844 |
|
|---|
| 845 |
|
|---|
| 846 | G4double G4PenelopeIonisation::CalculateDeltaFermi(G4double ene,G4int Z,
|
|---|
| 847 | G4double electronVolumeDensity)
|
|---|
| 848 | {
|
|---|
| 849 | G4double plasmaEnergyCoefficient = 1.377e-39; //(e*hbar)^2/(epsilon0*electron_mass)
|
|---|
| 850 | G4double plasmaEnergySquared = plasmaEnergyCoefficient*(electronVolumeDensity*m3);
|
|---|
| 851 | // std::sqrt(plasmaEnergySquared) is the plasma energy of the solid (MeV)
|
|---|
| 852 | G4double gam = 1.0+ene/electron_mass_c2;
|
|---|
| 853 | G4double gam2=gam*gam;
|
|---|
| 854 | G4double delta = 0.0;
|
|---|
| 855 |
|
|---|
| 856 | //Density effect
|
|---|
| 857 | G4double TST = ((G4double) Z)/(gam2*plasmaEnergySquared);
|
|---|
| 858 |
|
|---|
| 859 | G4double wl2 = 0.0;
|
|---|
| 860 | G4double fdel=0.0;
|
|---|
| 861 | G4double wr=0;
|
|---|
| 862 | G4double help1=0.0;
|
|---|
| 863 | size_t nbOsc = resonanceEnergy->find(Z)->second->size();
|
|---|
| 864 | for(size_t i=0;i<nbOsc;i++)
|
|---|
| 865 | {
|
|---|
| 866 | G4int occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i];
|
|---|
| 867 | wr = (*(resonanceEnergy->find(Z)->second))[i];
|
|---|
| 868 | fdel += occupNb/(wr*wr+wl2);
|
|---|
| 869 | }
|
|---|
| 870 | if (fdel < TST) return delta;
|
|---|
| 871 | help1 = (*(resonanceEnergy->find(Z)->second))[nbOsc-1];
|
|---|
| 872 | wl2 = help1*help1;
|
|---|
| 873 | do{
|
|---|
| 874 | wl2=wl2*2.0;
|
|---|
| 875 | fdel = 0.0;
|
|---|
| 876 | for (size_t ii=0;ii<nbOsc;ii++){
|
|---|
| 877 | G4int occupNb = (G4int) (*(occupationNumber->find(Z)->second))[ii];
|
|---|
| 878 | wr = (*(resonanceEnergy->find(Z)->second))[ii];
|
|---|
| 879 | fdel += occupNb/(wr*wr+wl2);
|
|---|
| 880 | }
|
|---|
| 881 | }while (fdel > TST);
|
|---|
| 882 | G4double wl2l=0.0;
|
|---|
| 883 | G4double wl2u = wl2;
|
|---|
| 884 | G4double control = 0.0;
|
|---|
| 885 | do{
|
|---|
| 886 | wl2=0.5*(wl2l+wl2u);
|
|---|
| 887 | fdel = 0.0;
|
|---|
| 888 | for (size_t jj=0;jj<nbOsc;jj++){
|
|---|
| 889 | G4int occupNb = (G4int) (*(occupationNumber->find(Z)->second))[jj];
|
|---|
| 890 | wr = (*(resonanceEnergy->find(Z)->second))[jj];
|
|---|
| 891 | fdel += occupNb/(wr*wr+wl2);
|
|---|
| 892 | }
|
|---|
| 893 | if (fdel > TST)
|
|---|
| 894 | {
|
|---|
| 895 | wl2l = wl2;
|
|---|
| 896 | }
|
|---|
| 897 | else
|
|---|
| 898 | {
|
|---|
| 899 | wl2u = wl2;
|
|---|
| 900 | }
|
|---|
| 901 | control = wl2u-wl2l-wl2*1e-12;
|
|---|
| 902 | }while(control>0);
|
|---|
| 903 |
|
|---|
| 904 | //Density correction effect
|
|---|
| 905 | for (size_t kk=0;kk<nbOsc;kk++){
|
|---|
| 906 | G4int occupNb = (G4int) (*(occupationNumber->find(Z)->second))[kk];
|
|---|
| 907 | wr = (*(resonanceEnergy->find(Z)->second))[kk];
|
|---|
| 908 | delta += occupNb*std::log(1.0+wl2/(wr*wr));
|
|---|
| 909 | }
|
|---|
| 910 | delta = (delta/((G4double) Z))-wl2/(gam2*plasmaEnergySquared);
|
|---|
| 911 | return delta;
|
|---|
| 912 | }
|
|---|
| 913 |
|
|---|
| 914 | G4double G4PenelopeIonisation::CalculateContinuous(G4double ene,G4double cutoff,
|
|---|
| 915 | G4int Z,G4double electronVolumeDensity,
|
|---|
| 916 | const G4ParticleDefinition& particle)
|
|---|
| 917 | {
|
|---|
| 918 | //Constants
|
|---|
| 919 | G4double gamma = 1.0+ene/electron_mass_c2;
|
|---|
| 920 | G4double gamma2 = gamma*gamma;
|
|---|
| 921 | G4double beta2 = (gamma2-1.0)/gamma2;
|
|---|
| 922 | G4double constant = pi*classic_electr_radius*classic_electr_radius
|
|---|
| 923 | *2.0*electron_mass_c2/beta2;
|
|---|
| 924 |
|
|---|
| 925 |
|
|---|
| 926 | G4double delta = CalculateDeltaFermi(ene,Z,electronVolumeDensity);
|
|---|
| 927 | G4int nbOsc = (G4int) resonanceEnergy->find(Z)->second->size();
|
|---|
| 928 | G4double S1 = 0.0;
|
|---|
| 929 | G4double stoppingPower = 0.0;
|
|---|
| 930 | for (G4int i=0;i<nbOsc;i++){
|
|---|
| 931 | G4double resEnergy = (*(resonanceEnergy->find(Z)->second))[i];
|
|---|
| 932 | if (&particle == G4Electron::Electron())
|
|---|
| 933 | {
|
|---|
| 934 | S1 = CalculateStoppingPowerForElectrons(ene,resEnergy,delta,cutoff);
|
|---|
| 935 | }
|
|---|
| 936 | else if (&particle == G4Positron::Positron())
|
|---|
| 937 | {
|
|---|
| 938 | S1 = CalculateStoppingPowerForPositrons(ene,resEnergy,delta,cutoff);
|
|---|
| 939 | }
|
|---|
| 940 | G4double occupNb = (*(occupationNumber->find(Z)->second))[i];
|
|---|
| 941 | stoppingPower += occupNb*constant*S1;
|
|---|
| 942 | }
|
|---|
| 943 |
|
|---|
| 944 | return stoppingPower;
|
|---|
| 945 | }
|
|---|
| 946 |
|
|---|
| 947 | G4double G4PenelopeIonisation::CalculateStoppingPowerForElectrons(G4double ene,G4double resEne,
|
|---|
| 948 | G4double delta,G4double cutoff)
|
|---|
| 949 | {
|
|---|
| 950 | //Calculate constants
|
|---|
| 951 | G4double gamma = 1.0+ene/electron_mass_c2;
|
|---|
| 952 | G4double gamma2 = gamma*gamma;
|
|---|
| 953 | G4double beta2 = (gamma2-1.0)/gamma2;
|
|---|
| 954 | G4double cps = ene*(ene+2.0*electron_mass_c2);
|
|---|
| 955 | G4double amol = (gamma-1.0)*(gamma-1.0)/gamma2;
|
|---|
| 956 | G4double sPower = 0.0;
|
|---|
| 957 | if (ene < resEne) return sPower;
|
|---|
| 958 |
|
|---|
| 959 | //Distant interactions
|
|---|
| 960 | G4double cp1s = (ene-resEne)*(ene-resEne+2.0*electron_mass_c2);
|
|---|
| 961 | G4double cp1 = std::sqrt(cp1s);
|
|---|
| 962 | G4double cp = std::sqrt(cps);
|
|---|
| 963 | G4double sdLong=0.0, sdTrans = 0.0, sdDist=0.0;
|
|---|
| 964 |
|
|---|
| 965 | //Distant longitudinal interactions
|
|---|
| 966 | G4double qm = 0.0;
|
|---|
| 967 |
|
|---|
| 968 | if (resEne > ene*(1e-6))
|
|---|
| 969 | {
|
|---|
| 970 | qm = std::sqrt((cp-cp1)*(cp-cp1)+(electron_mass_c2*electron_mass_c2))-electron_mass_c2;
|
|---|
| 971 | }
|
|---|
| 972 | else
|
|---|
| 973 | {
|
|---|
| 974 | qm = resEne*resEne/(beta2*2.0*electron_mass_c2);
|
|---|
| 975 | qm = qm*(1.0-0.5*qm/electron_mass_c2);
|
|---|
| 976 | }
|
|---|
| 977 |
|
|---|
| 978 | if (qm < resEne)
|
|---|
| 979 | {
|
|---|
| 980 | sdLong = std::log(resEne*(qm+2.0*electron_mass_c2)/(qm*(resEne+2.0*electron_mass_c2)));
|
|---|
| 981 | }
|
|---|
| 982 | else
|
|---|
| 983 | {
|
|---|
| 984 | sdLong = 0.0;
|
|---|
| 985 | }
|
|---|
| 986 |
|
|---|
| 987 | if (sdLong > 0) {
|
|---|
| 988 | sdTrans = std::max(std::log(gamma2)-beta2-delta,0.0);
|
|---|
| 989 | sdDist = sdTrans + sdLong;
|
|---|
| 990 | if (cutoff > resEne) sPower = sdDist;
|
|---|
| 991 | }
|
|---|
| 992 |
|
|---|
| 993 |
|
|---|
| 994 | // Close collisions (Moeller's cross section)
|
|---|
| 995 | G4double wl = std::max(cutoff,resEne);
|
|---|
| 996 | G4double wu = 0.5*ene;
|
|---|
| 997 |
|
|---|
| 998 | if (wl < (wu-1*eV)) wu=wl;
|
|---|
| 999 | wl = resEne;
|
|---|
| 1000 | if (wl > (wu-1*eV)) return sPower;
|
|---|
| 1001 | sPower += std::log(wu/wl)+(ene/(ene-wu))-(ene/(ene-wl))
|
|---|
| 1002 | + (2.0 - amol)*std::log((ene-wu)/(ene-wl))
|
|---|
| 1003 | + amol*((wu*wu)-(wl*wl))/(2.0*ene*ene);
|
|---|
| 1004 |
|
|---|
| 1005 | return sPower;
|
|---|
| 1006 | }
|
|---|
| 1007 |
|
|---|
| 1008 | G4double G4PenelopeIonisation::CalculateStoppingPowerForPositrons(G4double ene,G4double resEne,
|
|---|
| 1009 | G4double delta,G4double cutoff)
|
|---|
| 1010 | {
|
|---|
| 1011 | //Calculate constants
|
|---|
| 1012 | G4double gamma = 1.0+ene/electron_mass_c2;
|
|---|
| 1013 | G4double gamma2 = gamma*gamma;
|
|---|
| 1014 | G4double beta2 = (gamma2-1.0)/gamma2;
|
|---|
| 1015 | G4double cps = ene*(ene+2.0*electron_mass_c2);
|
|---|
| 1016 | G4double amol = (ene/(ene+electron_mass_c2)) * (ene/(ene+electron_mass_c2));
|
|---|
| 1017 | G4double help = (gamma+1.0)*(gamma+1.0);
|
|---|
| 1018 | G4double bha1 = amol*(2.0*help-1.0)/(gamma2-1.0);
|
|---|
| 1019 | G4double bha2 = amol*(3.0+1.0/help);
|
|---|
| 1020 | G4double bha3 = amol*2.0*gamma*(gamma-1.0)/help;
|
|---|
| 1021 | G4double bha4 = amol*(gamma-1.0)*(gamma-1.0)/help;
|
|---|
| 1022 |
|
|---|
| 1023 | G4double sPower = 0.0;
|
|---|
| 1024 | if (ene < resEne) return sPower;
|
|---|
| 1025 |
|
|---|
| 1026 | //Distant interactions
|
|---|
| 1027 | G4double cp1s = (ene-resEne)*(ene-resEne+2.0*electron_mass_c2);
|
|---|
| 1028 | G4double cp1 = std::sqrt(cp1s);
|
|---|
| 1029 | G4double cp = std::sqrt(cps);
|
|---|
| 1030 | G4double sdLong=0.0, sdTrans = 0.0, sdDist=0.0;
|
|---|
| 1031 |
|
|---|
| 1032 | //Distant longitudinal interactions
|
|---|
| 1033 | G4double qm = 0.0;
|
|---|
| 1034 |
|
|---|
| 1035 | if (resEne > ene*(1e-6))
|
|---|
| 1036 | {
|
|---|
| 1037 | qm = std::sqrt((cp-cp1)*(cp-cp1)+(electron_mass_c2*electron_mass_c2))-electron_mass_c2;
|
|---|
| 1038 | }
|
|---|
| 1039 | else
|
|---|
| 1040 | {
|
|---|
| 1041 | qm = resEne*resEne/(beta2*2.0*electron_mass_c2);
|
|---|
| 1042 | qm = qm*(1.0-0.5*qm/electron_mass_c2);
|
|---|
| 1043 | }
|
|---|
| 1044 |
|
|---|
| 1045 | if (qm < resEne)
|
|---|
| 1046 | {
|
|---|
| 1047 | sdLong = std::log(resEne*(qm+2.0*electron_mass_c2)/(qm*(resEne+2.0*electron_mass_c2)));
|
|---|
| 1048 | }
|
|---|
| 1049 | else
|
|---|
| 1050 | {
|
|---|
| 1051 | sdLong = 0.0;
|
|---|
| 1052 | }
|
|---|
| 1053 |
|
|---|
| 1054 | if (sdLong > 0) {
|
|---|
| 1055 | sdTrans = std::max(std::log(gamma2)-beta2-delta,0.0);
|
|---|
| 1056 | sdDist = sdTrans + sdLong;
|
|---|
| 1057 | if (cutoff > resEne) sPower = sdDist;
|
|---|
| 1058 | }
|
|---|
| 1059 |
|
|---|
| 1060 |
|
|---|
| 1061 | // Close collisions (Bhabha's cross section)
|
|---|
| 1062 | G4double wl = std::max(cutoff,resEne);
|
|---|
| 1063 | G4double wu = ene;
|
|---|
| 1064 |
|
|---|
| 1065 | if (wl < (wu-1*eV)) wu=wl;
|
|---|
| 1066 | wl = resEne;
|
|---|
| 1067 | if (wl > (wu-1*eV)) return sPower;
|
|---|
| 1068 | sPower += std::log(wu/wl)-bha1*(wu-wl)/ene
|
|---|
| 1069 | + bha2*((wu*wu)-(wl*wl))/(2.0*ene*ene)
|
|---|
| 1070 | - bha3*((wu*wu*wu)-(wl*wl*wl))/(3.0*ene*ene*ene)
|
|---|
| 1071 | + bha4*((wu*wu*wu*wu)-(wl*wl*wl*wl))/(4.0*ene*ene*ene*ene);
|
|---|
| 1072 |
|
|---|
| 1073 | return sPower;
|
|---|
| 1074 | }
|
|---|
| 1075 |
|
|---|
| 1076 | void G4PenelopeIonisation::CalculateDiscreteForPositrons(G4double ene,G4double cutoff,
|
|---|
| 1077 | G4int Z,G4double electronVolumeDensity)
|
|---|
| 1078 |
|
|---|
| 1079 | {
|
|---|
| 1080 | kineticEnergy1=ene;
|
|---|
| 1081 | cosThetaPrimary=1.0;
|
|---|
| 1082 | energySecondary=0.0;
|
|---|
| 1083 | cosThetaSecondary=1.0;
|
|---|
| 1084 | iOsc=-1;
|
|---|
| 1085 | //constants
|
|---|
| 1086 | G4double rb=ene+2.0*electron_mass_c2;
|
|---|
| 1087 | G4double gamma = 1.0+ene/electron_mass_c2;
|
|---|
| 1088 | G4double gamma2 = gamma*gamma;
|
|---|
| 1089 | G4double beta2 = (gamma2-1.0)/gamma2;
|
|---|
| 1090 | G4double amol = (gamma-1.0)*(gamma-1.0)/gamma2;
|
|---|
| 1091 | G4double cps = ene*rb;
|
|---|
| 1092 | G4double cp = std::sqrt(cps);
|
|---|
| 1093 | G4double help = (gamma+1.0)*(gamma+1.0);
|
|---|
| 1094 | G4double bha1 = amol*(2.0*help-1.0)/(gamma2-1.0);
|
|---|
| 1095 | G4double bha2 = amol*(3.0+1.0/help);
|
|---|
| 1096 | G4double bha3 = amol*2.0*gamma*(gamma-1.0)/help;
|
|---|
| 1097 | G4double bha4 = amol*(gamma-1.0)*(gamma-1.0)/help;
|
|---|
| 1098 |
|
|---|
| 1099 | G4double delta = CalculateDeltaFermi(ene,Z,electronVolumeDensity);
|
|---|
| 1100 | G4double distantTransvCS0 = std::max(std::log(gamma2)-beta2-delta,0.0);
|
|---|
| 1101 |
|
|---|
| 1102 | G4double rl,rl1;
|
|---|
| 1103 |
|
|---|
| 1104 | if (cutoff > ene) return; //delta rays are not generated
|
|---|
| 1105 |
|
|---|
| 1106 | G4DataVector* qm = new G4DataVector();
|
|---|
| 1107 | G4DataVector* cumulHardCS = new G4DataVector();
|
|---|
| 1108 | G4DataVector* typeOfInteraction = new G4DataVector();
|
|---|
| 1109 | G4DataVector* nbOfLevel = new G4DataVector();
|
|---|
| 1110 |
|
|---|
| 1111 |
|
|---|
| 1112 | //Hard close collisions with outer shells
|
|---|
| 1113 | G4double wmaxc = ene;
|
|---|
| 1114 | G4double closeCS0 = 0.0;
|
|---|
| 1115 | G4double closeCS = 0.0;
|
|---|
| 1116 | if (cutoff>0.1*eV)
|
|---|
| 1117 | {
|
|---|
| 1118 | rl=cutoff/ene;
|
|---|
| 1119 | rl1=1.0-rl;
|
|---|
| 1120 | if (rl < 1.0)
|
|---|
| 1121 | closeCS0 = (((1.0/rl)-1.0) + bha1*std::log(rl) + bha2*rl1
|
|---|
| 1122 | + (bha3/2.0)*((rl*rl)-1.0)
|
|---|
| 1123 | + (bha4/3.0)*(1.0-(rl*rl*rl)))/ene;
|
|---|
| 1124 | }
|
|---|
| 1125 |
|
|---|
| 1126 | // Cross sections for the different oscillators
|
|---|
| 1127 |
|
|---|
| 1128 | // totalHardCS contains the cumulative hard interaction cross section for the different
|
|---|
| 1129 | // excitable levels and the different interaction channels (close, distant, etc.),
|
|---|
| 1130 | // i.e.
|
|---|
| 1131 | // cumulHardCS[0] = 0.0
|
|---|
| 1132 | // cumulHardCS[1] = 1st excitable level (distant longitudinal only)
|
|---|
| 1133 | // cumulHardCS[2] = 1st excitable level (distant longitudinal + transverse)
|
|---|
| 1134 | // cumulHardCS[3] = 1st excitable level (distant longitudinal + transverse + close)
|
|---|
| 1135 | // cumulHardCS[4] = 1st excitable level (all channels) + 2nd excitable level (distant long only)
|
|---|
| 1136 | // etc.
|
|---|
| 1137 | // This is used for sampling the atomic level which is ionised and the channel of the
|
|---|
| 1138 | // interaction.
|
|---|
| 1139 | //
|
|---|
| 1140 | // For each index iFill of the cumulHardCS vector,
|
|---|
| 1141 | // nbOfLevel[iFill] contains the current excitable atomic level and
|
|---|
| 1142 | // typeOfInteraction[iFill] contains the current interaction channel, with the legenda:
|
|---|
| 1143 | // 1 = distant longitudinal interaction
|
|---|
| 1144 | // 2 = distant transverse interaction
|
|---|
| 1145 | // 3 = close collision
|
|---|
| 1146 | // 4 = close collision with outer shells (in this case nbOfLevel < 0 --> no binding energy)
|
|---|
| 1147 |
|
|---|
| 1148 |
|
|---|
| 1149 | G4int nOscil = ionizationEnergy->find(Z)->second->size();
|
|---|
| 1150 | G4double totalHardCS = 0.0;
|
|---|
| 1151 | G4double involvedElectrons = 0.0;
|
|---|
| 1152 | for (G4int i=0;i<nOscil;i++){
|
|---|
| 1153 | G4double wi = (*(resonanceEnergy->find(Z)->second))[i];
|
|---|
| 1154 | G4int occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i];
|
|---|
| 1155 | //Distant excitations
|
|---|
| 1156 | if (wi>cutoff && wi<ene)
|
|---|
| 1157 | {
|
|---|
| 1158 | if (wi>(1e-6*ene)){
|
|---|
| 1159 | G4double cpp=std::sqrt((ene-wi)*(ene-wi+2.0*electron_mass_c2));
|
|---|
| 1160 | qm->push_back(std::sqrt((cp-cpp)*(cp-cpp)+ electron_mass_c2 * electron_mass_c2)-electron_mass_c2);
|
|---|
| 1161 | }
|
|---|
| 1162 | else
|
|---|
| 1163 | {
|
|---|
| 1164 | qm->push_back(wi*wi/(beta2+2.0*electron_mass_c2));
|
|---|
| 1165 | }
|
|---|
| 1166 | //verificare che quando arriva qui il vettore ha SEMPRE l'i-esimo elemento
|
|---|
| 1167 | if ((*qm)[i] < wi)
|
|---|
| 1168 | {
|
|---|
| 1169 |
|
|---|
| 1170 | G4double distantLongitCS = occupNb*std::log(wi*((*qm)[i]+2.0*electron_mass_c2)/
|
|---|
| 1171 | ((*qm)[i]*(wi+2.0*electron_mass_c2)))/wi;
|
|---|
| 1172 | cumulHardCS->push_back(totalHardCS);
|
|---|
| 1173 | typeOfInteraction->push_back(1.0); //distant longitudinal
|
|---|
| 1174 | nbOfLevel->push_back((G4double) i); //only excitable level are counted
|
|---|
| 1175 | totalHardCS += distantLongitCS;
|
|---|
| 1176 |
|
|---|
| 1177 | G4double distantTransvCS = occupNb*distantTransvCS0/wi;
|
|---|
| 1178 |
|
|---|
| 1179 | cumulHardCS->push_back(totalHardCS);
|
|---|
| 1180 | typeOfInteraction->push_back(2.0); //distant tranverse
|
|---|
| 1181 | nbOfLevel->push_back((G4double) i);
|
|---|
| 1182 | totalHardCS += distantTransvCS;
|
|---|
| 1183 | }
|
|---|
| 1184 | }
|
|---|
| 1185 | else
|
|---|
| 1186 | {
|
|---|
| 1187 | qm->push_back(wi);
|
|---|
| 1188 | }
|
|---|
| 1189 | //close collisions
|
|---|
| 1190 | if(wi < wmaxc){
|
|---|
| 1191 | if (wi < cutoff) {
|
|---|
| 1192 | involvedElectrons += occupNb;
|
|---|
| 1193 | }
|
|---|
| 1194 | else
|
|---|
| 1195 | {
|
|---|
| 1196 | rl=wi/ene;
|
|---|
| 1197 | rl1=1.0-rl;
|
|---|
| 1198 | closeCS = occupNb*(((1.0/rl)-1.0)+bha1*std::log(rl)+bha2*rl1
|
|---|
| 1199 | + (bha3/2.0)*((rl*rl)-1.0)
|
|---|
| 1200 | + (bha4/3.0)*(1.0-(rl*rl*rl)))/ene;
|
|---|
| 1201 | cumulHardCS->push_back(totalHardCS);
|
|---|
| 1202 | typeOfInteraction->push_back(3.0); //close
|
|---|
| 1203 | nbOfLevel->push_back((G4double) i);
|
|---|
| 1204 | totalHardCS += closeCS;
|
|---|
| 1205 | }
|
|---|
| 1206 | }
|
|---|
| 1207 | } // loop on the levels
|
|---|
| 1208 |
|
|---|
| 1209 | cumulHardCS->push_back(totalHardCS);
|
|---|
| 1210 | typeOfInteraction->push_back(4.0); //close interaction with outer shells
|
|---|
| 1211 | nbOfLevel->push_back(-1.0);
|
|---|
| 1212 | totalHardCS += involvedElectrons*closeCS0;
|
|---|
| 1213 | cumulHardCS->push_back(totalHardCS); //this is the final value of the totalHardCS
|
|---|
| 1214 |
|
|---|
| 1215 | if (totalHardCS < 1e-30) {
|
|---|
| 1216 | kineticEnergy1=ene;
|
|---|
| 1217 | cosThetaPrimary=1.0;
|
|---|
| 1218 | energySecondary=0.0;
|
|---|
| 1219 | cosThetaSecondary=0.0;
|
|---|
| 1220 | iOsc=-1;
|
|---|
| 1221 | delete qm;
|
|---|
| 1222 | delete cumulHardCS;
|
|---|
| 1223 | delete typeOfInteraction;
|
|---|
| 1224 | delete nbOfLevel;
|
|---|
| 1225 | return;
|
|---|
| 1226 | }
|
|---|
| 1227 |
|
|---|
| 1228 |
|
|---|
| 1229 | //Selection of the active oscillator on the basis of the cumulative cross sections
|
|---|
| 1230 | G4double TST = totalHardCS*G4UniformRand();
|
|---|
| 1231 | G4int is=0;
|
|---|
| 1232 | G4int js= nbOfLevel->size();
|
|---|
| 1233 | do{
|
|---|
| 1234 | G4int it=(is+js)/2;
|
|---|
| 1235 | if (TST > (*cumulHardCS)[it]) is=it;
|
|---|
| 1236 | if (TST <= (*cumulHardCS)[it]) js=it;
|
|---|
| 1237 | }while((js-is) > 1);
|
|---|
| 1238 |
|
|---|
| 1239 | G4double UII=0.0;
|
|---|
| 1240 | G4double rkc=cutoff/ene;
|
|---|
| 1241 | G4double dde;
|
|---|
| 1242 | G4int kks;
|
|---|
| 1243 |
|
|---|
| 1244 | G4double sampledInteraction = (*typeOfInteraction)[is];
|
|---|
| 1245 | iOsc = (G4int) (*nbOfLevel)[is];
|
|---|
| 1246 |
|
|---|
| 1247 | //Generates the final state according to the sampled level and
|
|---|
| 1248 | //interaction channel
|
|---|
| 1249 |
|
|---|
| 1250 | if (sampledInteraction == 1.0) //Hard distant longitudinal collisions
|
|---|
| 1251 | {
|
|---|
| 1252 | dde= (*(resonanceEnergy->find(Z)->second))[iOsc];
|
|---|
| 1253 | kineticEnergy1=ene-dde;
|
|---|
| 1254 | G4double qs=(*qm)[iOsc]/(1.0+((*qm)[iOsc]/(2.0*electron_mass_c2)));
|
|---|
| 1255 | G4double q=qs/(std::pow((qs/dde)*(1.0+(0.5*dde/electron_mass_c2)),G4UniformRand())-(0.5*qs/electron_mass_c2));
|
|---|
| 1256 | G4double qtrev = q*(q+2.0*electron_mass_c2);
|
|---|
| 1257 | G4double cpps = kineticEnergy1*(kineticEnergy1+2.0*electron_mass_c2);
|
|---|
| 1258 | cosThetaPrimary = (cpps+cps-qtrev)/(2.0*cp*std::sqrt(cpps));
|
|---|
| 1259 | if (cosThetaPrimary>1.0) cosThetaPrimary=1.0;
|
|---|
| 1260 | //Energy and emission angle of the delta ray
|
|---|
| 1261 | kks = (G4int) (*(shellFlag->find(Z)->second))[iOsc];
|
|---|
| 1262 | if (kks>4)
|
|---|
| 1263 | {
|
|---|
| 1264 | energySecondary=dde;
|
|---|
| 1265 | }
|
|---|
| 1266 | else
|
|---|
| 1267 | {
|
|---|
| 1268 | energySecondary=dde-(*(ionizationEnergy->find(Z)->second))[iOsc];
|
|---|
| 1269 | }
|
|---|
| 1270 | cosThetaSecondary = 0.5*(dde*(ene+rb-dde)+qtrev)/std::sqrt(cps*qtrev);
|
|---|
| 1271 | if (cosThetaSecondary>1.0) cosThetaSecondary=1.0;
|
|---|
| 1272 | }
|
|---|
| 1273 |
|
|---|
| 1274 | else if (sampledInteraction == 2.0) //Hard distant transverse collisions
|
|---|
| 1275 | {
|
|---|
| 1276 | dde=(*(resonanceEnergy->find(Z)->second))[iOsc];
|
|---|
| 1277 | kineticEnergy1=ene-dde;
|
|---|
| 1278 | cosThetaPrimary=1.0;
|
|---|
| 1279 | //Energy and emission angle of the delta ray
|
|---|
| 1280 | kks = (G4int) (*(shellFlag->find(Z)->second))[iOsc];
|
|---|
| 1281 | if (kks>4)
|
|---|
| 1282 | {
|
|---|
| 1283 | energySecondary=dde;
|
|---|
| 1284 | }
|
|---|
| 1285 | else
|
|---|
| 1286 | {
|
|---|
| 1287 | energySecondary=dde-(*(ionizationEnergy->find(Z)->second))[iOsc];
|
|---|
| 1288 | }
|
|---|
| 1289 | cosThetaSecondary = 1.0;
|
|---|
| 1290 | }
|
|---|
| 1291 |
|
|---|
| 1292 | else if (sampledInteraction == 3.0 || sampledInteraction == 4.0) //Close interaction
|
|---|
| 1293 | {
|
|---|
| 1294 | if (sampledInteraction == 4.0) //interaction with inner shells
|
|---|
| 1295 | {
|
|---|
| 1296 | UII=0.0;
|
|---|
| 1297 | rkc = cutoff/ene;
|
|---|
| 1298 | iOsc = -1;
|
|---|
| 1299 | }
|
|---|
| 1300 | else
|
|---|
| 1301 | {
|
|---|
| 1302 | kks = (G4int) (*(shellFlag->find(Z)->second))[iOsc];
|
|---|
| 1303 | if (kks > 4) {
|
|---|
| 1304 | UII=0.0;
|
|---|
| 1305 | }
|
|---|
| 1306 | else
|
|---|
| 1307 | {
|
|---|
| 1308 | UII = (*(ionizationEnergy->find(Z)->second))[iOsc];
|
|---|
| 1309 | }
|
|---|
| 1310 | rkc = (*(resonanceEnergy->find(Z)->second))[iOsc]/ene;
|
|---|
| 1311 | }
|
|---|
| 1312 | G4double phi,rk;
|
|---|
| 1313 | do{
|
|---|
| 1314 | rk=rkc/(1.0-G4UniformRand()*(1.0-rkc));
|
|---|
| 1315 | phi = 1.0-rk*(bha1-rk*(bha2-rk*(bha3-bha4*rk)));
|
|---|
| 1316 | }while ( G4UniformRand() > phi);
|
|---|
| 1317 | //Energy and scattering angle (primary electron);
|
|---|
| 1318 | kineticEnergy1 = ene*(1.0-rk);
|
|---|
| 1319 | cosThetaPrimary = std::sqrt(kineticEnergy1*rb/(ene*(rb-(rk*ene))));
|
|---|
| 1320 | //Energy and scattering angle of the delta ray
|
|---|
| 1321 | energySecondary = ene-kineticEnergy1-UII;
|
|---|
| 1322 | cosThetaSecondary = std::sqrt(rk*ene*rb/(ene*(rk*ene+2.0*electron_mass_c2)));
|
|---|
| 1323 | }
|
|---|
| 1324 | else
|
|---|
| 1325 | {
|
|---|
| 1326 | G4String excep = "G4PenelopeIonisation - Error in the calculation of the final state";
|
|---|
| 1327 | G4Exception(excep);
|
|---|
| 1328 | }
|
|---|
| 1329 |
|
|---|
| 1330 | delete qm;
|
|---|
| 1331 | delete cumulHardCS;
|
|---|
| 1332 | delete typeOfInteraction;
|
|---|
| 1333 | delete nbOfLevel;
|
|---|
| 1334 |
|
|---|
| 1335 | return;
|
|---|
| 1336 | }
|
|---|
| 1337 |
|
|---|
| 1338 | // This stuff in needed in order to interface with the Cross Section Handler
|
|---|
| 1339 |
|
|---|
| 1340 | G4double G4PenelopeIonisation::CalculateCrossSectionsRatio(G4double ene,G4double cutoff,
|
|---|
| 1341 | G4int Z,G4double electronVolumeDensity,
|
|---|
| 1342 | const G4ParticleDefinition& particle)
|
|---|
| 1343 | {
|
|---|
| 1344 | //Constants
|
|---|
| 1345 | G4double gamma = 1.0+ene/electron_mass_c2;
|
|---|
| 1346 | G4double gamma2 = gamma*gamma;
|
|---|
| 1347 | G4double beta2 = (gamma2-1.0)/gamma2;
|
|---|
| 1348 | G4double constant = pi*classic_electr_radius*classic_electr_radius*2.0*electron_mass_c2/beta2;
|
|---|
| 1349 | G4double delta = CalculateDeltaFermi(ene,Z,electronVolumeDensity);
|
|---|
| 1350 | G4int nbOsc = (G4int) resonanceEnergy->find(Z)->second->size();
|
|---|
| 1351 | G4double S0 = 0.0, H0=0.0;
|
|---|
| 1352 | G4double softCS = 0.0;
|
|---|
| 1353 | G4double hardCS = 0.0;
|
|---|
| 1354 | for (G4int i=0;i<nbOsc;i++){
|
|---|
| 1355 | G4double resEnergy = (*(resonanceEnergy->find(Z)->second))[i];
|
|---|
| 1356 | if (&particle == G4Electron::Electron())
|
|---|
| 1357 | {
|
|---|
| 1358 | S0 = CrossSectionsRatioForElectrons(ene,resEnergy,delta,cutoff,1);
|
|---|
| 1359 | H0 = CrossSectionsRatioForElectrons(ene,resEnergy,delta,cutoff,2);
|
|---|
| 1360 | }
|
|---|
| 1361 | else if (&particle == G4Positron::Positron())
|
|---|
| 1362 | {
|
|---|
| 1363 | S0 = CrossSectionsRatioForPositrons(ene,resEnergy,delta,cutoff,1);
|
|---|
| 1364 | H0 = CrossSectionsRatioForPositrons(ene,resEnergy,delta,cutoff,2);
|
|---|
| 1365 | }
|
|---|
| 1366 | G4double occupNb = (*(occupationNumber->find(Z)->second))[i];
|
|---|
| 1367 | softCS += occupNb*constant*S0;
|
|---|
| 1368 | hardCS += occupNb*constant*H0;
|
|---|
| 1369 | }
|
|---|
| 1370 | G4double ratio = 0.0;
|
|---|
| 1371 | if (softCS+hardCS) ratio = (hardCS)/(softCS+hardCS);
|
|---|
| 1372 | return ratio;
|
|---|
| 1373 | }
|
|---|
| 1374 |
|
|---|
| 1375 |
|
|---|
| 1376 | G4double G4PenelopeIonisation::CrossSectionsRatioForElectrons(G4double ene,G4double resEne,
|
|---|
| 1377 | G4double delta,G4double cutoff,
|
|---|
| 1378 | G4int index)
|
|---|
| 1379 | {
|
|---|
| 1380 | //Calculate constants
|
|---|
| 1381 | G4double gamma = 1.0+ene/electron_mass_c2;
|
|---|
| 1382 | G4double gamma2 = gamma*gamma;
|
|---|
| 1383 | G4double beta2 = (gamma2-1.0)/gamma2;
|
|---|
| 1384 | G4double cps = ene*(ene+2.0*electron_mass_c2);
|
|---|
| 1385 | G4double amol = (ene/(ene+electron_mass_c2)) * (ene/(ene+electron_mass_c2)) ;
|
|---|
| 1386 | G4double hardCont = 0.0;
|
|---|
| 1387 | G4double softCont = 0.0;
|
|---|
| 1388 | if (ene < resEne) return 0.0;
|
|---|
| 1389 |
|
|---|
| 1390 | //Distant interactions
|
|---|
| 1391 | G4double cp1s = (ene-resEne)*(ene-resEne+2.0*electron_mass_c2);
|
|---|
| 1392 | G4double cp1 = std::sqrt(cp1s);
|
|---|
| 1393 | G4double cp = std::sqrt(cps);
|
|---|
| 1394 | G4double sdLong=0.0, sdTrans = 0.0, sdDist=0.0;
|
|---|
| 1395 |
|
|---|
| 1396 | //Distant longitudinal interactions
|
|---|
| 1397 | G4double qm = 0.0;
|
|---|
| 1398 |
|
|---|
| 1399 | if (resEne > ene*(1e-6))
|
|---|
| 1400 | {
|
|---|
| 1401 | qm = std::sqrt((cp-cp1)*(cp-cp1)+(electron_mass_c2*electron_mass_c2))-electron_mass_c2;
|
|---|
| 1402 | }
|
|---|
| 1403 | else
|
|---|
| 1404 | {
|
|---|
| 1405 | qm = resEne*resEne/(beta2*2.0*electron_mass_c2);
|
|---|
| 1406 | qm = qm*(1.0-0.5*qm/electron_mass_c2);
|
|---|
| 1407 | }
|
|---|
| 1408 |
|
|---|
| 1409 | if (qm < resEne)
|
|---|
| 1410 | {
|
|---|
| 1411 | sdLong = std::log(resEne*(qm+2.0*electron_mass_c2)/(qm*(resEne+2.0*electron_mass_c2)));
|
|---|
| 1412 | }
|
|---|
| 1413 | else
|
|---|
| 1414 | {
|
|---|
| 1415 | sdLong = 0.0;
|
|---|
| 1416 | }
|
|---|
| 1417 |
|
|---|
| 1418 | if (sdLong > 0) {
|
|---|
| 1419 | sdTrans = std::max(std::log(gamma2)-beta2-delta,0.0);
|
|---|
| 1420 | sdDist = sdTrans + sdLong;
|
|---|
| 1421 | if (cutoff > resEne)
|
|---|
| 1422 | {
|
|---|
| 1423 | softCont = sdDist/resEne;
|
|---|
| 1424 | }
|
|---|
| 1425 | else
|
|---|
| 1426 | {
|
|---|
| 1427 | hardCont = sdDist/resEne;
|
|---|
| 1428 | }
|
|---|
| 1429 | }
|
|---|
| 1430 |
|
|---|
| 1431 |
|
|---|
| 1432 | // Close collisions (Moeller's cross section)
|
|---|
| 1433 | G4double wl = std::max(cutoff,resEne);
|
|---|
| 1434 | G4double wu = 0.5*ene;
|
|---|
| 1435 |
|
|---|
| 1436 | if (wl < (wu-1*eV))
|
|---|
| 1437 | {
|
|---|
| 1438 | hardCont += (1.0/(ene-wu))-(1.0/(ene-wl))
|
|---|
| 1439 | - (1.0/wu)+(1.0/wl)
|
|---|
| 1440 | + (1.0-amol)*std::log(((ene-wu)*wl)/((ene-wl)*wu))/ene
|
|---|
| 1441 | + amol*(wu-wl)/(ene*ene);
|
|---|
| 1442 | wu=wl;
|
|---|
| 1443 | }
|
|---|
| 1444 |
|
|---|
| 1445 | wl = resEne;
|
|---|
| 1446 | if (wl > (wu-1*eV)) {
|
|---|
| 1447 | if (index == 1) return softCont;
|
|---|
| 1448 | if (index == 2) return hardCont;
|
|---|
| 1449 | }
|
|---|
| 1450 | softCont += (1.0/(ene-wu))-(1.0/(ene-wl))
|
|---|
| 1451 | - (1.0/wu)+(1.0/wl)
|
|---|
| 1452 | + (1.0-amol)*std::log(((ene-wu)*wl)/((ene-wl)*wu))/ene
|
|---|
| 1453 | + amol*(wu-wl)/(ene*ene);
|
|---|
| 1454 | if (index == 1) return softCont;
|
|---|
| 1455 | return hardCont;
|
|---|
| 1456 | }
|
|---|
| 1457 |
|
|---|
| 1458 | G4double G4PenelopeIonisation::CrossSectionsRatioForPositrons(G4double ene,G4double resEne,
|
|---|
| 1459 | G4double delta,G4double cutoff,G4int index)
|
|---|
| 1460 | {
|
|---|
| 1461 | //Calculate constants
|
|---|
| 1462 | G4double gamma = 1.0+ene/electron_mass_c2;
|
|---|
| 1463 | G4double gamma2 = gamma*gamma;
|
|---|
| 1464 | G4double beta2 = (gamma2-1.0)/gamma2;
|
|---|
| 1465 | G4double cps = ene*(ene+2.0*electron_mass_c2);
|
|---|
| 1466 | G4double amol = (ene/(ene+electron_mass_c2)) * (ene/(ene+electron_mass_c2)) ;
|
|---|
| 1467 | G4double help = (gamma+1.0)*(gamma+1.0);
|
|---|
| 1468 | G4double bha1 = amol*(2.0*help-1.0)/(gamma2-1.0);
|
|---|
| 1469 | G4double bha2 = amol*(3.0+1.0/help);
|
|---|
| 1470 | G4double bha3 = amol*2.0*gamma*(gamma-1.0)/help;
|
|---|
| 1471 | G4double bha4 = amol*(gamma-1.0)*(gamma-1.0)/help;
|
|---|
| 1472 | G4double hardCont = 0.0;
|
|---|
| 1473 | G4double softCont = 0.0;
|
|---|
| 1474 | if (ene < resEne) return 0.0;
|
|---|
| 1475 |
|
|---|
| 1476 |
|
|---|
| 1477 | //Distant interactions
|
|---|
| 1478 | G4double cp1s = (ene-resEne)*(ene-resEne+2.0*electron_mass_c2);
|
|---|
| 1479 | G4double cp1 = std::sqrt(cp1s);
|
|---|
| 1480 | G4double cp = std::sqrt(cps);
|
|---|
| 1481 | G4double sdLong=0.0, sdTrans = 0.0, sdDist=0.0;
|
|---|
| 1482 |
|
|---|
| 1483 | //Distant longitudinal interactions
|
|---|
| 1484 | G4double qm = 0.0;
|
|---|
| 1485 |
|
|---|
| 1486 | if (resEne > ene*(1e-6))
|
|---|
| 1487 | {
|
|---|
| 1488 | qm = std::sqrt((cp-cp1)*(cp-cp1)+(electron_mass_c2*electron_mass_c2))-electron_mass_c2;
|
|---|
| 1489 | }
|
|---|
| 1490 | else
|
|---|
| 1491 | {
|
|---|
| 1492 | qm = resEne*resEne/(beta2*2.0*electron_mass_c2);
|
|---|
| 1493 | qm = qm*(1.0-0.5*qm/electron_mass_c2);
|
|---|
| 1494 | }
|
|---|
| 1495 |
|
|---|
| 1496 | if (qm < resEne)
|
|---|
| 1497 | {
|
|---|
| 1498 | sdLong = std::log(resEne*(qm+2.0*electron_mass_c2)/(qm*(resEne+2.0*electron_mass_c2)));
|
|---|
| 1499 | }
|
|---|
| 1500 | else
|
|---|
| 1501 | {
|
|---|
| 1502 | sdLong = 0.0;
|
|---|
| 1503 | }
|
|---|
| 1504 |
|
|---|
| 1505 | if (sdLong > 0) {
|
|---|
| 1506 | sdTrans = std::max(std::log(gamma2)-beta2-delta,0.0);
|
|---|
| 1507 | sdDist = sdTrans + sdLong;
|
|---|
| 1508 | if (cutoff > resEne)
|
|---|
| 1509 | {
|
|---|
| 1510 | softCont = sdDist/resEne;
|
|---|
| 1511 | }
|
|---|
| 1512 | else
|
|---|
| 1513 | {
|
|---|
| 1514 | hardCont = sdDist/resEne;
|
|---|
| 1515 | }
|
|---|
| 1516 | }
|
|---|
| 1517 |
|
|---|
| 1518 |
|
|---|
| 1519 | // Close collisions (Bhabha's cross section)
|
|---|
| 1520 | G4double wl = std::max(cutoff,resEne);
|
|---|
| 1521 | G4double wu = ene;
|
|---|
| 1522 |
|
|---|
| 1523 | if (wl < (wu-1*eV)) {
|
|---|
| 1524 | hardCont += (1.0/wl)-(1.0/wu)-bha1*std::log(wu/wl)/ene
|
|---|
| 1525 | + bha2*(wu-wl)/(ene*ene) -bha3*((wu*wu)-(wl*wl))/(2.0*ene*ene*ene)
|
|---|
| 1526 | + bha4*((wu*wu*wu)-(wl*wl*wl))/(3.0*ene*ene*ene*ene);
|
|---|
| 1527 | wu=wl;
|
|---|
| 1528 | }
|
|---|
| 1529 | wl = resEne;
|
|---|
| 1530 | if (wl > (wu-1*eV))
|
|---|
| 1531 | {
|
|---|
| 1532 | if (index == 1) return softCont;
|
|---|
| 1533 | if (index == 2) return hardCont;
|
|---|
| 1534 | }
|
|---|
| 1535 | softCont += (1.0/wl)-(1.0/wu)-bha1*std::log(wu/wl)/ene
|
|---|
| 1536 | + bha2*(wu-wl)/(ene*ene) -bha3*((wu*wu)-(wl*wl))/(2.0*ene*ene*ene)
|
|---|
| 1537 | + bha4*((wu*wu*wu)-(wl*wl*wl))/(3.0*ene*ene*ene*ene);
|
|---|
| 1538 |
|
|---|
| 1539 | if (index == 1) return softCont;
|
|---|
| 1540 | return hardCont;
|
|---|
| 1541 | }
|
|---|