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

Last change on this file since 1143 was 1055, checked in by garnier, 17 years ago

maj sur la beta de geant 4.9.3

File size: 24.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// 090122 Be8 -> Alpha + Alpha
34// 090331 Change member shenXS and genspaXS object to pointer
35//
36#include "G4QMDReaction.hh"
37#include "G4QMDNucleus.hh"
38#include "G4QMDGroundStateNucleus.hh"
39
40#include "G4NistManager.hh"
41
42G4QMDReaction::G4QMDReaction()
43: system ( NULL )
44, deltaT ( 1 ) // in fsec
45, maxTime ( 100 ) // will have maxTime-th time step
46, gem ( true )
47, frag ( false )
48{
49
50 //090331
51 shenXS = new G4IonsShenCrossSection();
52 //genspaXS = new G4GeneralSpaceNNCrossSection();
53 meanField = new G4QMDMeanField();
54 collision = new G4QMDCollision();
55
56 excitationHandler = new G4ExcitationHandler;
57 evaporation = new G4Evaporation;
58 excitationHandler->SetEvaporation( evaporation );
59 setEvaporationCh();
60}
61
62
63
64G4QMDReaction::~G4QMDReaction()
65{
66 delete evaporation;
67 delete excitationHandler;
68 delete collision;
69 delete meanField;
70}
71
72
73
74G4HadFinalState* G4QMDReaction::ApplyYourself( const G4HadProjectile & projectile , G4Nucleus & target )
75{
76 //G4cout << "G4QMDReaction::ApplyYourself" << G4endl;
77
78 theParticleChange.Clear();
79
80 system = new G4QMDSystem;
81
82 G4int proj_Z = 0;
83 G4int proj_A = 0;
84 G4ParticleDefinition* proj_pd = ( G4ParticleDefinition* ) projectile.GetDefinition();
85 if ( proj_pd->GetParticleType() == "nucleus" )
86 {
87 proj_Z = proj_pd->GetAtomicNumber();
88 proj_A = proj_pd->GetAtomicMass();
89 }
90 else
91 {
92 proj_Z = (int)( proj_pd->GetPDGCharge()/eplus );
93 proj_A = 1;
94 }
95 G4int targ_Z = int ( target.GetZ() + 0.5 );
96 G4int targ_A = int ( target.GetN() + 0.5 );
97 G4ParticleDefinition* targ_pd = G4ParticleTable::GetParticleTable()->GetIon( targ_Z , targ_A , 0.0 );
98
99
100 G4NistManager* nistMan = G4NistManager::Instance();
101// G4Element* G4NistManager::FindOrBuildElement( targ_Z );
102
103 const G4DynamicParticle* proj_dp = new G4DynamicParticle ( proj_pd , projectile.Get4Momentum() );
104 const G4Element* targ_ele = nistMan->FindOrBuildElement( targ_Z );
105 G4double aTemp = projectile.GetMaterial()->GetTemperature();
106
107 //090331
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 );
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 * std::sqrt ( 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 G4QMDGroundStateNucleus* proj(NULL);
167 if ( projectile.GetDefinition()->GetParticleType() == "nucleus"
168 || projectile.GetDefinition()->GetParticleName() == "proton"
169 || projectile.GetDefinition()->GetParticleName() == "neutron" )
170 {
171
172 proj_Z = proj_pd->GetAtomicNumber();
173 proj_A = proj_pd->GetAtomicMass();
174
175 proj = new G4QMDGroundStateNucleus( proj_Z , proj_A );
176 //proj->ShowParticipants();
177
178 }
179
180 meanField->SetSystem ( proj );
181 proj->SetTotalPotential( meanField->GetTotalPotential() );
182 proj->CalEnergyAndAngularMomentumInCM();
183
184// Target
185 G4int iz = int ( target.GetZ() );
186 G4int ia = int ( target.GetN() );
187
188 G4QMDGroundStateNucleus* 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( deltaT );
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 if ( frag == true )
389 if ( elasticLike_energy == false ) elastic = false;
390// InelasticLike with Collision
391 if ( elasticLike_energy == false && withCollision == true ) elastic = false; // ielst = 3
392
393 }
394
395 }
396 else
397 {
398
399// numberOfSecondary != 2
400 elastic = false;
401
402 }
403
404//071115
405 //G4cout << elastic << G4endl;
406 // if elastic is true try again from sampling of impact parameter
407
408 if ( elastic == true )
409 {
410 // delete this nucleues
411 for ( std::vector< G4QMDNucleus* >::iterator
412 it = nucleuses.begin() ; it != nucleuses.end() ; it++ )
413 {
414 delete *it;
415 }
416 nucleuses.clear();
417 }
418 }
419
420
421// Statical Decay Phase
422
423 for ( std::vector< G4QMDNucleus* >::iterator it
424 = nucleuses.begin() ; it != nucleuses.end() ; it++ )
425 {
426
427/*
428 std::cout << "G4QMDRESULT "
429 << (*it)->GetAtomicNumber()
430 << " "
431 << (*it)->GetMassNumber()
432 << " "
433 << (*it)->Get4Momentum()
434 << " "
435 << (*it)->Get4Momentum().vect()
436 << " "
437 << (*it)->Get4Momentum().restMass()
438 << " "
439 << (*it)->GetNuclearMass()/GeV
440 << std::endl;
441*/
442
443 meanField->SetNucleus ( *it );
444
445 if ( (*it)->GetAtomicNumber() == 0 // neutron cluster
446 || (*it)->GetAtomicNumber() == (*it)->GetMassNumber() ) // proton cluster
447 {
448 // push back system
449 for ( G4int i = 0 ; i < (*it)->GetTotalNumberOfParticipant() ; i++ )
450 {
451 G4QMDParticipant* aP = new G4QMDParticipant( ( (*it)->GetParticipant( i ) )->GetDefinition() , ( (*it)->GetParticipant( i ) )->GetMomentum() , ( (*it)->GetParticipant( i ) )->GetPosition() );
452 system->SetParticipant ( aP );
453 }
454 continue;
455 }
456
457 G4double nucleus_e = std::sqrt ( std::pow ( (*it)->GetNuclearMass()/GeV , 2 ) + std::pow ( (*it)->Get4Momentum().vect().mag() , 2 ) );
458 G4LorentzVector nucleus_p4CM ( (*it)->Get4Momentum().vect() , nucleus_e );
459
460// std::cout << "G4QMDRESULT nucleus deltaQ " << deltaQ << std::endl;
461
462 G4int ia = (*it)->GetMassNumber();
463 G4int iz = (*it)->GetAtomicNumber();
464
465 G4LorentzVector lv ( G4ThreeVector( 0.0 ) , (*it)->GetExcitationEnergy()*GeV + G4ParticleTable::GetParticleTable()->GetIonTable()->GetIonMass( iz , ia ) );
466
467 G4Fragment* aFragment = new G4Fragment( ia , iz , lv );
468
469 G4ReactionProductVector* rv;
470 rv = excitationHandler->BreakItUp( *aFragment );
471 G4bool notBreak = true;
472 for ( G4ReactionProductVector::iterator itt
473 = rv->begin() ; itt != rv->end() ; itt++ )
474 {
475
476 notBreak = false;
477 // Secondary from this nucleus (*it)
478 G4ParticleDefinition* pd = (*itt)->GetDefinition();
479
480 G4LorentzVector p4 ( (*itt)->GetMomentum()/GeV , (*itt)->GetTotalEnergy()/GeV ); //in nucleus(*it) rest system
481 G4LorentzVector p4_CM = CLHEP::boostOf( p4 , -nucleus_p4CM.findBoostToCM() ); // Back to CM
482 G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB ); // Back to LAB
483
484
485//090122
486 //theParticleChange.AddSecondary( dp );
487 if ( !( pd->GetAtomicNumber() == 4 && pd->GetAtomicMass() == 8 ) )
488 {
489 G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV );
490 theParticleChange.AddSecondary( dp );
491 }
492 else
493 {
494 //Be8 -> Alpha + Alpha + Q
495 G4ThreeVector randomized_direction( G4UniformRand() , G4UniformRand() , G4UniformRand() );
496 randomized_direction = randomized_direction.unit();
497 G4double q_decay = (*itt)->GetMass() - 2*G4Alpha::Alpha()->GetPDGMass();
498 G4double p_decay = std::sqrt ( std::pow(G4Alpha::Alpha()->GetPDGMass()+q_decay/2,2) - std::pow(G4Alpha::Alpha()->GetPDGMass() , 2 ) );
499 G4LorentzVector p4_a1 ( p_decay*randomized_direction , G4Alpha::Alpha()->GetPDGMass()+q_decay/2 ); //in Be8 rest system
500
501 G4LorentzVector p4_a1_Be8 = CLHEP::boostOf ( p4_a1/GeV , -p4.findBoostToCM() );
502 G4LorentzVector p4_a1_CM = CLHEP::boostOf ( p4_a1_Be8 , -nucleus_p4CM.findBoostToCM() );
503 G4LorentzVector p4_a1_LAB = CLHEP::boostOf ( p4_a1_CM , boostBackToLAB );
504
505 G4LorentzVector p4_a2 ( -p_decay*randomized_direction , G4Alpha::Alpha()->GetPDGMass()+q_decay/2 ); //in Be8 rest system
506
507 G4LorentzVector p4_a2_Be8 = CLHEP::boostOf ( p4_a2/GeV , -p4.findBoostToCM() );
508 G4LorentzVector p4_a2_CM = CLHEP::boostOf ( p4_a2_Be8 , -nucleus_p4CM.findBoostToCM() );
509 G4LorentzVector p4_a2_LAB = CLHEP::boostOf ( p4_a2_CM , boostBackToLAB );
510
511 G4DynamicParticle* dp1 = new G4DynamicParticle( G4Alpha::Alpha() , p4_a1_LAB*GeV );
512 G4DynamicParticle* dp2 = new G4DynamicParticle( G4Alpha::Alpha() , p4_a2_LAB*GeV );
513 theParticleChange.AddSecondary( dp1 );
514 theParticleChange.AddSecondary( dp2 );
515 }
516//090122
517
518/*
519 std::cout
520 << "Regist Secondary "
521 << (*itt)->GetDefinition()->GetParticleName()
522 << " "
523 << (*itt)->GetMomentum()/GeV
524 << " "
525 << (*itt)->GetKineticEnergy()/GeV
526 << " "
527 << (*itt)->GetMass()/GeV
528 << " "
529 << (*itt)->GetTotalEnergy()/GeV
530 << " "
531 << (*itt)->GetTotalEnergy()/GeV * (*itt)->GetTotalEnergy()/GeV
532 - (*itt)->GetMomentum()/GeV * (*itt)->GetMomentum()/GeV
533 << " "
534 << nucleus_p4CM.findBoostToCM()
535 << " "
536 << p4
537 << " "
538 << p4_CM
539 << " "
540 << p4_LAB
541 << std::endl;
542*/
543
544 }
545 if ( notBreak == true )
546 {
547
548 G4ParticleDefinition* pd = G4ParticleTable::GetParticleTable()->GetIon( (*it)->GetAtomicNumber() , (*it)->GetMassNumber(), (*it)->GetExcitationEnergy()*GeV );
549 G4LorentzVector p4_CM = nucleus_p4CM;
550 G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB ); // Back to LAB
551 G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV );
552 theParticleChange.AddSecondary( dp );
553
554 }
555
556 for ( G4ReactionProductVector::iterator itt
557 = rv->begin() ; itt != rv->end() ; itt++ )
558 {
559 delete *itt;
560 }
561 delete rv;
562
563 delete aFragment;
564
565 }
566
567
568
569 for ( G4int i = 0 ; i < system->GetTotalNumberOfParticipant() ; i++ )
570 {
571
572 // Secondary particles
573
574 G4ParticleDefinition* pd = system->GetParticipant( i )->GetDefinition();
575 G4LorentzVector p4_CM = system->GetParticipant( i )->Get4Momentum();
576 G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB );
577 G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV );
578 theParticleChange.AddSecondary( dp );
579
580/*
581 G4cout << "G4QMDRESULT "
582 << "r" << i << " " << system->GetParticipant ( i ) -> GetPosition() << " "
583 << "p" << i << " " << system->GetParticipant ( i ) -> Get4Momentum()
584 << G4endl;
585*/
586
587 }
588
589 for ( std::vector< G4QMDNucleus* >::iterator it
590 = nucleuses.begin() ; it != nucleuses.end() ; it++ )
591 {
592 delete *it; // delete nulceuse
593 }
594 nucleuses.clear();
595
596 system->Clear();
597 delete system;
598
599 theParticleChange.SetStatusChange( stopAndKill );
600
601 return &theParticleChange;
602
603}
604
605
606
607void G4QMDReaction::calcOffSetOfCollision( G4double b ,
608G4ParticleDefinition* pd_proj ,
609G4ParticleDefinition* pd_targ ,
610G4double ptot , G4double etot , G4double bmax , G4ThreeVector boostToCM )
611{
612 G4double mass_proj = pd_proj->GetPDGMass()/GeV;
613 G4double mass_targ = pd_targ->GetPDGMass()/GeV;
614
615 G4double stot = std::sqrt ( etot*etot - ptot*ptot );
616
617 G4double pstt = std::sqrt ( ( stot*stot - ( mass_proj + mass_targ ) * ( mass_proj + mass_targ )
618 ) * ( stot*stot - ( mass_proj - mass_targ ) * ( mass_proj - mass_targ ) ) )
619 / ( 2.0 * stot );
620
621 G4double pzcc = pstt;
622 G4double eccm = stot - ( mass_proj + mass_targ );
623
624 G4int zp = pd_proj->GetAtomicNumber();
625 G4int ap = pd_proj->GetAtomicMass();
626 G4int zt = pd_targ->GetAtomicNumber();
627 G4int at = pd_targ->GetAtomicMass();
628
629 //G4double rmax0 = 8.0; // T.K dicide parameter value // for low energy
630 G4double rmax0 = bmax + 4.0;
631 G4double rmax = std::sqrt( rmax0*rmax0 + b*b );
632
633 G4double ccoul = 0.001439767;
634 G4double pcca = 1.0 - double ( zp * zt ) * ccoul / eccm / rmax - ( b / rmax )*( b / rmax );
635
636 G4double pccf = std::sqrt( pcca );
637
638 G4double aas = 2.0 * eccm * b / double ( zp * zt ) / ccoul;
639 G4double bbs = 1.0 / std::sqrt ( 1.0 + aas*aas );
640 G4double aas1 = ( 1.0 + aas * b / rmax ) * bbs;
641
642 G4double cost = 0.0;
643 G4double sint = 0.0;
644 G4double thet1 = 0.0;
645 G4double thet2 = 0.0;
646 if ( 1.0 - aas1*aas1 <= 0 || 1.0 - bbs*bbs <= 0.0 )
647 {
648 cost = 1.0;
649 sint = 0.0;
650 }
651 else
652 {
653 G4double aat1 = aas1 / std::sqrt ( 1.0 - aas1*aas1 );
654 G4double aat2 = bbs / std::sqrt ( 1.0 - bbs*bbs );
655
656 thet1 = std::atan ( aat1 );
657 thet2 = std::atan ( aat2 );
658
659// TK enter to else block
660 G4double theta = thet1 - thet2;
661 cost = std::cos( theta );
662 sint = std::sin( theta );
663 }
664
665 G4double rzpr = -rmax * cost * ( mass_targ ) / ( mass_proj + mass_targ );
666 G4double rzta = rmax * cost * ( mass_proj ) / ( mass_proj + mass_targ );
667
668 G4double rxpr = rmax / 2.0 * sint;
669
670 G4double rxta = -rxpr;
671
672
673 G4double pzpc = pzcc * ( cost * pccf + sint * b / rmax );
674 G4double pxpr = pzcc * ( -sint * pccf + cost * b / rmax );
675
676 G4double pztc = - pzpc;
677 G4double pxta = - pxpr;
678
679 G4double epc = std::sqrt ( pzpc*pzpc + pxpr*pxpr + mass_proj*mass_proj );
680 G4double etc = std::sqrt ( pztc*pztc + pxta*pxta + mass_targ*mass_targ );
681
682 G4double pzpr = pzpc;
683 G4double pzta = pztc;
684 G4double epr = epc;
685 G4double eta = etc;
686
687// CM -> NN
688 G4double gammacm = boostToCM.gamma();
689 //G4double betacm = -boostToCM.beta();
690 G4double betacm = boostToCM.z();
691 pzpr = pzpc + betacm * gammacm * ( gammacm / ( 1. + gammacm ) * pzpc * betacm + epc );
692 pzta = pztc + betacm * gammacm * ( gammacm / ( 1. + gammacm ) * pztc * betacm + etc );
693 epr = gammacm * ( epc + betacm * pzpc );
694 eta = gammacm * ( etc + betacm * pztc );
695
696 //G4double betpr = pzpr / epr;
697 //G4double betta = pzta / eta;
698
699 G4double gammpr = epr / ( mass_proj );
700 G4double gammta = eta / ( mass_targ );
701
702 pzta = pzta / double ( at );
703 pxta = pxta / double ( at );
704
705 pzpr = pzpr / double ( ap );
706 pxpr = pxpr / double ( ap );
707
708 G4double zeroz = 0.0;
709
710 rzpr = rzpr -zeroz;
711 rzta = rzta -zeroz;
712
713 // Set results
714 coulomb_collision_gamma_proj = gammpr;
715 coulomb_collision_rx_proj = rxpr;
716 coulomb_collision_rz_proj = rzpr;
717 coulomb_collision_px_proj = pxpr;
718 coulomb_collision_pz_proj = pzpr;
719
720 coulomb_collision_gamma_targ = gammta;
721 coulomb_collision_rx_targ = rxta;
722 coulomb_collision_rz_targ = rzta;
723 coulomb_collision_px_targ = pxta;
724 coulomb_collision_pz_targ = pzta;
725
726}
727
728
729
730void G4QMDReaction::setEvaporationCh()
731{
732
733 if ( gem == true )
734 evaporation->SetGEMChannel();
735 else
736 evaporation->SetDefaultChannel();
737
738}
Note: See TracBrowser for help on using the repository browser.