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 #include "G4MuElecElasticModel.hh"
00039 #include "G4PhysicalConstants.hh"
00040 #include "G4SystemOfUnits.hh"
00041
00042
00043
00044 using namespace std;
00045
00046
00047
00048 G4MuElecElasticModel::G4MuElecElasticModel(const G4ParticleDefinition*,
00049 const G4String& nam)
00050 :G4VEmModel(nam),isInitialised(false)
00051 {
00052 nistSi = G4NistManager::Instance()->FindOrBuildMaterial("G4_Si");
00053
00054 killBelowEnergy = 16.7 * eV;
00055 lowEnergyLimit = 0 * eV;
00056 lowEnergyLimitOfModel = 5 * eV;
00057 highEnergyLimit = 100. * MeV;
00058 SetLowEnergyLimit(lowEnergyLimit);
00059 SetHighEnergyLimit(highEnergyLimit);
00060
00061 verboseLevel= 0;
00062
00063
00064
00065
00066
00067
00068
00069 if( verboseLevel>0 )
00070 {
00071 G4cout << "MuElec Elastic model is constructed " << G4endl
00072 << "Energy range: "
00073 << lowEnergyLimit / eV << " eV - "
00074 << highEnergyLimit / keV << " keV"
00075 << G4endl;
00076 }
00077 fParticleChangeForGamma = 0;
00078 }
00079
00080
00081
00082 G4MuElecElasticModel::~G4MuElecElasticModel()
00083 {
00084
00085
00086 std::map< G4String,G4MuElecCrossSectionDataSet*,std::less<G4String> >::iterator pos;
00087 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
00088 {
00089 G4MuElecCrossSectionDataSet* table = pos->second;
00090 delete table;
00091 }
00092
00093
00094
00095 eVecm.clear();
00096
00097 }
00098
00099
00100
00101 void G4MuElecElasticModel::Initialise(const G4ParticleDefinition* ,
00102 const G4DataVector& )
00103 {
00104
00105 if (verboseLevel > 3)
00106 G4cout << "Calling G4MuElecElasticModel::Initialise()" << G4endl;
00107
00108
00109
00110 if (LowEnergyLimit() < lowEnergyLimit)
00111 {
00112 G4cout << "G4MuElecElasticModel: low energy limit increased from " <<
00113 LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
00114 SetLowEnergyLimit(lowEnergyLimit);
00115 }
00116
00117 if (HighEnergyLimit() > highEnergyLimit)
00118 {
00119 G4cout << "G4MuElecElasticModel: high energy limit decreased from " <<
00120 HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl;
00121 SetHighEnergyLimit(highEnergyLimit);
00122 }
00123
00124
00125
00126 G4double scaleFactor = 1e-18 * cm * cm;
00127
00128 G4String fileElectron("muelec/sigma_elastic_e_Si");
00129
00130 G4ParticleDefinition* electronDef = G4Electron::ElectronDefinition();
00131 G4String electron;
00132
00133
00134
00135 electron = electronDef->GetParticleName();
00136
00137 tableFile[electron] = fileElectron;
00138
00139 G4MuElecCrossSectionDataSet* tableE = new G4MuElecCrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
00140 tableE->LoadData(fileElectron);
00141 tableData[electron] = tableE;
00142
00143
00144
00145 char *path = getenv("G4LEDATA");
00146
00147 if (!path)
00148 {
00149 G4Exception("G4MuElecElasticModel::Initialise","em0006",FatalException,"G4LEDATA environment variable not set.");
00150 return;
00151 }
00152
00153 std::ostringstream eFullFileName;
00154 eFullFileName << path << "/muelec/sigmadiff_elastic_e_Si.dat";
00155 std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
00156
00157 if (!eDiffCrossSection)
00158 G4Exception("G4MuElecElasticModel::Initialise","em0003",FatalException,"Missing data file: /muelec/sigmadiff_elastic_e_Si.dat");
00159
00160 eTdummyVec.push_back(0.);
00161
00162 while(!eDiffCrossSection.eof())
00163 {
00164 double tDummy;
00165 double eDummy;
00166 eDiffCrossSection>>tDummy>>eDummy;
00167
00168
00169 if (tDummy != eTdummyVec.back())
00170 {
00171 eTdummyVec.push_back(tDummy);
00172 eVecm[tDummy].push_back(0.);
00173 }
00174
00175 eDiffCrossSection>>eDiffCrossSectionData[tDummy][eDummy];
00176
00177
00178 if (!eDiffCrossSection.eof()) eDiffCrossSectionData[tDummy][eDummy]*=scaleFactor;
00179
00180 if (eDummy != eVecm[tDummy].back()) eVecm[tDummy].push_back(eDummy);
00181
00182 }
00183
00184
00185
00186 if (verboseLevel > 2)
00187 G4cout << "Loaded cross section files for MuElec Elastic model" << G4endl;
00188
00189 if( verboseLevel>0 )
00190 {
00191 G4cout << "MuElec Elastic model is initialized " << G4endl
00192 << "Energy range: "
00193 << LowEnergyLimit() / eV << " eV - "
00194 << HighEnergyLimit() / keV << " keV"
00195 << G4endl;
00196 }
00197
00198 if (isInitialised) { return; }
00199 fParticleChangeForGamma = GetParticleChangeForGamma();
00200 isInitialised = true;
00201
00202
00203 }
00204
00205
00206
00207 G4double G4MuElecElasticModel::CrossSectionPerVolume(const G4Material* material,
00208 const G4ParticleDefinition* p,
00209 G4double ekin,
00210 G4double,
00211 G4double)
00212 {
00213 if (verboseLevel > 3)
00214 G4cout << "Calling CrossSectionPerVolume() of G4MuElecElasticModel" << G4endl;
00215
00216
00217
00218 G4double sigma=0;
00219
00220 G4double density = material->GetTotNbOfAtomsPerVolume();
00221
00222 if (material == nistSi || material->GetBaseMaterial() == nistSi)
00223 {
00224 const G4String& particleName = p->GetParticleName();
00225
00226 if (ekin < highEnergyLimit)
00227 {
00228
00229 if (ekin < lowEnergyLimitOfModel) ekin = lowEnergyLimitOfModel;
00230
00231
00232 std::map< G4String,G4MuElecCrossSectionDataSet*,std::less<G4String> >::iterator pos;
00233 pos = tableData.find(particleName);
00234
00235 if (pos != tableData.end())
00236 {
00237 G4MuElecCrossSectionDataSet* table = pos->second;
00238 if (table != 0)
00239 {
00240 sigma = table->FindValue(ekin);
00241 }
00242 }
00243 else
00244 {
00245 G4Exception("G4MuElecElasticModel::ComputeCrossSectionPerVolume","em0002",FatalException,"Model not applicable to particle type.");
00246 }
00247 }
00248
00249 if (verboseLevel > 3)
00250 {
00251 G4cout << "---> Kinetic energy(eV)=" << ekin/eV << G4endl;
00252 G4cout << " - Cross section per Si atom (cm^2)=" << sigma/cm/cm << G4endl;
00253 G4cout << " - Cross section per Si atom (cm^-1)=" << sigma*density/(1./cm) << G4endl;
00254 }
00255
00256 }
00257
00258 return sigma*density;
00259 }
00260
00261
00262
00263 void G4MuElecElasticModel::SampleSecondaries(std::vector<G4DynamicParticle*>* ,
00264 const G4MaterialCutsCouple* ,
00265 const G4DynamicParticle* aDynamicElectron,
00266 G4double,
00267 G4double)
00268 {
00269
00270 if (verboseLevel > 3)
00271 G4cout << "Calling SampleSecondaries() of G4MuElecElasticModel" << G4endl;
00272
00273 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
00274
00275 if (electronEnergy0 < killBelowEnergy)
00276 {
00277 fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
00278 fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
00279 return ;
00280 }
00281
00282 if (electronEnergy0>= killBelowEnergy && electronEnergy0 < highEnergyLimit)
00283 {
00284 G4double cosTheta = RandomizeCosTheta(electronEnergy0);
00285
00286 G4double phi = 2. * pi * G4UniformRand();
00287
00288 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
00289 G4ThreeVector xVers = zVers.orthogonal();
00290 G4ThreeVector yVers = zVers.cross(xVers);
00291
00292 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
00293 G4double yDir = xDir;
00294 xDir *= std::cos(phi);
00295 yDir *= std::sin(phi);
00296
00297 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
00298
00299 fParticleChangeForGamma->ProposeMomentumDirection(zPrimeVers.unit()) ;
00300
00301 fParticleChangeForGamma->SetProposedKineticEnergy(electronEnergy0);
00302 }
00303
00304 }
00305
00306
00307
00308 G4double G4MuElecElasticModel::Theta
00309 (G4ParticleDefinition * particleDefinition, G4double k, G4double integrDiff)
00310 {
00311
00312 G4double theta = 0.;
00313 G4double valueT1 = 0;
00314 G4double valueT2 = 0;
00315 G4double valueE21 = 0;
00316 G4double valueE22 = 0;
00317 G4double valueE12 = 0;
00318 G4double valueE11 = 0;
00319 G4double xs11 = 0;
00320 G4double xs12 = 0;
00321 G4double xs21 = 0;
00322 G4double xs22 = 0;
00323
00324
00325 if (particleDefinition == G4Electron::ElectronDefinition())
00326 {
00327 std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k);
00328 std::vector<double>::iterator t1 = t2-1;
00329
00330 std::vector<double>::iterator e12 = std::upper_bound(eVecm[(*t1)].begin(),eVecm[(*t1)].end(), integrDiff);
00331 std::vector<double>::iterator e11 = e12-1;
00332
00333 std::vector<double>::iterator e22 = std::upper_bound(eVecm[(*t2)].begin(),eVecm[(*t2)].end(), integrDiff);
00334 std::vector<double>::iterator e21 = e22-1;
00335
00336 valueT1 =*t1;
00337 valueT2 =*t2;
00338 valueE21 =*e21;
00339 valueE22 =*e22;
00340 valueE12 =*e12;
00341 valueE11 =*e11;
00342
00343 xs11 = eDiffCrossSectionData[valueT1][valueE11];
00344 xs12 = eDiffCrossSectionData[valueT1][valueE12];
00345 xs21 = eDiffCrossSectionData[valueT2][valueE21];
00346 xs22 = eDiffCrossSectionData[valueT2][valueE22];
00347
00348 }
00349
00350 if (xs11==0 || xs12==0 ||xs21==0 ||xs22==0) return (0.);
00351
00352 theta = QuadInterpolator( valueE11, valueE12,
00353 valueE21, valueE22,
00354 xs11, xs12,
00355 xs21, xs22,
00356 valueT1, valueT2,
00357 k, integrDiff );
00358
00359 return theta;
00360 }
00361
00362
00363
00364 G4double G4MuElecElasticModel::LinLogInterpolate(G4double e1,
00365 G4double e2,
00366 G4double e,
00367 G4double xs1,
00368 G4double xs2)
00369 {
00370 G4double d1 = std::log(xs1);
00371 G4double d2 = std::log(xs2);
00372 G4double value = std::exp(d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
00373 return value;
00374 }
00375
00376
00377
00378 G4double G4MuElecElasticModel::LogLogInterpolate(G4double e1,
00379 G4double e2,
00380 G4double e,
00381 G4double xs1,
00382 G4double xs2)
00383 {
00384 G4double a = (std::log10(xs2)-std::log10(xs1)) / (std::log10(e2)-std::log10(e1));
00385 G4double b = std::log10(xs2) - a*std::log10(e2);
00386 G4double sigma = a*std::log10(e) + b;
00387 G4double value = (std::pow(10.,sigma));
00388 return value;
00389 }
00390
00391
00392
00393 G4double G4MuElecElasticModel::QuadInterpolator(G4double e11, G4double e12,
00394 G4double e21, G4double e22,
00395 G4double xs11, G4double xs12,
00396 G4double xs21, G4double xs22,
00397 G4double t1, G4double t2,
00398 G4double t, G4double e)
00399 {
00400
00401 G4double interpolatedvalue1 = LinLogInterpolate(e11, e12, e, xs11, xs12);
00402 G4double interpolatedvalue2 = LinLogInterpolate(e21, e22, e, xs21, xs22);
00403 G4double value = LinLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
00404 return value;
00405 }
00406
00407
00408
00409 G4double G4MuElecElasticModel::RandomizeCosTheta(G4double k)
00410 {
00411 G4double integrdiff=0;
00412 G4double uniformRand=G4UniformRand();
00413 integrdiff = uniformRand;
00414
00415 G4double theta=0.;
00416 G4double cosTheta=0.;
00417 theta = Theta(G4Electron::ElectronDefinition(),k/eV,integrdiff);
00418
00419 cosTheta= std::cos(theta*pi/180);
00420
00421 return cosTheta;
00422 }