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

Last change on this file since 830 was 819, checked in by garnier, 17 years ago

import all except CVS

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