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

Last change on this file since 978 was 962, checked in by garnier, 15 years ago

update processes

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