00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 #include <map>
00030 #include <algorithm>
00031 #include <numeric>
00032
00033 #include <CLHEP/Random/Stat.h>
00034
00035 #include "G4QMDMeanField.hh"
00036 #include "G4QMDParameters.hh"
00037
00038 #include "G4PhysicalConstants.hh"
00039 #include "Randomize.hh"
00040
00041 G4QMDMeanField::G4QMDMeanField()
00042 : rclds ( 4.0 )
00043 , epsx ( -20.0 )
00044 , epscl ( 0.0001 )
00045 , irelcr ( 1 )
00046 {
00047
00048 G4QMDParameters* parameters = G4QMDParameters::GetInstance();
00049 wl = parameters->Get_wl();
00050 cl = parameters->Get_cl();
00051 rho0 = parameters->Get_rho0();
00052 hbc = parameters->Get_hbc();
00053 gamm = parameters->Get_gamm();
00054
00055 cpw = parameters->Get_cpw();
00056 cph = parameters->Get_cph();
00057 cpc = parameters->Get_cpc();
00058
00059 c0 = parameters->Get_c0();
00060 c3 = parameters->Get_c3();
00061 cs = parameters->Get_cs();
00062
00063
00064 c0w = 1.0/4.0/wl;
00065
00066 c0sw = std::sqrt( c0w );
00067 clw = 2.0 / std::sqrt ( 4.0 * pi * wl );
00068
00069
00070 c0g = - c0 / ( 2.0 * wl );
00071 c3g = - c3 / ( 4.0 * wl ) * gamm;
00072 csg = - cs / ( 2.0 * wl );
00073 pag = gamm - 1;
00074
00075 }
00076
00077
00078
00079 G4QMDMeanField::~G4QMDMeanField()
00080 {
00081 ;
00082 }
00083
00084
00085
00086 void G4QMDMeanField::SetSystem ( G4QMDSystem* aSystem )
00087 {
00088
00089
00090
00091 system = aSystem;
00092
00093 G4int n = system->GetTotalNumberOfParticipant();
00094
00095 pp2.clear();
00096 rr2.clear();
00097 rbij.clear();
00098 rha.clear();
00099 rhe.clear();
00100 rhc.clear();
00101
00102 rr2.resize( n );
00103 pp2.resize( n );
00104 rbij.resize( n );
00105 rha.resize( n );
00106 rhe.resize( n );
00107 rhc.resize( n );
00108
00109 for ( int i = 0 ; i < n ; i++ )
00110 {
00111 rr2[i].resize( n );
00112 pp2[i].resize( n );
00113 rbij[i].resize( n );
00114 rha[i].resize( n );
00115 rhe[i].resize( n );
00116 rhc[i].resize( n );
00117 }
00118
00119
00120 ffr.clear();
00121 ffp.clear();
00122 rh3d.clear();
00123
00124 ffr.resize( n );
00125 ffp.resize( n );
00126 rh3d.resize( n );
00127
00128 Cal2BodyQuantities();
00129
00130 }
00131
00132 void G4QMDMeanField::SetNucleus ( G4QMDNucleus* aNucleus )
00133 {
00134
00135
00136
00137 SetSystem( aNucleus );
00138
00139 G4double totalPotential = GetTotalPotential();
00140 aNucleus->SetTotalPotential( totalPotential );
00141
00142 aNucleus->CalEnergyAndAngularMomentumInCM();
00143
00144 }
00145
00146
00147
00148 void G4QMDMeanField::Cal2BodyQuantities()
00149 {
00150
00151 if ( system->GetTotalNumberOfParticipant() < 2 ) return;
00152
00153 for ( G4int j = 1 ; j < system->GetTotalNumberOfParticipant() ; j++ )
00154 {
00155
00156 G4ThreeVector rj = system->GetParticipant( j )->GetPosition();
00157 G4LorentzVector p4j = system->GetParticipant( j )->Get4Momentum();
00158
00159 for ( G4int i = 0 ; i < j ; i++ )
00160 {
00161
00162 G4ThreeVector ri = system->GetParticipant( i )->GetPosition();
00163 G4LorentzVector p4i = system->GetParticipant( i )->Get4Momentum();
00164
00165 G4ThreeVector rij = ri - rj;
00166 G4ThreeVector pij = (p4i - p4j).v();
00167 G4LorentzVector p4ij = p4i - p4j;
00168 G4ThreeVector bij = ( p4i + p4j ).boostVector();
00169 G4double gammaij = ( p4i + p4j ).gamma();
00170
00171 G4double eij = ( p4i + p4j ).e();
00172
00173 G4double rbrb = rij*bij;
00174
00175 G4double rij2 = rij*rij;
00176 G4double pij2 = pij*pij;
00177
00178 rbrb = irelcr * rbrb;
00179 G4double gamma2_ij = gammaij*gammaij;
00180
00181
00182 rr2[i][j] = rij2 + gamma2_ij * rbrb*rbrb;
00183 rr2[j][i] = rr2[i][j];
00184
00185 rbij[i][j] = gamma2_ij * rbrb;
00186 rbij[j][i] = - rbij[i][j];
00187
00188 pp2[i][j] = pij2
00189 + irelcr * ( - std::pow ( p4i.e() - p4j.e() , 2 )
00190 + gamma2_ij * std::pow ( ( ( p4i.m2() - p4j.m2() ) / eij ) , 2 ) );
00191
00192
00193 pp2[j][i] = pp2[i][j];
00194
00195
00196
00197 G4double expa1 = - rr2[i][j] * c0w;
00198
00199 G4double rh1;
00200 if ( expa1 > epsx )
00201 {
00202 rh1 = std::exp( expa1 );
00203 }
00204 else
00205 {
00206 rh1 = 0.0;
00207 }
00208
00209 G4int ibry = system->GetParticipant(i)->GetBaryonNumber();
00210 G4int jbry = system->GetParticipant(j)->GetBaryonNumber();
00211
00212
00213 rha[i][j] = ibry*jbry*rh1;
00214 rha[j][i] = rha[i][j];
00215
00216
00217
00218 G4double rrs2 = rr2[i][j] + epscl;
00219 G4double rrs = std::sqrt ( rrs2 );
00220
00221 G4int icharge = system->GetParticipant(i)->GetChargeInUnitOfEplus();
00222 G4int jcharge = system->GetParticipant(j)->GetChargeInUnitOfEplus();
00223
00224 G4double erf = 0.0;
00225
00226 if ( rrs*c0sw < 5.8 )
00227 erf = CLHEP::HepStat::erf ( rrs*c0sw );
00228 else
00229 erf = 1.0;
00230
00231 G4double erfij = erf/rrs;
00232
00233
00234 rhe[i][j] = icharge*jcharge * erfij;
00235
00236 rhe[j][i] = rhe[i][j];
00237
00238 rhc[i][j] = icharge*jcharge * ( - erfij + clw * rh1 ) / rrs2;
00239
00240 rhc[j][i] = rhc[i][j];
00241
00242 }
00243 }
00244 }
00245
00246
00247
00248 void G4QMDMeanField::Cal2BodyQuantities( G4int i )
00249 {
00250
00251
00252
00253 G4ThreeVector ri = system->GetParticipant( i )->GetPosition();
00254 G4LorentzVector p4i = system->GetParticipant( i )->Get4Momentum();
00255
00256
00257 for ( G4int j = 0 ; j < system->GetTotalNumberOfParticipant() ; j ++ )
00258 {
00259 if ( j == i ) continue;
00260
00261 G4ThreeVector rj = system->GetParticipant( j )->GetPosition();
00262 G4LorentzVector p4j = system->GetParticipant( j )->Get4Momentum();
00263
00264 G4ThreeVector rij = ri - rj;
00265 G4ThreeVector pij = (p4i - p4j).v();
00266 G4LorentzVector p4ij = p4i - p4j;
00267 G4ThreeVector bij = ( p4i + p4j ).boostVector();
00268 G4double gammaij = ( p4i + p4j ).gamma();
00269
00270 G4double eij = ( p4i + p4j ).e();
00271
00272 G4double rbrb = rij*bij;
00273
00274 G4double rij2 = rij*rij;
00275 G4double pij2 = pij*pij;
00276
00277 rbrb = irelcr * rbrb;
00278 G4double gamma2_ij = gammaij*gammaij;
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305
00306 rr2[i][j] = rij2 + gamma2_ij * rbrb*rbrb;
00307 rr2[j][i] = rr2[i][j];
00308
00309 rbij[i][j] = gamma2_ij * rbrb;
00310 rbij[j][i] = - rbij[i][j];
00311
00312 pp2[i][j] = pij2
00313 + irelcr * ( - std::pow ( p4i.e() - p4j.e() , 2 )
00314 + gamma2_ij * std::pow ( ( ( p4i.m2() - p4j.m2() ) / eij ) , 2 ) );
00315
00316 pp2[j][i] = pp2[i][j];
00317
00318
00319
00320 G4double expa1 = - rr2[i][j] * c0w;
00321
00322 G4double rh1;
00323 if ( expa1 > epsx )
00324 {
00325 rh1 = std::exp( expa1 );
00326 }
00327 else
00328 {
00329 rh1 = 0.0;
00330 }
00331
00332 G4int ibry = system->GetParticipant(i)->GetBaryonNumber();
00333 G4int jbry = system->GetParticipant(j)->GetBaryonNumber();
00334
00335
00336 rha[i][j] = ibry*jbry*rh1;
00337 rha[j][i] = rha[i][j];
00338
00339
00340
00341 G4double rrs2 = rr2[i][j] + epscl;
00342 G4double rrs = std::sqrt ( rrs2 );
00343
00344 G4int icharge = system->GetParticipant(i)->GetChargeInUnitOfEplus();
00345 G4int jcharge = system->GetParticipant(j)->GetChargeInUnitOfEplus();
00346
00347 G4double erf = 0.0;
00348
00349 if ( rrs*c0sw < 5.8 )
00350 erf = CLHEP::HepStat::erf ( rrs*c0sw );
00351 else
00352 erf = 1.0;
00353
00354 G4double erfij = erf/rrs;
00355
00356
00357 rhe[i][j] = icharge*jcharge * erfij;
00358
00359 rhe[j][i] = rhe[i][j];
00360
00361
00362
00363 rhc[i][j] = icharge*jcharge * ( - erfij + clw * rh1 ) / rrs2;
00364
00365 rhc[j][i] = rhc[i][j];
00366
00367 }
00368
00369 }
00370
00371
00372
00373 void G4QMDMeanField::CalGraduate()
00374 {
00375
00376 ffr.resize( system->GetTotalNumberOfParticipant() );
00377 ffp.resize( system->GetTotalNumberOfParticipant() );
00378 rh3d.resize( system->GetTotalNumberOfParticipant() );
00379
00380 for ( G4int i = 0 ; i < system->GetTotalNumberOfParticipant() ; i ++ )
00381 {
00382 G4double rho3 = 0.0;
00383 for ( G4int j = 0 ; j < system->GetTotalNumberOfParticipant() ; j ++ )
00384 {
00385 rho3 += rha[j][i];
00386 }
00387 rh3d[i] = std::pow ( rho3 , pag );
00388 }
00389
00390
00391 for ( G4int i = 0 ; i < system->GetTotalNumberOfParticipant() ; i ++ )
00392 {
00393
00394 G4ThreeVector ri = system->GetParticipant( i )->GetPosition();
00395 G4LorentzVector p4i = system->GetParticipant( i )->Get4Momentum();
00396
00397 G4ThreeVector betai = p4i.v()/p4i.e();
00398
00399
00400 G4double Vi = GetPotential( i );
00401 G4double p_zero = std::sqrt( p4i.e()*p4i.e() + 2*p4i.m()*Vi);
00402 G4ThreeVector betai_R = p4i.v()/p_zero;
00403 G4double mi_R = p4i.m()/p_zero;
00404
00405 ffr[i] = betai_R;
00406 ffp[i] = G4ThreeVector( 0.0 );
00407
00408 if ( false )
00409 {
00410 ffr[i] = betai;
00411 mi_R = 1.0;
00412 }
00413
00414 for ( G4int j = 0 ; j < system->GetTotalNumberOfParticipant() ; j ++ )
00415 {
00416
00417 G4ThreeVector rj = system->GetParticipant( j )->GetPosition();
00418 G4LorentzVector p4j = system->GetParticipant( j )->Get4Momentum();
00419
00420 G4double eij = p4i.e() + p4j.e();
00421
00422 G4int icharge = system->GetParticipant(i)->GetChargeInUnitOfEplus();
00423 G4int jcharge = system->GetParticipant(j)->GetChargeInUnitOfEplus();
00424
00425 G4int inuc = system->GetParticipant(i)->GetNuc();
00426 G4int jnuc = system->GetParticipant(j)->GetNuc();
00427
00428 G4double ccpp = c0g * rha[j][i]
00429 + c3g * rha[j][i] * ( rh3d[j] + rh3d[i] )
00430 + csg * rha[j][i] * jnuc * inuc
00431 * ( 1. - 2. * std::abs( jcharge - icharge ) )
00432 + cl * rhc[j][i];
00433 ccpp *= mi_R;
00434
00435
00436
00437
00438
00439
00440
00441
00442
00443 G4double grbb = - rbij[j][i];
00444 G4double ccrr = grbb * ccpp / eij;
00445
00446
00447
00448
00449
00450
00451
00452 G4ThreeVector rij = ri - rj;
00453 G4ThreeVector betaij = ( p4i + p4j ).v()/eij;
00454
00455 G4ThreeVector cij = betaij - betai;
00456
00457 ffr[i] = ffr[i] + 2*ccrr* ( rij + grbb*cij );
00458
00459 ffp[i] = ffp[i] - 2*ccpp* ( rij + grbb*betaij );
00460
00461 }
00462 }
00463
00464
00465
00466
00467 }
00468
00469
00470
00471 G4double G4QMDMeanField::GetPotential( G4int i )
00472 {
00473 G4int n = system->GetTotalNumberOfParticipant();
00474
00475 G4double rhoa = 0.0;
00476 G4double rho3 = 0.0;
00477 G4double rhos = 0.0;
00478 G4double rhoc = 0.0;
00479
00480
00481 G4int icharge = system->GetParticipant(i)->GetChargeInUnitOfEplus();
00482 G4int inuc = system->GetParticipant(i)->GetNuc();
00483
00484 for ( G4int j = 0 ; j < n ; j ++ )
00485 {
00486 G4int jcharge = system->GetParticipant(j)->GetChargeInUnitOfEplus();
00487 G4int jnuc = system->GetParticipant(j)->GetNuc();
00488
00489 rhoa += rha[j][i];
00490 rhoc += rhe[j][i];
00491 rhos += rha[j][i] * jnuc * inuc
00492 * ( 1 - 2 * std::abs ( jcharge - icharge ) );
00493 }
00494
00495 rho3 = std::pow ( rhoa , gamm );
00496
00497 G4double potential = c0 * rhoa
00498 + c3 * rho3
00499 + cs * rhos
00500 + cl * rhoc;
00501
00502 return potential;
00503 }
00504
00505
00506
00507 G4double G4QMDMeanField::GetTotalPotential()
00508 {
00509
00510 G4int n = system->GetTotalNumberOfParticipant();
00511
00512 std::vector < G4double > rhoa ( n , 0.0 );
00513 std::vector < G4double > rho3 ( n , 0.0 );
00514 std::vector < G4double > rhos ( n , 0.0 );
00515 std::vector < G4double > rhoc ( n , 0.0 );
00516
00517
00518 for ( G4int i = 0 ; i < n ; i ++ )
00519 {
00520 G4int icharge = system->GetParticipant(i)->GetChargeInUnitOfEplus();
00521 G4int inuc = system->GetParticipant(i)->GetNuc();
00522
00523 for ( G4int j = 0 ; j < n ; j ++ )
00524 {
00525 G4int jcharge = system->GetParticipant(j)->GetChargeInUnitOfEplus();
00526 G4int jnuc = system->GetParticipant(j)->GetNuc();
00527
00528 rhoa[i] += rha[j][i];
00529 rhoc[i] += rhe[j][i];
00530 rhos[i] += rha[j][i] * jnuc * inuc
00531 * ( 1 - 2 * std::abs ( jcharge - icharge ) );
00532 }
00533
00534 rho3[i] = std::pow ( rhoa[i] , gamm );
00535 }
00536
00537 G4double potential = c0 * std::accumulate( rhoa.begin() , rhoa.end() , 0.0 )
00538 + c3 * std::accumulate( rho3.begin() , rho3.end() , 0.0 )
00539 + cs * std::accumulate( rhos.begin() , rhos.end() , 0.0 )
00540 + cl * std::accumulate( rhoc.begin() , rhoc.end() , 0.0 );
00541
00542 return potential;
00543
00544 }
00545
00546
00547
00548 G4double G4QMDMeanField::calPauliBlockingFactor( G4int i )
00549 {
00550
00551 G4double pf = 0.0;
00552
00553 G4int icharge = system->GetParticipant(i)->GetChargeInUnitOfEplus();
00554
00555 for ( G4int j = 0 ; j < system->GetTotalNumberOfParticipant() ; j++ )
00556 {
00557
00558 G4int jcharge = system->GetParticipant(j)->GetChargeInUnitOfEplus();
00559 G4int jnuc = system->GetParticipant(j)->GetNuc();
00560
00561 if ( jcharge == icharge && jnuc == 1 )
00562 {
00563
00564
00565
00566
00567
00568
00569 G4double expa = -rr2[i][j]*cpw;
00570
00571
00572 if ( expa > epsx )
00573 {
00574 expa = expa - pp2[i][j]*cph;
00575
00576
00577
00578
00579
00580
00581 if ( expa > epsx )
00582 {
00583
00584 pf = pf + std::exp ( expa );
00585 }
00586 }
00587 }
00588
00589 }
00590
00591
00592 pf = ( pf - 1.0 ) * cpc;
00593
00594
00595
00596 return pf;
00597
00598 }
00599
00600
00601
00602 G4bool G4QMDMeanField::IsPauliBlocked( G4int i )
00603 {
00604 G4bool result = false;
00605
00606 if ( system->GetParticipant( i )->GetNuc() == 1 )
00607 {
00608 G4double pf = calPauliBlockingFactor( i );
00609 G4double rand = G4UniformRand();
00610 if ( pf > rand ) result = true;
00611 }
00612
00613 return result;
00614 }
00615
00616
00617
00618 void G4QMDMeanField::DoPropagation( G4double dt )
00619 {
00620
00621 G4double cc2 = 1.0;
00622 G4double cc1 = 1.0 - cc2;
00623 G4double cc3 = 1.0 / 2.0 / cc2;
00624
00625 G4double dt3 = dt * cc3;
00626 G4double dt1 = dt * ( cc1 - cc3 );
00627 G4double dt2 = dt * cc2;
00628
00629 CalGraduate();
00630
00631 G4int n = system->GetTotalNumberOfParticipant();
00632
00633
00634
00635 std::vector< G4ThreeVector > f0r, f0p;
00636 f0r.resize( n );
00637 f0p.resize( n );
00638
00639 for ( G4int i = 0 ; i < n ; i++ )
00640 {
00641 G4ThreeVector ri = system->GetParticipant( i )->GetPosition();
00642 G4ThreeVector p3i = system->GetParticipant( i )->GetMomentum();
00643
00644 ri += dt3* ffr[i];
00645 p3i += dt3* ffp[i];
00646
00647 f0r[i] = ffr[i];
00648 f0p[i] = ffp[i];
00649
00650 system->GetParticipant( i )->SetPosition( ri );
00651 system->GetParticipant( i )->SetMomentum( p3i );
00652
00653
00654 }
00655
00656
00657 Cal2BodyQuantities();
00658 CalGraduate();
00659
00660 for ( G4int i = 0 ; i < n ; i++ )
00661 {
00662 G4ThreeVector ri = system->GetParticipant( i )->GetPosition();
00663 G4ThreeVector p3i = system->GetParticipant( i )->GetMomentum();
00664
00665 ri += dt1* f0r[i] + dt2* ffr[i];
00666 p3i += dt1* f0p[i] + dt2* ffp[i];
00667
00668 system->GetParticipant( i )->SetPosition( ri );
00669 system->GetParticipant( i )->SetMomentum( p3i );
00670
00671
00672 }
00673
00674 Cal2BodyQuantities();
00675
00676 }
00677
00678
00679
00680 std::vector< G4QMDNucleus* > G4QMDMeanField::DoClusterJudgment()
00681 {
00682
00683
00684
00685 Cal2BodyQuantities();
00686
00687 G4double cpf2 = std::pow ( 1.5 * pi*pi * std::pow ( 4.0 * pi * wl , -1.5 )
00688 ,
00689 2./3. )
00690 * hbc * hbc;
00691 G4double rcc2 = rclds*rclds;
00692
00693 G4int n = system->GetTotalNumberOfParticipant();
00694 std::vector < G4double > rhoa;
00695 rhoa.resize ( n );
00696
00697 for ( G4int i = 0 ; i < n ; i++ )
00698 {
00699 rhoa[i] = 0.0;
00700
00701 if ( system->GetParticipant( i )->GetBaryonNumber() == 1 )
00702 {
00703 for ( G4int j = 0 ; j < n ; j++ )
00704 {
00705 if ( system->GetParticipant( j )->GetBaryonNumber() == 1 )
00706 rhoa[i] += rha[i][j];
00707 }
00708 }
00709
00710 rhoa[i] = std::pow ( rhoa[i] + 1 , 1.0/3.0 );
00711
00712 }
00713
00714
00715
00716 std::map < G4int , std::vector < G4int > > cluster_map;
00717 std::vector < G4bool > is_already_belong_some_cluster;
00718
00719
00720 std::multimap < G4int , G4int > comb_map;
00721 std::multimap < G4int , G4int > assign_map;
00722 assign_map.clear();
00723
00724 std::vector < G4int > mascl;
00725 std::vector < G4int > num;
00726 mascl.resize ( n );
00727 num.resize ( n );
00728 is_already_belong_some_cluster.resize ( n );
00729
00730 std::vector < G4int > is_assigned_to ( n , -1 );
00731 std::multimap < G4int , G4int > clusters;
00732
00733 for ( G4int i = 0 ; i < n ; i++ )
00734 {
00735 mascl[i] = 1;
00736 num[i] = 1;
00737
00738 is_already_belong_some_cluster[i] = false;
00739 }
00740
00741
00742 G4int nclst = 1;
00743 G4int ichek = 1;
00744
00745
00746 G4int id = 0;
00747 G4int cluster_id = -1;
00748 for ( G4int i = 0 ; i < n-1 ; i++ )
00749 {
00750
00751 G4bool hasThisCompany = false;
00752
00753
00754
00755 if ( system->GetParticipant( i )->GetBaryonNumber() == 1 )
00756 {
00757
00758
00759
00760
00761 G4int j1 = i + 1;
00762 for ( G4int j = j1 ; j < n ; j++ )
00763 {
00764
00765 std::vector < G4int > cluster_participants;
00766 if ( system->GetParticipant( j )->GetBaryonNumber() == 1 )
00767 {
00768 G4double rdist2 = rr2[ i ][ j ];
00769 G4double pdist2 = pp2[ i ][ j ];
00770
00771
00772 G4double pcc2 = cpf2
00773 * ( rhoa[ i ] + rhoa[ j ] )
00774 * ( rhoa[ i ] + rhoa[ j ] );
00775
00776
00777 if ( rdist2 < rcc2 && pdist2 < pcc2 )
00778 {
00779
00780
00781
00782
00783
00784
00785
00786
00787 if ( is_assigned_to [ j ] == -1 )
00788 {
00789 if ( is_assigned_to [ i ] == -1 )
00790 {
00791 if ( clusters.size() != 0 )
00792 {
00793 id = clusters.rbegin()->first + 1;
00794
00795 }
00796 else
00797 {
00798 id = 0;
00799 }
00800 clusters.insert ( std::multimap<G4int,G4int>::value_type ( id , i ) );
00801 is_assigned_to [ i ] = id;
00802 clusters.insert ( std::multimap<G4int,G4int>::value_type ( id , j ) );
00803 is_assigned_to [ j ] = id;
00804 }
00805 else
00806 {
00807 clusters.insert ( std::multimap<G4int,G4int>::value_type ( is_assigned_to [ i ] , j ) );
00808 is_assigned_to [ j ] = is_assigned_to [ i ];
00809 }
00810 }
00811 else
00812 {
00813
00814 if ( is_assigned_to [ i ] == -1 )
00815 {
00816 clusters.insert ( std::multimap<G4int,G4int>::value_type ( is_assigned_to [ j ] , i ) );
00817 is_assigned_to [ i ] = is_assigned_to [ j ];
00818 }
00819 else
00820 {
00821
00822 if ( is_assigned_to [ i ] != is_assigned_to [ j ] )
00823 {
00824
00825
00826
00827 std::multimap< G4int , G4int > clusters_tmp;
00828 G4int target_cluster_id;
00829 if ( is_assigned_to [ i ] > is_assigned_to [ j ] )
00830 target_cluster_id = is_assigned_to [ i ];
00831 else
00832 target_cluster_id = is_assigned_to [ j ];
00833
00834 for ( std::multimap< G4int , G4int >::iterator it
00835 = clusters.begin() ; it != clusters.end() ; it++ )
00836 {
00837
00838
00839 if ( it->first == target_cluster_id )
00840 {
00841
00842 is_assigned_to [ it->second ] = is_assigned_to [ j ];
00843 clusters_tmp.insert ( std::multimap<G4int,G4int>::value_type ( is_assigned_to [ j ] , it->second ) );
00844 }
00845 else
00846 {
00847 clusters_tmp.insert ( std::multimap<G4int,G4int>::value_type ( it->first , it->second ) );
00848 }
00849 }
00850
00851 clusters = clusters_tmp;
00852
00853
00854
00855 }
00856 }
00857 }
00858
00859
00860 comb_map.insert( std::multimap<G4int,G4int>::value_type ( i , j ) );
00861 cluster_participants.push_back ( j );
00862
00863
00864
00865 if ( assign_map.find( cluster_id ) == assign_map.end() )
00866 {
00867 is_already_belong_some_cluster[i] = true;
00868 assign_map.insert ( std::multimap<G4int,G4int>::value_type ( cluster_id , i ) );
00869 hasThisCompany = true;
00870 }
00871 assign_map.insert ( std::multimap<G4int,G4int>::value_type ( cluster_id , j ) );
00872 is_already_belong_some_cluster[j] = true;
00873
00874 }
00875
00876 if ( ichek == i )
00877 {
00878 nclst++;
00879 ichek++;
00880 }
00881 }
00882
00883 if ( cluster_participants.size() > 0 )
00884 {
00885
00886 cluster_map.insert ( std::pair < G4int , std::vector < G4int > > ( i , cluster_participants ) );
00887 }
00888 }
00889
00890 }
00891 if ( hasThisCompany == true ) cluster_id++;
00892 }
00893
00894
00895
00896
00897
00898
00899 std::multimap< G4int , G4int > sorted_cluster_map;
00900 for ( G4int i = 0 ; i <= id ; i++ )
00901 {
00902
00903
00904 sorted_cluster_map.insert ( std::multimap<G4int,G4int>::value_type ( (G4int) clusters.count( i ) , i ) );
00905
00906 }
00907
00908
00909
00910 std::vector < G4QMDNucleus* > result;
00911 for ( std::multimap < G4int , G4int >::reverse_iterator it
00912 = sorted_cluster_map.rbegin() ; it != sorted_cluster_map.rend() ; it ++)
00913 {
00914
00915
00916
00917 if ( it->first != 0 )
00918 {
00919 G4QMDNucleus* nucleus = new G4QMDNucleus();
00920 for ( std::multimap < G4int , G4int >::iterator itt
00921 = clusters.begin() ; itt != clusters.end() ; itt ++)
00922 {
00923
00924 if ( it->second == itt->first )
00925 {
00926 nucleus->SetParticipant( system->GetParticipant ( itt->second ) );
00927
00928 }
00929
00930 }
00931 result.push_back( nucleus );
00932 }
00933
00934 }
00935
00936
00937
00938 for ( std::vector < G4QMDNucleus* > ::iterator it
00939 = result.begin() ; it != result.end() ; it++ )
00940 {
00941 system->SubtractSystem ( *it );
00942 }
00943
00944 return result;
00945
00946 }
00947
00948
00949
00950 void G4QMDMeanField::Update()
00951 {
00952 SetSystem( system );
00953 }