Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4Scatterer.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// $Id: G4Scatterer.cc,v 1.16 2010-03-12 15:45:18 gunter Exp $ //
27//
28
29#include <vector>
30
31#include "globals.hh"
33#include "G4SystemOfUnits.hh"
34#include "G4ios.hh"
35#include "G4Scatterer.hh"
36#include "G4KineticTrack.hh"
37#include "G4ThreeVector.hh"
38#include "G4LorentzRotation.hh"
39#include "G4LorentzVector.hh"
40
41#include "G4CollisionNN.hh"
42#include "G4CollisionPN.hh"
44
46#include "G4HadTmpUtil.hh"
47#include "G4Pair.hh"
48
49// Declare the categories of collisions the Scatterer can handle
51
52
54{
55 Register aR;
56 G4ForEach<theChannels>::Apply(&aR, &collisions);
57}
58
59
61{
62 std::for_each(collisions.begin(), collisions.end(), G4Delete());
63 collisions.clear();
64}
65
67 const G4KineticTrack& trk2)
68{
69 G4double time = DBL_MAX;
70 G4double distance_fast;
72// G4cout << "zcomp=" << std::abs(mom1.vect().unit().z() -1 ) << G4endl;
73 G4double collisionTime;
74
75 if ( std::abs(mom1.vect().unit().z() -1 ) < 1e-6 )
76 {
78 G4double deltaz=position.z();
79 G4double velocity = mom1.z()/mom1.e() * c_light;
80
81 collisionTime=deltaz/velocity;
82 distance_fast=position.x()*position.x() + position.y()*position.y();
83 } else {
84
85 // The nucleons of the nucleus are FROZEN, ie. do not move..
86
88
89 G4ThreeVector velocity = mom1.vect()/mom1.e() * c_light; // mom1.boostVector() will exit on slightly negative mass
90 collisionTime = (position * velocity) / velocity.mag2(); // can't divide by /c_light;
91 position -= velocity * collisionTime;
92 distance_fast=position.mag2();
93
94// if ( collisionTime>0 ) G4cout << " dis1/2 square" << dis1 <<" "<< dis2 << G4endl;
95// collisionTime = GetTimeToClosestApproach(trk1,trk2);
96 }
97 if (collisionTime > 0)
98 {
99 static const G4double maxCrossSection = 500*millibarn;
100 if(0.7*pi*distance_fast>maxCrossSection) return time;
101
102
103 G4LorentzVector mom2(0,0,0,trk2.Get4Momentum().mag());
104
105// G4ThreeVector momLab = mom1.vect();// frozen Nucleus - mom2.vect();
106// G4ThreeVector posLab = trk1.GetPosition() - trk2.GetPosition();
107// G4double disLab=posLab * posLab - (posLab*momLab) * (posLab*momLab) /(momLab.mag2());
108
109 G4LorentzRotation toCMSFrame((-1)*(mom1 + mom2).boostVector());
110 mom1 = toCMSFrame * mom1;
111 mom2 = toCMSFrame * mom2;
112
113 G4LorentzVector coordinate1(trk1.GetPosition(), 100.);
114 G4LorentzVector coordinate2(trk2.GetPosition(), 100.);
115 G4ThreeVector pos = ((toCMSFrame * coordinate1).vect() -
116 (toCMSFrame * coordinate2).vect());
117
118 G4ThreeVector mom = mom1.vect() - mom2.vect();
119
120 // Calculate the impact parameter
121
122 G4double distance = pos * pos - (pos*mom) * (pos*mom) / (mom.mag2());
123
124// G4cout << " disDiff " << distance-disLab << " " << disLab
125// << " " << std::abs(distance-disLab)/distance << G4endl
126// << " mom/Lab " << mom << " " << momLab << G4endl
127// << " pos/Lab " << pos << " " << posLab
128// << G4endl;
129
130 if(pi*distance>maxCrossSection) return time;
131
132 // charged particles special
133 static const G4double maxChargedCrossSection = 200*millibarn;
134 if(std::abs(trk1.GetDefinition()->GetPDGCharge())>0.1 &&
135 std::abs(trk2.GetDefinition()->GetPDGCharge())>0.1 &&
136 pi*distance>maxChargedCrossSection) return time;
137
138 G4double sqrtS = (trk1.Get4Momentum() + trk2.Get4Momentum()).mag();
139 // neutrons special
140 if(( trk1.GetDefinition() == G4Neutron::Neutron() ||
141 trk1.GetDefinition() == G4Neutron::Neutron() ) &&
142 sqrtS>1.91*GeV && pi*distance>maxChargedCrossSection) return time;
143
144/*
145 * if(distance <= sqr(1.14*fermi))
146 * {
147 * time = collisionTime;
148 *
149 * *
150 * * G4cout << "Scatter distance/time: " << std::sqrt(distance)/fermi <<
151 * * " / "<< time/ns << G4endl;
152 * * G4ThreeVector pos1=trk1.GetPosition();
153 * * G4ThreeVector pos2=trk2.GetPosition();
154 * * G4LorentzVector xmom1 = trk1.Get4Momentum();
155 * * G4LorentzVector xmom2 = trk2.Get4Momentum();
156 * * G4cout << "position1: " << pos1.x() << " " << pos1.y() << " "
157 * * << pos1.z();
158 * * pos1+=(collisionTime*c_light/xmom1.e())*xmom1.vect();
159 * * G4cout << " straight line trprt: "
160 * * << pos1.x() << " " << pos1.y() << " "
161 * * << pos1.z() << G4endl;
162 * * G4cout << "position2: " << pos2.x() << " " << pos2.y() << " "
163 * * << pos2.z() << G4endl;
164 * * G4cout << "straight line distance 2 fixed:" << (pos1-pos2).mag()/fermi << G4endl;
165 * * pos2+= (collisionTime*c_light/xmom2.e())*xmom2.vect();
166 * * G4cout<< " straight line trprt: "
167 * * << pos2.x() << " " << pos2.y() << " "
168 * * << pos2.z() << G4endl;
169 * * G4cout << "straight line distance :" << (pos1-pos2).mag()/fermi << G4endl;
170 * *
171 * }
172 *
173 * if(1)
174 * return time;
175 */
176
177 if ((trk1.GetActualMass()+trk2.GetActualMass()) > sqrtS) return time;
178
179
180
181 G4VCollision* collision = FindCollision(trk1,trk2);
182 G4double totalCrossSection;
183 // The cross section is interpreted geometrically as an area
184 // Two particles are assumed to collide if their distance is < (totalCrossSection/pi)
185
186 if (collision != 0)
187 {
188 totalCrossSection = collision->CrossSection(trk1,trk2);
189 if ( totalCrossSection > 0 )
190 {
191/* G4cout << " totalCrossection = "<< totalCrossSection << ", trk1/2, s, e-m: "
192 * << trk1.GetDefinition()->GetParticleName()
193 * << " / "
194 * << trk2.GetDefinition()->GetParticleName()
195 * << ", "
196 * << (trk1.Get4Momentum()+trk2.Get4Momentum()).mag()
197 * << ", "
198 * << (trk1.Get4Momentum()+trk2.Get4Momentum()).mag()-
199 * trk1.Get4Momentum().mag() - trk2.Get4Momentum().mag()
200 * << G4endl;
201 */
202 if (distance <= totalCrossSection / pi)
203 {
204 time = collisionTime;
205 }
206 } else
207 {
208
209 // For debugging...
210 // G4cout << " totalCrossection = 0, trk1/2, s, e-m: "
211 // << trk1.GetDefinition()->GetParticleName()
212 // << " / "
213 // << trk2.GetDefinition()->GetParticleName()
214 // << ", "
215 // << (trk1.Get4Momentum()+trk2.Get4Momentum()).mag()
216 // << ", "
217 // << (trk1.Get4Momentum()+trk2.Get4Momentum()).mag()-
218 // trk1.Get4Momentum().mag() - trk2.Get4Momentum().mag()
219 // << G4endl;
220
221 }
222/*
223 * if(distance <= sqr(5.*fermi))
224 * {
225 * G4cout << " distance,xsect, std::sqrt(xsect/pi) : " << std::sqrt(distance)/fermi
226 * << " " << totalCrossSection/sqr(fermi)
227 * << " " << std::sqrt(totalCrossSection / pi)/fermi << G4endl;
228 * }
229 */
230
231 }
232 else
233 {
234 time = DBL_MAX;
235// /*
236 // For debugging
237//hpw G4cout << "G4Scatterer - collision not found: "
238//hpw << trk1.GetDefinition()->GetParticleName()
239//hpw << " - "
240//hpw << trk2.GetDefinition()->GetParticleName()
241//hpw << G4endl;
242 // End of debugging
243// */
244 }
245 }
246
247 else
248 {
249 /*
250 // For debugging
251 G4cout << "G4Scatterer - negative collisionTime"
252 << ": collisionTime = " << collisionTime
253 << ", position = " << position
254 << ", velocity = " << velocity
255 << G4endl;
256 // End of debugging
257 */
258 }
259
260 return time;
261}
262
264 const G4KineticTrack& trk2)
265{
266// G4double sqrtS = (trk1.Get4Momentum() + trk2.Get4Momentum()).mag();
267 G4LorentzVector pInitial=trk1.Get4Momentum() + trk2.Get4Momentum();
268 G4double energyBalance = pInitial.t();
269 G4double pxBalance = pInitial.vect().x();
270 G4double pyBalance = pInitial.vect().y();
271 G4double pzBalance = pInitial.vect().z();
272 G4int chargeBalance = G4lrint(trk1.GetDefinition()->GetPDGCharge()
273 + trk2.GetDefinition()->GetPDGCharge());
274 G4int baryonBalance = trk1.GetDefinition()->GetBaryonNumber()
275 + trk2.GetDefinition()->GetBaryonNumber();
276
277 G4VCollision* collision = FindCollision(trk1,trk2);
278 if (collision != 0)
279 {
280 G4double aCrossSection = collision->CrossSection(trk1,trk2);
281 if (aCrossSection > 0.0)
282 {
283
284
285 #ifdef debug_G4Scatterer
286 G4cout << "be4 FinalState 1(p,e,m): "
287 << trk1.Get4Momentum() << " "
288 << trk1.Get4Momentum().mag()
289 << ", 2: "
290 << trk2.Get4Momentum()<< " "
291 << trk2.Get4Momentum().mag() << " "
292 << G4endl;
293 #endif
294
295
296 G4KineticTrackVector* products = collision->FinalState(trk1,trk2);
297 if(!products || products->size() == 0) return products;
298
299 #ifdef debug_G4Scatterer
300 G4cout << "size of FS: "<<products->size()<<G4endl;
301 #endif
302
303 G4KineticTrack *final= products->operator[](0);
304
305
306 #ifdef debug_G4Scatterer
307 G4cout << " FinalState 1: "
308 << final->Get4Momentum()<< " "
309 << final->Get4Momentum().mag() ;
310 #endif
311
312 if(products->size() == 1) return products;
313 final=products->operator[](1);
314 #ifdef debug_G4Scatterer
315 G4cout << ", 2: "
316 << final->Get4Momentum() << " "
317 << final->Get4Momentum().mag() << " " << G4endl;
318 #endif
319
320 final= products->operator[](0);
321 G4LorentzVector pFinal=final->Get4Momentum();
322 if(products->size()==2)
323 {
324 final=products->operator[](1);
325 pFinal +=final->Get4Momentum();
326 }
327
328 #ifdef debug_G4Scatterer
329 if ( (pInitial-pFinal).mag() > 0.1*MeV )
330 {
331 G4cout << "G4Scatterer: momentum imbalance, pInitial= " <<pInitial << " pFinal= " <<pFinal<< G4endl;
332 }
333 G4cout << "Scatterer costh= " << trk1.Get4Momentum().vect().unit() *(products->operator[](0))->Get4Momentum().vect().unit()<< G4endl;
334 #endif
335
336 for(size_t hpw=0; hpw<products->size(); hpw++)
337 {
338 energyBalance-=products->operator[](hpw)->Get4Momentum().t();
339 pxBalance-=products->operator[](hpw)->Get4Momentum().vect().x();
340 pyBalance-=products->operator[](hpw)->Get4Momentum().vect().y();
341 pzBalance-=products->operator[](hpw)->Get4Momentum().vect().z();
342 chargeBalance-=G4lrint(products->operator[](hpw)->GetDefinition()->GetPDGCharge());
343 baryonBalance-=products->operator[](hpw)->GetDefinition()->GetBaryonNumber();
344 }
345 if(getenv("ScattererEnergyBalanceCheck"))
346 std::cout << "DEBUGGING energy balance A: "
347 <<energyBalance<<" "
348 <<pxBalance<<" "
349 <<pyBalance<<" "
350 <<pzBalance<<" "
351 <<chargeBalance<<" "
352 <<baryonBalance<<" "
353 <<G4endl;
354 if(chargeBalance !=0 )
355 {
356 G4cout << "track 1"<<trk1.GetDefinition()->GetParticleName()<<G4endl;
357 G4cout << "track 2"<<trk2.GetDefinition()->GetParticleName()<<G4endl;
358 for(size_t hpw=0; hpw<products->size(); hpw++)
359 {
360 G4cout << products->operator[](hpw)->GetDefinition()->GetParticleName()<<G4endl;
361 }
362 G4Exception("G4Scatterer", "im_r_matrix001", FatalException,
363 "Problem in ChargeBalance");
364 }
365 return products;
366 }
367 }
368
369 return NULL;
370}
371
372
373G4VCollision* G4Scatterer::FindCollision(const G4KineticTrack& trk1,
374 const G4KineticTrack& trk2)
375{
376 G4VCollision* collisionInCharge = 0;
377
378 size_t i;
379 for (i=0; i<collisions.size(); i++)
380 {
381 G4VCollision* component = collisions[i];
382 if (component->IsInCharge(trk1,trk2))
383 {
384 collisionInCharge = component;
385 break;
386 }
387 }
388// if(collisionInCharge)
389// {
390// G4cout << "found collision : "
391// << collisionInCharge->GetName()<< " "
392// << "for "
393// << trk1.GetDefinition()->GetParticleName()<<" + "
394// << trk2.GetDefinition()->GetParticleName()<<" "
395// << G4endl;;
396// }
397 return collisionInCharge;
398}
399
401 const G4KineticTrack& trk2)
402{
403 G4VCollision* collision = FindCollision(trk1,trk2);
404 G4double aCrossSection = 0;
405 if (collision != 0)
406 {
407 aCrossSection = collision->CrossSection(trk1,trk2);
408 }
409 return aCrossSection;
410}
411
412
413const std::vector<G4CollisionInitialState *> & G4Scatterer::
414GetCollisions(G4KineticTrack * aProjectile,
415 std::vector<G4KineticTrack *> & someCandidates,
416 G4double aCurrentTime)
417{
418 theCollisions.clear();
419 std::vector<G4KineticTrack *>::iterator j=someCandidates.begin();
420 for(; j != someCandidates.end(); ++j)
421 {
422 G4double collisionTime = GetTimeToInteraction(*aProjectile, **j);
423 if(collisionTime == DBL_MAX) // no collision
424 {
425 continue;
426 }
427 G4KineticTrackVector aTarget;
428 aTarget.push_back(*j);
429 theCollisions.push_back(
430 new G4CollisionInitialState(collisionTime+aCurrentTime, aProjectile, aTarget, this) );
431// G4cerr <<" !!!!!! debug collisions "<<collisionTime<<" "<<pkt->GetDefinition()->GetParticleName()<<G4endl;
432 }
433 return theCollisions;
434}
435
436
438GetFinalState(G4KineticTrack * aProjectile,
439 std::vector<G4KineticTrack *> & theTargets)
440{
441 G4KineticTrack target_reloc(*(theTargets[0]));
442 return Scatter(*aProjectile, target_reloc);
443}
@ FatalException
#define GROUP2(t1, t2)
Definition: G4Pair.hh:42
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
double z() const
Hep3Vector unit() const
double x() const
double mag2() const
double y() const
Hep3Vector vect() const
static void Apply()
Definition: G4Pair.hh:179
const G4ThreeVector & GetPosition() const
G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & GetTrackingMomentum() const
const G4LorentzVector & Get4Momentum() const
G4double GetActualMass() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4double GetPDGCharge() const
const G4String & GetParticleName() const
virtual G4KineticTrackVector * Scatter(const G4KineticTrack &trk1, const G4KineticTrack &trk2)
Definition: G4Scatterer.cc:263
virtual const std::vector< G4CollisionInitialState * > & GetCollisions(G4KineticTrack *aProjectile, std::vector< G4KineticTrack * > &someCandidates, G4double aCurrentTime)
Definition: G4Scatterer.cc:414
virtual G4double GetTimeToInteraction(const G4KineticTrack &trk1, const G4KineticTrack &trk2)
Definition: G4Scatterer.cc:66
virtual G4KineticTrackVector * GetFinalState(G4KineticTrack *aProjectile, std::vector< G4KineticTrack * > &theTargets)
Definition: G4Scatterer.cc:438
G4double GetCrossSection(const G4KineticTrack &trk1, const G4KineticTrack &trk2)
Definition: G4Scatterer.cc:400
virtual ~G4Scatterer()
Definition: G4Scatterer.cc:60
virtual G4double CrossSection(const G4KineticTrack &trk1, const G4KineticTrack &trk2) const
Definition: G4VCollision.cc:55
virtual G4KineticTrackVector * FinalState(const G4KineticTrack &trk1, const G4KineticTrack &trk2) const =0
virtual G4bool IsInCharge(const G4KineticTrack &trk1, const G4KineticTrack &trk2) const =0
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
int G4lrint(double ad)
Definition: templates.hh:163
#define DBL_MAX
Definition: templates.hh:83