// // ******************************************************************** // * License and Disclaimer * // * * // * The Geant4 software is copyright of the Copyright Holders of * // * the Geant4 Collaboration. It is provided under the terms and * // * conditions of the Geant4 Software License, included in the file * // * LICENSE and available at http://cern.ch/geant4/license . These * // * include a list of copyright holders. * // * * // * Neither the authors of this software system, nor their employing * // * institutes,nor the agencies providing financial support for this * // * work make any representation or warranty, express or implied, * // * regarding this software system or assume any liability for its * // * use. Please see the license in the file LICENSE and URL above * // * for the full disclaimer and the limitation of liability. * // * * // * This code implementation is the result of the scientific and * // * technical work of the GEANT4 collaboration. * // * By using, copying, modifying or distributing the software (or * // * any work based on the software) you agree to acknowledge its * // * use in resulting scientific publications, and indicate your * // * acceptance of all terms of the Geant4 Software license. * // ******************************************************************** // // // $Id: G4PAIxSection.cc,v 1.24 2008/05/30 16:04:40 grichine Exp $ // GEANT4 tag $Name: geant4-09-03-cand-01 $ // // // G4PAIxSection.cc -- class implementation file // // GEANT 4 class implementation file // // For information related to this code, please, contact // the Geant4 Collaboration. // // R&D: Vladimir.Grichine@cern.ch // // History: // // 13.05.03 V. Grichine, bug fixed for maxEnergyTransfer > max interval energy // 28.05.01 V.Ivanchenko minor changes to provide ANSI -wall compilation // 17.05.01 V. Grichine, low energy extension down to 10*keV of proton // 20.11.98 adapted to a new Material/SandiaTable interface, mma // 11.06.97 V. Grichine, 1st version // #include "G4PAIxSection.hh" #include "globals.hh" #include "G4ios.hh" #include "G4Poisson.hh" #include "G4Material.hh" #include "G4MaterialCutsCouple.hh" #include "G4SandiaTable.hh" using namespace std; /* ****************************************************************** // Init array of Lorentz factors const G4double G4PAIxSection::fLorentzFactor[22] = { 0.0 , 1.1 , 1.2 , 1.3 , 1.5 , 1.8 , 2.0 , 2.5 , 3.0 , 4.0 , 7.0 , 10.0 , 20.0 , 40.0 , 70.0 , 100.0 , 300.0 , 600.0 , 1000.0 , 3000.0 , 10000.0 , 50000.0 }; const G4int G4PAIxSection:: fRefGammaNumber = 29; // The number of gamma for creation of // spline (9) ***************************************************************** */ // Local class constants const G4double G4PAIxSection::fDelta = 0.005; // energy shift from interval border const G4double G4PAIxSection::fError = 0.005; // error in lin-log approximation const G4int G4PAIxSection::fMaxSplineSize = 500; // Max size of output spline // arrays ////////////////////////////////////////////////////////////////// // // Constructor // G4PAIxSection::G4PAIxSection(G4MaterialCutsCouple* matCC) { fDensity = matCC->GetMaterial()->GetDensity(); G4int matIndex = matCC->GetMaterial()->GetIndex(); fMaterialIndex = matIndex; fSandia = new G4SandiaTable(matIndex); G4int i, j; fMatSandiaMatrix = new G4OrderedTable(); for (i = 0; i < fSandia->GetMaxInterval()-1; i++) { fMatSandiaMatrix->push_back(new G4DataVector(5,0.)); } for (i = 0; i < fSandia->GetMaxInterval()-1; i++) { (*(*fMatSandiaMatrix)[i])[0] = fSandia->GetSandiaMatTable(i,0); for(j = 1; j < 5; j++) { (*(*fMatSandiaMatrix)[i])[j] = fSandia->GetSandiaMatTable(i,j)*fDensity; } } } //////////////////////////////////////////////////////////////// G4PAIxSection::G4PAIxSection(G4int materialIndex, G4double maxEnergyTransfer) { const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable(); G4int i, j; fMaterialIndex = materialIndex; fDensity = (*theMaterialTable)[materialIndex]->GetDensity(); fElectronDensity = (*theMaterialTable)[materialIndex]-> GetElectronDensity(); fIntervalNumber = (*theMaterialTable)[materialIndex]-> GetSandiaTable()->GetMatNbOfIntervals(); fIntervalNumber--; // G4cout< GetSandiaTable()->GetSandiaCofForMaterial(i-1,0) >= maxEnergyTransfer) || i > fIntervalNumber ) { fEnergyInterval[i] = maxEnergyTransfer; fIntervalNumber = i; break; } fEnergyInterval[i] = (*theMaterialTable)[materialIndex]-> GetSandiaTable()->GetSandiaCofForMaterial(i-1,0); fA1[i] = (*theMaterialTable)[materialIndex]-> GetSandiaTable()->GetSandiaCofForMaterial(i-1,1); fA2[i] = (*theMaterialTable)[materialIndex]-> GetSandiaTable()->GetSandiaCofForMaterial(i-1,2); fA3[i] = (*theMaterialTable)[materialIndex]-> GetSandiaTable()->GetSandiaCofForMaterial(i-1,3); fA4[i] = (*theMaterialTable)[materialIndex]-> GetSandiaTable()->GetSandiaCofForMaterial(i-1,4); // G4cout< 1.5*fDelta*(fEnergyInterval[i+1]+fEnergyInterval[i])) { continue; } else { for(j=i;jGetDensity(); fElectronDensity = (*theMaterialTable)[materialIndex]-> GetElectronDensity(); fIntervalNumber = intNumber; fIntervalNumber--; // G4cout<= maxEnergyTransfer ) || i > fIntervalNumber ) { fEnergyInterval[i] = maxEnergyTransfer; fIntervalNumber = i; break; } fEnergyInterval[i] = photoAbsCof[i-1][0]; fA1[i] = photoAbsCof[i-1][1]; fA2[i] = photoAbsCof[i-1][2]; fA3[i] = photoAbsCof[i-1][3]; fA4[i] = photoAbsCof[i-1][4]; // G4cout< 1.5*fDelta*(fEnergyInterval[i+1]+fEnergyInterval[i])) { continue; } else { for(j=i;j 1.5*fDelta*(fEnergyInterval[i+1]+fEnergyInterval[i])) { continue; } else { for( j = i; j < fIntervalNumber; j++ ) { fEnergyInterval[j] = fEnergyInterval[j+1]; fA1[j] = fA1[j+1]; fA2[j] = fA2[j+1]; fA3[j] = fA3[j+1]; fA4[j] = fA4[j+1]; } fIntervalNumber--; i--; } } /* ********************************* fSplineEnergy = new G4double[fMaxSplineSize]; fRePartDielectricConst = new G4double[fMaxSplineSize]; fImPartDielectricConst = new G4double[fMaxSplineSize]; fIntegralTerm = new G4double[fMaxSplineSize]; fDifPAIxSection = new G4double[fMaxSplineSize]; fIntegralPAIxSection = new G4double[fMaxSplineSize]; for(i=0;i= fIntegralPAIxSection[iTransfer] ) break; } if(iTransfer > fSplineNumber) iTransfer--; energyTransfer = fSplineEnergy[iTransfer]; if(iTransfer > 1) { energyTransfer -= (fSplineEnergy[iTransfer]-fSplineEnergy[iTransfer-1])*G4UniformRand(); } return energyTransfer; } ///////////////////////////////////////////////////////////////////////// // // Returns random Cerenkov energy loss over step G4double G4PAIxSection::GetStepCerenkovLoss( G4double step ) { G4long numOfCollisions; G4double meanNumber, loss = 0.0; // G4cout<<" G4PAIxSection::GetStepCerenkovLoss "<= fIntegralCerenkov[iTransfer] ) break; } if(iTransfer > fSplineNumber) iTransfer--; energyTransfer = fSplineEnergy[iTransfer]; if(iTransfer > 1) { energyTransfer -= (fSplineEnergy[iTransfer]-fSplineEnergy[iTransfer-1])*G4UniformRand(); } return energyTransfer; } ///////////////////////////////////////////////////////////////////////// // // Returns MM-Cerenkov energy transfer in one collision G4double G4PAIxSection::GetMMEnergyTransfer() { G4int iTransfer ; G4double energyTransfer, position; position = fIntegralMM[1]*G4UniformRand(); for( iTransfer = 1; iTransfer <= fSplineNumber; iTransfer++ ) { if( position >= fIntegralMM[iTransfer] ) break; } if(iTransfer > fSplineNumber) iTransfer--; energyTransfer = fSplineEnergy[iTransfer]; if(iTransfer > 1) { energyTransfer -= (fSplineEnergy[iTransfer]-fSplineEnergy[iTransfer-1])*G4UniformRand(); } return energyTransfer; } ///////////////////////////////////////////////////////////////////////// // // Returns random plasmon energy loss over step G4double G4PAIxSection::GetStepPlasmonLoss( G4double step ) { G4long numOfCollisions; G4double meanNumber, loss = 0.0; // G4cout<<" G4PAIxSection::GetStepPlasmonLoss "<= fIntegralPlasmon[iTransfer] ) break; } if(iTransfer > fSplineNumber) iTransfer--; energyTransfer = fSplineEnergy[iTransfer]; if(iTransfer > 1) { energyTransfer -= (fSplineEnergy[iTransfer]-fSplineEnergy[iTransfer-1])*G4UniformRand(); } return energyTransfer; } ///////////////////////////////////////////////////////////////////////// // // Returns random resonance energy loss over step G4double G4PAIxSection::GetStepResonanceLoss( G4double step ) { G4long numOfCollisions; G4double meanNumber, loss = 0.0; // G4cout<<" G4PAIxSection::GetStepCreLosnkovs "<= fIntegralResonance[iTransfer] ) break; } if(iTransfer > fSplineNumber) iTransfer--; energyTransfer = fSplineEnergy[iTransfer]; if(iTransfer > 1) { energyTransfer -= (fSplineEnergy[iTransfer]-fSplineEnergy[iTransfer-1])*G4UniformRand(); } return energyTransfer; } ///////////////////////////////////////////////////////////////////////// // // Returns Rutherford energy transfer in one collision G4double G4PAIxSection::GetRutherfordEnergyTransfer() { G4int iTransfer ; G4double energyTransfer, position; position = (fIntegralPlasmon[1]-fIntegralResonance[1])*G4UniformRand(); for( iTransfer = 1; iTransfer <= fSplineNumber; iTransfer++ ) { if( position >= (fIntegralPlasmon[iTransfer]-fIntegralResonance[iTransfer]) ) break; } if(iTransfer > fSplineNumber) iTransfer--; energyTransfer = fSplineEnergy[iTransfer]; if(iTransfer > 1) { energyTransfer -= (fSplineEnergy[iTransfer]-fSplineEnergy[iTransfer-1])*G4UniformRand(); } return energyTransfer; } ///////////////////////////////////////////////////////////////////////////// // // Init array of Lorentz factors // G4int G4PAIxSection::fNumberOfGammas = 111; const G4double G4PAIxSection::fLorentzFactor[112] = // fNumberOfGammas+1 { 0.0, 1.094989e+00, 1.107813e+00, 1.122369e+00, 1.138890e+00, 1.157642e+00, 1.178925e+00, 1.203082e+00, 1.230500e+00, 1.261620e+00, 1.296942e+00, // 10 1.337032e+00, 1.382535e+00, 1.434181e+00, 1.492800e+00, 1.559334e+00, 1.634850e+00, 1.720562e+00, 1.817845e+00, 1.928263e+00, 2.053589e+00, // 20 2.195835e+00, 2.357285e+00, 2.540533e+00, 2.748522e+00, 2.984591e+00, 3.252533e+00, 3.556649e+00, 3.901824e+00, 4.293602e+00, 4.738274e+00, // 30 5.242981e+00, 5.815829e+00, 6.466019e+00, 7.203990e+00, 8.041596e+00, 8.992288e+00, 1.007133e+01, 1.129606e+01, 1.268614e+01, 1.426390e+01, // 40 1.605467e+01, 1.808721e+01, 2.039417e+01, 2.301259e+01, 2.598453e+01, 2.935771e+01, 3.318630e+01, 3.753180e+01, 4.246399e+01, 4.806208e+01, // 50 5.441597e+01, 6.162770e+01, 6.981310e+01, 7.910361e+01, 8.964844e+01, 1.016169e+02, 1.152013e+02, 1.306197e+02, 1.481198e+02, 1.679826e+02, // 60 1.905270e+02, 2.161152e+02, 2.451581e+02, 2.781221e+02, 3.155365e+02, 3.580024e+02, 4.062016e+02, 4.609081e+02, 5.230007e+02, 5.934765e+02, // 70 6.734672e+02, 7.642575e+02, 8.673056e+02, 9.842662e+02, 1.117018e+03, 1.267692e+03, 1.438709e+03, 1.632816e+03, 1.853128e+03, 2.103186e+03, // 80 2.387004e+03, 2.709140e+03, 3.074768e+03, 3.489760e+03, 3.960780e+03, 4.495394e+03, 5.102185e+03, 5.790900e+03, 6.572600e+03, 7.459837e+03, // 90 8.466860e+03, 9.609843e+03, 1.090714e+04, 1.237959e+04, 1.405083e+04, 1.594771e+04, 1.810069e+04, 2.054434e+04, 2.331792e+04, 2.646595e+04, // 100 3.003901e+04, 3.409446e+04, 3.869745e+04, 4.392189e+04, 4.985168e+04, 5.658206e+04, 6.422112e+04, 7.289153e+04, 8.273254e+04, 9.390219e+04, // 110 1.065799e+05 }; /////////////////////////////////////////////////////////////////////// // // The number of gamma for creation of spline (near ion-min , G ~ 4 ) // const G4int G4PAIxSection::fRefGammaNumber = 29; // // end of G4PAIxSection implementation file // ////////////////////////////////////////////////////////////////////////////