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