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