G4CollisionComposite.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id: G4CollisionComposite.cc,v 1.9 2010-03-12 15:45:18 gunter Exp $ //
00028 
00029 #include "globals.hh"
00030 #include "G4SystemOfUnits.hh"
00031 #include "G4CollisionComposite.hh"
00032 #include "G4VCollision.hh"
00033 #include "G4CollisionVector.hh"
00034 #include "G4KineticTrack.hh"
00035 #include "G4KineticTrackVector.hh"
00036 #include "G4VCrossSectionSource.hh"
00037 #include "G4HadTmpUtil.hh"
00038 
00039 const G4int G4CollisionComposite::nPoints = 32;
00040 
00041 G4double G4CollisionComposite::theT[nPoints] = 
00042 {.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};
00043 
00044 G4CollisionComposite::G4CollisionComposite()
00045 { 
00046 }
00047 
00048 
00049 G4CollisionComposite::~G4CollisionComposite()
00050 { 
00051   std::for_each(components.begin(), components.end(), G4Delete());
00052 }
00053 
00054 
00055 G4double G4CollisionComposite::CrossSection(const G4KineticTrack& trk1, 
00056                                             const G4KineticTrack& trk2) const
00057 {
00058   G4double crossSect = 0.;
00059   const G4VCrossSectionSource* xSource = GetCrossSectionSource();
00060   if (xSource != 0)
00061   // There is a total cross section for this Collision
00062   {
00063     crossSect = xSource->CrossSection(trk1,trk2);
00064   }
00065   else
00066   {
00067     // waiting for mutable to enable buffering.
00068     const_cast<G4CollisionComposite *>(this)->BufferCrossSection(trk1.GetDefinition(), trk2.GetDefinition());
00069 //    G4cerr << "Buffer filled, reying with sqrts = "<< (trk1.Get4Momentum()+trk2.Get4Momentum()).mag() <<G4endl;
00070     crossSect = BufferedCrossSection(trk1,trk2);
00071   }
00072   return crossSect;
00073 }
00074 
00075 
00076 G4KineticTrackVector* G4CollisionComposite::FinalState(const G4KineticTrack& trk1, 
00077                                                           const G4KineticTrack& trk2) const
00078 {
00079   std::vector<G4double> cxCache;
00080   G4double partialCxSum = 0.0;
00081 
00082   size_t i;
00083   for (i=0; i<components.size(); i++) 
00084   {
00085     G4double partialCx;
00086 //    cout << "comp" << i << " " << components[i]()->GetName();
00087     if (components[i]->IsInCharge(trk1,trk2)) 
00088     {
00089       partialCx = components[i]->CrossSection(trk1,trk2);
00090     } 
00091     else 
00092     {
00093       partialCx = 0.0;
00094     }
00095 //    cout << "   cx=" << partialCx << endl;
00096     partialCxSum += partialCx;
00097     cxCache.push_back(partialCx);
00098   }
00099 
00100   G4double random = G4UniformRand()*partialCxSum;
00101   G4double running = 0;
00102   for (i=0; i<cxCache.size(); i++) 
00103   {
00104     running += cxCache[i];
00105     if (running > random) 
00106     {
00107       return components[i]->FinalState(trk1, trk2);
00108     }
00109   }
00110 //  G4cerr <<"in charge = "<<IsInCharge(trk1, trk2)<<G4endl;
00111 //  G4cerr <<"Cross-section = "<<CrossSection(trk1, trk2)/millibarn<<" "<<running<<" "<<cxCache.size()<<G4endl;
00112 //  G4cerr <<"Names = "<<trk1.GetDefinition()->GetParticleName()<<", "<<trk2.GetDefinition()->GetParticleName()<<G4endl;
00113 //  throw G4HadronicException(__FILE__, __LINE__, "G4CollisionComposite: no final state found!");
00114   return NULL;
00115 }
00116 
00117 
00118 G4bool G4CollisionComposite::IsInCharge(const G4KineticTrack& trk1, 
00119                                         const G4KineticTrack& trk2) const
00120 {
00121   G4bool isInCharge = false;
00122 
00123   // The composite is in charge if any of its components is in charge
00124 
00125   const G4CollisionVector* comps = GetComponents();
00126   if (comps)
00127     {
00128       G4CollisionVector::const_iterator iter;
00129       for (iter = comps->begin(); iter != comps->end(); ++iter)
00130         {
00131          if ( ((*iter))->IsInCharge(trk1,trk2) ) isInCharge = true;
00132         }
00133     }
00134 
00135   return isInCharge;
00136 }
00137 
00138 void G4CollisionComposite::
00139 BufferCrossSection(const G4ParticleDefinition * aP, const G4ParticleDefinition * bP)
00140 {
00141    // check if already buffered
00142    size_t i;
00143    for(i=0; i<theBuffer.size(); i++)
00144    {
00145      if(theBuffer[i].InCharge(aP, bP)) return;
00146    }
00147 //   G4cerr << "Buffering for "<<aP->GetParticleName()<<" "<<bP->GetParticleName()<<G4endl;
00148    
00149    // buffer the new one.
00150    G4CrossSectionBuffer aNewBuff(aP, bP);
00151    size_t maxE=nPoints;
00152    for(size_t tt=0; tt<maxE; tt++)
00153    {
00154      G4double aT = theT[tt]*GeV;
00155      G4double crossSect = 0;
00156      // The total cross-section is summed over all the component channels
00157      
00158      //A.R. 28-Sep-2012 Fix reproducibility problem
00159      //                 Assign the kinetic energy to the lightest of the
00160      //                 two particles, instead to the first one always.
00161      G4double atime = 0;
00162      G4double btime = 0;
00163      G4ThreeVector aPosition(0,0,0);
00164      G4ThreeVector bPosition(0,0,0);
00165      G4double aM = aP->GetPDGMass();
00166      G4double bM = bP->GetPDGMass();
00167      G4double aE = aM;
00168      G4double bE = bM;
00169      G4ThreeVector aMom(0,0,0);
00170      G4ThreeVector bMom(0,0,0);
00171      if ( aM <= bM ) {
00172       aE += aT;
00173       aMom = G4ThreeVector(0,0,std::sqrt(aE*aE-aM*aM));
00174      } else {
00175       bE += aT;
00176       bMom = G4ThreeVector(0,0,std::sqrt(bE*bE-bM*bM));
00177      }
00178      G4LorentzVector a4Momentum(aE, aMom);
00179      G4LorentzVector b4Momentum(bE, bMom);
00180      G4KineticTrack a(const_cast<G4ParticleDefinition *>(aP), atime, aPosition, a4Momentum);
00181      G4KineticTrack b(const_cast<G4ParticleDefinition *>(bP), btime, bPosition, b4Momentum);
00182      
00183      for (i=0; i<components.size(); i++)
00184      {
00185        if(components[i]->IsInCharge(a,b))
00186        {
00187          crossSect += components[i]->CrossSection(a,b);
00188        }
00189      }
00190      G4double sqrts = (a4Momentum+b4Momentum).mag();
00191      aNewBuff.push_back(sqrts, crossSect);
00192    }
00193    theBuffer.push_back(aNewBuff);
00194 //   theBuffer.back().Print();
00195 }
00196 
00197 
00198 G4double G4CollisionComposite::
00199 BufferedCrossSection(const G4KineticTrack& trk1, const G4KineticTrack& trk2) const
00200 {
00201    for(size_t i=0; i<theBuffer.size(); i++)
00202    {
00203      if(theBuffer[i].InCharge(trk1.GetDefinition(), trk2.GetDefinition())) 
00204      {
00205        return theBuffer[i].CrossSection(trk1, trk2);
00206      }
00207    }
00208    throw G4HadronicException(__FILE__, __LINE__, "G4CollisionComposite::BufferedCrossSection - Blitz !!");
00209    return 0;
00210 }
00211 

Generated on Mon May 27 17:47:54 2013 for Geant4 by  doxygen 1.4.7