| [1350] | 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: G4CrossSectionPairGG.cc,v 1.5 2010/11/18 11:01:01 gunter Exp $
|
|---|
| 27 | // $ GEANT4 tag $Name: geant4-09-04-ref-00 $
|
|---|
| 28 | //
|
|---|
| 29 | // Class G4CrossSectionPairGG
|
|---|
| 30 | //
|
|---|
| 31 | // smoothly join two cross section sets by scaling the second at a given
|
|---|
| 32 | // transition energy to match the first.
|
|---|
| 33 | //
|
|---|
| 34 | // Author: Gunter Folger
|
|---|
| 35 | // November 2009
|
|---|
| 36 | //
|
|---|
| 37 |
|
|---|
| 38 | #include "G4CrossSectionPairGG.hh"
|
|---|
| 39 |
|
|---|
| 40 | #include "globals.hh"
|
|---|
| 41 | #include "G4HadTmpUtil.hh"
|
|---|
| 42 | #include "G4NistManager.hh"
|
|---|
| 43 | #include "G4ThreeVector.hh"
|
|---|
| 44 |
|
|---|
| 45 | G4CrossSectionPairGG::G4CrossSectionPairGG(G4VCrossSectionDataSet * low,
|
|---|
| 46 | G4double Etransit):
|
|---|
| 47 | theLowX(low),
|
|---|
| 48 | ETransition(Etransit)
|
|---|
| 49 | {
|
|---|
| 50 | theHighX=new G4GlauberGribovCrossSection();
|
|---|
| 51 | verboseLevel=0;
|
|---|
| 52 | }
|
|---|
| 53 |
|
|---|
| 54 | G4CrossSectionPairGG::~G4CrossSectionPairGG()
|
|---|
| 55 | {
|
|---|
| 56 | // The cross section registry wil delete these
|
|---|
| 57 | // delete theLowX;
|
|---|
| 58 | // delete theHighX;
|
|---|
| 59 | }
|
|---|
| 60 |
|
|---|
| 61 | G4bool G4CrossSectionPairGG::IsIsoApplicable(const G4DynamicParticle* aParticle,
|
|---|
| 62 | G4int ZZ, G4int AA)
|
|---|
| 63 | {
|
|---|
| 64 | G4bool isApplicable(false);
|
|---|
| 65 | G4double Ekin=aParticle->GetKineticEnergy();
|
|---|
| 66 | if (Ekin < ETransition )
|
|---|
| 67 | {
|
|---|
| 68 | isApplicable = theLowX->IsIsoApplicable(aParticle,ZZ,AA);
|
|---|
| 69 | } else {
|
|---|
| 70 | isApplicable = theHighX->IsIsoApplicable(aParticle,ZZ,AA);
|
|---|
| 71 | }
|
|---|
| 72 |
|
|---|
| 73 | return isApplicable;
|
|---|
| 74 | }
|
|---|
| 75 |
|
|---|
| 76 | G4double G4CrossSectionPairGG::GetZandACrossSection(const G4DynamicParticle* aParticle,
|
|---|
| 77 | G4int ZZ, G4int AA,
|
|---|
| 78 | G4double aTemperature)
|
|---|
| 79 | {
|
|---|
| 80 | G4double Xsec(0.);
|
|---|
| 81 | std::vector<ParticleXScale>::iterator iter;
|
|---|
| 82 | iter=scale_factors.begin();
|
|---|
| 83 | G4ParticleDefinition * pDef=aParticle->GetDefinition();
|
|---|
| 84 | while ( iter !=scale_factors.end() && (*iter).first != pDef ) {++iter;}
|
|---|
| 85 |
|
|---|
| 86 |
|
|---|
| 87 | G4double Ekin=aParticle->GetKineticEnergy();
|
|---|
| 88 | if (Ekin < ETransition )
|
|---|
| 89 | {
|
|---|
| 90 | Xsec=theLowX->GetZandACrossSection(aParticle,ZZ,AA,aTemperature);
|
|---|
| 91 | } else {
|
|---|
| 92 | Xsec=theHighX->GetInelasticGlauberGribov(aParticle,ZZ,AA)
|
|---|
| 93 | * (*iter).second[ZZ];
|
|---|
| 94 | if ( verboseLevel > 2 )
|
|---|
| 95 | { G4cout << " scaling .." << ZZ << " " << AA << " " <<
|
|---|
| 96 | (*iter).second[ZZ]<< " " <<theHighX->GetInelasticGlauberGribov(aParticle,ZZ,AA) << " "
|
|---|
| 97 | << Xsec << G4endl;
|
|---|
| 98 | }
|
|---|
| 99 | }
|
|---|
| 100 |
|
|---|
| 101 | return Xsec;
|
|---|
| 102 | }
|
|---|
| 103 |
|
|---|
| 104 |
|
|---|
| 105 |
|
|---|
| 106 | void G4CrossSectionPairGG::BuildPhysicsTable(const G4ParticleDefinition& pDef)
|
|---|
| 107 | {
|
|---|
| 108 | theLowX->BuildPhysicsTable(pDef);
|
|---|
| 109 | theHighX->BuildPhysicsTable(pDef);
|
|---|
| 110 |
|
|---|
| 111 | G4NistManager* NistMan = G4NistManager::Instance();
|
|---|
| 112 | G4ParticleDefinition * myDef=const_cast<G4ParticleDefinition*>(&pDef);
|
|---|
| 113 | std::vector<ParticleXScale>::iterator iter;
|
|---|
| 114 | iter=scale_factors.begin();
|
|---|
| 115 | while ( iter !=scale_factors.end() && (*iter).first != myDef ) {++iter;}
|
|---|
| 116 |
|
|---|
| 117 | // new particle, initialise
|
|---|
| 118 |
|
|---|
| 119 | if ( iter == scale_factors.end() )
|
|---|
| 120 | {
|
|---|
| 121 | XS_factors factors (93);
|
|---|
| 122 | G4ThreeVector mom(0.0,0.0,1.0);
|
|---|
| 123 | G4DynamicParticle DynPart(myDef, mom, ETransition); // last is kinetic Energy
|
|---|
| 124 |
|
|---|
| 125 | if (verboseLevel > 0) {
|
|---|
| 126 | G4cout << "G4CrossSectionPairGG::BuildPhysicsTable for particle "
|
|---|
| 127 | << pDef.GetParticleName() << G4endl;
|
|---|
| 128 | }
|
|---|
| 129 | for (G4int aZ=1; aZ<93; ++aZ)
|
|---|
| 130 | {
|
|---|
| 131 | factors[aZ]=1.; // default, to give reasonable value if only high is applicable
|
|---|
| 132 | G4int AA=G4lrint(NistMan->GetAtomicMassAmu(aZ));
|
|---|
| 133 | G4bool isApplicable = theLowX->IsIsoApplicable(&DynPart, aZ, AA) &&
|
|---|
| 134 | theHighX->IsIsoApplicable(&DynPart, aZ, AA-aZ);
|
|---|
| 135 |
|
|---|
| 136 | if (isApplicable)
|
|---|
| 137 | {
|
|---|
| 138 | factors[aZ]=theLowX->GetZandACrossSection(&DynPart,aZ,AA,0) /
|
|---|
| 139 | theHighX->GetInelasticGlauberGribov(&DynPart,aZ,AA);
|
|---|
| 140 |
|
|---|
| 141 | }
|
|---|
| 142 | if (verboseLevel > 0) {
|
|---|
| 143 | G4cout << "Z=" << aZ<< ", A="<< AA << ", scale=" << factors[aZ];
|
|---|
| 144 | if ( verboseLevel == 1) { G4cout << G4endl; }
|
|---|
| 145 | else {
|
|---|
| 146 | if (isApplicable) {
|
|---|
| 147 | G4cout << ", low / high " << theLowX->GetZandACrossSection(&DynPart,aZ,AA,0)
|
|---|
| 148 | << " " << theHighX->GetInelasticGlauberGribov(&DynPart,aZ,AA) << G4endl;
|
|---|
| 149 | } else { G4cout << ", N/A" << G4endl; }
|
|---|
| 150 | }
|
|---|
| 151 | }
|
|---|
| 152 | }
|
|---|
| 153 | ParticleXScale forPart(myDef,factors);
|
|---|
| 154 | scale_factors.push_back(forPart);
|
|---|
| 155 | }
|
|---|
| 156 | }
|
|---|
| 157 | void G4CrossSectionPairGG::DumpPhysicsTable(const G4ParticleDefinition&)
|
|---|
| 158 | {
|
|---|
| 159 | }
|
|---|