Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QMDReaction.cc
Go to the documentation of this file.
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"
40#include "G4Pow.hh"
42#include "G4SystemOfUnits.hh"
43#include "G4NistManager.hh"
44
46#include "G4BGGPionElasticXS.hh"
52
53
55: G4HadronicInteraction("QMDModel")
56, system ( NULL )
57, deltaT ( 1 ) // in fsec (c=1)
58, maxTime ( 100 ) // will have maxTime-th time step
59, envelopF ( 1.05 ) // 10% for Peripheral reactions
60, gem ( true )
61, frag ( false )
62, secID( -1 )
63{
65 pipElNucXS = new G4BGGPionElasticXS(G4PionPlus::PionPlus() );
66 pipElNucXS->BuildPhysicsTable(*(G4PionPlus::PionPlus() ) );
67
68 pimElNucXS = new G4BGGPionElasticXS(G4PionMinus::PionMinus() );
69 pimElNucXS->BuildPhysicsTable(*(G4PionMinus::PionMinus() ) );
70
71 pipInelNucXS = new G4BGGPionInelasticXS(G4PionPlus::PionPlus() );
72 pipInelNucXS->BuildPhysicsTable(*(G4PionPlus::PionPlus() ) );
73
74 pimInelNucXS = new G4BGGPionInelasticXS(G4PionMinus::PionMinus() );
75 pimInelNucXS->BuildPhysicsTable(*(G4PionMinus::PionMinus() ) );
76
77 meanField = new G4QMDMeanField();
78 collision = new G4QMDCollision();
79
80 excitationHandler = new G4ExcitationHandler;
81 excitationHandler->SetDeexChannelsType( fCombined );
82 //fEvaporation - 8 default channels
83 //fCombined - 8 default + 60 GEM
84 //fGEM - 2 default + 66 GEM
85 evaporation = new G4Evaporation;
86 excitationHandler->SetEvaporation( evaporation );
87 setEvaporationCh();
88
89 coulomb_collision_gamma_proj = 0.0;
90 coulomb_collision_rx_proj = 0.0;
91 coulomb_collision_rz_proj = 0.0;
92 coulomb_collision_px_proj = 0.0;
93 coulomb_collision_pz_proj = 0.0;
94
95 coulomb_collision_gamma_targ = 0.0;
96 coulomb_collision_rx_targ = 0.0;
97 coulomb_collision_rz_targ = 0.0;
98 coulomb_collision_px_targ = 0.0;
99 coulomb_collision_pz_targ = 0.0;
100
101 secID = G4PhysicsModelCatalog::GetModelID( "model_QMDModel" );
102}
103
104
106{
107 delete evaporation;
108 delete excitationHandler;
109 delete collision;
110 delete meanField;
111}
112
113
115{
116 //G4cout << "G4QMDReaction::ApplyYourself" << G4endl;
117
119
120 system = new G4QMDSystem;
121
122 G4int proj_Z = 0;
123 G4int proj_A = 0;
124 const G4ParticleDefinition* proj_pd = ( const G4ParticleDefinition* ) projectile.GetDefinition();
125 if ( proj_pd->GetParticleType() == "nucleus" )
126 {
127 proj_Z = proj_pd->GetAtomicNumber();
128 proj_A = proj_pd->GetAtomicMass();
129 }
130 else
131 {
132 proj_Z = (int)( proj_pd->GetPDGCharge()/eplus );
133 proj_A = 1;
134 }
135 //G4int targ_Z = int ( target.GetZ() + 0.5 );
136 //G4int targ_A = int ( target.GetN() + 0.5 );
137 //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this)
138 G4int targ_Z = target.GetZ_asInt();
139 G4int targ_A = target.GetA_asInt();
140 const G4ParticleDefinition* targ_pd = G4IonTable::GetIonTable()->GetIon( targ_Z , targ_A , 0.0 );
141
142
143 //G4NistManager* nistMan = G4NistManager::Instance();
144// G4Element* G4NistManager::FindOrBuildElement( targ_Z );
145
146 const G4DynamicParticle* proj_dp = new G4DynamicParticle ( proj_pd , projectile.Get4Momentum() );
147 //const G4Element* targ_ele = nistMan->FindOrBuildElement( targ_Z );
148 //G4double aTemp = projectile.GetMaterial()->GetTemperature();
149
150 // Glauber-Gribov nucleus-nucleus cross section does not have GetIsoCrossSection,
151 // therefore call GetElementCrossSection instead.
152 //G4double xs_0 = theXS->GetIsoCrossSection ( proj_dp , targ_Z , targ_A );
153 G4double xs_0 = theXS->GetElementCrossSection( proj_dp , targ_Z , projectile.GetMaterial() );
154
155 // When the projectile is a pion
156 if (proj_pd == G4PionPlus::PionPlus() ) {
157 xs_0 = pipElNucXS->GetElementCrossSection(proj_dp, targ_Z, projectile.GetMaterial() ) +
158 pipInelNucXS->GetElementCrossSection(proj_dp, targ_Z, projectile.GetMaterial() );
159 } else if (proj_pd == G4PionMinus::PionMinus() ) {
160 xs_0 = pimElNucXS->GetElementCrossSection(proj_dp, targ_Z, projectile.GetMaterial() ) +
161 pimInelNucXS->GetElementCrossSection(proj_dp, targ_Z, projectile.GetMaterial() );
162 }
163
164 //G4double xs_0 = genspaXS->GetCrossSection ( proj_dp , targ_ele , aTemp );
165 //G4double xs_0 = theXS->GetCrossSection ( proj_dp , targ_ele , aTemp );
166 //110822
167
168 G4double bmax_0 = std::sqrt( xs_0 / pi );
169 //std::cout << "bmax_0 in fm (fermi) " << bmax_0/fermi << std::endl;
170
171 //delete proj_dp;
172
173 G4bool elastic = true;
174
175 std::vector< G4QMDNucleus* > nucleuses; // Secondary nuceluses
176 G4ThreeVector boostToReac; // ReactionSystem (CM or NN);
177 G4ThreeVector boostBackToLAB; // Reaction System to LAB;
178
179 G4LorentzVector targ4p( G4ThreeVector( 0.0 ) , targ_pd->GetPDGMass()/GeV );
180 G4ThreeVector boostLABtoCM = targ4p.findBoostToCM( proj_dp->Get4Momentum()/GeV ); // CM of target and proj;
181
182 G4double p1 = proj_dp->GetMomentum().mag()/GeV/proj_A;
183 G4double m1 = proj_dp->GetDefinition()->GetPDGMass()/GeV/proj_A;
184 G4double e1 = std::sqrt( p1*p1 + m1*m1 );
185 G4double e2 = targ_pd->GetPDGMass()/GeV/targ_A;
186 G4double beta_nn = -p1 / ( e1+e2 );
187
188 G4ThreeVector boostLABtoNN ( 0. , 0. , beta_nn ); // CM of NN;
189
190 G4double beta_nncm = ( - boostLABtoCM.beta() + boostLABtoNN.beta() ) / ( 1 - boostLABtoCM.beta() * boostLABtoNN.beta() ) ;
191
192 //std::cout << targ4p << std::endl;
193 //std::cout << proj_dp->Get4Momentum()<< std::endl;
194 //std::cout << beta_nncm << std::endl;
195 G4ThreeVector boostNNtoCM( 0. , 0. , beta_nncm ); //
196 G4ThreeVector boostCMtoNN( 0. , 0. , -beta_nncm ); //
197
198 boostToReac = boostLABtoNN;
199 boostBackToLAB = -boostLABtoNN;
200
201 delete proj_dp;
202
203 G4int icounter = 0;
204 G4int icounter_max = 1024;
205 while ( elastic ) // Loop checking, 11.03.2015, T. Koi
206 {
207 icounter++;
208 if ( icounter > icounter_max ) {
209 G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
210 break;
211 }
212
213// impact parameter
214 //G4double bmax = 1.05*(bmax_0/fermi); // 10% for Peripheral reactions
215 G4double bmax = envelopF*(bmax_0/fermi);
216 G4double b = bmax * std::sqrt ( G4UniformRand() );
217//071112
218 //G4double b = 0;
219 //G4double b = bmax;
220 //G4double b = bmax/1.05 * 0.7 * G4UniformRand();
221
222 //G4cout << "G4QMDRESULT bmax_0 = " << bmax_0/fermi << " fm, bmax = " << bmax << " fm , b = " << b << " fm " << G4endl;
223
224 G4double plab = projectile.GetTotalMomentum()/GeV;
225 G4double elab = ( projectile.GetKineticEnergy() + proj_pd->GetPDGMass() + targ_pd->GetPDGMass() )/GeV;
226
227 calcOffSetOfCollision( b , proj_pd , targ_pd , plab , elab , bmax , boostCMtoNN );
228
229// Projectile
230 G4LorentzVector proj4pLAB = projectile.Get4Momentum()/GeV;
231
232 G4QMDGroundStateNucleus* proj(NULL);
233 if ( projectile.GetDefinition()->GetParticleType() == "nucleus"
234 || projectile.GetDefinition()->GetParticleName() == "proton"
235 || projectile.GetDefinition()->GetParticleName() == "neutron" )
236 {
237
238 proj_Z = proj_pd->GetAtomicNumber();
239 proj_A = proj_pd->GetAtomicMass();
240
241 proj = new G4QMDGroundStateNucleus( proj_Z , proj_A );
242 //proj->ShowParticipants();
243
244
245 meanField->SetSystem ( proj );
246 proj->SetTotalPotential( meanField->GetTotalPotential() );
248
249 }
250
251// Target
252 //G4int iz = int ( target.GetZ() );
253 //G4int ia = int ( target.GetN() );
254 //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this)
255 G4int iz = int ( target.GetZ_asInt() );
256 G4int ia = int ( target.GetA_asInt() );
257
259
260 meanField->SetSystem (targ );
261 targ->SetTotalPotential( meanField->GetTotalPotential() );
263
264 //G4LorentzVector targ4p( G4ThreeVector( 0.0 ) , targ->GetNuclearMass()/GeV );
265// Boost Vector to CM
266 //boostToCM = targ4p.findBoostToCM( proj4pLAB );
267
268// Target
269 for ( G4int i = 0 ; i < targ->GetTotalNumberOfParticipant() ; i++ )
270 {
271
272 G4ThreeVector p0 = targ->GetParticipant( i )->GetMomentum();
273 G4ThreeVector r0 = targ->GetParticipant( i )->GetPosition();
274
275 G4ThreeVector p ( p0.x() + coulomb_collision_px_targ
276 , p0.y()
277 , p0.z() * coulomb_collision_gamma_targ + coulomb_collision_pz_targ );
278
279 G4ThreeVector r ( r0.x() + coulomb_collision_rx_targ
280 , r0.y()
281 , r0.z() / coulomb_collision_gamma_targ + coulomb_collision_rz_targ );
282
283 system->SetParticipant( new G4QMDParticipant( targ->GetParticipant( i )->GetDefinition() , p , r ) );
284 system->GetParticipant( i )->SetTarget();
285
286 }
287
288 G4LorentzVector proj4pCM = CLHEP::boostOf ( proj4pLAB , boostToReac );
289 G4LorentzVector targ4pCM = CLHEP::boostOf ( targ4p , boostToReac );
290
291// Projectile
292 //G4cout << "proj : " << proj << G4endl;
293 //if ( proj != NULL )
294 if ( proj_A != 1 )
295 {
296
297// projectile is nucleus
298
299 for ( G4int i = 0 ; i < proj->GetTotalNumberOfParticipant() ; i++ )
300 {
301
302 G4ThreeVector p0 = proj->GetParticipant( i )->GetMomentum();
303 G4ThreeVector r0 = proj->GetParticipant( i )->GetPosition();
304
305 G4ThreeVector p ( p0.x() + coulomb_collision_px_proj
306 , p0.y()
307 , p0.z() * coulomb_collision_gamma_proj + coulomb_collision_pz_proj );
308
309 G4ThreeVector r ( r0.x() + coulomb_collision_rx_proj
310 , r0.y()
311 , r0.z() / coulomb_collision_gamma_proj + coulomb_collision_rz_proj );
312
313 system->SetParticipant( new G4QMDParticipant( proj->GetParticipant( i )->GetDefinition() , p , r ) );
315 }
316
317 }
318 else
319 {
320
321// projectile is particle
322
323 // avoid multiple set in "elastic" loop
324 //G4cout << "system Total Participants : " << system->GetTotalNumberOfParticipant() << ", target : " << targ->GetTotalNumberOfParticipant() << G4endl;
326 {
327
329
330 G4ThreeVector p0( 0 );
331 G4ThreeVector r0( 0 );
332
333 G4ThreeVector p ( p0.x() + coulomb_collision_px_proj
334 , p0.y()
335 , p0.z() * coulomb_collision_gamma_proj + coulomb_collision_pz_proj );
336
337 G4ThreeVector r ( r0.x() + coulomb_collision_rx_proj
338 , r0.y()
339 , r0.z() / coulomb_collision_gamma_proj + coulomb_collision_rz_proj );
340
341 system->SetParticipant( new G4QMDParticipant( (G4ParticleDefinition*)projectile.GetDefinition() , p , r ) );
342 // This is not important becase only 1 projectile particle.
343 system->GetParticipant ( i )->SetProjectile();
344 }
345
346 }
347 //system->ShowParticipants();
348
349 delete targ;
350 delete proj;
351
352 meanField->SetSystem ( system );
353 collision->SetMeanField ( meanField );
354
355// Time Evolution
356 //std::cout << "Start time evolution " << std::endl;
357 //system->ShowParticipants();
358 for ( G4int i = 0 ; i < maxTime ; i++ )
359 {
360 //G4cout << " do Paropagate " << i << " th time step. " << G4endl;
361 meanField->DoPropagation( deltaT );
362 //system->ShowParticipants();
363 collision->CalKinematicsOfBinaryCollisions( deltaT );
364
365 if ( i / 10 * 10 == i )
366 {
367 //G4cout << i << " th time step. " << G4endl;
368 //system->ShowParticipants();
369 }
370 //system->ShowParticipants();
371 }
372 //system->ShowParticipants();
373
374
375 //std::cout << "Doing Cluster Judgment " << std::endl;
376
377 nucleuses = meanField->DoClusterJudgment();
378
379// Elastic Judgment
380
381 G4int numberOfSecondary = int ( nucleuses.size() ) + system->GetTotalNumberOfParticipant();
382
383 G4int sec_a_Z = 0;
384 G4int sec_a_A = 0;
385 const G4ParticleDefinition* sec_a_pd = NULL;
386 G4int sec_b_Z = 0;
387 G4int sec_b_A = 0;
388 const G4ParticleDefinition* sec_b_pd = NULL;
389
390 if ( numberOfSecondary == 2 )
391 {
392
393 G4bool elasticLike_system = false;
394 if ( nucleuses.size() == 2 )
395 {
396
397 sec_a_Z = nucleuses[0]->GetAtomicNumber();
398 sec_a_A = nucleuses[0]->GetMassNumber();
399 sec_b_Z = nucleuses[1]->GetAtomicNumber();
400 sec_b_A = nucleuses[1]->GetMassNumber();
401
402 if ( ( sec_a_Z == proj_Z && sec_a_A == proj_A && sec_b_Z == targ_Z && sec_b_A == targ_A )
403 || ( sec_a_Z == targ_Z && sec_a_A == targ_A && sec_b_Z == proj_Z && sec_b_A == proj_A ) )
404 {
405 elasticLike_system = true;
406 }
407
408 }
409 else if ( nucleuses.size() == 1 )
410 {
411
412 sec_a_Z = nucleuses[0]->GetAtomicNumber();
413 sec_a_A = nucleuses[0]->GetMassNumber();
414 sec_b_pd = system->GetParticipant( 0 )->GetDefinition();
415
416 if ( ( sec_a_Z == proj_Z && sec_a_A == proj_A && sec_b_pd == targ_pd )
417 || ( sec_a_Z == targ_Z && sec_a_A == targ_A && sec_b_pd == proj_pd ) )
418 {
419 elasticLike_system = true;
420 }
421
422 }
423 else
424 {
425
426 sec_a_pd = system->GetParticipant( 0 )->GetDefinition();
427 sec_b_pd = system->GetParticipant( 1 )->GetDefinition();
428
429 if ( ( sec_a_pd == proj_pd && sec_b_pd == targ_pd )
430 || ( sec_a_pd == targ_pd && sec_b_pd == proj_pd ) )
431 {
432 elasticLike_system = true;
433 }
434
435 }
436
437 if ( elasticLike_system == true )
438 {
439
440 G4bool elasticLike_energy = true;
441// Cal ExcitationEnergy
442 for ( G4int i = 0 ; i < int ( nucleuses.size() ) ; i++ )
443 {
444
445 //meanField->SetSystem( nucleuses[i] );
446 meanField->SetNucleus( nucleuses[i] );
447 //nucleuses[i]->SetTotalPotential( meanField->GetTotalPotential() );
448 //nucleuses[i]->CalEnergyAndAngularMomentumInCM();
449
450 if ( nucleuses[i]->GetExcitationEnergy()*GeV > 1.0*MeV ) elasticLike_energy = false;
451
452 }
453
454// Check Collision
455 G4bool withCollision = true;
456 if ( system->GetNOCollision() == 0 ) withCollision = false;
457
458// Final judegement for Inelasitc or Elastic;
459//
460// ElasticLike without Collision
461 //if ( elasticLike_energy == true && withCollision == false ) elastic = true; // ielst = 0
462// ElasticLike with Collision
463 //if ( elasticLike_energy == true && withCollision == true ) elastic = true; // ielst = 1
464// InelasticLike without Collision
465 //if ( elasticLike_energy == false ) elastic = false; // ielst = 2
466 if ( frag == true )
467 if ( elasticLike_energy == false ) elastic = false;
468// InelasticLike with Collision
469 if ( elasticLike_energy == false && withCollision == true ) elastic = false; // ielst = 3
470
471 }
472
473 }
474 else
475 {
476
477// numberOfSecondary != 2
478 elastic = false;
479
480 }
481
482//071115
483 //G4cout << elastic << G4endl;
484 // if elastic is true try again from sampling of impact parameter
485
486 if ( elastic == true )
487 {
488 // delete this nucleues
489 for ( std::vector< G4QMDNucleus* >::iterator
490 it = nucleuses.begin() ; it != nucleuses.end() ; it++ )
491 {
492 delete *it;
493 }
494 nucleuses.clear();
495 }
496
497 }
498
499
500// Statical Decay Phase
501
502 for ( std::vector< G4QMDNucleus* >::iterator it
503 = nucleuses.begin() ; it != nucleuses.end() ; it++ )
504 {
505
506/*
507 G4cout << "G4QMDRESULT "
508 << (*it)->GetAtomicNumber()
509 << " "
510 << (*it)->GetMassNumber()
511 << " "
512 << (*it)->Get4Momentum()
513 << " "
514 << (*it)->Get4Momentum().vect()
515 << " "
516 << (*it)->Get4Momentum().restMass()
517 << " "
518 << (*it)->GetNuclearMass()/GeV
519 << G4endl;
520*/
521
522 meanField->SetNucleus ( *it );
523
524 if ( (*it)->GetAtomicNumber() == 0 // neutron cluster
525 || (*it)->GetAtomicNumber() == (*it)->GetMassNumber() ) // proton cluster
526 {
527 // push back system
528 for ( G4int i = 0 ; i < (*it)->GetTotalNumberOfParticipant() ; i++ )
529 {
530 G4QMDParticipant* aP = new G4QMDParticipant( ( (*it)->GetParticipant( i ) )->GetDefinition() , ( (*it)->GetParticipant( i ) )->GetMomentum() , ( (*it)->GetParticipant( i ) )->GetPosition() );
531 system->SetParticipant ( aP );
532 }
533 continue;
534 }
535
536 G4double nucleus_e = std::sqrt ( G4Pow::GetInstance()->powN ( (*it)->GetNuclearMass()/GeV , 2 ) + G4Pow::GetInstance()->powN ( (*it)->Get4Momentum().vect().mag() , 2 ) );
537 G4LorentzVector nucleus_p4CM ( (*it)->Get4Momentum().vect() , nucleus_e );
538
539// std::cout << "G4QMDRESULT nucleus deltaQ " << deltaQ << std::endl;
540
541 G4int ia = (*it)->GetMassNumber();
542 G4int iz = (*it)->GetAtomicNumber();
543
544 G4LorentzVector lv ( G4ThreeVector( 0.0 ) , (*it)->GetExcitationEnergy()*GeV + G4IonTable::GetIonTable()->GetIonMass( iz , ia ) );
545
546 G4Fragment* aFragment = new G4Fragment( ia , iz , lv );
547
549 rv = excitationHandler->BreakItUp( *aFragment );
550 G4bool notBreak = true;
551 for ( G4ReactionProductVector::iterator itt
552 = rv->begin() ; itt != rv->end() ; itt++ )
553 {
554
555 notBreak = false;
556 // Secondary from this nucleus (*it)
557 const G4ParticleDefinition* pd = (*itt)->GetDefinition();
558
559 G4LorentzVector p4 ( (*itt)->GetMomentum()/GeV , (*itt)->GetTotalEnergy()/GeV ); //in nucleus(*it) rest system
560 G4LorentzVector p4_CM = CLHEP::boostOf( p4 , -nucleus_p4CM.findBoostToCM() ); // Back to CM
561 G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB ); // Back to LAB
562
563
564//090122
565 //theParticleChange.AddSecondary( dp );
566 if ( !( pd->GetAtomicNumber() == 4 && pd->GetAtomicMass() == 8 ) )
567 {
568 //G4cout << "pd out of notBreak loop : " << pd->GetParticleName() << G4endl;
569 G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV );
571 }
572 else
573 {
574 //Be8 -> Alpha + Alpha + Q
575 G4ThreeVector randomized_direction( G4UniformRand() , G4UniformRand() , G4UniformRand() );
576 randomized_direction = randomized_direction.unit();
577 G4double q_decay = (*itt)->GetMass() - 2*G4Alpha::Alpha()->GetPDGMass();
578 G4double p_decay = std::sqrt ( G4Pow::GetInstance()->powN(G4Alpha::Alpha()->GetPDGMass()+q_decay/2,2) - G4Pow::GetInstance()->powN(G4Alpha::Alpha()->GetPDGMass() , 2 ) );
579 G4LorentzVector p4_a1 ( p_decay*randomized_direction , G4Alpha::Alpha()->GetPDGMass()+q_decay/2 ); //in Be8 rest system
580
581 G4LorentzVector p4_a1_Be8 = CLHEP::boostOf ( p4_a1/GeV , -p4.findBoostToCM() );
582 G4LorentzVector p4_a1_CM = CLHEP::boostOf ( p4_a1_Be8 , -nucleus_p4CM.findBoostToCM() );
583 G4LorentzVector p4_a1_LAB = CLHEP::boostOf ( p4_a1_CM , boostBackToLAB );
584
585 G4LorentzVector p4_a2 ( -p_decay*randomized_direction , G4Alpha::Alpha()->GetPDGMass()+q_decay/2 ); //in Be8 rest system
586
587 G4LorentzVector p4_a2_Be8 = CLHEP::boostOf ( p4_a2/GeV , -p4.findBoostToCM() );
588 G4LorentzVector p4_a2_CM = CLHEP::boostOf ( p4_a2_Be8 , -nucleus_p4CM.findBoostToCM() );
589 G4LorentzVector p4_a2_LAB = CLHEP::boostOf ( p4_a2_CM , boostBackToLAB );
590
591 G4DynamicParticle* dp1 = new G4DynamicParticle( G4Alpha::Alpha() , p4_a1_LAB*GeV );
592 G4DynamicParticle* dp2 = new G4DynamicParticle( G4Alpha::Alpha() , p4_a2_LAB*GeV );
595 }
596//090122
597
598/*
599 G4cout
600 << "Regist Secondary "
601 << (*itt)->GetDefinition()->GetParticleName()
602 << " "
603 << (*itt)->GetMomentum()/GeV
604 << " "
605 << (*itt)->GetKineticEnergy()/GeV
606 << " "
607 << (*itt)->GetMass()/GeV
608 << " "
609 << (*itt)->GetTotalEnergy()/GeV
610 << " "
611 << (*itt)->GetTotalEnergy()/GeV * (*itt)->GetTotalEnergy()/GeV
612 - (*itt)->GetMomentum()/GeV * (*itt)->GetMomentum()/GeV
613 << " "
614 << nucleus_p4CM.findBoostToCM()
615 << " "
616 << p4
617 << " "
618 << p4_CM
619 << " "
620 << p4_LAB
621 << G4endl;
622*/
623
624 }
625 if ( notBreak == true )
626 {
627
628 const G4ParticleDefinition* pd = G4IonTable::GetIonTable()->GetIon( (*it)->GetAtomicNumber() , (*it)->GetMassNumber(), (*it)->GetExcitationEnergy()*GeV );
629 //G4cout << "pd in notBreak loop : " << pd->GetParticleName() << G4endl;
630 G4LorentzVector p4_CM = nucleus_p4CM;
631 G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB ); // Back to LAB
632 G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV );
634
635 }
636
637 for ( G4ReactionProductVector::iterator itt
638 = rv->begin() ; itt != rv->end() ; itt++ )
639 {
640 delete *itt;
641 }
642 delete rv;
643
644 delete aFragment;
645
646 }
647
648
649
650 for ( G4int i = 0 ; i < system->GetTotalNumberOfParticipant() ; i++ )
651 {
652 // Secondary particles
653
654 const G4ParticleDefinition* pd = system->GetParticipant( i )->GetDefinition();
655 G4LorentzVector p4_CM = system->GetParticipant( i )->Get4Momentum();
656 G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB );
657 G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV );
659 //G4cout << "In the last theParticleChange loop : " << pd->GetParticleName() << G4endl;
660
661/*
662 G4cout << "G4QMDRESULT "
663 << "r" << i << " " << system->GetParticipant ( i ) -> GetPosition() << " "
664 << "p" << i << " " << system->GetParticipant ( i ) -> Get4Momentum()
665 << G4endl;
666*/
667
668 }
669
670 for ( std::vector< G4QMDNucleus* >::iterator it
671 = nucleuses.begin() ; it != nucleuses.end() ; it++ )
672 {
673 delete *it; // delete nulceuse
674 }
675 nucleuses.clear();
676
677 system->Clear();
678 delete system;
679
681
682 for (G4int i = 0; i < G4int(theParticleChange.GetNumberOfSecondaries() ); i++)
683 {
684 //G4cout << "Particle : " << theParticleChange.GetSecondary(i)->GetParticle()->GetParticleDefinition()->GetParticleName() << G4endl;
685 //G4cout << "KEnergy : " << theParticleChange.GetSecondary(i)->GetParticle()->GetKineticEnergy() << G4endl;
686 //G4cout << "modelID : " << theParticleChange.GetSecondary(i)->GetCreatorModelID() << G4endl;
688 }
689
690 return &theParticleChange;
691
692}
693
694
695
696void G4QMDReaction::calcOffSetOfCollision( G4double b ,
697const G4ParticleDefinition* pd_proj ,
698const G4ParticleDefinition* pd_targ ,
699G4double ptot , G4double etot , G4double bmax , G4ThreeVector boostToCM )
700{
701
702 G4double mass_proj = pd_proj->GetPDGMass()/GeV;
703 G4double mass_targ = pd_targ->GetPDGMass()/GeV;
704
705 G4double stot = std::sqrt ( etot*etot - ptot*ptot );
706
707 G4double pstt = std::sqrt ( ( stot*stot - ( mass_proj + mass_targ ) * ( mass_proj + mass_targ )
708 ) * ( stot*stot - ( mass_proj - mass_targ ) * ( mass_proj - mass_targ ) ) )
709 / ( 2.0 * stot );
710
711 G4double pzcc = pstt;
712 G4double eccm = stot - ( mass_proj + mass_targ );
713
714 G4int zp = 1;
715 G4int ap = 1;
716 if ( pd_proj->GetParticleType() == "nucleus" )
717 {
718 zp = pd_proj->GetAtomicNumber();
719 ap = pd_proj->GetAtomicMass();
720 }
721 else
722 {
723 // proton, neutron, mesons
724 zp = int ( pd_proj->GetPDGCharge()/eplus + 0.5 );
725 // ap = 1;
726 }
727
728
729 G4int zt = pd_targ->GetAtomicNumber();
730 G4int at = pd_targ->GetAtomicMass();
731
732
733 // Check the ramx0 value
734 //G4double rmax0 = 8.0; // T.K dicide parameter value // for low energy
735 G4double rmax0 = bmax + 4.0;
736 G4double rmax = std::sqrt( rmax0*rmax0 + b*b );
737
738 G4double ccoul = 0.001439767;
739 G4double pcca = 1.0 - double ( zp * zt ) * ccoul / eccm / rmax - ( b / rmax )*( b / rmax );
740
741 G4double pccf = std::sqrt( pcca );
742
743 //Fix for neutral particles
744 G4double aas1 = 0.0;
745 G4double bbs = 0.0;
746
747 if ( zp != 0 )
748 {
749 G4double aas = 2.0 * eccm * b / double ( zp * zt ) / ccoul;
750 bbs = 1.0 / std::sqrt ( 1.0 + aas*aas );
751 aas1 = ( 1.0 + aas * b / rmax ) * bbs;
752 }
753
754 G4double cost = 0.0;
755 G4double sint = 0.0;
756 G4double thet1 = 0.0;
757 G4double thet2 = 0.0;
758 if ( 1.0 - aas1*aas1 <= 0 || 1.0 - bbs*bbs <= 0.0 )
759 {
760 cost = 1.0;
761 sint = 0.0;
762 }
763 else
764 {
765 G4double aat1 = aas1 / std::sqrt ( 1.0 - aas1*aas1 );
766 G4double aat2 = bbs / std::sqrt ( 1.0 - bbs*bbs );
767
768 thet1 = std::atan ( aat1 );
769 thet2 = std::atan ( aat2 );
770
771// TK enter to else block
772 G4double theta = thet1 - thet2;
773 cost = std::cos( theta );
774 sint = std::sin( theta );
775 }
776
777 G4double rzpr = -rmax * cost * ( mass_targ ) / ( mass_proj + mass_targ );
778 G4double rzta = rmax * cost * ( mass_proj ) / ( mass_proj + mass_targ );
779
780 G4double rxpr = rmax / 2.0 * sint;
781
782 G4double rxta = -rxpr;
783
784
785 G4double pzpc = pzcc * ( cost * pccf + sint * b / rmax );
786 G4double pxpr = pzcc * ( -sint * pccf + cost * b / rmax );
787
788 G4double pztc = - pzpc;
789 G4double pxta = - pxpr;
790
791 G4double epc = std::sqrt ( pzpc*pzpc + pxpr*pxpr + mass_proj*mass_proj );
792 G4double etc = std::sqrt ( pztc*pztc + pxta*pxta + mass_targ*mass_targ );
793
794 G4double pzpr = pzpc;
795 G4double pzta = pztc;
796 G4double epr = epc;
797 G4double eta = etc;
798
799// CM -> NN
800 G4double gammacm = boostToCM.gamma();
801 //G4double betacm = -boostToCM.beta();
802 G4double betacm = boostToCM.z();
803 pzpr = pzpc + betacm * gammacm * ( gammacm / ( 1. + gammacm ) * pzpc * betacm + epc );
804 pzta = pztc + betacm * gammacm * ( gammacm / ( 1. + gammacm ) * pztc * betacm + etc );
805 epr = gammacm * ( epc + betacm * pzpc );
806 eta = gammacm * ( etc + betacm * pztc );
807
808 //G4double betpr = pzpr / epr;
809 //G4double betta = pzta / eta;
810
811 G4double gammpr = epr / ( mass_proj );
812 G4double gammta = eta / ( mass_targ );
813
814 pzta = pzta / double ( at );
815 pxta = pxta / double ( at );
816
817 pzpr = pzpr / double ( ap );
818 pxpr = pxpr / double ( ap );
819
820 G4double zeroz = 0.0;
821
822 rzpr = rzpr -zeroz;
823 rzta = rzta -zeroz;
824
825 // Set results
826 coulomb_collision_gamma_proj = gammpr;
827 coulomb_collision_rx_proj = rxpr;
828 coulomb_collision_rz_proj = rzpr;
829 coulomb_collision_px_proj = pxpr;
830 coulomb_collision_pz_proj = pzpr;
831
832 coulomb_collision_gamma_targ = gammta;
833 coulomb_collision_rx_targ = rxta;
834 coulomb_collision_rz_targ = rzta;
835 coulomb_collision_px_targ = pxta;
836 coulomb_collision_pz_targ = pzta;
837
838}
839
840
841
842void G4QMDReaction::setEvaporationCh()
843{
844
845 if ( gem == true )
846 evaporation->SetGEMChannel();
847 else
848 evaporation->SetDefaultChannel();
849
850}
851
852void G4QMDReaction::ModelDescription(std::ostream& outFile) const
853{
854 outFile << "Lorentz covarianted Quantum Molecular Dynamics model for nucleus (particle) vs nucleus reactions\n";
855}
@ stopAndKill
std::vector< G4ReactionProduct * > G4ReactionProductVector
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
double beta() const
Definition: SpaceVectorP.cc:26
double z() const
Hep3Vector unit() const
double x() const
double y() const
double mag() const
double gamma() const
Definition: SpaceVectorP.cc:35
Hep3Vector findBoostToCM() const
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat) final
void BuildPhysicsTable(const G4ParticleDefinition &) final
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat) final
G4ParticleDefinition * GetDefinition() const
G4LorentzVector Get4Momentum() const
G4ThreeVector GetMomentum() const
void SetGEMChannel()
void SetDefaultChannel()
void SetEvaporation(G4VEvaporation *ptr, G4bool isLocal=false)
G4ReactionProductVector * BreakItUp(const G4Fragment &theInitialState)
void SetDeexChannelsType(G4DeexChannelType val)
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
std::size_t GetNumberOfSecondaries() const
G4HadSecondary * GetSecondary(size_t i)
const G4Material * GetMaterial() const
G4double GetTotalMomentum() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
void SetCreatorModelID(G4int id)
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
static G4IonTable * GetIonTable()
Definition: G4IonTable.cc:170
G4double GetIonMass(G4int Z, G4int A, G4int nL=0, G4int lvl=0) const
Definition: G4IonTable.cc:1517
G4int GetA_asInt() const
Definition: G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105
G4int GetAtomicNumber() const
const G4String & GetParticleType() const
G4int GetAtomicMass() const
G4double GetPDGCharge() const
const G4String & GetParticleName() const
static G4int GetModelID(const G4int modelIndex)
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:97
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:97
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powN(G4double x, G4int n) const
Definition: G4Pow.cc:162
void SetMeanField(G4QMDMeanField *meanfield)
void CalKinematicsOfBinaryCollisions(G4double)
G4double GetTotalPotential()
void SetNucleus(G4QMDNucleus *aSystem)
void DoPropagation(G4double)
std::vector< G4QMDNucleus * > DoClusterJudgment()
void SetSystem(G4QMDSystem *aSystem)
void SetTotalPotential(G4double x)
Definition: G4QMDNucleus.hh:62
void CalEnergyAndAngularMomentumInCM()
G4ThreeVector GetPosition()
const G4ParticleDefinition * GetDefinition()
G4LorentzVector Get4Momentum()
G4ThreeVector GetMomentum()
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
virtual void ModelDescription(std::ostream &outFile) const
G4QMDParticipant * GetParticipant(G4int i)
Definition: G4QMDSystem.hh:83
G4int GetTotalNumberOfParticipant()
Definition: G4QMDSystem.hh:78
void Clear()
Definition: G4QMDSystem.cc:68
void SetParticipant(G4QMDParticipant *particle)
Definition: G4QMDSystem.hh:52
G4int GetNOCollision()
Definition: G4QMDSystem.hh:93
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
HepLorentzVector boostOf(const HepLorentzVector &vec, const Hep3Vector &betaVector)