source: trunk/source/processes/hadronic/models/qmd/src/G4QMDMeanField.cc @ 962

Last change on this file since 962 was 962, checked in by garnier, 15 years ago

update processes

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