| 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: G4BremsstrahlungCrossSectionHandler.cc,v 1.12 2009/09/27 10:47:42 sincerti Exp $
|
|---|
| 27 | // GEANT4 tag $Name: geant4-09-03 $
|
|---|
| 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 | //
|
|---|
| 42 | // 10.10.2001 MGP Revision to improve code quality and consistency with design
|
|---|
| 43 | // 21.01.2003 VI cut per region
|
|---|
| 44 | // 03.03.2009 LP Added public method to make a easier migration of
|
|---|
| 45 | // G4LowEnergyBremsstrahlung to G4LivermoreBremsstrahlungModel
|
|---|
| 46 | //
|
|---|
| 47 | // 15 Jul 2009 Nicolas A. Karakatsanis
|
|---|
| 48 | //
|
|---|
| 49 | // - BuildCrossSectionForMaterials method was revised in order to calculate the
|
|---|
| 50 | // logarithmic values of the loaded data.
|
|---|
| 51 | // It retrieves the data values from the G4EMLOW data files but, then, calculates the
|
|---|
| 52 | // respective log values and loads them to seperate data structures.
|
|---|
| 53 | // The EM data sets, initialized this way, contain both non-log and log values.
|
|---|
| 54 | // These initialized data sets can enhance the computing performance of data interpolation
|
|---|
| 55 | // operations
|
|---|
| 56 | //
|
|---|
| 57 | //
|
|---|
| 58 | //
|
|---|
| 59 | //
|
|---|
| 60 | // -------------------------------------------------------------------
|
|---|
| 61 |
|
|---|
| 62 | #include "G4BremsstrahlungCrossSectionHandler.hh"
|
|---|
| 63 | #include "G4eBremsstrahlungSpectrum.hh"
|
|---|
| 64 | #include "G4DataVector.hh"
|
|---|
| 65 | #include "G4CompositeEMDataSet.hh"
|
|---|
| 66 | #include "G4VDataSetAlgorithm.hh"
|
|---|
| 67 | #include "G4SemiLogInterpolation.hh"
|
|---|
| 68 | #include "G4VEMDataSet.hh"
|
|---|
| 69 | #include "G4EMDataSet.hh"
|
|---|
| 70 | #include "G4Material.hh"
|
|---|
| 71 | #include "G4ProductionCutsTable.hh"
|
|---|
| 72 |
|
|---|
| 73 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 74 |
|
|---|
| 75 | G4BremsstrahlungCrossSectionHandler::G4BremsstrahlungCrossSectionHandler(const G4VEnergySpectrum* spec,
|
|---|
| 76 | G4VDataSetAlgorithm* )
|
|---|
| 77 | : theBR(spec)
|
|---|
| 78 | {
|
|---|
| 79 | interp = new G4SemiLogInterpolation();
|
|---|
| 80 | }
|
|---|
| 81 |
|
|---|
| 82 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 83 |
|
|---|
| 84 | G4BremsstrahlungCrossSectionHandler::~G4BremsstrahlungCrossSectionHandler()
|
|---|
| 85 | {
|
|---|
| 86 | delete interp;
|
|---|
| 87 | }
|
|---|
| 88 |
|
|---|
| 89 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 90 |
|
|---|
| 91 | std::vector<G4VEMDataSet*>*
|
|---|
| 92 | G4BremsstrahlungCrossSectionHandler::BuildCrossSectionsForMaterials(const G4DataVector& energyVector,
|
|---|
| 93 | const G4DataVector* energyCuts)
|
|---|
| 94 | {
|
|---|
| 95 | std::vector<G4VEMDataSet*>* set = new std::vector<G4VEMDataSet*>;
|
|---|
| 96 |
|
|---|
| 97 | G4DataVector* energies;
|
|---|
| 98 | G4DataVector* cs;
|
|---|
| 99 |
|
|---|
| 100 | G4DataVector* log_energies;
|
|---|
| 101 | G4DataVector* log_cs;
|
|---|
| 102 |
|
|---|
| 103 | G4int nOfBins = energyVector.size();
|
|---|
| 104 |
|
|---|
| 105 | const G4ProductionCutsTable* theCoupleTable=
|
|---|
| 106 | G4ProductionCutsTable::GetProductionCutsTable();
|
|---|
| 107 | size_t numOfCouples = theCoupleTable->GetTableSize();
|
|---|
| 108 |
|
|---|
| 109 | for (size_t m=0; m<numOfCouples; m++) {
|
|---|
| 110 |
|
|---|
| 111 | const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(m);
|
|---|
| 112 | const G4Material* material= couple->GetMaterial();
|
|---|
| 113 | const G4ElementVector* elementVector = material->GetElementVector();
|
|---|
| 114 | const G4double* nAtomsPerVolume = material->GetVecNbOfAtomsPerVolume();
|
|---|
| 115 | G4int nElements = material->GetNumberOfElements();
|
|---|
| 116 |
|
|---|
| 117 | G4double tcut = (*energyCuts)[m];
|
|---|
| 118 |
|
|---|
| 119 | G4VDataSetAlgorithm* algo = interp->Clone();
|
|---|
| 120 | G4VEMDataSet* setForMat = new G4CompositeEMDataSet(algo,1.,1.);
|
|---|
| 121 |
|
|---|
| 122 | for (G4int i=0; i<nElements; i++) {
|
|---|
| 123 |
|
|---|
| 124 | G4int Z = (G4int) ((*elementVector)[i]->GetZ());
|
|---|
| 125 |
|
|---|
| 126 | energies = new G4DataVector;
|
|---|
| 127 | cs = new G4DataVector;
|
|---|
| 128 |
|
|---|
| 129 | log_energies = new G4DataVector;
|
|---|
| 130 | log_cs = new G4DataVector;
|
|---|
| 131 |
|
|---|
| 132 | G4double density = nAtomsPerVolume[i];
|
|---|
| 133 |
|
|---|
| 134 | for (G4int bin=0; bin<nOfBins; bin++) {
|
|---|
| 135 |
|
|---|
| 136 | G4double e = energyVector[bin];
|
|---|
| 137 | energies->push_back(e);
|
|---|
| 138 | if (e==0.) e=1e-300;
|
|---|
| 139 | log_energies->push_back(std::log10(e));
|
|---|
| 140 | G4double value = 0.0;
|
|---|
| 141 |
|
|---|
| 142 | if(e > tcut) {
|
|---|
| 143 | G4double elemCs = FindValue(Z, e);
|
|---|
| 144 |
|
|---|
| 145 | value = theBR->Probability(Z, tcut, e, e);
|
|---|
| 146 |
|
|---|
| 147 | value *= elemCs*density;
|
|---|
| 148 | }
|
|---|
| 149 | cs->push_back(value);
|
|---|
| 150 |
|
|---|
| 151 | if (value==0.) value=1e-300;
|
|---|
| 152 | log_cs->push_back(std::log10(value));
|
|---|
| 153 | }
|
|---|
| 154 | G4VDataSetAlgorithm* algol = interp->Clone();
|
|---|
| 155 |
|
|---|
| 156 | //G4VEMDataSet* elSet = new G4EMDataSet(i,energies,cs,algol,1.,1.);
|
|---|
| 157 |
|
|---|
| 158 | G4VEMDataSet* elSet = new G4EMDataSet(i,energies,cs,log_energies,log_cs,algol,1.,1.);
|
|---|
| 159 |
|
|---|
| 160 | setForMat->AddComponent(elSet);
|
|---|
| 161 | }
|
|---|
| 162 | set->push_back(setForMat);
|
|---|
| 163 | }
|
|---|
| 164 |
|
|---|
| 165 | return set;
|
|---|
| 166 | }
|
|---|
| 167 |
|
|---|
| 168 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 169 |
|
|---|
| 170 | G4double G4BremsstrahlungCrossSectionHandler::GetCrossSectionAboveThresholdForElement(G4double energy,
|
|---|
| 171 | G4double cutEnergy,
|
|---|
| 172 | G4int Z)
|
|---|
| 173 | {
|
|---|
| 174 | G4double value = 0.;
|
|---|
| 175 | if(energy > cutEnergy)
|
|---|
| 176 | {
|
|---|
| 177 | G4double elemCs = FindValue(Z, energy);
|
|---|
| 178 | value = theBR->Probability(Z,cutEnergy, energy, energy);
|
|---|
| 179 | value *= elemCs;
|
|---|
| 180 | }
|
|---|
| 181 | return value;
|
|---|
| 182 | }
|
|---|