// // ******************************************************************** // * 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: G4CollisionComposite.cc,v 1.7 2007/07/05 21:43:24 dennis Exp $ // #include "globals.hh" #include "G4CollisionComposite.hh" #include "G4VCollision.hh" #include "G4CollisionVector.hh" #include "G4KineticTrack.hh" #include "G4KineticTrackVector.hh" #include "G4VCrossSectionSource.hh" #include "G4HadTmpUtil.hh" const G4int G4CollisionComposite::nPoints = 32; G4double G4CollisionComposite::theT[nPoints] = {.01, .03, .05, .1, .15, .2, .3, .4, .5, .6, .7, .8, .9, 1.0, 1.2, 1.4, 1.6, 1.8, 2.0, 2.2, 2.5, 3.0, 3.5, 4.0, 5.0, 6.0, 8.0, 10., 15, 20, 50, 100}; G4CollisionComposite::G4CollisionComposite() { } G4CollisionComposite::~G4CollisionComposite() { std::for_each(components.begin(), components.end(), G4Delete()); } G4double G4CollisionComposite::CrossSection(const G4KineticTrack& trk1, const G4KineticTrack& trk2) const { G4double crossSect = 0.; const G4VCrossSectionSource* xSource = GetCrossSectionSource(); if (xSource != 0) // There is a total cross section for this Collision { crossSect = xSource->CrossSection(trk1,trk2); } else { // waiting for mutable to enable buffering. const_cast(this)->BufferCrossSection(trk1.GetDefinition(), trk2.GetDefinition()); // G4cerr << "Buffer filled, reying with sqrts = "<< (trk1.Get4Momentum()+trk2.Get4Momentum()).mag() < cxCache; G4double partialCxSum = 0.0; size_t i; for (i=0; iGetName(); if (components[i]->IsInCharge(trk1,trk2)) { partialCx = components[i]->CrossSection(trk1,trk2); } else { partialCx = 0.0; } // cout << " cx=" << partialCx << endl; partialCxSum += partialCx; cxCache.push_back(partialCx); } G4double random = G4UniformRand()*partialCxSum; G4double running = 0; for (i=0; i random) { return components[i]->FinalState(trk1, trk2); } } // G4cerr <<"in charge = "<GetParticleName()<begin(); iter != comps->end(); ++iter) { if ( ((*iter))->IsInCharge(trk1,trk2) ) isInCharge = true; } } return isInCharge; } void G4CollisionComposite:: BufferCrossSection(const G4ParticleDefinition * aP, const G4ParticleDefinition * bP) { // check if already buffered size_t i; for(i=0; iGetParticleName()<<" "<GetParticleName()<GetPDGMass(); G4double aE = aM+aT; G4ThreeVector aMom(0,0,std::sqrt(aE*aE-aM*aM)); G4LorentzVector a4Momentum(aE, aMom); G4KineticTrack a(const_cast(aP), atime, aPosition, a4Momentum); G4double btime = 0; G4ThreeVector bPosition(0,0,0); G4ThreeVector bMom(0,0,0*MeV); G4double bE = bP->GetPDGMass(); G4LorentzVector b4Momentum(bE, bMom); G4KineticTrack b(const_cast(bP), btime, bPosition, b4Momentum); for (i=0; iIsInCharge(a,b)) { crossSect += components[i]->CrossSection(a,b); } } G4double sqrts = (a4Momentum+b4Momentum).mag(); aNewBuff.push_back(sqrts, crossSect); } theBuffer.push_back(aNewBuff); // theBuffer.back().Print(); } G4double G4CollisionComposite:: BufferedCrossSection(const G4KineticTrack& trk1, const G4KineticTrack& trk2) const { for(size_t i=0; i