Geant4_10
G4QMDMeanField.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 // 081120 Add Update by T. Koi
27 //
28 
29 #include <map>
30 #include <algorithm>
31 #include <numeric>
32 
33 #include <cmath>
34 #include <CLHEP/Random/Stat.h>
35 
36 #include "G4QMDMeanField.hh"
37 #include "G4QMDParameters.hh"
38 
39 #include "G4PhysicalConstants.hh"
40 #include "Randomize.hh"
41 
43 : rclds ( 4.0 ) // distance for cluster judgement
44 , epsx ( -20.0 ) // gauss term
45 , epscl ( 0.0001 ) // coulomb term
46 , irelcr ( 1 )
47 {
48 
50  wl = parameters->Get_wl();
51  cl = parameters->Get_cl();
52  rho0 = parameters->Get_rho0();
53  hbc = parameters->Get_hbc();
54  gamm = parameters->Get_gamm();
55 
56  cpw = parameters->Get_cpw();
57  cph = parameters->Get_cph();
58  cpc = parameters->Get_cpc();
59 
60  c0 = parameters->Get_c0();
61  c3 = parameters->Get_c3();
62  cs = parameters->Get_cs();
63 
64 // distance
65  c0w = 1.0/4.0/wl;
66  //c3w = 1.0/4.0/wl; //no need
67  c0sw = std::sqrt( c0w );
68  clw = 2.0 / std::sqrt ( 4.0 * pi * wl );
69 
70 // graduate
71  c0g = - c0 / ( 2.0 * wl );
72  c3g = - c3 / ( 4.0 * wl ) * gamm;
73  csg = - cs / ( 2.0 * wl );
74  pag = gamm - 1;
75 
76  system = NULL; // will be set through SetSystem method
77 }
78 
79 
80 
82 {
83  ;
84 }
85 
86 
87 
89 {
90 
91  //std::cout << "QMDMeanField SetSystem" << std::endl;
92 
93  system = aSystem;
94 
96 
97  pp2.clear();
98  rr2.clear();
99  rbij.clear();
100  rha.clear();
101  rhe.clear();
102  rhc.clear();
103 
104  rr2.resize( n );
105  pp2.resize( n );
106  rbij.resize( n );
107  rha.resize( n );
108  rhe.resize( n );
109  rhc.resize( n );
110 
111  for ( int i = 0 ; i < n ; i++ )
112  {
113  rr2[i].resize( n );
114  pp2[i].resize( n );
115  rbij[i].resize( n );
116  rha[i].resize( n );
117  rhe[i].resize( n );
118  rhc[i].resize( n );
119  }
120 
121 
122  ffr.clear();
123  ffp.clear();
124  rh3d.clear();
125 
126  ffr.resize( n );
127  ffp.resize( n );
128  rh3d.resize( n );
129 
131 
132 }
133 
135 {
136 
137  //std::cout << "QMDMeanField SetNucleus" << std::endl;
138 
139  SetSystem( aNucleus );
140 
141  G4double totalPotential = GetTotalPotential();
142  aNucleus->SetTotalPotential( totalPotential );
143 
145 
146 }
147 
148 
149 
151 {
152 
153  if ( system->GetTotalNumberOfParticipant() < 2 ) return;
154 
155  for ( G4int j = 1 ; j < system->GetTotalNumberOfParticipant() ; j++ )
156  {
157 
158  G4ThreeVector rj = system->GetParticipant( j )->GetPosition();
159  G4LorentzVector p4j = system->GetParticipant( j )->Get4Momentum();
160 
161  for ( G4int i = 0 ; i < j ; i++ )
162  {
163 
164  G4ThreeVector ri = system->GetParticipant( i )->GetPosition();
165  G4LorentzVector p4i = system->GetParticipant( i )->Get4Momentum();
166 
167  G4ThreeVector rij = ri - rj;
168  G4ThreeVector pij = (p4i - p4j).v();
169  G4LorentzVector p4ij = p4i - p4j;
170  G4ThreeVector bij = ( p4i + p4j ).boostVector();
171  G4double gammaij = ( p4i + p4j ).gamma();
172 
173  G4double eij = ( p4i + p4j ).e();
174 
175  G4double rbrb = rij*bij;
176 // G4double bij2 = bij*bij;
177  G4double rij2 = rij*rij;
178  G4double pij2 = pij*pij;
179 
180  rbrb = irelcr * rbrb;
181  G4double gamma2_ij = gammaij*gammaij;
182 
183 
184  rr2[i][j] = rij2 + gamma2_ij * rbrb*rbrb;
185  rr2[j][i] = rr2[i][j];
186 
187  rbij[i][j] = gamma2_ij * rbrb;
188  rbij[j][i] = - rbij[i][j];
189 
190  pp2[i][j] = pij2
191  + irelcr * ( - std::pow ( p4i.e() - p4j.e() , 2 )
192  + gamma2_ij * std::pow ( ( ( p4i.m2() - p4j.m2() ) / eij ) , 2 ) );
193 
194 
195  pp2[j][i] = pp2[i][j];
196 
197 // Gauss term
198 
199  G4double expa1 = - rr2[i][j] * c0w;
200 
201  G4double rh1;
202  if ( expa1 > epsx )
203  {
204  rh1 = std::exp( expa1 );
205  }
206  else
207  {
208  rh1 = 0.0;
209  }
210 
211  G4int ibry = system->GetParticipant(i)->GetBaryonNumber();
212  G4int jbry = system->GetParticipant(j)->GetBaryonNumber();
213 
214 
215  rha[i][j] = ibry*jbry*rh1;
216  rha[j][i] = rha[i][j];
217 
218 // Coulomb terms
219 
220  G4double rrs2 = rr2[i][j] + epscl;
221  G4double rrs = std::sqrt ( rrs2 );
222 
223  G4int icharge = system->GetParticipant(i)->GetChargeInUnitOfEplus();
224  G4int jcharge = system->GetParticipant(j)->GetChargeInUnitOfEplus();
225 
226  G4double xerf = 0.0;
227  // T. K. add this protection. 5.8 is good enough for double
228  if ( rrs*c0sw < 5.8 ) {
229  //erf = G4RandStat::erf ( rrs*c0sw );
230  //Restore to CLHEP for avoiding compilation error in MT
231  //erf = CLHEP::HepStat::erf ( rrs*c0sw );
232  //Use cmath
233 #if defined WIN32-VC
234  xerf = CLHEP::HepStat::erf ( rrs*c0sw );
235 #else
236  xerf = erf ( rrs*c0sw );
237 #endif
238  } else {
239  xerf = 1.0;
240  }
241 
242  G4double erfij = xerf/rrs;
243 
244 
245  rhe[i][j] = icharge*jcharge * erfij;
246 
247  rhe[j][i] = rhe[i][j];
248 
249  rhc[i][j] = icharge*jcharge * ( - erfij + clw * rh1 ) / rrs2;
250 
251  rhc[j][i] = rhc[i][j];
252 
253  } // i
254  } // j
255 }
256 
257 
258 
260 {
261 
262  //std::cout << "Cal2BodyQuantities " << i << std::endl;
263 
264  G4ThreeVector ri = system->GetParticipant( i )->GetPosition();
265  G4LorentzVector p4i = system->GetParticipant( i )->Get4Momentum();
266 
267 
268  for ( G4int j = 0 ; j < system->GetTotalNumberOfParticipant() ; j ++ )
269  {
270  if ( j == i ) continue;
271 
272  G4ThreeVector rj = system->GetParticipant( j )->GetPosition();
273  G4LorentzVector p4j = system->GetParticipant( j )->Get4Momentum();
274 
275  G4ThreeVector rij = ri - rj;
276  G4ThreeVector pij = (p4i - p4j).v();
277  G4LorentzVector p4ij = p4i - p4j;
278  G4ThreeVector bij = ( p4i + p4j ).boostVector();
279  G4double gammaij = ( p4i + p4j ).gamma();
280 
281  G4double eij = ( p4i + p4j ).e();
282 
283  G4double rbrb = rij*bij;
284 // G4double bij2 = bij*bij;
285  G4double rij2 = rij*rij;
286  G4double pij2 = pij*pij;
287 
288  rbrb = irelcr * rbrb;
289  G4double gamma2_ij = gammaij*gammaij;
290 
291 /*
292  G4double rbrb = 0.0;
293  G4double beta2_ij = 0.0;
294  G4double rij2 = 0.0;
295  G4double pij2 = 0.0;
296 
297 //
298  G4LorentzVector p4ip4j = p4i + p4j;
299  G4double eij = p4ip4j.e();
300 
301  G4ThreeVector r = ri - rj;
302  G4LorentzVector p4 = p4i - p4j;
303 
304  rbrb = r.x()*p4ip4j.x()/eij
305  + r.y()*p4ip4j.y()/eij
306  + r.z()*p4ip4j.z()/eij;
307 
308  beta2_ij = ( p4ip4j.x()*p4ip4j.x() + p4ip4j.y()*p4ip4j.y() + p4ip4j.z()*p4ip4j.z() ) / ( eij*eij );
309  rij2 = r*r;
310  pij2 = p4.v()*p4.v();
311 
312  rbrb = irelcr * rbrb;
313 
314  G4double gamma2_ij = 1 / ( 1 - beta2_ij );
315 */
316 
317  rr2[i][j] = rij2 + gamma2_ij * rbrb*rbrb;
318  rr2[j][i] = rr2[i][j];
319 
320  rbij[i][j] = gamma2_ij * rbrb;
321  rbij[j][i] = - rbij[i][j];
322 
323  pp2[i][j] = pij2
324  + irelcr * ( - std::pow ( p4i.e() - p4j.e() , 2 )
325  + gamma2_ij * std::pow ( ( ( p4i.m2() - p4j.m2() ) / eij ) , 2 ) );
326 
327  pp2[j][i] = pp2[i][j];
328 
329 // Gauss term
330 
331  G4double expa1 = - rr2[i][j] * c0w;
332 
333  G4double rh1;
334  if ( expa1 > epsx )
335  {
336  rh1 = std::exp( expa1 );
337  }
338  else
339  {
340  rh1 = 0.0;
341  }
342 
343  G4int ibry = system->GetParticipant(i)->GetBaryonNumber();
344  G4int jbry = system->GetParticipant(j)->GetBaryonNumber();
345 
346 
347  rha[i][j] = ibry*jbry*rh1;
348  rha[j][i] = rha[i][j];
349 
350 // Coulomb terms
351 
352  G4double rrs2 = rr2[i][j] + epscl;
353  G4double rrs = std::sqrt ( rrs2 );
354 
355  G4int icharge = system->GetParticipant(i)->GetChargeInUnitOfEplus();
356  G4int jcharge = system->GetParticipant(j)->GetChargeInUnitOfEplus();
357 
358  G4double xerf = 0.0;
359  // T. K. add this protection. 5.8 is good enough for double
360  if ( rrs*c0sw < 5.8 ) {
361  //xerf = G4RandStat::erf ( rrs*c0sw );
362  //Use cmath
363 #if defined WIN32-VC
364  xerf = CLHEP::HepStat::erf ( rrs*c0sw );
365 #else
366  xerf = erf ( rrs*c0sw );
367 #endif
368  } else {
369  xerf = 1.0;
370  }
371 
372  G4double erfij = xerf/rrs;
373 
374 
375  rhe[i][j] = icharge*jcharge * erfij;
376 
377  rhe[j][i] = rhe[i][j];
378 
379 // G4double clw;
380 
381  rhc[i][j] = icharge*jcharge * ( - erfij + clw * rh1 ) / rrs2;
382 
383  rhc[j][i] = rhc[i][j];
384 
385  }
386 
387 }
388 
389 
390 
392 {
393 
394  ffr.resize( system->GetTotalNumberOfParticipant() );
395  ffp.resize( system->GetTotalNumberOfParticipant() );
396  rh3d.resize( system->GetTotalNumberOfParticipant() );
397 
398  for ( G4int i = 0 ; i < system->GetTotalNumberOfParticipant() ; i ++ )
399  {
400  G4double rho3 = 0.0;
401  for ( G4int j = 0 ; j < system->GetTotalNumberOfParticipant() ; j ++ )
402  {
403  rho3 += rha[j][i];
404  }
405  rh3d[i] = std::pow ( rho3 , pag );
406  }
407 
408 
409  for ( G4int i = 0 ; i < system->GetTotalNumberOfParticipant() ; i ++ )
410  {
411 
412  G4ThreeVector ri = system->GetParticipant( i )->GetPosition();
413  G4LorentzVector p4i = system->GetParticipant( i )->Get4Momentum();
414 
415  G4ThreeVector betai = p4i.v()/p4i.e();
416 
417 // R-JQMD
418  G4double Vi = GetPotential( i );
419  G4double p_zero = std::sqrt( p4i.e()*p4i.e() + 2*p4i.m()*Vi);
420  G4ThreeVector betai_R = p4i.v()/p_zero;
421  G4double mi_R = p4i.m()/p_zero;
422 //
423  ffr[i] = betai_R;
424  ffp[i] = G4ThreeVector( 0.0 );
425 
426  if ( false )
427  {
428  ffr[i] = betai;
429  mi_R = 1.0;
430  }
431 
432  for ( G4int j = 0 ; j < system->GetTotalNumberOfParticipant() ; j ++ )
433  {
434 
435  G4ThreeVector rj = system->GetParticipant( j )->GetPosition();
436  G4LorentzVector p4j = system->GetParticipant( j )->Get4Momentum();
437 
438  G4double eij = p4i.e() + p4j.e();
439 
440  G4int icharge = system->GetParticipant(i)->GetChargeInUnitOfEplus();
441  G4int jcharge = system->GetParticipant(j)->GetChargeInUnitOfEplus();
442 
443  G4int inuc = system->GetParticipant(i)->GetNuc();
444  G4int jnuc = system->GetParticipant(j)->GetNuc();
445 
446  G4double ccpp = c0g * rha[j][i]
447  + c3g * rha[j][i] * ( rh3d[j] + rh3d[i] )
448  + csg * rha[j][i] * jnuc * inuc
449  * ( 1. - 2. * std::abs( jcharge - icharge ) )
450  + cl * rhc[j][i];
451  ccpp *= mi_R;
452 
453 /*
454  G4cout << c0g << " " << c3g << " " << csg << " " << cl << G4endl;
455  G4cout << "ccpp " << i << " " << j << " " << ccpp << G4endl;
456  G4cout << "rha[j][i] " << rha[j][i] << G4endl;
457  G4cout << "rh3d " << rh3d[j] << " " << rh3d[i] << G4endl;
458  G4cout << "rhc[j][i] " << rhc[j][i] << G4endl;
459 */
460 
461  G4double grbb = - rbij[j][i];
462  G4double ccrr = grbb * ccpp / eij;
463 
464 /*
465  G4cout << "ccrr " << ccrr << G4endl;
466  G4cout << "grbb " << grbb << G4endl;
467 */
468 
469 
470  G4ThreeVector rij = ri - rj;
471  G4ThreeVector betaij = ( p4i + p4j ).v()/eij;
472 
473  G4ThreeVector cij = betaij - betai;
474 
475  ffr[i] = ffr[i] + 2*ccrr* ( rij + grbb*cij );
476 
477  ffp[i] = ffp[i] - 2*ccpp* ( rij + grbb*betaij );
478 
479  }
480  }
481 
482  //std::cout << "gradu 0 " << ffr[0] << " " << ffp[0] << std::endl;
483  //std::cout << "gradu 1 " << ffr[1] << " " << ffp[1] << std::endl;
484 
485 }
486 
487 
488 
490 {
491  G4int n = system->GetTotalNumberOfParticipant();
492 
493  G4double rhoa = 0.0;
494  G4double rho3 = 0.0;
495  G4double rhos = 0.0;
496  G4double rhoc = 0.0;
497 
498 
499  G4int icharge = system->GetParticipant(i)->GetChargeInUnitOfEplus();
500  G4int inuc = system->GetParticipant(i)->GetNuc();
501 
502  for ( G4int j = 0 ; j < n ; j ++ )
503  {
504  G4int jcharge = system->GetParticipant(j)->GetChargeInUnitOfEplus();
505  G4int jnuc = system->GetParticipant(j)->GetNuc();
506 
507  rhoa += rha[j][i];
508  rhoc += rhe[j][i];
509  rhos += rha[j][i] * jnuc * inuc
510  * ( 1 - 2 * std::abs ( jcharge - icharge ) );
511  }
512 
513  rho3 = std::pow ( rhoa , gamm );
514 
515  G4double potential = c0 * rhoa
516  + c3 * rho3
517  + cs * rhos
518  + cl * rhoc;
519 
520  return potential;
521 }
522 
523 
524 
526 {
527 
528  G4int n = system->GetTotalNumberOfParticipant();
529 
530  std::vector < G4double > rhoa ( n , 0.0 );
531  std::vector < G4double > rho3 ( n , 0.0 );
532  std::vector < G4double > rhos ( n , 0.0 );
533  std::vector < G4double > rhoc ( n , 0.0 );
534 
535 
536  for ( G4int i = 0 ; i < n ; i ++ )
537  {
538  G4int icharge = system->GetParticipant(i)->GetChargeInUnitOfEplus();
539  G4int inuc = system->GetParticipant(i)->GetNuc();
540 
541  for ( G4int j = 0 ; j < n ; j ++ )
542  {
543  G4int jcharge = system->GetParticipant(j)->GetChargeInUnitOfEplus();
544  G4int jnuc = system->GetParticipant(j)->GetNuc();
545 
546  rhoa[i] += rha[j][i];
547  rhoc[i] += rhe[j][i];
548  rhos[i] += rha[j][i] * jnuc * inuc
549  * ( 1 - 2 * std::abs ( jcharge - icharge ) );
550  }
551 
552  rho3[i] = std::pow ( rhoa[i] , gamm );
553  }
554 
555  G4double potential = c0 * std::accumulate( rhoa.begin() , rhoa.end() , 0.0 )
556  + c3 * std::accumulate( rho3.begin() , rho3.end() , 0.0 )
557  + cs * std::accumulate( rhos.begin() , rhos.end() , 0.0 )
558  + cl * std::accumulate( rhoc.begin() , rhoc.end() , 0.0 );
559 
560  return potential;
561 
562 }
563 
564 
565 
566 G4double G4QMDMeanField::calPauliBlockingFactor( G4int i )
567 {
568 
569  G4double pf = 0.0;
570 // i is supposed beyond total number of Participant()
571  G4int icharge = system->GetParticipant(i)->GetChargeInUnitOfEplus();
572 
573  for ( G4int j = 0 ; j < system->GetTotalNumberOfParticipant() ; j++ )
574  {
575 
576  G4int jcharge = system->GetParticipant(j)->GetChargeInUnitOfEplus();
577  G4int jnuc = system->GetParticipant(j)->GetNuc();
578 
579  if ( jcharge == icharge && jnuc == 1 )
580  {
581 
582 /*
583  G4cout << "Pauli i j " << i << " " << j << G4endl;
584  G4cout << "Pauli icharge " << icharge << G4endl;
585  G4cout << "Pauli jcharge " << jcharge << G4endl;
586 */
587  G4double expa = -rr2[i][j]*cpw;
588 
589 
590  if ( expa > epsx )
591  {
592  expa = expa - pp2[i][j]*cph;
593 /*
594  G4cout << "Pauli cph " << cph << G4endl;
595  G4cout << "Pauli pp2 " << pp2[i][j] << G4endl;
596  G4cout << "Pauli expa " << expa << G4endl;
597  G4cout << "Pauli epsx " << epsx << G4endl;
598 */
599  if ( expa > epsx )
600  {
601 // std::cout << "Pauli phase " << pf << std::endl;
602  pf = pf + std::exp ( expa );
603  }
604  }
605  }
606 
607  }
608 
609 
610  pf = ( pf - 1.0 ) * cpc;
611 
612  //std::cout << "Pauli pf " << pf << std::endl;
613 
614  return pf;
615 
616 }
617 
618 
619 
621 {
622  G4bool result = false;
623 
624  if ( system->GetParticipant( i )->GetNuc() == 1 )
625  {
626  G4double pf = calPauliBlockingFactor( i );
627  G4double rand = G4UniformRand();
628  if ( pf > rand ) result = true;
629  }
630 
631  return result;
632 }
633 
634 
635 
637 {
638 
639  G4double cc2 = 1.0;
640  G4double cc1 = 1.0 - cc2;
641  G4double cc3 = 1.0 / 2.0 / cc2;
642 
643  G4double dt3 = dt * cc3;
644  G4double dt1 = dt * ( cc1 - cc3 );
645  G4double dt2 = dt * cc2;
646 
647  CalGraduate();
648 
649  G4int n = system->GetTotalNumberOfParticipant();
650 
651 // 1st Step
652 
653  std::vector< G4ThreeVector > f0r, f0p;
654  f0r.resize( n );
655  f0p.resize( n );
656 
657  for ( G4int i = 0 ; i < n ; i++ )
658  {
659  G4ThreeVector ri = system->GetParticipant( i )->GetPosition();
660  G4ThreeVector p3i = system->GetParticipant( i )->GetMomentum();
661 
662  ri += dt3* ffr[i];
663  p3i += dt3* ffp[i];
664 
665  f0r[i] = ffr[i];
666  f0p[i] = ffp[i];
667 
668  system->GetParticipant( i )->SetPosition( ri );
669  system->GetParticipant( i )->SetMomentum( p3i );
670 
671 // we do not need set total momentum by ourselvs
672  }
673 
674 // 2nd Step
676  CalGraduate();
677 
678  for ( G4int i = 0 ; i < n ; i++ )
679  {
680  G4ThreeVector ri = system->GetParticipant( i )->GetPosition();
681  G4ThreeVector p3i = system->GetParticipant( i )->GetMomentum();
682 
683  ri += dt1* f0r[i] + dt2* ffr[i];
684  p3i += dt1* f0p[i] + dt2* ffp[i];
685 
686  system->GetParticipant( i )->SetPosition( ri );
687  system->GetParticipant( i )->SetMomentum( p3i );
688 
689 // we do not need set total momentum by ourselvs
690  }
691 
693 
694 }
695 
696 
697 
698 std::vector< G4QMDNucleus* > G4QMDMeanField::DoClusterJudgment()
699 {
700 
701  //std::cout << "MeanField DoClusterJudgemnt" << std::endl;
702 
704 
705  G4double cpf2 = std::pow ( 1.5 * pi*pi * std::pow ( 4.0 * pi * wl , -1.5 )
706  ,
707  2./3. )
708  * hbc * hbc;
709  G4double rcc2 = rclds*rclds;
710 
711  G4int n = system->GetTotalNumberOfParticipant();
712  std::vector < G4double > rhoa;
713  rhoa.resize ( n );
714 
715  for ( G4int i = 0 ; i < n ; i++ )
716  {
717  rhoa[i] = 0.0;
718 
719  if ( system->GetParticipant( i )->GetBaryonNumber() == 1 )
720  {
721  for ( G4int j = 0 ; j < n ; j++ )
722  {
723  if ( system->GetParticipant( j )->GetBaryonNumber() == 1 )
724  rhoa[i] += rha[i][j];
725  }
726  }
727 
728  rhoa[i] = std::pow ( rhoa[i] + 1 , 1.0/3.0 );
729 
730  }
731 
732 // identification of the cluster
733 
734  std::map < G4int , std::vector < G4int > > cluster_map;
735  std::vector < G4bool > is_already_belong_some_cluster;
736 
737  // cluster_id participant_id
738  std::multimap < G4int , G4int > comb_map;
739  std::multimap < G4int , G4int > assign_map;
740  assign_map.clear();
741 
742  std::vector < G4int > mascl;
743  std::vector < G4int > num;
744  mascl.resize ( n );
745  num.resize ( n );
746  is_already_belong_some_cluster.resize ( n );
747 
748  std::vector < G4int > is_assigned_to ( n , -1 );
749  std::multimap < G4int , G4int > clusters;
750 
751  for ( G4int i = 0 ; i < n ; i++ )
752  {
753  mascl[i] = 1;
754  num[i] = 1;
755 
756  is_already_belong_some_cluster[i] = false;
757  }
758 
759 
760  G4int nclst = 1;
761  G4int ichek = 1;
762 
763 
764  G4int id = 0;
765  G4int cluster_id = -1;
766  for ( G4int i = 0 ; i < n-1 ; i++ )
767  {
768 
769  G4bool hasThisCompany = false;
770 // Check only for bryons?
771 // std::cout << "Check Baryon " << i << std::endl;
772 
773  if ( system->GetParticipant( i )->GetBaryonNumber() == 1 )
774  {
775 
776 // if ( is_already_belong_some_cluster[i] != true )
777 // {
778  //G4int j1 = ichek + 1;
779  G4int j1 = i + 1;
780  for ( G4int j = j1 ; j < n ; j++ )
781  {
782 
783  std::vector < G4int > cluster_participants;
784  if ( system->GetParticipant( j )->GetBaryonNumber() == 1 )
785  {
786  G4double rdist2 = rr2[ i ][ j ];
787  G4double pdist2 = pp2[ i ][ j ];
788  //G4double rdist2 = rr2[ num[i] ][ num[j] ];
789  //G4double pdist2 = pp2[ num[i] ][ num[j] ];
790  G4double pcc2 = cpf2
791  * ( rhoa[ i ] + rhoa[ j ] )
792  * ( rhoa[ i ] + rhoa[ j ] );
793 
794 // Check phase space: close enough?
795  if ( rdist2 < rcc2 && pdist2 < pcc2 )
796  {
797 
798 /*
799  G4cout << "G4QMDRESULT "
800  << i << " " << j << " " << id << " "
801  << is_assigned_to [ i ] << " " << is_assigned_to [ j ]
802  << G4endl;
803 */
804 
805  if ( is_assigned_to [ j ] == -1 )
806  {
807  if ( is_assigned_to [ i ] == -1 )
808  {
809  if ( clusters.size() != 0 )
810  {
811  id = clusters.rbegin()->first + 1;
812  //std::cout << "id is increare " << id << std::endl;
813  }
814  else
815  {
816  id = 0;
817  }
818  clusters.insert ( std::multimap<G4int,G4int>::value_type ( id , i ) );
819  is_assigned_to [ i ] = id;
820  clusters.insert ( std::multimap<G4int,G4int>::value_type ( id , j ) );
821  is_assigned_to [ j ] = id;
822  }
823  else
824  {
825  clusters.insert ( std::multimap<G4int,G4int>::value_type ( is_assigned_to [ i ] , j ) );
826  is_assigned_to [ j ] = is_assigned_to [ i ];
827  }
828  }
829  else
830  {
831 // j is already belong to some cluester
832  if ( is_assigned_to [ i ] == -1 )
833  {
834  clusters.insert ( std::multimap<G4int,G4int>::value_type ( is_assigned_to [ j ] , i ) );
835  is_assigned_to [ i ] = is_assigned_to [ j ];
836  }
837  else
838  {
839 // i has companion
840  if ( is_assigned_to [ i ] != is_assigned_to [ j ] )
841  {
842 // move companions to the cluster
843 //
844  //std::cout << "combine " << is_assigned_to [ i ] << " to " << is_assigned_to [ j ] << std::endl;
845  std::multimap< G4int , G4int > clusters_tmp;
846  G4int target_cluster_id;
847  if ( is_assigned_to [ i ] > is_assigned_to [ j ] )
848  target_cluster_id = is_assigned_to [ i ];
849  else
850  target_cluster_id = is_assigned_to [ j ];
851 
852  for ( std::multimap< G4int , G4int >::iterator it
853  = clusters.begin() ; it != clusters.end() ; it++ )
854  {
855 
856  //std::cout << it->first << " " << it->second << " " << target_cluster_id << std::endl;
857  if ( it->first == target_cluster_id )
858  {
859  //std::cout << "move " << it->first << " " << it->second << std::endl;
860  is_assigned_to [ it->second ] = is_assigned_to [ j ];
861  clusters_tmp.insert ( std::multimap<G4int,G4int>::value_type ( is_assigned_to [ j ] , it->second ) );
862  }
863  else
864  {
865  clusters_tmp.insert ( std::multimap<G4int,G4int>::value_type ( it->first , it->second ) );
866  }
867  }
868 
869  clusters = clusters_tmp;
870  //id = clusters.rbegin()->first;
871  //id = target_cluster_id;
872  //std::cout << "id " << id << std::endl;
873  }
874  }
875  }
876 
877  //std::cout << "combination " << i << " " << j << std::endl;
878  comb_map.insert( std::multimap<G4int,G4int>::value_type ( i , j ) );
879  cluster_participants.push_back ( j );
880 
881 
882 
883  if ( assign_map.find( cluster_id ) == assign_map.end() )
884  {
885  is_already_belong_some_cluster[i] = true;
886  assign_map.insert ( std::multimap<G4int,G4int>::value_type ( cluster_id , i ) );
887  hasThisCompany = true;
888  }
889  assign_map.insert ( std::multimap<G4int,G4int>::value_type ( cluster_id , j ) );
890  is_already_belong_some_cluster[j] = true;
891 
892  }
893 
894  if ( ichek == i )
895  {
896  nclst++;
897  ichek++;
898  }
899  }
900 
901  if ( cluster_participants.size() > 0 )
902  {
903 // cluster , participant
904  cluster_map.insert ( std::pair < G4int , std::vector < G4int > > ( i , cluster_participants ) );
905  }
906  }
907 // }
908  }
909  if ( hasThisCompany == true ) cluster_id++;
910  }
911 
912  //std::cout << " id " << id << std::endl;
913 
914 // sort
915 // Heavy cluster comes first
916 // size cluster_id
917  std::multimap< G4int , G4int > sorted_cluster_map;
918  for ( G4int i = 0 ; i <= id ; i++ ) // << "<=" because id is highest cluster nubmer.
919  {
920 
921  //std::cout << i << " cluster has " << clusters.count( i ) << " nucleons." << std::endl;
922  sorted_cluster_map.insert ( std::multimap<G4int,G4int>::value_type ( (G4int) clusters.count( i ) , i ) );
923 
924  }
925 
926 
927 // create nucleus from devided clusters
928  std::vector < G4QMDNucleus* > result;
929  for ( std::multimap < G4int , G4int >::reverse_iterator it
930  = sorted_cluster_map.rbegin() ; it != sorted_cluster_map.rend() ; it ++)
931  {
932 
933  //G4cout << "Add Participants to cluseter " << it->second << G4endl;
934 
935  if ( it->first != 0 )
936  {
937  G4QMDNucleus* nucleus = new G4QMDNucleus();
938  for ( std::multimap < G4int , G4int >::iterator itt
939  = clusters.begin() ; itt != clusters.end() ; itt ++)
940  {
941 
942  if ( it->second == itt->first )
943  {
944  nucleus->SetParticipant( system->GetParticipant ( itt->second ) );
945  //G4cout << "Add Participants " << itt->second << " " << system->GetParticipant ( itt->second )->GetPosition() << G4endl;
946  }
947 
948  }
949  result.push_back( nucleus );
950  }
951 
952  }
953 
954 // delete participants from current system
955 
956  for ( std::vector < G4QMDNucleus* > ::iterator it
957  = result.begin() ; it != result.end() ; it++ )
958  {
959  system->SubtractSystem ( *it );
960  }
961 
962  return result;
963 
964 }
965 
966 
967 
969 {
970  SetSystem( system );
971 }
G4ThreeVector GetPosition()
G4double Get_rho0()
void SetParticipant(G4QMDParticipant *particle)
Definition: G4QMDSystem.hh:51
CLHEP::Hep3Vector G4ThreeVector
void SetPosition(G4ThreeVector r)
G4double Get_hbc()
void SubtractSystem(G4QMDSystem *)
Definition: G4QMDSystem.cc:59
G4int GetChargeInUnitOfEplus()
Hep3Vector v() const
G4double G4NeutronHPJENDLHEData::G4double result
void SetNucleus(G4QMDNucleus *aSystem)
int G4int
Definition: G4Types.hh:78
void Cal2BodyQuantities()
G4bool IsPauliBlocked(G4int)
G4ThreeVector GetMomentum()
G4double Get_cpc()
Char_t n[5]
#define G4UniformRand()
Definition: Randomize.hh:87
static G4QMDParameters * GetInstance()
bool G4bool
Definition: G4Types.hh:79
G4double GetPotential(G4int)
G4int GetTotalNumberOfParticipant()
Definition: G4QMDSystem.hh:60
void SetMomentum(G4ThreeVector p)
G4double Get_gamm()
G4double GetTotalPotential()
tuple v
Definition: test.py:18
std::vector< G4QMDNucleus * > DoClusterJudgment()
G4QMDParticipant * GetParticipant(G4int i)
Definition: G4QMDSystem.hh:62
G4LorentzVector Get4Momentum()
double m2() const
void DoPropagation(G4double)
void SetSystem(G4QMDSystem *aSystem)
void SetTotalPotential(G4double x)
Definition: G4QMDNucleus.hh:62
static double erf(double x)
G4double Get_cpw()
G4double Get_cph()
double G4double
Definition: G4Types.hh:76
void CalEnergyAndAngularMomentumInCM()