source: trunk/source/processes/hadronic/models/im_r_matrix/src/G4CollisionComposite.cc@ 1201

Last change on this file since 1201 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 6.7 KB
Line 
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.7 2007/07/05 21:43:24 dennis Exp $ //
28
29#include "globals.hh"
30#include "G4CollisionComposite.hh"
31#include "G4VCollision.hh"
32#include "G4CollisionVector.hh"
33#include "G4KineticTrack.hh"
34#include "G4KineticTrackVector.hh"
35#include "G4VCrossSectionSource.hh"
36#include "G4HadTmpUtil.hh"
37
38const G4int G4CollisionComposite::nPoints = 32;
39
40G4double G4CollisionComposite::theT[nPoints] =
41{.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};
42
43G4CollisionComposite::G4CollisionComposite()
44{
45}
46
47
48G4CollisionComposite::~G4CollisionComposite()
49{
50 std::for_each(components.begin(), components.end(), G4Delete());
51}
52
53
54G4double G4CollisionComposite::CrossSection(const G4KineticTrack& trk1,
55 const G4KineticTrack& trk2) const
56{
57 G4double crossSect = 0.;
58 const G4VCrossSectionSource* xSource = GetCrossSectionSource();
59 if (xSource != 0)
60 // There is a total cross section for this Collision
61 {
62 crossSect = xSource->CrossSection(trk1,trk2);
63 }
64 else
65 {
66 // waiting for mutable to enable buffering.
67 const_cast<G4CollisionComposite *>(this)->BufferCrossSection(trk1.GetDefinition(), trk2.GetDefinition());
68// G4cerr << "Buffer filled, reying with sqrts = "<< (trk1.Get4Momentum()+trk2.Get4Momentum()).mag() <<G4endl;
69 crossSect = BufferedCrossSection(trk1,trk2);
70 }
71 return crossSect;
72}
73
74
75G4KineticTrackVector* G4CollisionComposite::FinalState(const G4KineticTrack& trk1,
76 const G4KineticTrack& trk2) const
77{
78 std::vector<G4double> cxCache;
79 G4double partialCxSum = 0.0;
80
81 size_t i;
82 for (i=0; i<components.size(); i++)
83 {
84 G4double partialCx;
85// cout << "comp" << i << " " << components[i]()->GetName();
86 if (components[i]->IsInCharge(trk1,trk2))
87 {
88 partialCx = components[i]->CrossSection(trk1,trk2);
89 }
90 else
91 {
92 partialCx = 0.0;
93 }
94// cout << " cx=" << partialCx << endl;
95 partialCxSum += partialCx;
96 cxCache.push_back(partialCx);
97 }
98
99 G4double random = G4UniformRand()*partialCxSum;
100 G4double running = 0;
101 for (i=0; i<cxCache.size(); i++)
102 {
103 running += cxCache[i];
104 if (running > random)
105 {
106 return components[i]->FinalState(trk1, trk2);
107 }
108 }
109// G4cerr <<"in charge = "<<IsInCharge(trk1, trk2)<<G4endl;
110// G4cerr <<"Cross-section = "<<CrossSection(trk1, trk2)/millibarn<<" "<<running<<" "<<cxCache.size()<<G4endl;
111// G4cerr <<"Names = "<<trk1.GetDefinition()->GetParticleName()<<", "<<trk2.GetDefinition()->GetParticleName()<<G4endl;
112// throw G4HadronicException(__FILE__, __LINE__, "G4CollisionComposite: no final state found!");
113 return NULL;
114}
115
116
117G4bool G4CollisionComposite::IsInCharge(const G4KineticTrack& trk1,
118 const G4KineticTrack& trk2) const
119{
120 G4bool isInCharge = false;
121
122 // The composite is in charge if any of its components is in charge
123
124 const G4CollisionVector* comps = GetComponents();
125 if (comps)
126 {
127 G4CollisionVector::const_iterator iter;
128 for (iter = comps->begin(); iter != comps->end(); ++iter)
129 {
130 if ( ((*iter))->IsInCharge(trk1,trk2) ) isInCharge = true;
131 }
132 }
133
134 return isInCharge;
135}
136
137void G4CollisionComposite::
138BufferCrossSection(const G4ParticleDefinition * aP, const G4ParticleDefinition * bP)
139{
140 // check if already buffered
141 size_t i;
142 for(i=0; i<theBuffer.size(); i++)
143 {
144 if(theBuffer[i].InCharge(aP, bP)) return;
145 }
146// G4cerr << "Buffering for "<<aP->GetParticleName()<<" "<<bP->GetParticleName()<<G4endl;
147
148 // buffer the new one.
149 G4CrossSectionBuffer aNewBuff(aP, bP);
150 size_t maxE=nPoints;
151 for(size_t tt=0; tt<maxE; tt++)
152 {
153 G4double aT = theT[tt]*GeV;
154 G4double crossSect = 0;
155 // The total cross-section is summed over all the component channels
156
157 G4double atime = 0;
158 G4ThreeVector aPosition(0,0,0);
159 G4double aM = aP->GetPDGMass();
160 G4double aE = aM+aT;
161 G4ThreeVector aMom(0,0,std::sqrt(aE*aE-aM*aM));
162 G4LorentzVector a4Momentum(aE, aMom);
163 G4KineticTrack a(const_cast<G4ParticleDefinition *>(aP), atime, aPosition, a4Momentum);
164
165 G4double btime = 0;
166 G4ThreeVector bPosition(0,0,0);
167 G4ThreeVector bMom(0,0,0*MeV);
168 G4double bE = bP->GetPDGMass();
169 G4LorentzVector b4Momentum(bE, bMom);
170 G4KineticTrack b(const_cast<G4ParticleDefinition *>(bP), btime, bPosition, b4Momentum);
171
172 for (i=0; i<components.size(); i++)
173 {
174 if(components[i]->IsInCharge(a,b))
175 {
176 crossSect += components[i]->CrossSection(a,b);
177 }
178 }
179 G4double sqrts = (a4Momentum+b4Momentum).mag();
180 aNewBuff.push_back(sqrts, crossSect);
181 }
182 theBuffer.push_back(aNewBuff);
183// theBuffer.back().Print();
184}
185
186
187G4double G4CollisionComposite::
188BufferedCrossSection(const G4KineticTrack& trk1, const G4KineticTrack& trk2) const
189{
190 for(size_t i=0; i<theBuffer.size(); i++)
191 {
192 if(theBuffer[i].InCharge(trk1.GetDefinition(), trk2.GetDefinition()))
193 {
194 return theBuffer[i].CrossSection(trk1, trk2);
195 }
196 }
197 throw G4HadronicException(__FILE__, __LINE__, "G4CollisionComposite::BufferedCrossSection - Blitz !!");
198 return 0;
199}
200
Note: See TracBrowser for help on using the repository browser.