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

Last change on this file since 1347 was 1347, checked in by garnier, 13 years ago

geant4 tag 9.4

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