00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048 #include "G4VAtomDeexcitation.hh"
00049 #include "G4SystemOfUnits.hh"
00050 #include "G4ParticleDefinition.hh"
00051 #include "G4DynamicParticle.hh"
00052 #include "G4Step.hh"
00053 #include "G4Region.hh"
00054 #include "G4RegionStore.hh"
00055 #include "G4MaterialCutsCouple.hh"
00056 #include "G4MaterialCutsCouple.hh"
00057 #include "G4Material.hh"
00058 #include "G4Element.hh"
00059 #include "G4ElementVector.hh"
00060 #include "Randomize.hh"
00061 #include "G4VParticleChange.hh"
00062
00063
00064
00065 G4VAtomDeexcitation::G4VAtomDeexcitation(const G4String& modname,
00066 const G4String& pname)
00067 : lowestKinEnergy(keV), verbose(1), name(modname), namePIXE(pname),
00068 nameElectronPIXE(""), isActive(false), flagAuger(false), flagPIXE(false)
00069 {
00070 vdyn.reserve(5);
00071 theCoupleTable = 0;
00072 }
00073
00074
00075
00076 G4VAtomDeexcitation::~G4VAtomDeexcitation()
00077 {}
00078
00079
00080
00081 void G4VAtomDeexcitation::InitialiseAtomicDeexcitation()
00082 {
00083
00084 theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
00085 size_t numOfCouples = theCoupleTable->GetTableSize();
00086
00087
00088 if(0 == numOfCouples) { numOfCouples = 1; }
00089
00090 activeDeexcitationMedia.resize(numOfCouples, false);
00091 activeAugerMedia.resize(numOfCouples, false);
00092 activePIXEMedia.resize(numOfCouples, false);
00093 activeZ.resize(93, false);
00094
00095
00096 if( !isActive ) { return; }
00097
00098
00099 size_t nRegions = deRegions.size();
00100
00101 if(0 == nRegions) {
00102 SetDeexcitationActiveRegion("World",isActive,flagAuger,flagPIXE);
00103 nRegions = 1;
00104 }
00105
00106 if(0 < verbose) {
00107 G4cout << G4endl;
00108 G4cout << "### === Deexcitation model " << name
00109 << " is activated for " << nRegions;
00110 if(1 == nRegions) { G4cout << " region:" << G4endl; }
00111 else { G4cout << " regions:" << G4endl;}
00112 }
00113
00114
00115 G4RegionStore* regionStore = G4RegionStore::GetInstance();
00116 for(size_t j=0; j<nRegions; ++j) {
00117 const G4Region* reg = regionStore->GetRegion(activeRegions[j], false);
00118 const G4ProductionCuts* rpcuts = reg->GetProductionCuts();
00119 if(0 < verbose) {
00120 G4cout << " " << activeRegions[j] << G4endl;
00121 }
00122
00123 for(size_t i=0; i<numOfCouples; ++i) {
00124 const G4MaterialCutsCouple* couple =
00125 theCoupleTable->GetMaterialCutsCouple(i);
00126 if (couple->GetProductionCuts() == rpcuts) {
00127 activeDeexcitationMedia[i] = deRegions[j];
00128 activeAugerMedia[i] = AugerRegions[j];
00129 activePIXEMedia[i] = PIXERegions[j];
00130 const G4Material* mat = couple->GetMaterial();
00131 const G4ElementVector* theElementVector =
00132 mat->GetElementVector();
00133 G4int nelm = mat->GetNumberOfElements();
00134 if(deRegions[j]) {
00135 for(G4int k=0; k<nelm; ++k) {
00136 G4int Z = G4lrint(((*theElementVector)[k])->GetZ());
00137 if(Z > 5 && Z < 93) {
00138 activeZ[Z] = true;
00139
00140 }
00141 }
00142 }
00143 }
00144 }
00145 }
00146
00147
00148 InitialiseForNewRun();
00149
00150 if(0 < verbose && flagPIXE) {
00151 G4cout << "### === PIXE model for hadrons: " << namePIXE
00152 << " " << IsPIXEActive()
00153 << G4endl;
00154 G4cout << "### === PIXE model for e+-: " << nameElectronPIXE
00155 << " " << IsPIXEActive()
00156 << G4endl;
00157 }
00158 }
00159
00160
00161
00162 void
00163 G4VAtomDeexcitation::SetDeexcitationActiveRegion(const G4String& rname,
00164 G4bool valDeexcitation,
00165 G4bool valAuger,
00166 G4bool valPIXE)
00167 {
00168 G4String ss = rname;
00169
00170
00171
00172 if(ss == "world" || ss == "World" || ss == "WORLD") {
00173 ss = "DefaultRegionForTheWorld";
00174 }
00175 size_t n = deRegions.size();
00176 if(n > 0) {
00177 for(size_t i=0; i<n; ++i) {
00178
00179
00180 if(ss == activeRegions[i]) {
00181 deRegions[i] = valDeexcitation;
00182 AugerRegions[i] = valAuger;
00183 PIXERegions[i] = valPIXE;
00184 return;
00185 }
00186 }
00187 }
00188
00189 activeRegions.push_back(ss);
00190 deRegions.push_back(valDeexcitation);
00191 AugerRegions.push_back(valAuger);
00192 PIXERegions.push_back(valPIXE);
00193 }
00194
00195
00196
00197 void
00198 G4VAtomDeexcitation::AlongStepDeexcitation(std::vector<G4Track*>& tracks,
00199 const G4Step& step,
00200 G4double& eLossMax,
00201 G4int coupleIndex)
00202 {
00203 G4double truelength = step.GetStepLength();
00204 if(!flagPIXE && !activePIXEMedia[coupleIndex]) { return; }
00205 if(eLossMax <= 0.0 || truelength <= 0.0) { return; }
00206
00207
00208 const G4StepPoint* preStep = step.GetPreStepPoint();
00209 G4ThreeVector prePos = preStep->GetPosition();
00210 G4ThreeVector delta = step.GetPostStepPoint()->GetPosition() - prePos;
00211 G4double preTime = preStep->GetGlobalTime();
00212 G4double dt = step.GetPostStepPoint()->GetGlobalTime() - preTime;
00213
00214
00215 const G4Track* track = step.GetTrack();
00216 const G4ParticleDefinition* part = track->GetDefinition();
00217 G4double ekin = preStep->GetKineticEnergy();
00218
00219
00220 G4double gCut = (*theCoupleTable->GetEnergyCutsVector(0))[coupleIndex];
00221 G4double eCut = DBL_MAX;
00222 if(CheckAugerActiveRegion(coupleIndex)) {
00223 eCut = (*theCoupleTable->GetEnergyCutsVector(1))[coupleIndex];
00224 }
00225
00226
00227
00228
00229 const G4Material* material = preStep->GetMaterial();
00230 const G4ElementVector* theElementVector = material->GetElementVector();
00231 const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
00232 G4int nelm = material->GetNumberOfElements();
00233
00234
00235 for(G4int i=0; i<nelm; ++i) {
00236 G4int Z = G4lrint((*theElementVector)[i]->GetZ());
00237 if(activeZ[Z] && Z < 93) {
00238 G4int nshells = std::min(9,(*theElementVector)[i]->GetNbOfAtomicShells());
00239 G4double rho = truelength*theAtomNumDensityVector[i];
00240
00241 for(G4int ii=0; ii<nshells; ++ii) {
00242 G4AtomicShellEnumerator as = G4AtomicShellEnumerator(ii);
00243 const G4AtomicShell* shell = GetAtomicShell(Z, as);
00244 G4double bindingEnergy = shell->BindingEnergy();
00245
00246 if(gCut > bindingEnergy) { break; }
00247
00248 if(eLossMax > bindingEnergy) {
00249 G4double sig = rho*
00250 GetShellIonisationCrossSectionPerAtom(part, Z, as, ekin, material);
00251
00252
00253 if(sig > 0.0) {
00254 G4double mfp = 1.0/sig;
00255 G4double stot = 0.0;
00256
00257
00258 do {
00259 stot -= mfp*std::log(G4UniformRand());
00260 if( stot > 1.0 || eLossMax < bindingEnergy) { break; }
00261
00262 vdyn.clear();
00263 GenerateParticles(&vdyn, shell, Z, gCut, eCut);
00264 G4int nsec = vdyn.size();
00265 if(nsec > 0) {
00266 G4ThreeVector r = prePos + stot*delta;
00267 G4double time = preTime + stot*dt;
00268 for(G4int j=0; j<nsec; ++j) {
00269 G4DynamicParticle* dp = vdyn[j];
00270 G4double e = dp->GetKineticEnergy();
00271
00272
00273 if(eLossMax >= e) {
00274 eLossMax -= e;
00275 G4Track* t = new G4Track(dp, time, r);
00276 tracks.push_back(t);
00277 } else {
00278 delete dp;
00279 }
00280 }
00281 }
00282 } while (stot < 1.0);
00283 }
00284 }
00285 }
00286 }
00287 }
00288 return;
00289 }
00290
00291