source: trunk/source/processes/hadronic/models/qmd/src/G4QMDReaction.cc @ 1198

Last change on this file since 1198 was 1196, checked in by garnier, 15 years ago

update CVS release candidate geant4.9.3.01

File size: 25.0 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// 080505 Fixed and changed sampling method of impact parameter by T. Koi
27// 080602 Fix memory leaks by T. Koi
28// 080612 Delete unnecessary dependency and unused functions
29//        Change criterion of reaction by T. Koi
30// 081107 Add UnUseGEM (then use the default channel of G4Evaporation)
31//            UseFrag (chage criterion of a inelastic reaction)
32//        Fix bug in nucleon projectiles  by T. Koi   
33// 090122 Be8 -> Alpha + Alpha
34// 090331 Change member shenXS and genspaXS object to pointer
35// 091119 Fix for incidence of neutral particles
36//
37#include "G4QMDReaction.hh"
38#include "G4QMDNucleus.hh"
39#include "G4QMDGroundStateNucleus.hh"
40
41#include "G4NistManager.hh"
42
43G4QMDReaction::G4QMDReaction()
44: system ( NULL )
45, deltaT ( 1 ) // in fsec
46, maxTime ( 100 ) // will have maxTime-th time step
47, gem ( true )
48, frag ( false )
49{
50
51   //090331
52   shenXS = new G4IonsShenCrossSection();
53   //genspaXS = new G4GeneralSpaceNNCrossSection();
54   meanField = new G4QMDMeanField();
55   collision = new G4QMDCollision();
56
57   excitationHandler = new G4ExcitationHandler;
58   evaporation = new G4Evaporation;
59   excitationHandler->SetEvaporation( evaporation );
60   setEvaporationCh();
61}
62
63
64
65G4QMDReaction::~G4QMDReaction()
66{
67   delete evaporation; 
68   delete excitationHandler;
69   delete collision;
70   delete meanField;
71}
72
73
74
75G4HadFinalState* G4QMDReaction::ApplyYourself( const G4HadProjectile & projectile , G4Nucleus & target )
76{
77   //G4cout << "G4QMDReaction::ApplyYourself" << G4endl;
78
79   theParticleChange.Clear();
80
81   system = new G4QMDSystem;
82
83   G4int proj_Z = 0;
84   G4int proj_A = 0;
85   G4ParticleDefinition* proj_pd = ( G4ParticleDefinition* ) projectile.GetDefinition();
86   if ( proj_pd->GetParticleType() == "nucleus" )
87   {
88      proj_Z = proj_pd->GetAtomicNumber();
89      proj_A = proj_pd->GetAtomicMass();
90   }
91   else
92   {
93      proj_Z = (int)( proj_pd->GetPDGCharge()/eplus );
94      proj_A = 1;
95   }
96   G4int targ_Z = int ( target.GetZ() + 0.5 );
97   G4int targ_A = int ( target.GetN() + 0.5 );
98   G4ParticleDefinition* targ_pd = G4ParticleTable::GetParticleTable()->GetIon( targ_Z , targ_A , 0.0 );
99
100
101   G4NistManager* nistMan = G4NistManager::Instance();
102//   G4Element* G4NistManager::FindOrBuildElement( targ_Z );
103
104     const G4DynamicParticle* proj_dp = new G4DynamicParticle ( proj_pd , projectile.Get4Momentum() );
105     const G4Element* targ_ele =  nistMan->FindOrBuildElement( targ_Z ); 
106     G4double aTemp = projectile.GetMaterial()->GetTemperature();
107
108     //090331
109     G4double xs_0 = shenXS->GetCrossSection ( proj_dp , targ_ele , aTemp );
110     //G4double xs_0 = genspaXS->GetCrossSection ( proj_dp , targ_ele , aTemp );
111     G4double bmax_0 = std::sqrt( xs_0 / pi );
112     //std::cout << "bmax_0 in fm (fermi) " <<  bmax_0/fermi << std::endl;
113
114     //delete proj_dp;
115
116   G4bool elastic = true;
117   
118   std::vector< G4QMDNucleus* > nucleuses; // Secondary nuceluses
119   G4ThreeVector boostToReac; // ReactionSystem (CM or NN);
120   G4ThreeVector boostBackToLAB; // Reaction System to LAB;
121
122   G4LorentzVector targ4p( G4ThreeVector( 0.0 ) , targ_pd->GetPDGMass()/GeV );
123   G4ThreeVector boostLABtoCM = targ4p.findBoostToCM( proj_dp->Get4Momentum()/GeV ); // CM of target and proj;
124
125   G4double p1 = proj_dp->GetMomentum().mag()/GeV/proj_A; 
126   G4double m1 = proj_dp->GetDefinition()->GetPDGMass()/GeV/proj_A;
127   G4double e1 = std::sqrt( p1*p1 + m1*m1 ); 
128   G4double e2 = targ_pd->GetPDGMass()/GeV/targ_A;
129   G4double beta_nn = -p1 / ( e1+e2 );
130
131   G4ThreeVector boostLABtoNN ( 0. , 0. , beta_nn ); // CM of NN;
132
133   G4double beta_nncm = ( - boostLABtoCM.beta() + boostLABtoNN.beta() ) / ( 1 - boostLABtoCM.beta() * boostLABtoNN.beta() ) ; 
134
135   //std::cout << targ4p << std::endl;
136   //std::cout << proj_dp->Get4Momentum()<< std::endl;
137   //std::cout << beta_nncm << std::endl;
138   G4ThreeVector boostNNtoCM( 0. , 0. , beta_nncm ); //
139   G4ThreeVector boostCMtoNN( 0. , 0. , -beta_nncm ); //
140
141   boostToReac = boostLABtoNN; 
142   boostBackToLAB = -boostLABtoNN; 
143
144   delete proj_dp; 
145
146   while ( elastic ) 
147   {
148
149// impact parameter
150      G4double bmax = 1.05*(bmax_0/fermi);  // 10% for Peripheral reactions
151      G4double b = bmax * std::sqrt ( G4UniformRand() );
152//071112
153      //G4double b = 0;
154      //G4double b = bmax;
155      //G4double b = bmax/1.05 * 0.7 * G4UniformRand();
156
157      //G4cout << "G4QMDRESULT bmax_0 = " << bmax_0/fermi << " fm, bmax = " << bmax << " fm , b = " << b  << " fm " << G4endl;
158
159      G4double plab = projectile.GetTotalMomentum()/GeV;
160      G4double elab = (projectile.GetKineticEnergy() + proj_pd->GetPDGMass() + targ_pd->GetPDGMass() )/GeV;
161
162      calcOffSetOfCollision( b , proj_pd , targ_pd , plab , elab , bmax , boostCMtoNN );
163
164// Projectile
165      G4LorentzVector proj4pLAB = projectile.Get4Momentum()/GeV;
166
167      G4QMDGroundStateNucleus* proj(NULL); 
168      if ( projectile.GetDefinition()->GetParticleType() == "nucleus" 
169        || projectile.GetDefinition()->GetParticleName() == "proton"
170        || projectile.GetDefinition()->GetParticleName() == "neutron" )
171      {
172
173         proj_Z = proj_pd->GetAtomicNumber();
174         proj_A = proj_pd->GetAtomicMass();
175
176         proj = new G4QMDGroundStateNucleus( proj_Z , proj_A );
177         //proj->ShowParticipants();
178
179      }
180
181      meanField->SetSystem ( proj );
182      proj->SetTotalPotential( meanField->GetTotalPotential() );
183      proj->CalEnergyAndAngularMomentumInCM();
184
185// Target
186      G4int iz = int ( target.GetZ() );
187      G4int ia = int ( target.GetN() );
188
189      G4QMDGroundStateNucleus* targ = new G4QMDGroundStateNucleus( iz , ia );
190
191      meanField->SetSystem (targ );
192      targ->SetTotalPotential( meanField->GetTotalPotential() );
193      targ->CalEnergyAndAngularMomentumInCM();
194   
195      //G4LorentzVector targ4p( G4ThreeVector( 0.0 ) , targ->GetNuclearMass()/GeV );
196// Boost Vector to CM
197      //boostToCM = targ4p.findBoostToCM( proj4pLAB );
198
199//    Target
200      for ( G4int i = 0 ; i < targ->GetTotalNumberOfParticipant() ; i++ )
201      {
202
203         G4ThreeVector p0 = targ->GetParticipant( i )->GetMomentum();
204         G4ThreeVector r0 = targ->GetParticipant( i )->GetPosition();
205
206         G4ThreeVector p ( p0.x() + coulomb_collision_px_targ
207                         , p0.y() 
208                         , p0.z() * coulomb_collision_gamma_targ + coulomb_collision_pz_targ ); 
209
210         G4ThreeVector r ( r0.x() + coulomb_collision_rx_targ
211                         , r0.y() 
212                         , r0.z() / coulomb_collision_gamma_targ + coulomb_collision_rz_targ ); 
213     
214         system->SetParticipant( new G4QMDParticipant( targ->GetParticipant( i )->GetDefinition() , p , r ) );
215         system->GetParticipant( i )->SetTarget();
216
217      }
218
219      G4LorentzVector proj4pCM = CLHEP::boostOf ( proj4pLAB , boostToReac );
220      G4LorentzVector targ4pCM = CLHEP::boostOf ( targ4p , boostToReac );
221
222
223//    Projectile
224      if ( proj != NULL )
225      {
226
227//    projectile is nucleus
228
229         for ( G4int i = 0 ; i < proj->GetTotalNumberOfParticipant() ; i++ )
230         {
231
232            G4ThreeVector p0 = proj->GetParticipant( i )->GetMomentum();
233            G4ThreeVector r0 = proj->GetParticipant( i )->GetPosition();
234
235            G4ThreeVector p ( p0.x() + coulomb_collision_px_proj
236                            , p0.y() 
237                            , p0.z() * coulomb_collision_gamma_proj + coulomb_collision_pz_proj ); 
238
239            G4ThreeVector r ( r0.x() + coulomb_collision_rx_proj
240                            , r0.y() 
241                            , r0.z() / coulomb_collision_gamma_proj + coulomb_collision_rz_proj ); 
242     
243            system->SetParticipant( new G4QMDParticipant( proj->GetParticipant( i )->GetDefinition() , p  , r ) );
244            system->GetParticipant ( i + targ->GetTotalNumberOfParticipant() )->SetProjectile();
245         }
246
247      }
248      else
249      {
250
251//       projectile is particle
252
253         G4int i = targ->GetTotalNumberOfParticipant(); 
254     
255         G4ThreeVector p0( 0 ); 
256         G4ThreeVector r0( 0 );
257
258
259         G4ThreeVector p ( p0.x() + coulomb_collision_px_proj
260                         , p0.y() 
261                         , p0.z() * coulomb_collision_gamma_proj + coulomb_collision_pz_proj ); 
262
263         G4ThreeVector r ( r0.x() + coulomb_collision_rx_proj
264                         , r0.y() 
265                         , r0.z() / coulomb_collision_gamma_proj + coulomb_collision_rz_proj ); 
266
267         system->SetParticipant( new G4QMDParticipant( (G4ParticleDefinition*)projectile.GetDefinition() , p , r ) );
268         system->GetParticipant ( i )->SetProjectile();
269      }
270
271      delete targ;
272      delete proj;
273
274
275   meanField->SetSystem ( system );
276   collision->SetMeanField ( meanField );
277
278// Time Evolution
279   //std::cout << "Start time evolution " << std::endl;
280   //system->ShowParticipants();
281   for ( G4int i = 0 ; i < maxTime ; i++ )
282   {
283      //G4cout << " do Paropagate " << i << " th time step. " << G4endl;
284      meanField->DoPropagation( deltaT );
285      //system->ShowParticipants();
286      collision->CalKinematicsOfBinaryCollisions( deltaT );
287
288      if ( i / 10 * 10 == i ) 
289      {
290         //G4cout << i << " th time step. " << G4endl;
291         //system->ShowParticipants();
292      } 
293      //system->ShowParticipants();
294   }
295   //system->ShowParticipants();
296
297
298   //std::cout << "Doing Cluster Judgment " << std::endl;
299
300   nucleuses = meanField->DoClusterJudgment();
301
302// Elastic Judgment 
303
304   G4int numberOfSecondary = int ( nucleuses.size() ) + system->GetTotalNumberOfParticipant(); 
305
306   G4int sec_a_Z = 0;
307   G4int sec_a_A = 0;
308   G4ParticleDefinition* sec_a_pd = NULL;
309   G4int sec_b_Z = 0;
310   G4int sec_b_A = 0;
311   G4ParticleDefinition* sec_b_pd = NULL;
312
313   if ( numberOfSecondary == 2 )
314   {
315
316      G4bool elasticLike_system = false;
317      if ( nucleuses.size() == 2 ) 
318      {
319
320         sec_a_Z = nucleuses[0]->GetAtomicNumber();
321         sec_a_A = nucleuses[0]->GetMassNumber();
322         sec_b_Z = nucleuses[1]->GetAtomicNumber();
323         sec_b_A = nucleuses[1]->GetMassNumber();
324
325         if ( ( sec_a_Z == proj_Z && sec_a_A == proj_A && sec_b_Z == targ_Z && sec_b_A == targ_A )
326           || ( sec_a_Z == targ_Z && sec_a_A == targ_A && sec_b_Z == proj_Z && sec_b_A == proj_A ) )
327         {
328            elasticLike_system = true;
329         } 
330
331      }
332      else if ( nucleuses.size() == 1 ) 
333      {
334
335         sec_a_Z = nucleuses[0]->GetAtomicNumber();
336         sec_a_A = nucleuses[0]->GetMassNumber();
337         sec_b_pd = system->GetParticipant( 0 )->GetDefinition();
338
339         if ( ( sec_a_Z == proj_Z && sec_a_A == proj_A && sec_b_pd == targ_pd )
340           || ( sec_a_Z == targ_Z && sec_a_A == targ_A && sec_b_pd == proj_pd ) )
341         {
342            elasticLike_system = true;
343         } 
344
345      } 
346      else
347      {
348
349         sec_a_pd = system->GetParticipant( 0 )->GetDefinition();
350         sec_b_pd = system->GetParticipant( 1 )->GetDefinition();
351 
352         if ( ( sec_a_pd == proj_pd && sec_b_pd == targ_pd ) 
353           || ( sec_a_pd == targ_pd && sec_b_pd == proj_pd ) ) 
354         {
355            elasticLike_system = true;
356         } 
357
358      } 
359
360      if ( elasticLike_system == true )
361      {
362
363         G4bool elasticLike_energy = true;
364//    Cal ExcitationEnergy
365         for ( G4int i = 0 ; i < int ( nucleuses.size() ) ; i++ )
366         { 
367
368            //meanField->SetSystem( nucleuses[i] );
369            meanField->SetNucleus( nucleuses[i] );
370            //nucleuses[i]->SetTotalPotential( meanField->GetTotalPotential() );
371            //nucleuses[i]->CalEnergyAndAngularMomentumInCM();
372
373            if ( nucleuses[i]->GetExcitationEnergy()*GeV > 1.0*MeV ) elasticLike_energy = false; 
374
375         } 
376
377//    Check Collision
378         G4bool withCollision = true;
379         if ( system->GetNOCollision() == 0 ) withCollision = false;
380
381//    Final judegement for Inelasitc or Elastic;
382//
383//       ElasticLike without Collision
384         //if ( elasticLike_energy == true && withCollision == false ) elastic = true;  // ielst = 0
385//       ElasticLike with Collision
386         //if ( elasticLike_energy == true && withCollision == true ) elastic = true;   // ielst = 1
387//       InelasticLike without Collision
388         //if ( elasticLike_energy == false ) elastic = false;                          // ielst = 2               
389         if ( frag == true )
390            if ( elasticLike_energy == false ) elastic = false;
391//       InelasticLike with Collision
392         if ( elasticLike_energy == false && withCollision == true ) elastic = false; // ielst = 3
393
394      }
395
396      }
397      else
398      {
399
400//       numberOfSecondary != 2
401         elastic = false;
402
403      }
404
405//071115
406      //G4cout << elastic << G4endl;
407      // if elastic is true try again from sampling of impact parameter
408
409      if ( elastic == true )
410      {
411         // delete this nucleues
412         for ( std::vector< G4QMDNucleus* >::iterator
413               it = nucleuses.begin() ; it != nucleuses.end() ; it++ )
414         {
415            delete *it;
416         }
417         nucleuses.clear();
418      }
419   } 
420
421
422// Statical Decay Phase
423
424   for ( std::vector< G4QMDNucleus* >::iterator it
425       = nucleuses.begin() ; it != nucleuses.end() ; it++ )
426   {
427
428/*
429      std::cout << "G4QMDRESULT "
430                << (*it)->GetAtomicNumber()
431                << " "
432                << (*it)->GetMassNumber()
433                << " "
434                << (*it)->Get4Momentum()
435                << " "
436                << (*it)->Get4Momentum().vect()
437                << " "
438                << (*it)->Get4Momentum().restMass()
439                << " "
440                << (*it)->GetNuclearMass()/GeV
441                << std::endl;
442*/
443
444      meanField->SetNucleus ( *it );
445
446      if ( (*it)->GetAtomicNumber() == 0  // neutron cluster
447        || (*it)->GetAtomicNumber() == (*it)->GetMassNumber() ) // proton cluster
448      {
449         // push back system
450         for ( G4int i = 0 ; i < (*it)->GetTotalNumberOfParticipant() ; i++ )
451         {
452            G4QMDParticipant* aP = new G4QMDParticipant( ( (*it)->GetParticipant( i ) )->GetDefinition() , ( (*it)->GetParticipant( i ) )->GetMomentum() , ( (*it)->GetParticipant( i ) )->GetPosition() ); 
453            system->SetParticipant ( aP ); 
454         } 
455         continue; 
456      }
457
458      G4double nucleus_e = std::sqrt ( std::pow ( (*it)->GetNuclearMass()/GeV , 2 ) + std::pow ( (*it)->Get4Momentum().vect().mag() , 2 ) );
459      G4LorentzVector nucleus_p4CM ( (*it)->Get4Momentum().vect() , nucleus_e ); 
460
461//      std::cout << "G4QMDRESULT nucleus deltaQ " << deltaQ << std::endl;
462
463      G4int ia = (*it)->GetMassNumber();
464      G4int iz = (*it)->GetAtomicNumber();
465
466      G4LorentzVector lv ( G4ThreeVector( 0.0 ) , (*it)->GetExcitationEnergy()*GeV + G4ParticleTable::GetParticleTable()->GetIonTable()->GetIonMass( iz , ia ) );
467
468      G4Fragment* aFragment = new G4Fragment( ia , iz , lv );
469
470      G4ReactionProductVector* rv;
471      rv = excitationHandler->BreakItUp( *aFragment );
472      G4bool notBreak = true;
473      for ( G4ReactionProductVector::iterator itt
474          = rv->begin() ; itt != rv->end() ; itt++ )
475      {
476
477          notBreak = false;
478          // Secondary from this nucleus (*it)
479          G4ParticleDefinition* pd = (*itt)->GetDefinition();
480
481          G4LorentzVector p4 ( (*itt)->GetMomentum()/GeV , (*itt)->GetTotalEnergy()/GeV );  //in nucleus(*it) rest system
482          G4LorentzVector p4_CM = CLHEP::boostOf( p4 , -nucleus_p4CM.findBoostToCM() );  // Back to CM
483          G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB ); // Back to LAB 
484
485
486//090122
487          //theParticleChange.AddSecondary( dp );
488          if ( !( pd->GetAtomicNumber() == 4 && pd->GetAtomicMass() == 8 ) )
489          {
490             G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV ); 
491             theParticleChange.AddSecondary( dp ); 
492          }
493          else
494          {
495             //Be8 -> Alpha + Alpha + Q
496             G4ThreeVector randomized_direction( G4UniformRand() , G4UniformRand() , G4UniformRand() );
497             randomized_direction = randomized_direction.unit();
498             G4double q_decay = (*itt)->GetMass() - 2*G4Alpha::Alpha()->GetPDGMass();
499             G4double p_decay = std::sqrt ( std::pow(G4Alpha::Alpha()->GetPDGMass()+q_decay/2,2) - std::pow(G4Alpha::Alpha()->GetPDGMass() , 2 ) ); 
500             G4LorentzVector p4_a1 ( p_decay*randomized_direction , G4Alpha::Alpha()->GetPDGMass()+q_decay/2 );  //in Be8 rest system
501             
502             G4LorentzVector p4_a1_Be8 = CLHEP::boostOf ( p4_a1/GeV , -p4.findBoostToCM() );
503             G4LorentzVector p4_a1_CM = CLHEP::boostOf ( p4_a1_Be8 , -nucleus_p4CM.findBoostToCM() );
504             G4LorentzVector p4_a1_LAB = CLHEP::boostOf ( p4_a1_CM , boostBackToLAB );
505
506             G4LorentzVector p4_a2 ( -p_decay*randomized_direction , G4Alpha::Alpha()->GetPDGMass()+q_decay/2 );  //in Be8 rest system
507             
508             G4LorentzVector p4_a2_Be8 = CLHEP::boostOf ( p4_a2/GeV , -p4.findBoostToCM() );
509             G4LorentzVector p4_a2_CM = CLHEP::boostOf ( p4_a2_Be8 , -nucleus_p4CM.findBoostToCM() );
510             G4LorentzVector p4_a2_LAB = CLHEP::boostOf ( p4_a2_CM , boostBackToLAB );
511             
512             G4DynamicParticle* dp1 = new G4DynamicParticle( G4Alpha::Alpha() , p4_a1_LAB*GeV ); 
513             G4DynamicParticle* dp2 = new G4DynamicParticle( G4Alpha::Alpha() , p4_a2_LAB*GeV ); 
514             theParticleChange.AddSecondary( dp1 ); 
515             theParticleChange.AddSecondary( dp2 ); 
516          }
517//090122
518
519/*
520          std::cout
521                << "Regist Secondary "
522                << (*itt)->GetDefinition()->GetParticleName()
523                << " "
524                << (*itt)->GetMomentum()/GeV
525                << " "
526                << (*itt)->GetKineticEnergy()/GeV
527                << " "
528                << (*itt)->GetMass()/GeV
529                << " "
530                << (*itt)->GetTotalEnergy()/GeV
531                << " "
532                << (*itt)->GetTotalEnergy()/GeV * (*itt)->GetTotalEnergy()/GeV
533                 - (*itt)->GetMomentum()/GeV * (*itt)->GetMomentum()/GeV
534                << " "
535                << nucleus_p4CM.findBoostToCM()
536                << " "
537                << p4
538                << " "
539                << p4_CM
540                << " "
541                << p4_LAB
542                << std::endl;
543*/
544
545      }
546      if ( notBreak == true )
547      {
548
549         G4ParticleDefinition* pd = G4ParticleTable::GetParticleTable()->GetIon( (*it)->GetAtomicNumber() , (*it)->GetMassNumber(), (*it)->GetExcitationEnergy()*GeV );
550         G4LorentzVector p4_CM = nucleus_p4CM;
551         G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB ); // Back to LAB 
552         G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV ); 
553         theParticleChange.AddSecondary( dp ); 
554
555      }
556
557      for ( G4ReactionProductVector::iterator itt
558          = rv->begin() ; itt != rv->end() ; itt++ )
559      {
560          delete *itt;
561      }
562      delete rv;
563
564      delete aFragment;
565
566   }
567
568
569
570   for ( G4int i = 0 ; i < system->GetTotalNumberOfParticipant() ; i++ )
571   {
572
573      // Secondary particles
574
575      G4ParticleDefinition* pd = system->GetParticipant( i )->GetDefinition();
576      G4LorentzVector p4_CM = system->GetParticipant( i )->Get4Momentum();
577      G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB );
578      G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV ); 
579      theParticleChange.AddSecondary( dp ); 
580
581/*
582      G4cout << "G4QMDRESULT "
583      << "r" << i << " " << system->GetParticipant ( i ) -> GetPosition() << " "
584      << "p" << i << " " << system->GetParticipant ( i ) -> Get4Momentum()
585      << G4endl;
586*/
587
588   }
589
590   for ( std::vector< G4QMDNucleus* >::iterator it
591       = nucleuses.begin() ; it != nucleuses.end() ; it++ )
592   {
593      delete *it;  // delete nulceuse
594   }
595   nucleuses.clear();
596
597   system->Clear();
598   delete system; 
599
600   theParticleChange.SetStatusChange( stopAndKill );
601
602   return &theParticleChange;
603
604}
605
606
607
608void G4QMDReaction::calcOffSetOfCollision( G4double b , 
609G4ParticleDefinition* pd_proj , 
610G4ParticleDefinition* pd_targ , 
611G4double ptot , G4double etot , G4double bmax , G4ThreeVector boostToCM )
612{
613   G4double mass_proj = pd_proj->GetPDGMass()/GeV;
614   G4double mass_targ = pd_targ->GetPDGMass()/GeV;
615 
616   G4double stot = std::sqrt ( etot*etot - ptot*ptot );
617
618   G4double pstt = std::sqrt ( ( stot*stot - ( mass_proj + mass_targ ) * ( mass_proj + mass_targ ) 
619                  ) * ( stot*stot - ( mass_proj - mass_targ ) * ( mass_proj - mass_targ ) ) ) 
620                 / ( 2.0 * stot );
621
622   G4double pzcc = pstt;
623   G4double eccm = stot - ( mass_proj + mass_targ );
624   
625   G4int zp = pd_proj->GetAtomicNumber();
626   G4int ap = pd_proj->GetAtomicMass();
627   G4int zt = pd_targ->GetAtomicNumber();
628   G4int at = pd_targ->GetAtomicMass();
629
630   //G4double rmax0 = 8.0;  // T.K dicide parameter value  // for low energy
631   G4double rmax0 = bmax + 4.0;
632   G4double rmax = std::sqrt( rmax0*rmax0 + b*b );
633
634   G4double ccoul = 0.001439767;
635   G4double pcca = 1.0 - double ( zp * zt ) * ccoul / eccm / rmax - ( b / rmax )*( b / rmax );
636
637   G4double pccf = std::sqrt( pcca );
638
639   //Fix for neutral particles
640   G4double aas1 = 0.0;
641   G4double bbs = 0.0;
642
643   if ( zp != 0 )
644   {
645      G4double aas = 2.0 * eccm * b / double ( zp * zt ) / ccoul;
646      bbs = 1.0 / std::sqrt ( 1.0 + aas*aas );
647      aas1 = ( 1.0 + aas * b / rmax ) * bbs;
648   }
649
650   G4double cost = 0.0;
651   G4double sint = 0.0;
652   G4double thet1 = 0.0;
653   G4double thet2 = 0.0;
654   if ( 1.0 - aas1*aas1 <= 0 || 1.0 - bbs*bbs <= 0.0 )   
655   {
656      cost = 1.0;
657      sint = 0.0;
658   } 
659   else 
660   {
661      G4double aat1 = aas1 / std::sqrt ( 1.0 - aas1*aas1 );
662      G4double aat2 = bbs / std::sqrt ( 1.0 - bbs*bbs );
663
664      thet1 = std::atan ( aat1 );
665      thet2 = std::atan ( aat2 );
666
667//    TK enter to else block 
668      G4double theta = thet1 - thet2;
669      cost = std::cos( theta );
670      sint = std::sin( theta );
671   }
672
673   G4double rzpr = -rmax * cost * ( mass_targ ) / ( mass_proj + mass_targ );
674   G4double rzta =  rmax * cost * ( mass_proj ) / ( mass_proj + mass_targ );
675
676   G4double rxpr = rmax / 2.0 * sint;
677
678   G4double rxta = -rxpr;
679
680
681   G4double pzpc = pzcc * (  cost * pccf + sint * b / rmax ); 
682   G4double pxpr = pzcc * ( -sint * pccf + cost * b / rmax ); 
683
684   G4double pztc = - pzpc;
685   G4double pxta = - pxpr;
686
687   G4double epc = std::sqrt ( pzpc*pzpc + pxpr*pxpr + mass_proj*mass_proj );
688   G4double etc = std::sqrt ( pztc*pztc + pxta*pxta + mass_targ*mass_targ );
689
690   G4double pzpr = pzpc;
691   G4double pzta = pztc;
692   G4double epr = epc;
693   G4double eta = etc;
694
695// CM -> NN
696   G4double gammacm = boostToCM.gamma();
697   //G4double betacm = -boostToCM.beta();
698   G4double betacm = boostToCM.z();
699   pzpr = pzpc + betacm * gammacm * ( gammacm / ( 1. + gammacm ) * pzpc * betacm + epc );
700   pzta = pztc + betacm * gammacm * ( gammacm / ( 1. + gammacm ) * pztc * betacm + etc );
701   epr = gammacm * ( epc + betacm * pzpc );
702   eta = gammacm * ( etc + betacm * pztc );
703
704   //G4double betpr = pzpr / epr;
705   //G4double betta = pzta / eta;
706
707   G4double gammpr = epr / ( mass_proj );
708   G4double gammta = eta / ( mass_targ );
709     
710   pzta = pzta / double ( at );
711   pxta = pxta / double ( at );
712     
713   pzpr = pzpr / double ( ap );
714   pxpr = pxpr / double ( ap );
715
716   G4double zeroz = 0.0; 
717
718   rzpr = rzpr -zeroz;
719   rzta = rzta -zeroz;
720
721   // Set results
722   coulomb_collision_gamma_proj = gammpr;
723   coulomb_collision_rx_proj = rxpr;
724   coulomb_collision_rz_proj = rzpr;
725   coulomb_collision_px_proj = pxpr;
726   coulomb_collision_pz_proj = pzpr;
727
728   coulomb_collision_gamma_targ = gammta;
729   coulomb_collision_rx_targ = rxta;
730   coulomb_collision_rz_targ = rzta;
731   coulomb_collision_px_targ = pxta;
732   coulomb_collision_pz_targ = pzta;
733
734}
735
736
737
738void G4QMDReaction::setEvaporationCh()
739{
740
741   if ( gem == true ) 
742      evaporation->SetGEMChannel();
743   else
744      evaporation->SetDefaultChannel();
745
746}
Note: See TracBrowser for help on using the repository browser.