Geant4  10.01.p02
G4CollisionComposite.cc
Go to the documentation of this file.
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 //
27 // $Id: G4CollisionComposite.cc,v 1.9 2010-03-12 15:45:18 gunter Exp $ //
28 
29 #include "globals.hh"
30 #include "G4SystemOfUnits.hh"
31 #include "G4CollisionComposite.hh"
32 #include "G4VCollision.hh"
33 #include "G4CollisionVector.hh"
34 #include "G4KineticTrack.hh"
35 #include "G4KineticTrackVector.hh"
36 #include "G4VCrossSectionSource.hh"
37 #include "G4HadTmpUtil.hh"
38 #include "G4AutoLock.hh"
39 
41 
42 const G4double G4CollisionComposite::theT[nPoints] =
43 {.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};
44 
46 {
48 }
49 
50 
52 {
54  std::for_each(components.begin(), components.end(), G4Delete());
55 }
56 
57 
59  const G4KineticTrack& trk2) const
60 {
61  G4double crossSect = 0.;
63  if (xSource != 0)
64  // There is a total cross section for this Collision
65  {
66  crossSect = xSource->CrossSection(trk1,trk2);
67  }
68  else
69  {
71  // waiting for mutable to enable buffering.
72  const_cast<G4CollisionComposite *>(this)->BufferCrossSection(trk1.GetDefinition(), trk2.GetDefinition());
73 // G4cerr << "Buffer filled, reying with sqrts = "<< (trk1.Get4Momentum()+trk2.Get4Momentum()).mag() <<G4endl;
74  crossSect = BufferedCrossSection(trk1,trk2);
75  }
76  return crossSect;
77 }
78 
79 
81  const G4KineticTrack& trk2) const
82 {
83  std::vector<G4double> cxCache;
84  G4double partialCxSum = 0.0;
85 
86  size_t i;
87  for (i=0; i<components.size(); i++)
88  {
89  G4double partialCx;
90 // cout << "comp" << i << " " << components[i]()->GetName();
91  if (components[i]->IsInCharge(trk1,trk2))
92  {
93  partialCx = components[i]->CrossSection(trk1,trk2);
94  }
95  else
96  {
97  partialCx = 0.0;
98  }
99 // cout << " cx=" << partialCx << endl;
100  partialCxSum += partialCx;
101  cxCache.push_back(partialCx);
102  }
103 
104  G4double random = G4UniformRand()*partialCxSum;
105  G4double running = 0;
106  for (i=0; i<cxCache.size(); i++)
107  {
108  running += cxCache[i];
109  if (running > random)
110  {
111  return components[i]->FinalState(trk1, trk2);
112  }
113  }
114 // G4cerr <<"in charge = "<<IsInCharge(trk1, trk2)<<G4endl;
115 // G4cerr <<"Cross-section = "<<CrossSection(trk1, trk2)/millibarn<<" "<<running<<" "<<cxCache.size()<<G4endl;
116 // G4cerr <<"Names = "<<trk1.GetDefinition()->GetParticleName()<<", "<<trk2.GetDefinition()->GetParticleName()<<G4endl;
117 // throw G4HadronicException(__FILE__, __LINE__, "G4CollisionComposite: no final state found!");
118  return NULL;
119 }
120 
121 
123  const G4KineticTrack& trk2) const
124 {
125  G4bool isInCharge = false;
126 
127  // The composite is in charge if any of its components is in charge
128 
129  const G4CollisionVector* comps = GetComponents();
130  if (comps)
131  {
132  G4CollisionVector::const_iterator iter;
133  for (iter = comps->begin(); iter != comps->end(); ++iter)
134  {
135  if ( ((*iter))->IsInCharge(trk1,trk2) ) isInCharge = true;
136  }
137  }
138 
139  return isInCharge;
140 }
141 
144 {
145  // check if already buffered
146  size_t i;
147  for(i=0; i<theBuffer.size(); i++)
148  {
149  if(theBuffer[i].InCharge(aP, bP)) return;
150  }
151 // G4cerr << "Buffering for "<<aP->GetParticleName()<<" "<<bP->GetParticleName()<<G4endl;
152 
153  // buffer the new one.
154  G4CrossSectionBuffer aNewBuff(aP, bP);
155  size_t maxE=nPoints;
156  for(size_t tt=0; tt<maxE; tt++)
157  {
158  G4double aT = theT[tt]*GeV;
159  G4double crossSect = 0;
160  // The total cross-section is summed over all the component channels
161 
162  //A.R. 28-Sep-2012 Fix reproducibility problem
163  // Assign the kinetic energy to the lightest of the
164  // two particles, instead to the first one always.
165  G4double atime = 0;
166  G4double btime = 0;
167  G4ThreeVector aPosition(0,0,0);
168  G4ThreeVector bPosition(0,0,0);
169  G4double aM = aP->GetPDGMass();
170  G4double bM = bP->GetPDGMass();
171  G4double aE = aM;
172  G4double bE = bM;
173  G4ThreeVector aMom(0,0,0);
174  G4ThreeVector bMom(0,0,0);
175  if ( aM <= bM ) {
176  aE += aT;
177  aMom = G4ThreeVector(0,0,std::sqrt(aE*aE-aM*aM));
178  } else {
179  bE += aT;
180  bMom = G4ThreeVector(0,0,std::sqrt(bE*bE-bM*bM));
181  }
182  G4LorentzVector a4Momentum(aE, aMom);
183  G4LorentzVector b4Momentum(bE, bMom);
184  G4KineticTrack a(aP, atime, aPosition, a4Momentum);
185  G4KineticTrack b(bP, btime, bPosition, b4Momentum);
186 
187  for (i=0; i<components.size(); i++)
188  {
189  if(components[i]->IsInCharge(a,b))
190  {
191  crossSect += components[i]->CrossSection(a,b);
192  }
193  }
194  G4double sqrts = (a4Momentum+b4Momentum).mag();
195  aNewBuff.push_back(sqrts, crossSect);
196  }
197  theBuffer.push_back(aNewBuff);
198 // theBuffer.back().Print();
199 }
200 
201 
203 BufferedCrossSection(const G4KineticTrack& trk1, const G4KineticTrack& trk2) const
204 {
205  for(size_t i=0; i<theBuffer.size(); i++)
206  {
207  if(theBuffer[i].InCharge(trk1.GetDefinition(), trk2.GetDefinition()))
208  {
209  return theBuffer[i].CrossSection(trk1, trk2);
210  }
211  }
212  throw G4HadronicException(__FILE__, __LINE__, "G4CollisionComposite::BufferedCrossSection - Blitz !!");
213  return 0;
214 }
215 
virtual const G4CollisionVector * GetComponents() const
std::vector< G4CrossSectionBuffer > theBuffer
CLHEP::Hep3Vector G4ThreeVector
virtual const G4VCrossSectionSource * GetCrossSectionSource() const
static const G4int nPoints
G4CollisionVector components
#define G4MUTEXINIT(mutex)
Definition: G4Threading.hh:177
G4double BufferedCrossSection(const G4KineticTrack &trk1, const G4KineticTrack &trk2) const
G4double a
Definition: TRTMaterials.hh:39
virtual G4double CrossSection(const G4KineticTrack &trk1, const G4KineticTrack &trk2) const
int G4int
Definition: G4Types.hh:78
void BufferCrossSection(const G4ParticleDefinition *aP, const G4ParticleDefinition *bP)
#define G4UniformRand()
Definition: Randomize.hh:93
bool G4bool
Definition: G4Types.hh:79
std::vector< G4VCollision * > G4CollisionVector
static const double GeV
Definition: G4SIunits.hh:196
virtual G4KineticTrackVector * FinalState(const G4KineticTrack &trk1, const G4KineticTrack &trk2) const
G4double GetPDGMass() const
virtual G4bool IsInCharge(const G4KineticTrack &trk1, const G4KineticTrack &trk2) const
void push_back(G4double S, G4double x)
double G4double
Definition: G4Types.hh:76
static const G4double theT[]
const G4ParticleDefinition * GetDefinition() const
#define G4MUTEXDESTROY(mutex)
Definition: G4Threading.hh:178
virtual G4double CrossSection(const G4KineticTrack &trk1, const G4KineticTrack &trk2) const =0
CLHEP::HepLorentzVector G4LorentzVector