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

Last change on this file since 1330 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

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