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

Last change on this file since 847 was 819, checked in by garnier, 16 years ago

import all except CVS

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