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

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

update processes

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