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

Last change on this file was 1347, checked in by garnier, 14 years ago

geant4 tag 9.4

File size: 25.9 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//    R-JQMD
398      G4double Vi = GetPotential( i );
399      G4double p_zero = std::sqrt( p4i.e()*p4i.e() + 2*p4i.m()*Vi);
400      G4ThreeVector betai_R = p4i.v()/p_zero;
401      G4double mi_R = p4i.m()/p_zero;
402//
403      ffr[i] = betai_R;
404      ffp[i] = G4ThreeVector( 0.0 );
405
406      if ( false ) 
407      {
408         ffr[i] = betai;
409         mi_R = 1.0;
410      }
411
412      for ( G4int j = 0 ; j < system->GetTotalNumberOfParticipant() ; j ++ )
413      {
414
415         G4ThreeVector rj = system->GetParticipant( j )->GetPosition(); 
416         G4LorentzVector p4j = system->GetParticipant( j )->Get4Momentum(); 
417
418         G4double eij = p4i.e() + p4j.e(); 
419
420         G4int icharge = system->GetParticipant(i)->GetChargeInUnitOfEplus();
421         G4int jcharge = system->GetParticipant(j)->GetChargeInUnitOfEplus();
422
423         G4int inuc = system->GetParticipant(i)->GetNuc();
424         G4int jnuc = system->GetParticipant(j)->GetNuc();
425
426         G4double ccpp = c0g * rha[j][i]
427                       + c3g * rha[j][i] * ( rh3d[j] + rh3d[i] )
428                       + csg * rha[j][i] * jnuc * inuc
429                           * ( 1. - 2. * std::abs( jcharge - icharge ) )
430                       + cl * rhc[j][i];
431         ccpp *= mi_R;
432
433/*
434           std::cout << c0g << " " << c3g << " " << csg << " " << cl << std::endl;
435           std::cout << "ccpp " << i << " " << j << " " << ccpp << std::endl;
436           std::cout << "rha[j][i] " << rha[j][i] << std::endl;
437           std::cout << "rh3d " << rh3d[j] << " " << rh3d[i] << std::endl;
438           std::cout << "rhc[j][i] " << rhc[j][i] << std::endl;
439*/
440
441         G4double grbb = - rbij[j][i];
442         G4double ccrr = grbb * ccpp / eij;
443
444/*
445           std::cout << "ccrr " << ccrr << std::endl;
446           std::cout << "grbb " << grbb << std::endl;
447*/
448
449
450         G4ThreeVector rij = ri - rj;   
451         G4ThreeVector betaij =  ( p4i + p4j ).v()/eij;   
452
453         G4ThreeVector cij = betaij - betai;   
454
455         ffr[i] = ffr[i] + 2*ccrr* ( rij + grbb*cij );
456           
457         ffp[i] = ffp[i] - 2*ccpp* ( rij + grbb*betaij );
458
459      }
460   }
461
462   //std::cout << "gradu 0 " << ffr[0] << " " << ffp[0] << std::endl;
463   //std::cout << "gradu 1 " << ffr[1] << " " << ffp[1] << std::endl;
464
465}
466
467
468
469G4double G4QMDMeanField::GetPotential( G4int i )
470{
471   G4int n = system->GetTotalNumberOfParticipant();
472
473   G4double rhoa = 0.0;
474   G4double rho3 = 0.0;
475   G4double rhos = 0.0;
476   G4double rhoc = 0.0;
477
478
479   G4int icharge = system->GetParticipant(i)->GetChargeInUnitOfEplus();
480   G4int inuc = system->GetParticipant(i)->GetNuc();
481
482   for ( G4int j = 0 ; j < n ; j ++ )
483   {
484      G4int jcharge = system->GetParticipant(j)->GetChargeInUnitOfEplus();
485      G4int jnuc = system->GetParticipant(j)->GetNuc();
486
487      rhoa += rha[j][i];
488      rhoc += rhe[j][i];
489      rhos += rha[j][i] * jnuc * inuc
490                * ( 1 - 2 * std::abs ( jcharge - icharge ) );
491   }
492
493   rho3 = std::pow ( rhoa , gamm );
494
495   G4double potential = c0 * rhoa
496                      + c3 * rho3
497                      + cs * rhos
498                      + cl * rhoc;
499
500   return potential;
501}
502
503
504
505G4double G4QMDMeanField::GetTotalPotential()
506{
507
508   G4int n = system->GetTotalNumberOfParticipant();
509
510   std::vector < G4double > rhoa ( n , 0.0 ); 
511   std::vector < G4double > rho3 ( n , 0.0 ); 
512   std::vector < G4double > rhos ( n , 0.0 ); 
513   std::vector < G4double > rhoc ( n , 0.0 ); 
514   
515
516   for ( G4int i = 0 ; i < n ; i ++ )
517   {
518      G4int icharge = system->GetParticipant(i)->GetChargeInUnitOfEplus();
519      G4int inuc = system->GetParticipant(i)->GetNuc();
520
521      for ( G4int j = 0 ; j < n ; j ++ )
522      {
523         G4int jcharge = system->GetParticipant(j)->GetChargeInUnitOfEplus();
524         G4int jnuc = system->GetParticipant(j)->GetNuc();
525
526         rhoa[i] += rha[j][i];
527         rhoc[i] += rhe[j][i];
528         rhos[i] += rha[j][i] * jnuc * inuc
529                   * ( 1 - 2 * std::abs ( jcharge - icharge ) );
530      }
531
532      rho3[i] = std::pow ( rhoa[i] , gamm );
533   }
534
535   G4double potential = c0 * std::accumulate( rhoa.begin() , rhoa.end() , 0.0 ) 
536                      + c3 * std::accumulate( rho3.begin() , rho3.end() , 0.0 ) 
537                      + cs * std::accumulate( rhos.begin() , rhos.end() , 0.0 ) 
538                      + cl * std::accumulate( rhoc.begin() , rhoc.end() , 0.0 );
539
540   return potential;
541
542}
543
544
545
546G4double G4QMDMeanField::calPauliBlockingFactor( G4int i )
547{
548
549   G4double pf = 0.0;
550// i is supposed beyond total number of Participant()
551   G4int icharge = system->GetParticipant(i)->GetChargeInUnitOfEplus();
552
553   for ( G4int j = 0 ; j < system->GetTotalNumberOfParticipant() ; j++ )
554   {
555
556      G4int jcharge = system->GetParticipant(j)->GetChargeInUnitOfEplus();
557      G4int jnuc = system->GetParticipant(j)->GetNuc();
558
559      if ( jcharge == icharge && jnuc == 1 )
560      {
561
562/*
563   std::cout << "Pauli i j " << i << " " << j << std::endl;
564   std::cout << "Pauli icharge " << icharge << std::endl;
565   std::cout << "Pauli jcharge " << jcharge << std::endl;
566*/
567         G4double expa = -rr2[i][j]*cpw;   
568
569
570         if ( expa > epsx ) 
571         {
572            expa = expa - pp2[i][j]*cph;
573/*
574   std::cout << "Pauli cph " << cph << std::endl;
575   std::cout << "Pauli pp2 " << pp2[i][j] << std::endl;
576   std::cout << "Pauli expa " <<  expa  << std::endl;
577   std::cout << "Pauli epsx " <<  epsx  << std::endl;
578*/
579            if ( expa > epsx ) 
580            { 
581//   std::cout << "Pauli phase " <<  pf  << std::endl;
582               pf = pf + std::exp ( expa );
583            }
584         }
585      }
586     
587   }
588
589
590   pf  = ( pf - 1.0 ) * cpc;
591
592   //std::cout << "Pauli pf " << pf << std::endl;
593
594   return pf;
595   
596}
597
598
599
600G4bool G4QMDMeanField::IsPauliBlocked( G4int i )
601{
602    G4bool result = false; 
603   
604    if ( system->GetParticipant( i )->GetNuc() == 1 )
605    {
606       G4double pf = calPauliBlockingFactor( i );
607       G4double rand = G4UniformRand(); 
608       if ( pf > rand ) result = true;
609    }
610
611    return result; 
612}
613
614
615
616void G4QMDMeanField::DoPropagation( G4double dt )
617{
618 
619   G4double c2 = 1.0; 
620   G4double c1 = 1.0 - c2; 
621   G4double c3 = 1.0 / 2.0 / c2; 
622
623   G4double dt3 = dt * c3;
624   G4double dt1 = dt * ( c1 - c3 );
625   G4double dt2 = dt * c2;
626
627   CalGraduate(); 
628
629   G4int n = system->GetTotalNumberOfParticipant();
630
631// 1st Step
632
633   std::vector< G4ThreeVector > f0r, f0p;
634   f0r.resize( n );
635   f0p.resize( n );
636
637   for ( G4int i = 0 ; i < n ; i++ )
638   {
639      G4ThreeVector ri = system->GetParticipant( i )->GetPosition(); 
640      G4ThreeVector p3i = system->GetParticipant( i )->GetMomentum(); 
641
642      ri += dt3* ffr[i];
643      p3i += dt3* ffp[i];
644
645      f0r[i] = ffr[i];
646      f0p[i] = ffp[i];
647     
648      system->GetParticipant( i )->SetPosition( ri ); 
649      system->GetParticipant( i )->SetMomentum( p3i ); 
650
651//    we do not need set total momentum by ourselvs
652   }
653
654// 2nd Step
655   Cal2BodyQuantities(); 
656   CalGraduate(); 
657
658   for ( G4int i = 0 ; i < n ; i++ )
659   {
660      G4ThreeVector ri = system->GetParticipant( i )->GetPosition(); 
661      G4ThreeVector p3i = system->GetParticipant( i )->GetMomentum(); 
662
663      ri += dt1* f0r[i] + dt2* ffr[i];
664      p3i += dt1* f0p[i] + dt2* ffp[i];
665
666      system->GetParticipant( i )->SetPosition( ri ); 
667      system->GetParticipant( i )->SetMomentum( p3i ); 
668
669//    we do not need set total momentum by ourselvs
670   }
671
672   Cal2BodyQuantities(); 
673
674
675}
676
677
678
679std::vector< G4QMDNucleus* > G4QMDMeanField::DoClusterJudgment()
680{
681
682   //std::cout << "MeanField DoClusterJudgemnt" << std::endl;
683
684   Cal2BodyQuantities(); 
685
686   G4double cpf2 = std::pow ( 1.5 * pi*pi * std::pow ( 4.0 * pi * wl , -1.5 ) 
687                              , 
688                              2./3. )
689                 * hbc * hbc;
690   G4double rcc2 = rclds*rclds;
691
692   G4int n = system->GetTotalNumberOfParticipant();
693   std::vector < G4double > rhoa;
694   rhoa.resize ( n );
695
696   for ( G4int i = 0 ; i < n ; i++ )
697   {
698      rhoa[i] = 0.0;
699
700      if ( system->GetParticipant( i )->GetBaryonNumber() == 1 ) 
701      {
702      for ( G4int j = 0 ; j < n ; j++ )
703      {
704         if ( system->GetParticipant( j )->GetBaryonNumber() == 1 ) 
705         rhoa[i] += rha[i][j];
706      }
707      }
708
709      rhoa[i] = std::pow ( rhoa[i] + 1 , 1.0/3.0 );
710
711   }
712
713// identification of the cluster
714
715   std::map < G4int , std::vector < G4int > > cluster_map;
716   std::vector < G4bool > is_already_belong_some_cluster;
717
718   //         cluster_id   participant_id
719   std::multimap < G4int , G4int > comb_map; 
720   std::multimap < G4int , G4int > assign_map; 
721   assign_map.clear(); 
722
723   std::vector < G4int > mascl;
724   std::vector < G4int > num;
725   mascl.resize ( n );
726   num.resize ( n );
727   is_already_belong_some_cluster.resize ( n );
728
729   std::vector < G4int > is_assigned_to ( n , -1 );
730   std::multimap < G4int , G4int > clusters;
731
732   for ( G4int i = 0 ; i < n ; i++ )
733   {
734      mascl[i] = 1;
735      num[i] = 1;
736
737      is_already_belong_some_cluster[i] = false;
738   }
739
740
741   G4int nclst = 1;
742   G4int ichek = 1;
743
744
745   G4int id = 0;
746   G4int cluster_id = -1; 
747   for ( G4int i = 0 ; i < n-1 ; i++ )
748   { 
749
750      G4bool hasThisCompany = false;
751//    Check only for bryons?
752//      std::cout << "Check Baryon " << i << std::endl;
753
754      if ( system->GetParticipant( i )->GetBaryonNumber() == 1 ) 
755      {
756
757//      if ( is_already_belong_some_cluster[i] != true ) 
758//      {
759      //G4int j1 = ichek + 1;
760      G4int j1 = i + 1;
761      for ( G4int j = j1 ; j < n ; j++ )
762      {
763
764         std::vector < G4int > cluster_participants;
765         if ( system->GetParticipant( j )->GetBaryonNumber() == 1 ) 
766         {
767         G4double rdist2 = rr2[ i ][ j ];
768         G4double pdist2 = pp2[ i ][ j ];
769         //G4double rdist2 = rr2[ num[i] ][ num[j] ];
770         //G4double pdist2 = pp2[ num[i] ][ num[j] ];
771         G4double pcc2 = cpf2
772                       * ( rhoa[ i ] + rhoa[ j ] )
773                       * ( rhoa[ i ] + rhoa[ j ] );
774
775//       Check phase space: close enough?
776         if ( rdist2 < rcc2 && pdist2 < pcc2 ) 
777         {
778
779/*
780            std::cout << "G4QMDRESULT "
781                      << i << " " << j << " " << id << " "
782                      << is_assigned_to [ i ] << " " << is_assigned_to [ j ]
783                      << std::endl;
784*/
785
786            if ( is_assigned_to [ j ] == -1 )
787            {
788               if ( is_assigned_to [ i ] == -1 )
789               {
790                  if ( clusters.size() != 0 )
791                  {
792                     id = clusters.rbegin()->first + 1; 
793                     //std::cout << "id is increare " << id << std::endl;
794                  }
795                  else
796                  {
797                     id = 0;
798                  }
799                  clusters.insert ( std::multimap<G4int,G4int>::value_type ( id , i ) );
800                  is_assigned_to [ i ] = id; 
801                  clusters.insert ( std::multimap<G4int,G4int>::value_type ( id , j ) );
802                  is_assigned_to [ j ] = id; 
803               }
804               else
805               {
806                  clusters.insert ( std::multimap<G4int,G4int>::value_type ( is_assigned_to [ i ] , j ) );
807                  is_assigned_to [ j ] = is_assigned_to [ i ]; 
808               }
809            } 
810            else
811            {
812//             j is already belong to some cluester
813               if ( is_assigned_to [ i ] == -1 )
814               {
815                  clusters.insert ( std::multimap<G4int,G4int>::value_type ( is_assigned_to [ j ] , i ) );
816                  is_assigned_to [ i ] = is_assigned_to [ j ]; 
817               }
818               else
819               {
820//             i has companion 
821                  if ( is_assigned_to [ i ] != is_assigned_to [ j ] )
822                  {
823//                   move companions to the cluster   
824//                 
825                     //std::cout <<  "combine " << is_assigned_to [ i ] << " to " << is_assigned_to [ j ] << std::endl;
826                     std::multimap< G4int , G4int > clusters_tmp;
827                     G4int target_cluster_id; 
828                     if ( is_assigned_to [ i ] > is_assigned_to [ j ] )
829                        target_cluster_id = is_assigned_to [ i ];
830                     else
831                        target_cluster_id = is_assigned_to [ j ];
832                     
833                     for ( std::multimap< G4int , G4int >::iterator it
834                         = clusters.begin() ; it != clusters.end() ; it++ ) 
835                     {
836
837                        //std::cout << it->first << " " << it->second  << " " << target_cluster_id << std::endl;
838                        if ( it->first == target_cluster_id )
839                        { 
840                           //std::cout << "move " <<  it->first << " " << it->second << std::endl;
841                           is_assigned_to [ it->second ] = is_assigned_to [ j ]; 
842                           clusters_tmp.insert ( std::multimap<G4int,G4int>::value_type (  is_assigned_to [ j ] , it->second ) );
843                        }
844                        else
845                        {
846                           clusters_tmp.insert ( std::multimap<G4int,G4int>::value_type ( it->first , it->second ) );
847                        }
848                     }
849
850                     clusters = clusters_tmp;
851                     //id = clusters.rbegin()->first;
852                     //id = target_cluster_id;
853                     //std::cout <<  "id " << id << std::endl;
854                  }
855               }
856            } 
857
858            //std::cout << "combination " << i << " " << j << std::endl;
859            comb_map.insert( std::multimap<G4int,G4int>::value_type ( i , j ) ); 
860            cluster_participants.push_back ( j );
861
862
863
864            if ( assign_map.find( cluster_id ) == assign_map.end() ) 
865            { 
866               is_already_belong_some_cluster[i] = true;
867               assign_map.insert ( std::multimap<G4int,G4int>::value_type ( cluster_id , i ) );
868               hasThisCompany = true; 
869            } 
870            assign_map.insert ( std::multimap<G4int,G4int>::value_type ( cluster_id , j ) ); 
871            is_already_belong_some_cluster[j] = true;
872
873         } 
874
875         if ( ichek == i ) 
876         {
877            nclst++;
878            ichek++;
879         }
880         }
881
882         if ( cluster_participants.size() > 0 ) 
883         {
884//                                         cluster , participant
885            cluster_map.insert ( std::pair < G4int , std::vector < G4int > > ( i , cluster_participants ) );
886         }
887      }
888//      }
889      } 
890      if ( hasThisCompany == true ) cluster_id++;
891   } 
892
893   //std::cout << " id " << id << std::endl;
894
895// sort
896// Heavy cluster comes first
897//                size    cluster_id
898   std::multimap< G4int , G4int > sorted_cluster_map; 
899   for ( G4int i = 0 ; i <= id ; i++ )  // << "<=" because id is highest cluster nubmer.
900   {
901 
902      //std::cout << i << " cluster has " << clusters.count( i )  << " nucleons." << std::endl;
903      sorted_cluster_map.insert ( std::multimap<G4int,G4int>::value_type ( clusters.count( i ) , i ) );
904
905   }
906
907
908// create nucleus from devided clusters
909   std::vector < G4QMDNucleus* > result;
910   for ( std::multimap < G4int , G4int >::reverse_iterator it
911       = sorted_cluster_map.rbegin() ; it != sorted_cluster_map.rend()  ; it ++) 
912   {
913
914      //G4cout << "Add Participants to cluseter " << it->second << G4endl;
915
916      if ( it->first != 0 ) 
917      {
918         G4QMDNucleus* nucleus = new G4QMDNucleus();
919         for ( std::multimap < G4int , G4int >::iterator itt
920             = clusters.begin() ; itt != clusters.end()  ; itt ++)
921         {
922
923            if ( it->second == itt->first )
924            {
925               nucleus->SetParticipant( system->GetParticipant ( itt->second ) );
926               //G4cout << "Add Participants " << itt->second << " "  << system->GetParticipant ( itt->second )->GetPosition() << G4endl;
927             }
928
929         }
930         result.push_back( nucleus );
931      } 
932
933   }
934
935// delete participants from current system
936
937   for ( std::vector < G4QMDNucleus* > ::iterator it
938       = result.begin() ; it != result.end() ; it++ ) 
939   {
940      system->SubtractSystem ( *it ); 
941   }
942   
943   return result;
944   
945}
946
947
948
949void G4QMDMeanField::Update() 
950{ 
951   SetSystem( system ); 
952}
Note: See TracBrowser for help on using the repository browser.