source: trunk/source/processes/hadronic/models/qmd/src/G4QMDCollision.cc @ 1253

Last change on this file since 1253 was 1055, checked in by garnier, 15 years ago

maj sur la beta de geant 4.9.3

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