Geant4  10.01.p01
G4QMDGroundStateNucleus.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 // 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
27 //
29 
30 #include "G4NucleiProperties.hh"
31 
32 #include "G4Proton.hh"
33 #include "G4Neutron.hh"
34 
35 #include "G4PhysicalConstants.hh"
36 #include "Randomize.hh"
37 
39 : maxTrial ( 1000 )
40 , r00 ( 1.124 ) // radius parameter for Woods-Saxon [fm]
41 , r01 ( 0.5 ) // radius parameter for Woods-Saxon
42 , saa ( 0.2 ) // diffuse parameter for initial Woods-Saxon shape
43 , rada ( 0.9 ) // cutoff parameter
44 , radb ( 0.3 ) // cutoff parameter
45 , dsam ( 1.5 ) // minimum distance for same particle [fm]
46 , ddif ( 1.0 ) // minimum distance for different particle
47 , edepth ( 0.0 )
48 , epse ( 0.000001 ) // torelance for energy in [GeV]
49 , meanfield ( NULL )
50 {
51 
52  //std::cout << " G4QMDGroundStateNucleus( G4int z , G4int a ) Begin " << z << " " << a << std::endl;
53 
54  dsam2 = dsam*dsam;
55  ddif2 = ddif*ddif;
56 
58 
59  hbc = parameters->Get_hbc();
60  gamm = parameters->Get_gamm();
61  cpw = parameters->Get_cpw();
62  cph = parameters->Get_cph();
63  epsx = parameters->Get_epsx();
64  cpc = parameters->Get_cpc();
65 
66  cdp = parameters->Get_cdp();
67  c0p = parameters->Get_c0p();
68  c3p = parameters->Get_c3p();
69  csp = parameters->Get_csp();
70  clp = parameters->Get_clp();
71 
72  //edepth = 0.0;
73 
74  for ( int i = 0 ; i < a ; i++ )
75  {
76 
78 
79  if ( i < z )
80  {
81  pd = G4Proton::Proton();
82  }
83  else
84  {
85  pd = G4Neutron::Neutron();
86  }
87 
88  G4ThreeVector p( 0.0 );
89  G4ThreeVector r( 0.0 );
90  G4QMDParticipant* aParticipant = new G4QMDParticipant( pd , p , r );
91  SetParticipant( aParticipant );
92 
93  }
94 
95  G4double radious = r00 * std::pow ( double ( GetMassNumber() ) , 1.0/3.0 );
96 
97  rt00 = radious - r01;
98  radm = radious - rada * ( gamm - 1.0 ) + radb;
99  rmax = 1.0 / ( 1.0 + std::exp ( -rt00/saa ) );
100 
101  //maxTrial = 1000;
102 
103  //Nucleon primary or target case;
104  if ( z == 1 && a == 1 ) { // Hydrogen Case or proton primary
106  ebini = 0.0;
107  return;
108  }
109  else if ( z == 0 && a == 1 ) { // Neutron primary
111  ebini = 0.0;
112  return;
113  }
114 
115 
116  meanfield = new G4QMDMeanField();
117  meanfield->SetSystem( this );
118 
119  //std::cout << "G4QMDGroundStateNucleus( G4int z , G4int a ) packNucleons Begin ( z , a ) ( " << z << ", " << a << " )" << std::endl;
120  packNucleons();
121  //std::cout << "G4QMDGroundStateNucleus( G4int z , G4int a ) packNucleons End" << std::endl;
122 
123  delete meanfield;
124 
125 }
126 
127 
128 
130 {
131 
132  //std::cout << "G4QMDGroundStateNucleus::packNucleons" << std::endl;
133 
135 
136  G4double ebin00 = ebini * 0.001;
137 
138  G4double ebin0 = 0.0;
139  G4double ebin1 = 0.0;
140 
141  if ( GetMassNumber() != 4 )
142  {
143  ebin0 = ( ebini - 0.5 ) * 0.001;
144  ebin1 = ( ebini + 0.5 ) * 0.001;
145  }
146  else
147  {
148  ebin0 = ( ebini - 1.5 ) * 0.001;
149  ebin1 = ( ebini + 1.5 ) * 0.001;
150  }
151 
152 {
153  G4int n0Try = 0;
154  G4bool isThisOK = false;
155  while ( n0Try < maxTrial )
156  {
157  n0Try++;
158  //std::cout << "TKDB packNucleons n0Try " << n0Try << std::endl;
159 
160 // Sampling Position
161 
162  //std::cout << "TKDB Sampling Position " << std::endl;
163 
164  G4bool areThesePsOK = false;
165  G4int npTry = 0;
166  while ( npTry < maxTrial )
167  {
168  //std::cout << "TKDB Sampling Position npTry " << npTry << std::endl;
169  npTry++;
170  G4int i = 0;
171  if ( samplingPosition( i ) )
172  {
173  //std::cout << "packNucleons samplingPosition 0 succeed " << std::endl;
174  for ( i = 1 ; i < GetMassNumber() ; i++ )
175  {
176  //std::cout << "packNucleons samplingPosition " << i << " trying " << std::endl;
177  if ( !( samplingPosition( i ) ) )
178  {
179  //std::cout << "packNucleons samplingPosition " << i << " failed" << std::endl;
180  break;
181  }
182  }
183  if ( i == GetMassNumber() )
184  {
185  //std::cout << "packNucleons samplingPosition all scucceed " << std::endl;
186  areThesePsOK = true;
187  break;
188  }
189  }
190  }
191  if ( areThesePsOK == false ) continue;
192 
193  //std::cout << "TKDB Sampling Position End" << std::endl;
194 
195 // Calculate Two-body quantities
196 
198  std::vector< G4double > rho_a ( GetMassNumber() , 0.0 );
199  std::vector< G4double > rho_s ( GetMassNumber() , 0.0 );
200  std::vector< G4double > rho_c ( GetMassNumber() , 0.0 );
201 
202  rho_l.resize ( GetMassNumber() , 0.0 );
203  d_pot.resize ( GetMassNumber() , 0.0 );
204 
205  for ( G4int i = 0 ; i < GetMassNumber() ; i++ )
206  {
207  for ( G4int j = 0 ; j < GetMassNumber() ; j++ )
208  {
209 
210  rho_a[ i ] += meanfield->GetRHA( j , i );
211  G4int k = 0;
212  if ( participants[j]->GetDefinition() != participants[i]->GetDefinition() )
213  {
214  k = 1;
215  }
216 
217  rho_s[ i ] += meanfield->GetRHA( j , i )*( 1.0 - 2.0 * k ); // OK?
218 
219  rho_c[ i ] += meanfield->GetRHE( j , i );
220  }
221 
222  }
223 
224  for ( G4int i = 0 ; i < GetMassNumber() ; i++ )
225  {
226  rho_l[i] = cdp * ( rho_a[i] + 1.0 );
227  d_pot[i] = c0p * rho_a[i]
228  + c3p * std::pow ( rho_a[i] , gamm )
229  + csp * rho_s[i]
230  + clp * rho_c[i];
231 
232  //std::cout << "d_pot[i] " << i << " " << d_pot[i] << std::endl;
233  }
234 
235 
236 // Sampling Momentum
237 
238  //std::cout << "TKDB Sampling Momentum " << std::endl;
239 
240  phase_g.clear();
241  phase_g.resize( GetMassNumber() , 0.0 );
242 
243  //std::cout << "TKDB Sampling Momentum 1st " << std::endl;
244  G4bool isThis1stMOK = false;
245  G4int nmTry = 0;
246  while ( nmTry < maxTrial )
247  {
248  nmTry++;
249  G4int i = 0;
250  if ( samplingMomentum( i ) )
251  {
252  isThis1stMOK = true;
253  break;
254  }
255  }
256  if ( isThis1stMOK == false ) continue;
257 
258  //std::cout << "TKDB Sampling Momentum 2nd so on" << std::endl;
259 
260  G4bool areTheseMsOK = true;
261  nmTry = 0;
262  while ( nmTry < maxTrial )
263  {
264  nmTry++;
265  G4int i = 0;
266  for ( i = 1 ; i < GetMassNumber() ; i++ )
267  {
268  //std::cout << "TKDB packNucleons samplingMomentum try " << i << std::endl;
269  if ( !( samplingMomentum( i ) ) )
270  {
271  //std::cout << "TKDB packNucleons samplingMomentum " << i << " failed" << std::endl;
272  areTheseMsOK = false;
273  break;
274  }
275  }
276  if ( i == GetMassNumber() )
277  {
278  areTheseMsOK = true;
279  }
280 
281  break;
282  }
283  if ( areTheseMsOK == false ) continue;
284 
285 // Kill Angluar Momentum
286 
287  //std::cout << "TKDB Sampling Kill Angluar Momentum " << std::endl;
288 
290 
291 
292 // Check Binding Energy
293 
294  //std::cout << "packNucleons Check Binding Energy Begin " << std::endl;
295 
296  G4double ekinal = 0.0;
297  for ( int i = 0 ; i < GetMassNumber() ; i++ )
298  {
299  ekinal += participants[i]->GetKineticEnergy();
300  }
301 
303 
304  G4double totalPotentialE = meanfield->GetTotalPotential();
305  G4double ebinal = ( totalPotentialE + ekinal ) / double ( GetMassNumber() );
306 
307  //std::cout << "packNucleons totalPotentialE " << totalPotentialE << std::endl;
308  //std::cout << "packNucleons ebinal " << ebinal << std::endl;
309  //std::cout << "packNucleons ekinal " << ekinal << std::endl;
310 
311  if ( ebinal < ebin0 || ebinal > ebin1 )
312  {
313  //std::cout << "packNucleons ebin0 " << ebin0 << std::endl;
314  //std::cout << "packNucleons ebin1 " << ebin1 << std::endl;
315  //std::cout << "packNucleons ebinal " << ebinal << std::endl;
316  //std::cout << "packNucleons Check Binding Energy Failed " << std::endl;
317  continue;
318  }
319 
320  //std::cout << "packNucleons Check Binding Energy End = OK " << std::endl;
321 
322 
323 // Energy Adujstment
324 
325  G4double dtc = 1.0;
326  G4double frg = -0.1;
327  G4double rdf0 = 0.5;
328 
329  G4double edif0 = ebinal - ebin00;
330 
331  G4double cfrc = 0.0;
332  if ( 0 < edif0 )
333  {
334  cfrc = frg;
335  }
336  else
337  {
338  cfrc = -frg;
339  }
340 
341  G4int ifrc = 1;
342 
343  G4int neaTry = 0;
344 
345  G4bool isThisEAOK = false;
346  while ( neaTry < maxTrial )
347  {
348  neaTry++;
349 
350  G4double edif = ebinal - ebin00;
351 
352  //std::cout << "TKDB edif " << edif << std::endl;
353  if ( std::abs ( edif ) < epse )
354  {
355 
356  isThisEAOK = true;
357  //std::cout << "isThisEAOK " << isThisEAOK << std::endl;
358  break;
359  }
360 
361  G4int jfrc = 0;
362  if ( edif < 0.0 )
363  {
364  jfrc = 1;
365  }
366  else
367  {
368  jfrc = -1;
369  }
370 
371  if ( jfrc != ifrc )
372  {
373  cfrc = -rdf0 * cfrc;
374  dtc = rdf0 * dtc;
375  }
376 
377  if ( jfrc == ifrc && std::abs( edif0 ) < std::abs( edif ) )
378  {
379  cfrc = -rdf0 * cfrc;
380  dtc = rdf0 * dtc;
381  }
382 
383  ifrc = jfrc;
384  edif0 = edif;
385 
387 
388  for ( int i = 0 ; i < GetMassNumber() ; i++ )
389  {
390  G4ThreeVector ri = participants[i]->GetPosition();
391  G4ThreeVector p_i = participants[i]->GetMomentum();
392 
393  ri += dtc * ( meanfield->GetFFr(i) - cfrc * ( meanfield->GetFFp(i) ) );
394  p_i += dtc * ( meanfield->GetFFp(i) + cfrc * ( meanfield->GetFFr(i) ) );
395 
396  participants[i]->SetPosition( ri );
397  participants[i]->SetMomentum( p_i );
398  }
399 
400  ekinal = 0.0;
401 
402  for ( int i = 0 ; i < GetMassNumber() ; i++ )
403  {
404  ekinal += participants[i]->GetKineticEnergy();
405  }
406 
408  totalPotentialE = meanfield->GetTotalPotential();
409 
410  ebinal = ( totalPotentialE + ekinal ) / double ( GetMassNumber() );
411 
412  }
413  //std::cout << "isThisEAOK " << isThisEAOK << std::endl;
414  if ( isThisEAOK == false ) continue;
415 
416  isThisOK = true;
417  //std::cout << "isThisOK " << isThisOK << std::endl;
418  break;
419 
420  }
421 
422  if ( isThisOK == false )
423  {
424  G4cout << "GroundStateNucleus state cannot be created. Try again with another parameters." << G4endl;
425  }
426 
427  //std::cout << "packNucleons End" << std::endl;
428  return;
429 
430 }
431 
432 
433 
434 
435 
436 // Start packing
437 // position
438 
439  G4int n0Try = 0;
440 
441  while ( n0Try < maxTrial )
442  {
443  if ( samplingPosition( 0 ) )
444  {
445  G4int i = 0;
446  for ( i = 1 ; i < GetMassNumber() ; i++ )
447  {
448  if ( !( samplingPosition( i ) ) )
449  {
450  break;
451  }
452  }
453  if ( i == GetMassNumber() ) break;
454  }
455  n0Try++;
456  }
457 
458  if ( n0Try > maxTrial )
459  {
460  G4cout << "GroundStateNucleus state cannot be created. Try again with another parameters." << G4endl;
461  return;
462  }
463 
465  std::vector< G4double > rho_a ( GetMassNumber() , 0.0 );
466  std::vector< G4double > rho_s ( GetMassNumber() , 0.0 );
467  std::vector< G4double > rho_c ( GetMassNumber() , 0.0 );
468 
469  rho_l.resize ( GetMassNumber() , 0.0 );
470  d_pot.resize ( GetMassNumber() , 0.0 );
471 
472  for ( int i = 0 ; i < GetMassNumber() ; i++ )
473  {
474  for ( int j = 0 ; j < GetMassNumber() ; j++ )
475  {
476 
477  rho_a[ i ] += meanfield->GetRHA( j , i );
478  G4int k = 0;
479  if ( participants[i]->GetDefinition() != participants[j]->GetDefinition() )
480  {
481  k = 1;
482  }
483 
484  rho_s[ i ] += meanfield->GetRHA( j , i )*( 1.0 - 2.0 * k ); // OK?
485 
486  rho_c[ i ] += meanfield->GetRHE( j , i );
487  }
488  }
489 
490  for ( int i = 0 ; i < GetMassNumber() ; i++ )
491  {
492  rho_l[i] = cdp * ( rho_a[i] + 1.0 );
493  d_pot[i] = c0p * rho_a[i]
494  + c3p * std::pow ( rho_a[i] , gamm )
495  + csp * rho_s[i]
496  + clp * rho_c[i];
497  }
498 
499 
500 
501 
502 
503 
504 // momentum
505 
506  phase_g.resize( GetMassNumber() , 0.0 );
507 
508  //G4int i = 0;
509  samplingMomentum( 0 );
510 
511  G4int n1Try = 0;
512  //G4int maxTry = 1000;
513 
514  while ( n1Try < maxTrial )
515  {
516  if ( samplingPosition( 0 ) )
517  {
518  G4int i = 0;
519  G4bool isThisOK = true;
520  for ( i = 1 ; i < GetMassNumber() ; i++ )
521  {
522  if ( !( samplingMomentum( i ) ) )
523  {
524  isThisOK = false;
525  break;
526  }
527  }
528  if ( isThisOK == true ) break;
529  //if ( i == GetMassNumber() ) break;
530  }
531  n1Try++;
532  }
533 
534  if ( n1Try > maxTrial )
535  {
536  G4cout << "GroundStateNucleus state cannot be created. Try again with another parameters." << G4endl;
537  return;
538  }
539 
540 //
541 
542 // Shift nucleus to thier CM frame and kill angular momentum
544 
545 // Check binding energy
546 
547  G4double ekinal = 0.0;
548  for ( int i = 0 ; i < GetMassNumber() ; i++ )
549  {
550  ekinal += participants[i]->GetKineticEnergy();
551  }
552 
554  G4double totalPotentialE = meanfield->GetTotalPotential();
555  G4double ebinal = ( totalPotentialE + ekinal ) / double ( GetMassNumber() );
556 
557  if ( ebinal < ebin0 || ebinal > ebin1 )
558  {
559  // Retry From Position
560  }
561 
562 
563 // Adjust by frictional cooling or heating
564 
565  G4double dtc = 1.0;
566  G4double frg = -0.1;
567  G4double rdf0 = 0.5;
568 
569  G4double edif0 = ebinal - ebin00;
570 
571  G4double cfrc = 0.0;
572  if ( 0 < edif0 )
573  {
574  cfrc = frg;
575  }
576  else
577  {
578  cfrc = -frg;
579  }
580 
581  G4int ifrc = 1;
582 
583  G4int ntryACH = 0;
584 
585  G4bool isThisOK = false;
586  while ( ntryACH < maxTrial )
587  {
588 
589  G4double edif = ebinal - ebin00;
590  if ( std::abs ( edif ) < epse )
591  {
592  isThisOK = true;
593  break;
594  }
595 
596  G4int jfrc = 0;
597  if ( edif < 0.0 )
598  {
599  jfrc = 1;
600  }
601  else
602  {
603  jfrc = -1;
604  }
605 
606  if ( jfrc != ifrc )
607  {
608  cfrc = -rdf0 * cfrc;
609  dtc = rdf0 * dtc;
610  }
611 
612  if ( jfrc == ifrc && std::abs( edif0 ) < std::abs( edif ) )
613  {
614  cfrc = -rdf0 * cfrc;
615  dtc = rdf0 * dtc;
616  }
617 
618  ifrc = jfrc;
619  edif0 = edif;
620 
622 
623  for ( int i = 0 ; i < GetMassNumber() ; i++ )
624  {
625  G4ThreeVector ri = participants[i]->GetPosition();
626  G4ThreeVector p_i = participants[i]->GetMomentum();
627 
628  ri += dtc * ( meanfield->GetFFr(i) - cfrc * ( meanfield->GetFFp(i) ) );
629  p_i += dtc * ( meanfield->GetFFp(i) + cfrc * ( meanfield->GetFFr(i) ) );
630 
631  participants[i]->SetPosition( ri );
632  participants[i]->SetMomentum( p_i );
633  }
634 
635  ekinal = 0.0;
636 
637  for ( int i = 0 ; i < GetMassNumber() ; i++ )
638  {
639  ekinal += participants[i]->GetKineticEnergy();
640  }
641 
643  totalPotentialE = meanfield->GetTotalPotential();
644 
645  ebinal = ( totalPotentialE + ekinal ) / double ( GetMassNumber() );
646 
647 
648  ntryACH++;
649  }
650 
651  if ( isThisOK )
652  {
653  return;
654  }
655 
656 }
657 
658 
660 {
661 
662  G4bool result = false;
663 
664  G4int nTry = 0;
665 
666  while ( nTry < maxTrial )
667  {
668 
669  //std::cout << "samplingPosition i th particle, nTtry " << i << " " << nTry << std::endl;
670 
671  G4double rwod = -1.0;
672  G4double rrr = 0.0;
673 
674  G4double rx = 0.0;
675  G4double ry = 0.0;
676  G4double rz = 0.0;
677 
678  while ( G4UniformRand() * rmax > rwod )
679  {
680 
681  G4double rsqr = 10.0;
682  while ( rsqr > 1.0 )
683  {
684  rx = 1.0 - 2.0 * G4UniformRand();
685  ry = 1.0 - 2.0 * G4UniformRand();
686  rz = 1.0 - 2.0 * G4UniformRand();
687  rsqr = rx*rx + ry*ry + rz*rz;
688  }
689  rrr = radm * std::sqrt ( rsqr );
690  rwod = 1.0 / ( 1.0 + std::exp ( ( rrr - rt00 ) / saa ) );
691 
692  }
693 
694  participants[i]->SetPosition( G4ThreeVector( rx , ry , rz )*radm );
695 
696  if ( i == 0 )
697  {
698  result = true;
699  return result;
700  }
701 
702 // i > 1 ( Second Particle or later )
703 // Check Distance to others
704 
705  G4bool isThisOK = true;
706  for ( G4int j = 0 ; j < i ; j++ )
707  {
708 
709  G4double r2 = participants[j]->GetPosition().diff2( participants[i]->GetPosition() );
710  G4double dmin2 = 0.0;
711 
712  if ( participants[j]->GetDefinition() == participants[i]->GetDefinition() )
713  {
714  dmin2 = dsam2;
715  }
716  else
717  {
718  dmin2 = ddif2;
719  }
720 
721  //std::cout << "distance between j and i " << j << " " << i << " " << r2 << " " << dmin2 << std::endl;
722  if ( r2 < dmin2 )
723  {
724  isThisOK = false;
725  break;
726  }
727 
728  }
729 
730  if ( isThisOK == true )
731  {
732  result = true;
733  return result;
734  }
735 
736  nTry++;
737 
738  }
739 
740 // Here return "false"
741  return result;
742 }
743 
744 
745 
747 {
748 
749  //std::cout << "TKDB samplingMomentum for " << i << std::endl;
750 
751  G4bool result = false;
752 
753  G4double pfm = hbc * std::pow ( ( 3.0 / 2.0 * pi*pi * rho_l[i] ) , 1./3. );
754 
755  if ( 10 < GetMassNumber() && -5.5 < ebini )
756  {
757  pfm = pfm * ( 1.0 + 0.2 * std::sqrt( std::abs( 8.0 + ebini ) / 8.0 ) );
758  }
759 
760  //std::cout << "TKDB samplingMomentum pfm " << pfm << std::endl;
761 
762  std::vector< G4double > phase;
763  phase.resize( i+1 ); // i start from 0
764 
765  G4int ntry = 0;
766 // 710
767  while ( ntry < maxTrial )
768  {
769 
770  //std::cout << " TKDB ntry " << ntry << std::endl;
771  ntry++;
772 
773  G4double ke = DBL_MAX;
774 
775  G4int tkdb_i =0;
776 // 700
777  while ( ke + d_pot [i] > edepth )
778  {
779 
780  G4double psqr = 10.0;
781  G4double px = 0.0;
782  G4double py = 0.0;
783  G4double pz = 0.0;
784 
785  while ( psqr > 1.0 )
786  {
787  px = 1.0 - 2.0*G4UniformRand();
788  py = 1.0 - 2.0*G4UniformRand();
789  pz = 1.0 - 2.0*G4UniformRand();
790 
791  psqr = px*px + py*py + pz*pz;
792  }
793 
794  G4ThreeVector p ( px , py , pz );
795  p = pfm * p;
796  participants[i]->SetMomentum( p );
797  G4LorentzVector p4 = participants[i]->Get4Momentum();
798  //ke = p4.e() - p4.restMass();
799  ke = participants[i]->GetKineticEnergy();
800 
801 
802  tkdb_i++;
803  if ( tkdb_i > maxTrial ) return result; // return false
804 
805  }
806 
807  //std::cout << "TKDB ke d_pot[i] " << ke << " " << d_pot[i] << std::endl;
808 
809 
810  if ( i == 0 )
811  {
812  result = true;
813  return result;
814  }
815 
816  G4bool isThisOK = true;
817 
818  // Check Pauli principle
819 
820  phase[ i ] = 0.0;
821 
822  //std::cout << "TKDB Check Pauli Principle " << i << std::endl;
823 
824  for ( G4int j = 0 ; j < i ; j++ )
825  {
826  phase[ j ] = 0.0;
827  //std::cout << "TKDB Check Pauli Principle i , j " << i << " , " << j << std::endl;
828  G4double expa = 0.0;
829  if ( participants[j]->GetDefinition() == participants[i]->GetDefinition() )
830  {
831 
832  expa = - meanfield->GetRR2(i,j) * cpw;
833 
834  if ( expa > epsx )
835  {
836  G4ThreeVector p_i = participants[i]->GetMomentum();
837  G4ThreeVector pj = participants[j]->GetMomentum();
838  G4double dist2_p = p_i.diff2( pj );
839 
840  dist2_p = dist2_p*cph;
841  expa = expa - dist2_p;
842 
843  if ( expa > epsx )
844  {
845 
846  phase[j] = std::exp ( expa );
847 
848  if ( phase[j] * cpc > 0.2 )
849  {
850 /*
851  G4cout << "TKDB Check Pauli Principle A i , j " << i << " , " << j << G4endl;
852  G4cout << "TKDB Check Pauli Principle phase[j] " << phase[j] << G4endl;
853  G4cout << "TKDB Check Pauli Principle phase[j]*cpc > 0.2 " << phase[j]*cpc << G4endl;
854 */
855  isThisOK = false;
856  break;
857  }
858  if ( ( phase_g[j] + phase[j] ) * cpc > 0.5 )
859  {
860 /*
861  G4cout << "TKDB Check Pauli Principle B i , j " << i << " , " << j << G4endl;
862  G4cout << "TKDB Check Pauli Principle B phase_g[j] " << phase_g[j] << G4endl;
863  G4cout << "TKDB Check Pauli Principle B phase[j] " << phase[j] << G4endl;
864  G4cout << "TKDB Check Pauli Principle B phase_g[j] + phase[j] ) * cpc > 0.5 " << ( phase_g[j] + phase[j] ) * cpc << G4endl;
865 */
866  isThisOK = false;
867  break;
868  }
869 
870  phase[i] += phase[j];
871  if ( phase[i] * cpc > 0.3 )
872  {
873 /*
874  G4cout << "TKDB Check Pauli Principle C i , j " << i << " , " << j << G4endl;
875  G4cout << "TKDB Check Pauli Principle C phase[i] " << phase[i] << G4endl;
876  G4cout << "TKDB Check Pauli Principle C phase[i] * cpc > 0.3 " << phase[i] * cpc << G4endl;
877 */
878  isThisOK = false;
879  break;
880  }
881 
882  //std::cout << "TKDB Check Pauli Principle OK i , j " << i << " , " << j << std::endl;
883 
884  }
885  else
886  {
887  //std::cout << "TKDB Check Pauli Principle OK i , j " << i << " , " << j << std::endl;
888  }
889 
890  }
891  else
892  {
893  //std::cout << "TKDB Check Pauli Principle OK i , j " << i << " , " << j << std::endl;
894  }
895 
896  }
897  else
898  {
899  //std::cout << "TKDB Check Pauli Principle OK i , j " << i << " , " << j << std::endl;
900  }
901 
902  }
903 
904  if ( isThisOK == true )
905  {
906 
907  phase_g[i] = phase[i];
908 
909  for ( int j = 0 ; j < i ; j++ )
910  {
911  phase_g[j] += phase[j];
912  }
913 
914  result = true;
915  return result;
916  }
917 
918  }
919 
920  return result;
921 
922 }
923 
924 
925 
927 {
928 
929 // CalEnergyAndAngularMomentumInCM();
930 
931  //std::vector< G4ThreeVector > p ( GetMassNumber() , 0.0 );
932  //std::vector< G4ThreeVector > r ( GetMassNumber() , 0.0 );
933 
934 // Move to cm system
935 
936  G4ThreeVector pcm_tmp ( 0.0 );
937  G4ThreeVector rcm_tmp ( 0.0 );
938  G4double sumMass = 0.0;
939 
940  for ( G4int i = 0 ; i < GetMassNumber() ; i++ )
941  {
942  pcm_tmp += participants[i]->GetMomentum();
943  rcm_tmp += participants[i]->GetPosition() * participants[i]->GetMass();
944  sumMass += participants[i]->GetMass();
945  }
946 
947  pcm_tmp = pcm_tmp/GetMassNumber();
948  rcm_tmp = rcm_tmp/sumMass;
949 
950  for ( G4int i = 0 ; i < GetMassNumber() ; i++ )
951  {
952  participants[i]->SetMomentum( participants[i]->GetMomentum() - pcm_tmp );
953  participants[i]->SetPosition( participants[i]->GetPosition() - rcm_tmp );
954  }
955 
956 // kill the angular momentum
957 
958  G4ThreeVector ll ( 0.0 );
959  for ( G4int i = 0 ; i < GetMassNumber() ; i++ )
960  {
961  ll += participants[i]->GetPosition().cross ( participants[i]->GetMomentum() );
962  }
963 
964  G4double rr[3][3];
965  G4double ss[3][3];
966  for ( G4int i = 0 ; i < 3 ; i++ )
967  {
968  for ( G4int j = 0 ; j < 3 ; j++ )
969  {
970  rr [i][j] = 0.0;
971 
972  if ( i == j )
973  {
974  ss [i][j] = 1.0;
975  }
976  else
977  {
978  ss [i][j] = 0.0;
979  }
980  }
981  }
982 
983  for ( G4int i = 0 ; i < GetMassNumber() ; i++ )
984  {
985  G4ThreeVector r = participants[i]->GetPosition();
986  rr[0][0] += r.y() * r.y() + r.z() * r.z();
987  rr[1][0] += - r.y() * r.x();
988  rr[2][0] += - r.z() * r.x();
989  rr[0][1] += - r.x() * r.y();
990  rr[1][1] += r.z() * r.z() + r.x() * r.x();
991  rr[2][1] += - r.z() * r.y();
992  rr[2][0] += - r.x() * r.z();
993  rr[2][1] += - r.y() * r.z();
994  rr[2][2] += r.x() * r.x() + r.y() * r.y();
995  }
996 
997  for ( G4int i = 0 ; i < 3 ; i++ )
998  {
999  G4double x = rr [i][i];
1000  for ( G4int j = 0 ; j < 3 ; j++ )
1001  {
1002  rr[i][j] = rr[i][j] / x;
1003  ss[i][j] = ss[i][j] / x;
1004  }
1005  for ( G4int j = 0 ; j < 3 ; j++ )
1006  {
1007  if ( j != i )
1008  {
1009  G4double y = rr [j][i];
1010  for ( G4int k = 0 ; k < 3 ; k++ )
1011  {
1012  rr[j][k] += -y * rr[i][k];
1013  ss[j][k] += -y * ss[i][k];
1014  }
1015  }
1016  }
1017  }
1018 
1019  G4double opl[3];
1020  G4double rll[3];
1021 
1022  rll[0] = ll.x();
1023  rll[1] = ll.y();
1024  rll[2] = ll.z();
1025 
1026  for ( G4int i = 0 ; i < 3 ; i++ )
1027  {
1028  opl[i] = 0.0;
1029 
1030  for ( G4int j = 0 ; j < 3 ; j++ )
1031  {
1032  opl[i] += ss[i][j]*rll[j];
1033  }
1034  }
1035 
1036  for ( G4int i = 0 ; i < GetMassNumber() ; i++ )
1037  {
1038  G4ThreeVector p_i = participants[i]->GetMomentum() ;
1039  G4ThreeVector ri = participants[i]->GetPosition() ;
1040  G4ThreeVector opl_v ( opl[0] , opl[1] , opl[2] );
1041 
1042  p_i += ri.cross(opl_v);
1043 
1044  participants[i]->SetMomentum( p_i );
1045  }
1046 
1047 }
G4double GetRHA(G4int i, G4int j)
G4int GetAtomicNumber()
Definition: G4QMDNucleus.cc:88
std::vector< G4double > rho_l
G4int GetMassNumber()
Definition: G4QMDNucleus.cc:67
void SetParticipant(G4QMDParticipant *particle)
Definition: G4QMDSystem.hh:51
CLHEP::Hep3Vector G4ThreeVector
G4ThreeVector GetFFr(G4int i)
G4double Get_hbc()
G4double z
Definition: TRTMaterials.hh:39
G4double Get_c0p()
const G4double pi
G4double a
Definition: TRTMaterials.hh:39
int G4int
Definition: G4Types.hh:78
void Cal2BodyQuantities()
G4double Get_cpc()
#define G4UniformRand()
Definition: Randomize.hh:95
G4GLOB_DLL std::ostream G4cout
static G4QMDParameters * GetInstance()
bool G4bool
Definition: G4Types.hh:79
G4double Get_csp()
G4double Get_cdp()
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double Get_epsx()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4double Get_gamm()
G4double GetTotalPotential()
std::vector< G4QMDParticipant * > participants
Definition: G4QMDSystem.hh:72
G4double GetRHE(G4int i, G4int j)
static G4double GetBindingEnergy(const G4int A, const G4int Z)
G4double GetRR2(G4int i, G4int j)
G4double Get_clp()
void SetSystem(G4QMDSystem *aSystem)
G4double Get_c3p()
G4ThreeVector GetFFp(G4int i)
#define G4endl
Definition: G4ios.hh:61
G4double Get_cpw()
std::vector< G4double > d_pot
G4double Get_cph()
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
std::vector< G4double > phase_g
CLHEP::HepLorentzVector G4LorentzVector