[819] | 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 | // |
---|
[1055] | 26 | // $Id: G4BremsstrahlungCrossSectionHandler.cc,v 1.10 2009/03/03 11:19:18 pandola Exp $ |
---|
| 27 | // GEANT4 tag $Name: geant4-09-03-beta-cand-01 $ |
---|
[819] | 28 | // |
---|
| 29 | // ------------------------------------------------------------------- |
---|
| 30 | // |
---|
| 31 | // GEANT4 Class file |
---|
| 32 | // |
---|
| 33 | // |
---|
| 34 | // File name: G4BremsstrahlungCrossSectionHandler |
---|
| 35 | // |
---|
| 36 | // Author: V.Ivanchenko (Vladimir.Ivanchenko@cern.ch) |
---|
| 37 | // |
---|
| 38 | // Creation date: 25 September 2001 |
---|
| 39 | // |
---|
| 40 | // Modifications: |
---|
| 41 | // 10.10.2001 MGP Revision to improve code quality and consistency with design |
---|
| 42 | // 21.01.2003 VI cut per region |
---|
[1055] | 43 | // 03.03.2009 LP Added public method to make a easier migration of |
---|
| 44 | // G4LowEnergyBremsstrahlung to G4LivermoreBremsstrahlungModel |
---|
[819] | 45 | // |
---|
| 46 | // ------------------------------------------------------------------- |
---|
| 47 | |
---|
| 48 | #include "G4BremsstrahlungCrossSectionHandler.hh" |
---|
| 49 | #include "G4eBremsstrahlungSpectrum.hh" |
---|
| 50 | #include "G4DataVector.hh" |
---|
| 51 | #include "G4CompositeEMDataSet.hh" |
---|
| 52 | #include "G4VDataSetAlgorithm.hh" |
---|
| 53 | #include "G4SemiLogInterpolation.hh" |
---|
| 54 | #include "G4VEMDataSet.hh" |
---|
| 55 | #include "G4EMDataSet.hh" |
---|
| 56 | #include "G4Material.hh" |
---|
| 57 | #include "G4ProductionCutsTable.hh" |
---|
| 58 | |
---|
[1055] | 59 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 60 | |
---|
[819] | 61 | G4BremsstrahlungCrossSectionHandler::G4BremsstrahlungCrossSectionHandler(const G4VEnergySpectrum* spec, |
---|
| 62 | G4VDataSetAlgorithm* ) |
---|
| 63 | : theBR(spec) |
---|
| 64 | { |
---|
| 65 | interp = new G4SemiLogInterpolation(); |
---|
| 66 | } |
---|
| 67 | |
---|
[1055] | 68 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
[819] | 69 | |
---|
| 70 | G4BremsstrahlungCrossSectionHandler::~G4BremsstrahlungCrossSectionHandler() |
---|
| 71 | { |
---|
| 72 | delete interp; |
---|
| 73 | } |
---|
| 74 | |
---|
[1055] | 75 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
[819] | 76 | |
---|
| 77 | std::vector<G4VEMDataSet*>* |
---|
| 78 | G4BremsstrahlungCrossSectionHandler::BuildCrossSectionsForMaterials(const G4DataVector& energyVector, |
---|
| 79 | const G4DataVector* energyCuts) |
---|
| 80 | { |
---|
| 81 | std::vector<G4VEMDataSet*>* set = new std::vector<G4VEMDataSet*>; |
---|
| 82 | |
---|
| 83 | G4DataVector* energies; |
---|
| 84 | G4DataVector* cs; |
---|
| 85 | G4int nOfBins = energyVector.size(); |
---|
| 86 | |
---|
| 87 | const G4ProductionCutsTable* theCoupleTable= |
---|
| 88 | G4ProductionCutsTable::GetProductionCutsTable(); |
---|
| 89 | size_t numOfCouples = theCoupleTable->GetTableSize(); |
---|
| 90 | |
---|
| 91 | for (size_t m=0; m<numOfCouples; m++) { |
---|
| 92 | |
---|
| 93 | const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(m); |
---|
| 94 | const G4Material* material= couple->GetMaterial(); |
---|
| 95 | const G4ElementVector* elementVector = material->GetElementVector(); |
---|
| 96 | const G4double* nAtomsPerVolume = material->GetVecNbOfAtomsPerVolume(); |
---|
| 97 | G4int nElements = material->GetNumberOfElements(); |
---|
| 98 | |
---|
| 99 | G4double tcut = (*energyCuts)[m]; |
---|
| 100 | |
---|
| 101 | G4VDataSetAlgorithm* algo = interp->Clone(); |
---|
| 102 | G4VEMDataSet* setForMat = new G4CompositeEMDataSet(algo,1.,1.); |
---|
| 103 | |
---|
| 104 | for (G4int i=0; i<nElements; i++) { |
---|
| 105 | |
---|
| 106 | G4int Z = (G4int) ((*elementVector)[i]->GetZ()); |
---|
| 107 | energies = new G4DataVector; |
---|
| 108 | cs = new G4DataVector; |
---|
| 109 | G4double density = nAtomsPerVolume[i]; |
---|
| 110 | |
---|
| 111 | for (G4int bin=0; bin<nOfBins; bin++) { |
---|
| 112 | |
---|
| 113 | G4double e = energyVector[bin]; |
---|
| 114 | energies->push_back(e); |
---|
| 115 | G4double value = 0.0; |
---|
| 116 | |
---|
| 117 | if(e > tcut) { |
---|
| 118 | G4double elemCs = FindValue(Z, e); |
---|
| 119 | |
---|
| 120 | value = theBR->Probability(Z, tcut, e, e); |
---|
| 121 | |
---|
| 122 | value *= elemCs*density; |
---|
| 123 | } |
---|
| 124 | cs->push_back(value); |
---|
| 125 | } |
---|
| 126 | G4VDataSetAlgorithm* algol = interp->Clone(); |
---|
| 127 | G4VEMDataSet* elSet = new G4EMDataSet(i,energies,cs,algol,1.,1.); |
---|
| 128 | setForMat->AddComponent(elSet); |
---|
| 129 | } |
---|
| 130 | set->push_back(setForMat); |
---|
| 131 | } |
---|
| 132 | |
---|
| 133 | return set; |
---|
| 134 | } |
---|
[1055] | 135 | |
---|
| 136 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 137 | |
---|
| 138 | G4double G4BremsstrahlungCrossSectionHandler::GetCrossSectionAboveThresholdForElement(G4double energy, |
---|
| 139 | G4double cutEnergy, |
---|
| 140 | G4int Z) |
---|
| 141 | { |
---|
| 142 | G4double value = 0.; |
---|
| 143 | if(energy > cutEnergy) |
---|
| 144 | { |
---|
| 145 | G4double elemCs = FindValue(Z, energy); |
---|
| 146 | value = theBR->Probability(Z,cutEnergy, energy, energy); |
---|
| 147 | value *= elemCs; |
---|
| 148 | } |
---|
| 149 | return value; |
---|
| 150 | } |
---|