Geant4  10.02
G4QMDCollision.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 // 080602 Fix memory leaks by T. Koi
27 // 081120 Add deltaT in signature of CalKinematicsOfBinaryCollisions
28 // Add several required updating of Mean Filed
29 // Modified handling of absorption case by T. Koi
30 // 090126 Fix in absorption case by T. Koi
31 // 090331 Fix for gamma participant by T. Koi
32 //
33 #include "G4QMDCollision.hh"
34 #include "G4Scatterer.hh"
35 #include "G4Pow.hh"
36 #include "G4Exp.hh"
37 #include "G4Log.hh"
38 #include "G4PhysicalConstants.hh"
39 #include "G4SystemOfUnits.hh"
40 #include "Randomize.hh"
41 
43 : fdeltar ( 4.0 )
44 , fbcmax0 ( 1.323142 ) // NN maximum impact parameter
45 , fbcmax1 ( 2.523 ) // others maximum impact parameter
46 // , sig0 ( 55 ) // NN cross section
47 //110617 fix for gcc 4.6 compilation warnings
48 //, sig1 ( 200 ) // others cross section
49 , fepse ( 0.0001 )
50 {
51  //These two pointers will be set through SetMeanField method
52  theSystem=NULL;
53  theMeanField=NULL;
54  theScatterer = new G4Scatterer();
55 }
56 
57 /*
58 G4QMDCollision::G4QMDCollision( const G4QMDCollision& obj )
59 : fdeltar ( obj.fdeltar )
60 , fbcmax0 ( obj.fbcmax0 ) // NN maximum impact parameter
61 , fbcmax1 ( obj.fbcmax1 ) // others maximum impact parameter
62 , fepse ( obj.fepse )
63 {
64 
65  if ( obj.theSystem != NULL ) {
66  theSystem = new G4QMDSystem;
67  *theSystem = *obj.theSystem;
68  } else {
69  theSystem = NULL;
70  }
71  if ( obj.theMeanField != NULL ) {
72  theMeanField = new G4QMDMeanField;
73  *theMeanField = *obj.theMeanField;
74  } else {
75  theMeanField = NULL;
76  }
77  theScatterer = new G4Scatterer();
78  *theScatterer = *obj.theScatterer;
79 }
80 
81 G4QMDCollision & G4QMDCollision::operator= ( const G4QMDCollision& obj)
82 {
83  fdeltar = obj.fdeltar;
84  fbcmax0 = obj.fbcmax1;
85  fepse = obj.fepse;
86 
87  if ( obj.theSystem != NULL ) {
88  delete theSystem;
89  theSystem = new G4QMDSystem;
90  *theSystem = *obj.theSystem;
91  } else {
92  theSystem = NULL;
93  }
94  if ( obj.theMeanField != NULL ) {
95  delete theMeanField;
96  theMeanField = new G4QMDMeanField;
97  *theMeanField = *obj.theMeanField;
98  } else {
99  theMeanField = NULL;
100  }
101  delete theScatterer;
102  theScatterer = new G4Scatterer();
103  *theScatterer = *obj.theScatterer;
104 
105  return *this;
106 }
107 */
108 
109 
111 {
112  //if ( theSystem != NULL ) delete theSystem;
113  //if ( theMeanField != NULL ) delete theMeanField;
114  delete theScatterer;
115 }
116 
117 
119 {
120  G4double deltaT = dt;
121 
123 //081118
124  //G4int nb = 0;
125  for ( G4int i = 0 ; i < n ; i++ )
126  {
129  //nb += theSystem->GetParticipant( i )->GetBaryonNumber();
130  }
131  //G4cout << "nb = " << nb << " n = " << n << G4endl;
132 
133 
134 //071101
135  for ( G4int i = 0 ; i < n ; i++ )
136  {
137 
138  //std::cout << i << " " << theSystem->GetParticipant( i )->GetDefinition()->GetParticleName() << " " << theSystem->GetParticipant( i )->GetPosition() << std::endl;
139 
141  {
142 
143  G4bool decayed = false;
144 
148 
150 
151  G4double eini = theMeanField->GetTotalPotential() + p40.e();
152 
154  G4int i0 = 0;
155 
156  G4bool isThisEnergyOK = false;
157 
158  G4int maximumNumberOfTrial=4;
159  for ( G4int ii = 0 ; ii < maximumNumberOfTrial ; ii++ )
160  {
161 
162  //G4LorentzVector p4 = theSystem->GetParticipant( i )->Get4Momentum();
163  G4LorentzVector p400 = p40;
164 
165  p400 *= GeV;
166  //G4KineticTrack kt( theSystem->GetParticipant( i )->GetDefinition() , 0.0 , (theSystem->GetParticipant( i )->GetPosition())*fermi , p4 );
167  G4KineticTrack kt( pd0 , 0.0 , r0*fermi , p400 );
168  //std::cout << "G4KineticTrack " << i << " " << kt.GetDefinition()->GetParticleName() << kt.GetPosition() << std::endl;
169  G4KineticTrackVector* secs = NULL;
170  secs = kt.Decay();
171  G4int id = 0;
172  G4double et = 0;
173  if ( secs )
174  {
175  for ( G4KineticTrackVector::iterator it
176  = secs->begin() ; it != secs->end() ; it++ )
177  {
178 /*
179  G4cout << "G4KineticTrack"
180  << " " << (*it)->GetDefinition()->GetParticleName()
181  << " " << (*it)->Get4Momentum()
182  << " " << (*it)->GetPosition()/fermi
183  << G4endl;
184 */
185  if ( id == 0 )
186  {
187  theSystem->GetParticipant( i )->SetDefinition( (*it)->GetDefinition() );
188  theSystem->GetParticipant( i )->SetMomentum( (*it)->Get4Momentum().v()/GeV );
189  theSystem->GetParticipant( i )->SetPosition( (*it)->GetPosition()/fermi );
190  //theMeanField->Cal2BodyQuantities( i );
191  et += (*it)->Get4Momentum().e()/GeV;
192  }
193  if ( id > 0 )
194  {
195  // Append end;
196  theSystem->SetParticipant ( new G4QMDParticipant ( (*it)->GetDefinition() , (*it)->Get4Momentum().v()/GeV , (*it)->GetPosition()/fermi ) );
197  et += (*it)->Get4Momentum().e()/GeV;
198  if ( id > 1 )
199  {
200  //081118
201  //G4cout << "G4QMDCollision id >2; id= " << id << G4endl;
202  }
203  }
204  id++; // number of daughter particles
205 
206  delete *it;
207  }
208 
209  theMeanField->Update();
210  i0 = id-1; // 0 enter to i
211 
212  delete secs;
213  }
214 
215 // EnergyCheck
216 
217  G4double efin = theMeanField->GetTotalPotential() + et;
218  //std::cout << std::abs ( eini - efin ) - fepse << std::endl;
219 // std::cout << std::abs ( eini - efin ) - fepse*10 << std::endl;
220 // *10 TK
221  if ( std::abs ( eini - efin ) < fepse*10 )
222  {
223  // Energy OK
224  isThisEnergyOK = true;
225  break;
226  }
227  else
228  {
229 
230  theSystem->GetParticipant( i )->SetDefinition( pd0 );
231  theSystem->GetParticipant( i )->SetPosition( r0 );
232  theSystem->GetParticipant( i )->SetMomentum( p0 );
233 
234  for ( G4int i0i = 0 ; i0i < id-1 ; i0i++ )
235  {
236  //081118
237  //std::cout << "Decay Energitically Blocked deleteing " << i0i+n0 << std::endl;
238  theSystem->DeleteParticipant( i0i+n0 );
239  }
240  //081103
241  theMeanField->Update();
242  }
243 
244  }
245 
246 
247 // Pauli Check
248  if ( isThisEnergyOK == true )
249  {
250  if ( theMeanField->IsPauliBlocked ( i ) != true )
251  {
252 
253  G4bool allOK = true;
254  for ( G4int i0i = 0 ; i0i < i0 ; i0i++ )
255  {
256  if ( theMeanField->IsPauliBlocked ( i0i+n0 ) == true )
257  {
258  allOK = false;
259  break;
260  }
261  }
262 
263  if ( allOK )
264  {
265  decayed = true; //Decay Succeeded
266  }
267  }
268 
269  }
270 //
271 
272  if ( decayed )
273  {
274  //081119
275  //G4cout << "Decay Suceeded! " << std::endl;
277  for ( G4int i0i = 0 ; i0i < i0 ; i0i++ )
278  {
279  theSystem->GetParticipant( i0i+n0 )->SetHitMark();
280  }
281 
282  }
283  else
284  {
285 
286 // Decay Blocked and re-enter orginal participant;
287 
288  if ( isThisEnergyOK == true ) // for false case already done
289  {
290 
291  theSystem->GetParticipant( i )->SetDefinition( pd0 );
292  theSystem->GetParticipant( i )->SetPosition( r0 );
293  theSystem->GetParticipant( i )->SetMomentum( p0 );
294 
295  for ( G4int i0i = 0 ; i0i < i0 ; i0i++ )
296  {
297  //081118
298  //std::cout << "Decay Blocked deleteing " << i0i+n0 << std::endl;
299  theSystem->DeleteParticipant( i0+n0-i0i-1 );
300  }
301  //081103
302  theMeanField->Update();
303  }
304 
305  }
306 
307  } //shortlive
308  } // go next participant
309 //071101
310 
311 
313 
314 //081118
315  //for ( G4int i = 1 ; i < n ; i++ )
316  for ( G4int i = 1 ; i < theSystem->GetTotalNumberOfParticipant() ; i++ )
317  {
318 
319  //std::cout << "Collision i " << i << std::endl;
320  if ( theSystem->GetParticipant( i )->IsThisHit() ) continue;
321 
324  G4double rmi = theSystem->GetParticipant( i )->GetMass();
326 //090331 gamma
327  if ( pdi->GetPDGMass() == 0.0 ) continue;
328 
329  //std::cout << " p4i00 " << p4i << std::endl;
330  for ( G4int j = 0 ; j < i ; j++ )
331  {
332 
333 
334 /*
335  G4cout << "Collision " << i << " " << theSystem->GetParticipant( i )->IsThisProjectile() << G4endl;
336  G4cout << "Collision " << j << " " << theSystem->GetParticipant( j )->IsThisProjectile() << G4endl;
337  G4cout << "Collision " << i << " " << theSystem->GetParticipant( i )->IsThisTarget() << G4endl;
338  G4cout << "Collision " << j << " " << theSystem->GetParticipant( j )->IsThisTarget() << G4endl;
339 */
340 
341  // Only 1 Collision allowed for each particle in a time step.
342  //081119
343  if ( theSystem->GetParticipant( i )->IsThisHit() ) continue;
344  if ( theSystem->GetParticipant( j )->IsThisHit() ) continue;
345 
346  //std::cout << "Collision " << i << " " << j << std::endl;
347 
348  // Do not allow collision between nucleons in target/projectile til its first collision.
350  {
351  if ( theSystem->GetParticipant( j )->IsThisProjectile() ) continue;
352  }
353  else if ( theSystem->GetParticipant( i )->IsThisTarget() )
354  {
355  if ( theSystem->GetParticipant( j )->IsThisTarget() ) continue;
356  }
357 
358 
361  G4double rmj = theSystem->GetParticipant( j )->GetMass();
363 //090331 gamma
364  if ( pdj->GetPDGMass() == 0.0 ) continue;
365 
366  G4double rr2 = theMeanField->GetRR2( i , j );
367 
368 // Here we assume elab (beam momentum less than 5 GeV/n )
369  if ( rr2 > fdeltar*fdeltar ) continue;
370 
371  //G4double s = (p4i+p4j)*(p4i+p4j);
372  //G4double srt = std::sqrt ( s );
373 
374  G4double srt = std::sqrt( (p4i+p4j)*(p4i+p4j) );
375 
376  G4double cutoff = 0.0;
377  G4double fbcmax = 0.0;
378  //110617 fix for gcc 4.6 compilation warnings
379  //G4double sig = 0.0;
380 
381  if ( rmi < 0.94 && rmj < 0.94 )
382  {
383 // nucleon or pion case
384  cutoff = rmi + rmj + 0.02;
385  fbcmax = fbcmax0;
386  //110617 fix for gcc 4.6 compilation warnings
387  //sig = sig0;
388  }
389  else
390  {
391  cutoff = rmi + rmj;
392  fbcmax = fbcmax1;
393  //110617 fix for gcc compilation warnings
394  //sig = sig1;
395  }
396 
397  //std::cout << "Collision cutoff " << i << " " << j << " " << cutoff << std::endl;
398  if ( srt < cutoff ) continue;
399 
400  G4ThreeVector dr = ri - rj;
401  G4double rsq = dr*dr;
402 
403  G4double pij = p4i*p4j;
404  G4double pidr = p4i.vect()*dr;
405  G4double pjdr = p4j.vect()*dr;
406 
407  G4double aij = 1.0 - ( rmi*rmj /pij ) * ( rmi*rmj /pij );
408  G4double bij = pidr / rmi - pjdr*rmi/pij;
409  G4double cij = rsq + ( pidr / rmi ) * ( pidr / rmi );
410  G4double brel = std::sqrt ( std::abs ( cij - bij*bij/aij ) );
411 
412  if ( brel > fbcmax ) continue;
413  //std::cout << "collisions3 " << std::endl;
414 
415  G4double bji = -pjdr/rmj + pidr * rmj /pij;
416 
417  G4double ti = ( pidr/rmi - bij / aij ) * p4i.e() / rmi;
418  G4double tj = (-pjdr/rmj - bji / aij ) * p4j.e() / rmj;
419 
420 
421 /*
422  G4cout << "collisions4 p4i " << p4i << G4endl;
423  G4cout << "collisions4 ri " << ri << G4endl;
424  G4cout << "collisions4 p4j " << p4j << G4endl;
425  G4cout << "collisions4 rj " << rj << G4endl;
426  G4cout << "collisions4 dr " << dr << G4endl;
427  G4cout << "collisions4 pij " << pij << G4endl;
428  G4cout << "collisions4 aij " << aij << G4endl;
429  G4cout << "collisions4 bij bji " << bij << " " << bji << G4endl;
430  G4cout << "collisions4 pidr pjdr " << pidr << " " << pjdr << G4endl;
431  G4cout << "collisions4 p4i.e() p4j.e() " << p4i.e() << " " << p4j.e() << G4endl;
432  G4cout << "collisions4 rmi rmj " << rmi << " " << rmj << G4endl;
433  G4cout << "collisions4 " << ti << " " << tj << G4endl;
434 */
435  if ( std::abs ( ti + tj ) > deltaT ) continue;
436  //std::cout << "collisions4 " << std::endl;
437 
438  G4ThreeVector beta = ( p4i + p4j ).boostVector();
439 
440  G4LorentzVector p = p4i;
441  G4LorentzVector p4icm = p.boost( p.findBoostToCM ( p4j ) );
442  G4ThreeVector pcm = p4icm.vect();
443 
444  G4double prcm = pcm.mag();
445 
446  if ( prcm <= 0.00001 ) continue;
447  //std::cout << "collisions5 " << std::endl;
448 
449  G4bool energetically_forbidden = !( CalFinalStateOfTheBinaryCollision ( i , j ) ); // Use Geant4 Collision Library
450  //G4bool energetically_forbidden = !( CalFinalStateOfTheBinaryCollisionJQMD ( sig , cutoff , pcm , prcm , srt, beta , gamma , i , j ) ); // JQMD Elastic
451 
452 /*
453  G4bool pauli_blocked = false;
454  if ( energetically_forbidden == false ) // result true
455  {
456  if ( theMeanField->IsPauliBlocked ( i ) == true || theMeanField->IsPauliBlocked ( j ) == true )
457  {
458  pauli_blocked = true;
459  //std::cout << "G4QMDRESULT Collsion Pauli Blocked " << std::endl;
460  }
461  }
462  else
463  {
464  if ( theMeanField->IsPauliBlocked ( i ) == true || theMeanField->IsPauliBlocked ( j ) == true )
465  pauli_blocked = false;
466  //std::cout << "G4QMDRESULT Collsion Blocked " << std::endl;
467  }
468 */
469 
470 /*
471  G4cout << "G4QMDRESULT Collsion initial p4 i and j "
472  << p4i << " " << p4j
473  << G4endl;
474 */
475 // 081118
476  //if ( energetically_forbidden == true || pauli_blocked == true )
477  if ( energetically_forbidden == true )
478  {
479 
480  //G4cout << " energetically_forbidden " << G4endl;
481 // Collsion not allowed then re enter orginal participants
482 // Now only momentum, becasuse we only consider elastic scattering of nucleons
483 
484  theSystem->GetParticipant( i )->SetMomentum( p4i.vect() );
485  theSystem->GetParticipant( i )->SetDefinition( pdi );
486  theSystem->GetParticipant( i )->SetPosition( ri );
487 
488  theSystem->GetParticipant( j )->SetMomentum( p4j.vect() );
489  theSystem->GetParticipant( j )->SetDefinition( pdj );
490  theSystem->GetParticipant( j )->SetPosition( rj );
491 
494 
495  }
496  else
497  {
498 
499 
500  G4bool absorption = false;
501  if ( n == theSystem->GetTotalNumberOfParticipant()+1 ) absorption = true;
502  if ( absorption )
503  {
504  //G4cout << "Absorption happend " << G4endl;
505  i = i-1;
506  n = n-1;
507  }
508 
509 // Collsion allowed (really happened)
510 
511  // Unset Projectile/Target flag
513  if ( !absorption ) theSystem->GetParticipant( j )->UnsetInitialMark();
514 
516  if ( !absorption ) theSystem->GetParticipant( j )->SetHitMark();
517 
519 
520 /*
521  G4cout << "G4QMDRESULT Collsion Really Happened between "
522  << i << " and " << j
523  << G4endl;
524  G4cout << "G4QMDRESULT Collsion initial p4 i and j "
525  << p4i << " " << p4j
526  << G4endl;
527  G4cout << "G4QMDRESULT Collsion after p4 i and j "
528  << theSystem->GetParticipant( i )->Get4Momentum()
529  << " "
530  << theSystem->GetParticipant( j )->Get4Momentum()
531  << G4endl;
532  G4cout << "G4QMDRESULT Collsion Diff "
533  << p4i + p4j - theSystem->GetParticipant( i )->Get4Momentum() - theSystem->GetParticipant( j )->Get4Momentum()
534  << G4endl;
535  G4cout << "G4QMDRESULT Collsion initial r i and j "
536  << ri << " " << rj
537  << G4endl;
538  G4cout << "G4QMDRESULT Collsion after r i and j "
539  << theSystem->GetParticipant( i )->GetPosition()
540  << " "
541  << theSystem->GetParticipant( j )->GetPosition()
542  << G4endl;
543 */
544 
545 
546  }
547 
548  }
549 
550  }
551 
552 
553 }
554 
555 
556 
558 {
559 
560 //081103
561  //G4cout << "CalFinalStateOfTheBinaryCollision " << i << " " << j << " " << theSystem->GetTotalNumberOfParticipant() << G4endl;
562 
563  G4bool result = false;
564  G4bool energyOK = false;
565  G4bool pauliOK = false;
566  G4bool abs = false;
567  G4QMDParticipant* absorbed = NULL;
568 
571 
572 //071031
573 
575 
576  G4double eini = epot + p4i.e() + p4j.e();
577 
578 //071031
579  // will use KineticTrack
582  G4LorentzVector p4i0 = p4i*GeV;
583  G4LorentzVector p4j0 = p4j*GeV;
586 
587  for ( G4int iitry = 0 ; iitry < 4 ; iitry++ )
588  {
589 
590  abs = false;
591 
592  G4KineticTrack kt1( pdi0 , 0.0 , ri0 , p4i0 );
593  G4KineticTrack kt2( pdj0 , 0.0 , rj0 , p4j0 );
594 
595  G4LorentzVector p4ix_new;
596  G4LorentzVector p4jx_new;
597  G4KineticTrackVector* secs = NULL;
598  secs = theScatterer->Scatter( kt1 , kt2 );
599 
600  //std::cout << "G4QMDSCATTERER BEFORE " << kt1.GetDefinition()->GetParticleName() << " " << kt1.Get4Momentum()/GeV << " " << kt1.GetPosition()/fermi << std::endl;
601  //std::cout << "G4QMDSCATTERER BEFORE " << kt2.GetDefinition()->GetParticleName() << " " << kt2.Get4Momentum()/GeV << " " << kt2.GetPosition()/fermi << std::endl;
602  //std::cout << "THESCATTERER " << theScatterer->GetCrossSection ( kt1 , kt2 )/millibarn << " " << elastic << " " << sig << std::endl;
603 
604 
605  if ( secs )
606  {
607  G4int iti = 0;
608  if ( secs->size() == 2 )
609  {
610  for ( G4KineticTrackVector::iterator it
611  = secs->begin() ; it != secs->end() ; it++ )
612  {
613  if ( iti == 0 )
614  {
615  theSystem->GetParticipant( i )->SetDefinition( (*it)->GetDefinition() );
616  p4ix_new = (*it)->Get4Momentum()/GeV;
617  //std::cout << "THESCATTERER " << (*it)->GetDefinition()->GetParticleName() << std::endl;
618  theSystem->GetParticipant( i )->SetMomentum( p4ix_new.v() );
619  }
620  if ( iti == 1 )
621  {
622  theSystem->GetParticipant( j )->SetDefinition( (*it)->GetDefinition() );
623  p4jx_new = (*it)->Get4Momentum()/GeV;
624  //std::cout << "THESCATTERER " << p4jx_new.e()-p4jx_new.m() << std::endl;
625  theSystem->GetParticipant( j )->SetMomentum( p4jx_new.v() );
626  }
627  //std::cout << "G4QMDSCATTERER AFTER " << (*it)->GetDefinition()->GetParticleName() << " " << (*it)->Get4Momentum()/GeV << std::endl;
628  iti++;
629  }
630  }
631  else if ( secs->size() == 1 )
632  {
633 //081118
634  abs = true;
635  //G4cout << "G4QMDCollision pion absrorption " << secs->front()->GetDefinition()->GetParticleName() << G4endl;
636  //secs->front()->Decay();
637  theSystem->GetParticipant( i )->SetDefinition( secs->front()->GetDefinition() );
638  p4ix_new = secs->front()->Get4Momentum()/GeV;
639  theSystem->GetParticipant( i )->SetMomentum( p4ix_new.v() );
640 
641  }
642 
643 //081118
644  if ( secs->size() > 2 )
645  {
646 
647  G4cout << "G4QMDCollision secs size > 2; " << secs->size() << G4endl;
648 
649  for ( G4KineticTrackVector::iterator it
650  = secs->begin() ; it != secs->end() ; it++ )
651  {
652  G4cout << "G4QMDSCATTERER AFTER " << (*it)->GetDefinition()->GetParticleName() << " " << (*it)->Get4Momentum()/GeV << G4endl;
653  }
654 
655  }
656 
657  // deleteing KineticTrack
658  for ( G4KineticTrackVector::iterator it
659  = secs->begin() ; it != secs->end() ; it++ )
660  {
661  delete *it;
662  }
663 
664  delete secs;
665  }
666 //071031
667 
668  if ( !abs )
669  {
672  }
673  else
674  {
675  absorbed = theSystem->EraseParticipant( j );
676  theMeanField->Update();
677  }
678 
680 
681  G4double efin = epot + p4ix_new.e() + p4jx_new.e();
682 
683  //std::cout << "Collision NEW epot " << i << " " << j << " " << epot << " " << std::abs ( eini - efin ) - fepse << std::endl;
684 
685 /*
686  G4cout << "Collision efin " << i << " " << j << " " << efin << G4endl;
687  G4cout << "Collision " << i << " " << j << " " << std::abs ( eini - efin ) << " " << fepse << G4endl;
688  G4cout << "Collision " << std::abs ( eini - efin ) << " " << fepse << G4endl;
689 */
690 
691 //071031
692  if ( std::abs ( eini - efin ) < fepse )
693  {
694  // Collison OK
695  //std::cout << "collisions6" << std::endl;
696  //std::cout << "collisions before " << p4i << " " << p4j << std::endl;
697  //std::cout << "collisions after " << theSystem->GetParticipant( i )->Get4Momentum() << " " << theSystem->GetParticipant( j )->Get4Momentum() << std::endl;
698  //std::cout << "collisions dif " << ( p4i + p4j ) - ( theSystem->GetParticipant( i )->Get4Momentum() + theSystem->GetParticipant( j )->Get4Momentum() ) << std::endl;
699  //std::cout << "collisions before " << ri0/fermi << " " << rj0/fermi << std::endl;
700  //std::cout << "collisions after " << theSystem->GetParticipant( i )->GetPosition() << " " << theSystem->GetParticipant( j )->GetPosition() << std::endl;
701  energyOK = true;
702  break;
703  }
704  else
705  {
706  //G4cout << "Energy Not OK " << G4endl;
707  if ( abs )
708  {
709  //G4cout << "TKDB reinsert j " << G4endl;
710  theSystem->InsertParticipant( absorbed , j );
711  theMeanField->Update();
712  }
713  // do not need reinsert in no absroption case
714  }
715 //071031
716  }
717 
718 // Energetically forbidden collision
719 
720  if ( energyOK )
721  {
722  // Pauli Check
723  //G4cout << "Pauli Checking " << theSystem->GetTotalNumberOfParticipant() << G4endl;
724  if ( !abs )
725  {
726  if ( !( theMeanField->IsPauliBlocked ( i ) == true || theMeanField->IsPauliBlocked ( j ) == true ) )
727  {
728  //G4cout << "Binary Collision Happen " << theSystem->GetTotalNumberOfParticipant() << G4endl;
729  pauliOK = true;
730  }
731  }
732  else
733  {
734  //if ( theMeanField->IsPauliBlocked ( i ) == false )
735  //090126 i-1 cause jth is erased
736  if ( theMeanField->IsPauliBlocked ( i-1 ) == false )
737  {
738  //G4cout << "Absorption Happen " << theSystem->GetTotalNumberOfParticipant() << G4endl;
739  delete absorbed;
740  pauliOK = true;
741  }
742  }
743 
744 
745  if ( pauliOK )
746  {
747  result = true;
748  }
749  else
750  {
751  //G4cout << "Pauli Blocked" << G4endl;
752  if ( abs )
753  {
754  //G4cout << "TKDB reinsert j pauli block" << G4endl;
755  theSystem->InsertParticipant( absorbed , j );
756  theMeanField->Update();
757  }
758  }
759  }
760 
761  return result;
762 
763 }
764 
765 
766 
768 {
769 
770  //G4cout << "CalFinalStateOfTheBinaryCollisionJQMD" << G4endl;
771 
772  G4bool result = true;
773 
775  G4double rmi = theSystem->GetParticipant( i )->GetMass();
777 
779  G4double rmj = theSystem->GetParticipant( j )->GetMass();
781 
782  G4double pr = prcm;
783 
784  G4double c2 = pcm.z()/pr;
785 
786  G4double csrt = srt - cutoff;
787 
788  //G4double pri = prcm;
789  //G4double prf = sqrt ( 0.25 * srt*srt -rm2 );
790 
791  G4double asrt = srt - rmi - rmj;
792  G4double pra = prcm;
793 
794 
795 
796  G4double elastic = 0.0;
797 
798  if ( zi == zj )
799  {
800  if ( csrt < 0.4286 )
801  {
802  elastic = 35.0 / ( 1. + csrt * 100.0 ) + 20.0;
803  }
804  else
805  {
806  elastic = ( - std::atan( ( csrt - 0.4286 ) * 1.5 - 0.8 )
807  * 2. / pi + 1.0 ) * 9.65 + 7.0;
808  }
809  }
810  else
811  {
812  if ( csrt < 0.4286 )
813  {
814  elastic = 28.0 / ( 1. + csrt * 100.0 ) + 27.0;
815  }
816  else
817  {
818  elastic = ( - std::atan( ( csrt - 0.4286 ) * 1.5 - 0.8 )
819  * 2. / pi + 1.0 ) * 12.34 + 10.0;
820  }
821  }
822 
823 // std::cout << "Collision csrt " << i << " " << j << " " << csrt << std::endl;
824 // std::cout << "Collision elstic " << i << " " << j << " " << elastic << std::endl;
825 
826 
827 // std::cout << "Collision sig " << i << " " << j << " " << sig << std::endl;
828  if ( G4UniformRand() > elastic / sig )
829  {
830  //std::cout << "Inelastic " << std::endl;
831  //std::cout << "elastic/sig " << elastic/sig << std::endl;
832  return result;
833  }
834  else
835  {
836  //std::cout << "Elastic " << std::endl;
837  }
838 // std::cout << "Collision ELSTIC " << i << " " << j << std::endl;
839 
840 
841  G4double as = G4Pow::GetInstance()->powN ( 3.65 * asrt , 6 );
842  G4double a = 6.0 * as / (1.0 + as);
843  G4double ta = -2.0 * pra*pra;
844  G4double x = G4UniformRand();
845  G4double t1 = G4Log( (1-x) * G4Exp(2.*a*ta) + x ) / a;
846  G4double c1 = 1.0 - t1/ta;
847 
848  if( std::abs(c1) > 1.0 ) c1 = 2.0 * x - 1.0;
849 
850 /*
851  G4cout << "Collision as " << i << " " << j << " " << as << G4endl;
852  G4cout << "Collision a " << i << " " << j << " " << a << G4endl;
853  G4cout << "Collision ta " << i << " " << j << " " << ta << G4endl;
854  G4cout << "Collision x " << i << " " << j << " " << x << G4endl;
855  G4cout << "Collision t1 " << i << " " << j << " " << t1 << G4endl;
856  G4cout << "Collision c1 " << i << " " << j << " " << c1 << G4endl;
857 */
858  t1 = 2.0*pi*G4UniformRand();
859 // std::cout << "Collision t1 " << i << " " << j << " " << t1 << std::endl;
860  G4double t2 = 0.0;
861  if ( pcm.x() == 0.0 && pcm.y() == 0 )
862  {
863  t2 = 0.0;
864  }
865  else
866  {
867  t2 = std::atan2( pcm.y() , pcm.x() );
868  }
869 // std::cout << "Collision t2 " << i << " " << j << " " << t2 << std::endl;
870 
871  G4double s1 = std::sqrt ( 1.0 - c1*c1 );
872  G4double s2 = std::sqrt ( 1.0 - c2*c2 );
873 
874  G4double ct1 = std::cos(t1);
875  G4double st1 = std::sin(t1);
876 
877  G4double ct2 = std::cos(t2);
878  G4double st2 = std::sin(t2);
879 
880  G4double ss = c2*s1*ct1 + s2*c1;
881 
882  pcm.setX( pr * ( ss*ct2 - s1*st1*st2) );
883  pcm.setY( pr * ( ss*st2 + s1*st1*ct2) );
884  pcm.setZ( pr * ( c1*c2 - s1*s2*ct1) );
885 
886 // std::cout << "Collision pcm " << i << " " << j << " " << pcm << std::endl;
887 
889 
890  G4double eini = epot + p4i.e() + p4j.e();
891  G4double etwo = p4i.e() + p4j.e();
892 
893 /*
894  G4cout << "Collision epot " << i << " " << j << " " << epot << G4endl;
895  G4cout << "Collision eini " << i << " " << j << " " << eini << G4endl;
896  G4cout << "Collision etwo " << i << " " << j << " " << etwo << G4endl;
897 */
898 
899 
900  for ( G4int itry = 0 ; itry < 4 ; itry++ )
901  {
902 
903  G4double eicm = std::sqrt ( rmi*rmi + pcm*pcm );
904  G4double pibeta = pcm*beta;
905 
906  G4double trans = gamma * ( gamma * pibeta / ( gamma + 1 ) + eicm );
907 
908  G4ThreeVector pi_new = beta*trans + pcm;
909 
910  G4double ejcm = std::sqrt ( rmj*rmj + pcm*pcm );
911  trans = gamma * ( gamma * pibeta / ( gamma + 1 ) + ejcm );
912 
913  G4ThreeVector pj_new = beta*trans - pcm;
914 
915 //
916 // Delete old
917 // Add new Particitipants
918 //
919 // Now only change momentum ( Beacuse we only have elastic sctter of nucleon
920 // In future Definition also will be change
921 //
922 
923  theSystem->GetParticipant( i )->SetMomentum( pi_new );
924  theSystem->GetParticipant( j )->SetMomentum( pj_new );
925 
926  G4double pi_new_e = (theSystem->GetParticipant( i )->Get4Momentum()).e();
927  G4double pj_new_e = (theSystem->GetParticipant( j )->Get4Momentum()).e();
928 
931 
933 
934  G4double efin = epot + pi_new_e + pj_new_e ;
935 
936  //std::cout << "Collision NEW epot " << i << " " << j << " " << epot << " " << std::abs ( eini - efin ) - fepse << std::endl;
937 /*
938  G4cout << "Collision efin " << i << " " << j << " " << efin << G4endl;
939  G4cout << "Collision " << i << " " << j << " " << std::abs ( eini - efin ) << " " << fepse << G4endl;
940  G4cout << "Collision " << std::abs ( eini - efin ) << " " << fepse << G4endl;
941 */
942 
943 //071031
944  if ( std::abs ( eini - efin ) < fepse )
945  {
946  // Collison OK
947  //std::cout << "collisions6" << std::endl;
948  //std::cout << "collisions before " << p4i << " " << p4j << std::endl;
949  //std::cout << "collisions after " << theSystem->GetParticipant( i )->Get4Momentum() << " " << theSystem->GetParticipant( j )->Get4Momentum() << std::endl;
950  //std::cout << "collisions dif " << ( p4i + p4j ) - ( theSystem->GetParticipant( i )->Get4Momentum() + theSystem->GetParticipant( j )->Get4Momentum() ) << std::endl;
951  //std::cout << "collisions before " << rix/fermi << " " << rjx/fermi << std::endl;
952  //std::cout << "collisions after " << theSystem->GetParticipant( i )->GetPosition() << " " << theSystem->GetParticipant( j )->GetPosition() << std::endl;
953  }
954 //071031
955 
956  if ( std::abs ( eini - efin ) < fepse ) return result; // Collison OK
957 
958  G4double cona = ( eini - efin + etwo ) / gamma;
959  G4double fac2 = 1.0 / ( 4.0 * cona*cona * pr*pr ) *
960  ( ( cona*cona - ( rmi*rmi + rmj*rmj ) )*( cona*cona - ( rmi*rmi + rmj*rmj ) )
961  - 4.0 * rmi*rmi * rmj*rmj );
962 
963  if ( fac2 > 0 )
964  {
965  G4double fact = std::sqrt ( fac2 );
966  pcm = fact*pcm;
967  }
968 
969 
970  }
971 
972 // Energetically forbidden collision
973  result = false;
974 
975  return result;
976 
977 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4ThreeVector GetPosition()
static c2_factory< G4double > c2
G4QMDMeanField * theMeanField
G4double powN(G4double x, G4int n) const
Definition: G4Pow.cc:128
void SetParticipant(G4QMDParticipant *particle)
Definition: G4QMDSystem.hh:51
const G4double fact
CLHEP::Hep3Vector G4ThreeVector
G4QMDSystem * theSystem
G4KineticTrackVector * Decay()
void SetPosition(G4ThreeVector r)
G4int GetChargeInUnitOfEplus()
G4double a
Definition: TRTMaterials.hh:39
const G4ParticleDefinition * GetDefinition()
G4QMDParticipant * EraseParticipant(G4int i)
Definition: G4QMDSystem.hh:56
int G4int
Definition: G4Types.hh:78
void Cal2BodyQuantities()
G4bool IsPauliBlocked(G4int)
G4ThreeVector GetMomentum()
void SetDefinition(const G4ParticleDefinition *pd)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4Scatterer * theScatterer
bool G4bool
Definition: G4Types.hh:79
static const double GeV
Definition: G4SIunits.hh:214
void CalKinematicsOfBinaryCollisions(G4double)
G4int GetTotalNumberOfParticipant()
Definition: G4QMDSystem.hh:60
const G4int n
void SetMomentum(G4ThreeVector p)
static const G4double c1
G4double GetTotalPotential()
G4QMDParticipant * GetParticipant(G4int i)
Definition: G4QMDSystem.hh:62
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double GetPDGMass() const
static const double pi
Definition: G4SIunits.hh:74
void DeleteParticipant(G4int i)
Definition: G4QMDSystem.hh:57
G4LorentzVector Get4Momentum()
G4double GetRR2(G4int i, G4int j)
virtual G4KineticTrackVector * Scatter(const G4KineticTrack &trk1, const G4KineticTrack &trk2) const
Definition: G4Scatterer.cc:284
G4bool CalFinalStateOfTheBinaryCollisionJQMD(G4double, G4double, G4ThreeVector, G4double, G4double, G4ThreeVector, G4double, G4int, G4int)
const G4double x[NPOINTSGL]
G4bool CalFinalStateOfTheBinaryCollision(G4int, G4int)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void InsertParticipant(G4QMDParticipant *particle, G4int j)
Definition: G4QMDSystem.cc:110
void IncrementCollisionCounter()
Definition: G4QMDSystem.hh:64
static const double fermi
Definition: G4SIunits.hh:102
const G4double r0
G4double elastic(Particle const *const p1, Particle const *const p2)
CLHEP::HepLorentzVector G4LorentzVector