DSQSS  1.1
dla.cc
説明を見る。
00001 //
00002 // バグ・今後の改良項目:
00003 //
00004 // サイトタイプが複数存在する場合 SiteProperty::EBASE が一般には
00005 // まちまちになるが,そのとき,平均ループ長が帯磁率に比例しなくなる.
00006 //
00007 // 粒子数が保存しない場合.つまりワームがもとの位置に戻ってくる前に
00008 // 通常のバーテックスで消滅するような場合もサポートすること.
00009 //
00010 
00011 //######################################################################
00012 //####
00013 //####  World-Line Monte Carlo simulation
00014 //####                       by the Directed-Loop Algorithm 
00015 //####
00016 //####                                 Mar.03 / 2005, Naoki Kawashima
00017 //####
00018 //######################################################################
00019 
00020 //######################################################################
00021 //####
00022 //####  World-Line Monte Carlo simulation
00023 //####                       by the non-Vertex Directed-Loop Algorithm 
00024 //####
00025 //####                                 Nov.11 / 2007, Yasuyuki Kato
00026 //####
00027 //######################################################################
00028 
00049 #define Version 1.1
00050 //#include <mpi.h>
00051 #include "dla.hpp"
00052 
00053 #ifdef GRAPHIC
00054 #include "graphic.cc"
00055 #endif
00056 
00057 #ifdef CHAINJOB
00058 #define killtime 3600 //(sec.)
00059 //#define killtime 32400
00060 #include "chainjob.hpp"
00061 #endif
00062 
00063 #ifdef PAPI
00064 #include <papi.h>
00065 #endif
00066 
00067 //######################################################################
00068 
00069 //int main(int argc, char** argv) {
00070 int main(int argc, char* argv[]) {
00071 
00072   // MPI  Initialization & Getting Size and Rank
00073 #ifdef MULTI
00074 #ifdef KASHIWA
00075     MPI_Init( &argc , &argv );
00076     MPI_Comm_size( MPI_COMM_WORLD , &N_PROC );
00077     MPI_Comm_rank( MPI_COMM_WORLD , &I_PROC );
00078 #else
00079     MPI::Init(argc,argv);
00080     I_PROC = MPI::COMM_WORLD.Get_rank();
00081     N_PROC = MPI::COMM_WORLD.Get_size();
00082 #endif
00083 
00084   //  int itag = 0;
00085   printf(
00086          ">>> The program is being run with MPI mode.( IP= %3d / NP= %3d ) \n\n", 
00087          I_PROC, N_PROC
00088          );
00089 #endif
00090 
00091 #ifdef DEB
00092   printf("\n\n>>> The program is being run in DEBUG mode.\n\n\n");
00093 #endif
00094 
00095 #ifdef GRAPHIC
00096   printf(">>> EGGX library is used.\n");
00097 #endif
00098 
00099 #ifdef CHAINJOB
00100   printf(">>CHAIN JOB\n");
00101   gettimeofday(&s_time, NULL);
00102 #endif
00103 
00104   Parameter P(argc,argv);
00105 
00106   TheSegmentPool.init(P.NSEGMAX);
00107   TheVertexPool.init(P.NVERMAX);
00108 
00109   RND.setSeed(P.SEED,32);
00110 
00111 #ifdef PAPI
00112   int ierr;
00113   long_long nflops;
00114   float real_time,proc_time,mflops;
00115 #endif
00116   
00117 #ifdef PAPI
00118   ierr = PAPI_flips(&real_time,&proc_time,&nflops,&mflops);
00119 #endif
00120 
00121   Simulation Sim(P);
00122 
00123 #ifdef PAPI
00124   ierr = PAPI_flips(&real_time,&proc_time,&nflops,&mflops);
00125 
00126   cout << "real_time = "<< real_time << endl;
00127   cout << "proc_time = "<< proc_time << endl;
00128   cout << "MFLOPS    = "<< mflops << endl;
00129 #endif 
00130 
00131 #ifdef MULTI
00132 #ifdef KASHIWA
00133   MPI_Finalize();
00134 #else
00135   MPI::Finalize();
00136 #endif
00137 #endif
00138 
00139 }
00140 //
00141 //######################################################################
00142  
00143 
00144 inline Simulation::Simulation(Parameter& P0) :
00145   P(P0), ALG(P.ALGFILE), LAT(P.LATFILE, ALG), MSR(P,LAT,ALG){
00146 
00147 #ifdef GRAPHIC
00148   g_init();
00149 #endif
00150   reset_counters();
00151 
00152 #ifdef CHAINJOB
00153   BinaryIO();
00154   if(!isChainjob){
00155     set_NCYC();
00156   }else{
00157     gettimeofday(&e_time, NULL);
00158     cout<<"Reading bin file takes "<<dtime(s_time,e_time)<<" sec."<<endl;
00159   }
00160   isEnd = false;
00161 #else
00162   if (P.RUNTYPE==0)set_NCYC();
00163 #endif
00164 
00165   int index;
00166 
00167   if (P.RUNTYPE >> 0 ){
00168  
00169     // =========================  edit  sakakura  =========================
00170     // =========================  edit  sakakura  =========================
00171 
00172     int nkink,nkink0,nkink1;
00173     double ediag,ediag0,ediag1;
00174 
00175     double dee;
00176 //    bool judge;
00177     int judge;
00178     int index0,index1,indexx;
00179     double hh[P.NREP];
00180 
00181     index  = I_PROC;
00182 
00184     Algorithm* alg = new Algorithm[P.NREP];
00185 
00186     for (int i=0; i<P.NREP; i++ ){
00187       alg[i].set_i(i);
00188       alg[i].set_alg(i);
00189     }
00190 
00191     if (P.RUNTYPE == 1){
00192       for (int i=0; i<P.NREP; i++) {
00193         hh[i]=P.VF + P.DF*i;
00194       }
00196     }
00197 
00198     if (P.RUNTYPE == 2){
00200       for (int i=0; i<P.NREP; i++) {
00201         hh[i]=P.VB + P.DB*i;
00202       }
00204     }  
00205 
00206     //set_NCYC();
00207     int is,ie,ikey,i;
00208 
00209     ikey = 1;
00210     is = 1;
00211     ie = P.NREP-2;
00212     //  judge = false;
00213     judge = 1;//true;
00214     index = I_PROC;
00215 
00216     MPI_Status recv_status;
00217 
00218     for (ISET=ISETstart; ISET<P.NSET; ISET++) {
00219       for (i=is; i<ie; i+=2) {
00221         if (judge == 1 && ISET !=ISETstart){
00222 
00223 #ifdef KASHIWA
00224           if (I_PROC == i+1) MPI_Send(&index , 1, MPI_INT, i  , 0, MPI_COMM_WORLD);
00225           if (I_PROC == i  ) MPI_Recv(&index0 , 1, MPI_INT, i+1  , 0, MPI_COMM_WORLD, &recv_status);
00226           if (I_PROC == i  ) MPI_Send(&index , 1, MPI_INT, i+1, 0, MPI_COMM_WORLD);
00227           if (I_PROC == i+1) MPI_Recv(&index1, 1, MPI_INT, i  , 0, MPI_COMM_WORLD, &recv_status);
00228 #else
00229           if (I_PROC == i+1) MPI::COMM_WORLD.Send(&index , 1, MPI::INT, i  , 0);
00230           if (I_PROC == i  ) MPI::COMM_WORLD.Recv(&index0, 1, MPI::INT, i+1, 0);
00231           if (I_PROC == i  ) MPI::COMM_WORLD.Send(&index , 1, MPI::INT, i+1, 0);
00232           if (I_PROC == i+1) MPI::COMM_WORLD.Recv(&index1, 1, MPI::INT, i  , 0);
00233 #endif
00234 
00235           if (I_PROC == i   ) index=index0;
00236           if (I_PROC == i+1 ) index=index1;
00237         }
00238       }
00239 
00240     //if (judge == true) {
00241       if (judge == 1) {
00242         if (P.RUNTYPE==1) RepChangeF(alg[index]);
00243         if (P.RUNTYPE==2) RepChangeB(hh[index]);
00244         set_NCYC();
00245       }
00246 
00247       if (ikey == 0) {
00248         ikey = 1;
00249         is = 1; ie=P.NREP-2;
00250       }
00251       else{
00252         ikey = 0;
00253         is = 0; ie=P.NREP-1;
00254       }
00255 
00256       Set(index);
00257       MSR.summary();
00258       //  ene = MSR.get_data();
00259 
00260       ediag = MSR.get_ediag();
00261 
00262       if (P.RUNTYPE==1){
00263         ediag = MSR.get_ediag();
00264       }   
00265 
00266       if (P.RUNTYPE==2){
00267         ediag = MSR.get_ediag();
00268         nkink = MSR.get_nkink();
00269       }   
00270 
00271 //    judge = false;
00272       judge = 0;
00273 
00274       if (ISET == P.NSET-1) break;
00275 
00276       for (i=is; i<ie; i+=2) {
00277 
00278         if (P.RUNTYPE==2){
00279 #ifdef KASHIWA
00280           if (I_PROC == i   ) MPI_Send(&nkink , 1, MPI_INT, i+1,0, MPI_COMM_WORLD);
00281           if (I_PROC == i+1 ) MPI_Recv(&nkink0, 1, MPI_INT, i  ,0, MPI_COMM_WORLD, &recv_status);
00282 #else
00283           if (I_PROC == i   ) MPI::COMM_WORLD.Send(&nkink , 1, MPI::INT, i+1,0);
00284           if (I_PROC == i+1 ) MPI::COMM_WORLD.Recv(&nkink0, 1, MPI::INT, i  ,0);
00285 #endif
00286         }
00287 
00288 #ifdef KASHIWA
00289         if (I_PROC == i  ) MPI_Send(&ediag , 1, MPI_DOUBLE, i+1,0, MPI_COMM_WORLD);
00290         if (I_PROC == i+1) MPI_Recv(&ediag0, 1, MPI_DOUBLE, i  ,0, MPI_COMM_WORLD, &recv_status);
00291         if (I_PROC == i  ) MPI_Send(&index , 1, MPI_INT, i+1, 0, MPI_COMM_WORLD);
00292         if (I_PROC == i+1) MPI_Recv(&index0, 1, MPI_INT, i  , 0, MPI_COMM_WORLD, &recv_status);
00293 #else
00294         if (I_PROC == i  ) MPI::COMM_WORLD.Send(&ediag , 1, MPI::DOUBLE, i+1,0);
00295         if (I_PROC == i+1) MPI::COMM_WORLD.Recv(&ediag0, 1, MPI::DOUBLE, i  ,0);
00296         if (I_PROC == i  ) MPI::COMM_WORLD.Send(&index , 1, MPI::INT, i+1, 0);
00297         if (I_PROC == i+1) MPI::COMM_WORLD.Recv(&index0, 1, MPI::INT, i  , 0);
00298 #endif
00299 
00300         if (I_PROC == i+1 ) {
00301 
00302           if (P.RUNTYPE == 1)dee = exp(  LAT.BETA * (ediag - ediag0) \
00303                                          * ( hh[index] - hh[index0]) ) ;
00304 
00305           if (P.RUNTYPE == 2)dee = exp( (hh[index]-hh[index0])*(ediag-ediag0) )\
00306                                * ( pow( hh[index0]/hh[index] ,nkink-nkink0 ) );
00307  
00309           if (dee > RND.Uniform()) judge = 1;
00310           // for debug//       judge=false;
00311 
00313           //  printf("rank=%d index(i,i+1)=%d %d judge=%d\n",\
00314           //  I_PROC,index0,index,judge);
00316 
00317         }
00318 
00319 #ifdef KASHIWA
00320         if (I_PROC == i+1 ) MPI_Send(&judge , 1, MPI_INT, i,   0, MPI_COMM_WORLD);
00321         if (I_PROC == i   ) MPI_Recv(&judge,  1, MPI_INT, i+1, 0, MPI_COMM_WORLD, &recv_status);
00322 #else 
00323         if (I_PROC == i+1 ) MPI::COMM_WORLD.Send(&judge , 1, MPI::INT, i,   0);
00324         if (I_PROC == i   ) MPI::COMM_WORLD.Recv(&judge,  1, MPI::INT, i+1, 0);
00325 //      if (I_PROC == i+1 ) MPI::COMM_WORLD.Send(&judge , 1, MPI::BOOL, i,   0);
00326 //      if (I_PROC == i   ) MPI::COMM_WORLD.Recv(&judge,  1, MPI::BOOL, i+1, 0);
00327 #endif
00328       }
00329 
00330     }
00331 
00332     MSR.summary_REP(index);
00333 
00334     ISETstart = 0;
00335 
00336     // =========================  edit  sakakura  =========================
00337  
00338   }else{
00340     int idum=0;
00341     for (ISET=ISETstart; ISET<P.NSET; ISET++) Set(idum);
00342     ISETstart=0;
00343     MSR.summary();
00344   }
00346 
00347   P.openfile();
00348   fprintf(P.FOUT,"C This is DSQSS ver.%.1f\n\n",Version);
00349   LAT.show_param(P.FOUT);
00350   P.dump(P.FOUT);
00351   MSR.show(P.FOUT);
00352 
00353   int ns_used = TheSegmentPool.number_of_used_elements();
00354   int nv_used = TheVertexPool.number_of_used_elements();
00355   int nrvi_used = TheRVIPool.number_of_used_elements();
00356   fprintf(P.FOUT,"I [the maximum number of segments]          = %d\n", ns_used);
00357   fprintf(P.FOUT,"I [the maximum number of vertices]          = %d\n", nv_used);
00358   fprintf(P.FOUT,"I [the maximum number of reg. vertex info.] = %d\n", nrvi_used);
00359   // === edit sakakura ===
00360   if (P.RUNTYPE > 0) {
00361     fprintf(P.FOUT,"I [the index of replica exchange ]         = %d\n", index);}
00362 
00363   P.closefile();
00364 #ifdef GRAPHIC
00365   g_clear();
00366 #endif
00367 
00368 #ifdef CHAINJOB
00369   isEnd = true;
00370   cjobout = new ofstream();
00371   (*cjobout).open(CJOBFILE, ios::out | ios::binary);//reset
00372   end_cjob();
00373 #endif
00374 
00375 }
00376 
00377 //======================================================================
00378 
00379 void Simulation::reset_counters() {
00380   ISET = -1;
00381   IMCSE = -1;
00382   IMCSD = -1;
00383   IMCS = -1;
00384   ICYC = -1;
00385 
00386   ISETstart  = 0;
00387   IMCSDstart = 0;
00388   IMCSstart  = 0;
00389 
00390 }
00391 
00392 //======================================================================
00393 
00394 inline Simulation::~Simulation() {
00395   //  printf("*** Destroying Simulation\n");
00396 }
00397 
00398 //======================================================================
00399 
00400 void Simulation::set_NCYC() {
00401   // edit Suzuki
00402   // 逆温度に対する規格化
00403 #ifdef NORM
00404   double vol = (double)LAT.NSITE;
00405 #else
00406   double vol = LAT.BETA * (double)LAT.NSITE;
00407 #endif
00408   // end edit Suzuki
00409 
00410   double path;
00411   int ncyc=1;
00412 
00413   int NSAMP = P.NMCSE/10;
00414   int*  ncycSAMP = new int[NSAMP]; //edit sakakura
00415   //  int ncycSAMP[NSAMP];
00416   
00417   for (IMCSE=0; IMCSE<P.NMCSE; IMCSE++) {
00418 
00419     path = 0.0;
00420     ncyc = 0;
00421     // 次の行の LAT.NSITE == 2 && ncyc%2 == 0 の指定は
00422     // S=1/2 反強磁性ハイゼンベルクモデルの特殊な非エルゴード性
00423     // を排除するため
00424     while ( path < vol || ( LAT.NSITE==2 && ncyc%2==0 ) ) {
00425       double p = Cycle();
00426       path += p;
00427       ncyc++;
00428     }
00429 
00430     ncycSAMP[IMCSE%NSAMP]= ncyc;
00431   };
00432 
00433   int ncycSUM = 0;
00434   for(int isamp = 0; isamp < NSAMP; isamp++){
00435     ncycSUM += ncycSAMP[isamp];
00436   }
00437 
00438   P.NCYC = ncycSUM / NSAMP ;
00439 #ifdef GRAPHIC
00440   g_draw();
00441 #endif
00442 }
00443 
00444 //======================================================================
00445 
00446 // === edit rep sakakrua ===
00447 //void Simulation::Set() {
00448 void Simulation::Set(int ix) {
00449   // === edit rep sakakrua ===
00450 
00451   ICYC=-1;
00452 
00453   for (IMCSD=IMCSDstart; IMCSD<P.NMCSD; IMCSD++) {
00454 #ifdef CHAINJOB
00455     gettimeofday(&e_time, NULL);
00456     if(dtime(s_time,e_time) > killtime ){IMCS=0;write_cjob();end_job();}
00457     isChainjob=false;
00458 #endif
00459     Sweep();
00460   }
00461   IMCSDstart = 0;
00462 #ifdef CHAINJOB
00463   if(!isChainjob){MSR.setinit();}else{isChainjob=false;}
00464 #else
00465   MSR.setinit();//kota for (int i=0; i<NACC; i++) ACC[i].reset();
00466 #endif
00467 
00468   for (IMCS=IMCSstart; IMCS<P.NMCS;  IMCS++) {
00469 #ifdef CHAINJOB
00470     gettimeofday(&e_time, NULL);
00471     if(dtime(s_time,e_time) > killtime ){write_cjob();end_job();}
00472 #endif
00473     Sweep();
00474     MSR.measure();
00475     MSR.accumulate_length(AMP);
00476   }
00477   IMCSstart = 0;
00478 
00479   // === edit rep sakakrua ===
00480   //MSR.setsummary();
00481   MSR.setsummary(ix);
00482   // === edit rep sakakrua ===
00483 
00484 }
00485 
00486 //======================================================================
00487 
00488 void Simulation::Sweep() {
00489 
00490   double len = 0.0;
00491   for (ICYC=0; ICYC<P.NCYC; ICYC++) {
00492     len += Cycle();
00493   }
00494   AMP = len/(double)P.NCYC;
00495 #ifdef GRAPHIC
00496   g_draw();
00497 #endif
00498 
00499 }
00500 
00501 //======================================================================
00502 
00503 double Simulation::Cycle() {
00504   if ( ! PlaceWorm() ) return 0.0;
00505   return MoveHead(); 
00506 
00507 }
00508 
00509 //======================================================================
00510 
00511 bool Simulation::PlaceWorm() {
00512 
00513 #ifdef MONITOR
00514   printf("\n#####################################################\n");
00515   printf("Simulation::PlaceWorm> Start.\n");
00516 #endif
00517 
00518   // ワーム初期バーテックスの設定
00519 
00520   Vertex& Vorg = W.origin();
00521   W.setCurrentVertex(Vorg);
00522     
00523   // ワーム初期位置の決定
00524 
00525   int s = RND.Int(LAT.NSITE);
00526   //  edit Suzuki
00527   // 逆温度に対する補正(beta -> 1で規格化)
00528 #ifdef NORM
00529   double t =  RND.Uniform();
00530 #else
00531   double t = (double)LAT.BETA * RND.Uniform();
00532 #endif
00533  
00534   // end edit Suzuki
00535   //kota//  乱数を用いて初期位置(サイト番号、位置)を決定
00536 
00537 #ifdef MONITOR
00538   printf("  The worm is placed at (s= %d, t= %8.3f)\n", s+1, t);
00539 #endif
00540 
00541   // ワームテールの初期化
00542   Site& ST = LAT.S(s);
00543 
00544   SiteProperty& SP = ST.Property();
00545 
00546   VertexProperty& VP = SP.getVertexProperty();
00547   //サイト、サイトプロパティ、バーテックスプロパティの各オブジェクトの設定
00548   Vorg.BareVertex::init( t , VP );
00549 
00550 
00551   //バーテックス初期化
00552   Vorg.setTime(t);
00553   //tをセット
00554 
00555   // セグメントを探す
00556   Segment& S0 = ST.findS(t);
00557   // 初期方向
00558   int x = S0.X();
00559 
00560   SiteInitialConfiguration& IC = SP.getInitialConfiguration(x);
00561 
00562   int NCH = IC.NCH; //Mに依存 when M=1,NCH=2
00563   int c;
00564   double r = RND.Uniform();
00565   for (c=0; c<NCH; c++) {
00566     if ( r < IC.CH[c].PROB ) break;
00567   } //初期状態の選択 alg.xmlの<Site>から
00568 
00569   ScatteringChannel& CH = IC.CH[c];
00570 
00571   int out  = CH.OUT;
00572   int xout = CH.XOUT;
00573 
00574   if ( out == DIR::UNDEF ) return false; // ワームは生成されない
00575   // セグメントをワームテールで切断 "UNDEF=-1"
00576   Segment& S1 = S0.cut( Vorg , 0 );
00577 
00578   // 方向に応じて初期セグメントを選ぶ
00579 
00580   if ( out == UORD::DOWN ) { 
00581     W.setCurrentSegment(S0);
00582   } else {       // moving up initially
00583     W.setCurrentSegment(S1);
00584   }
00585 
00586   // 初期タイプの設定
00587 
00588   W.setXBEHIND(xout);
00589 
00590 #ifdef MONITOR
00591   printf("\nSimulation::PlaceWorm> End.\n");
00592 #endif
00593   return true; // ワームが生成された
00594 }
00595 
00596 //======================================================================
00597 
00598 double Simulation::MoveHead() {
00599 
00600 #ifdef MONITOR
00601   printf("\n#####################################################\n");
00602   printf("Simulation::MoveHead> Start.\n");
00603   dump();
00604   printf("--------------------------------\n");
00605 #endif
00606 
00607   double len=0.0;
00608   double len0=0.0;
00609   double len1=0.0;
00610   EndOfCycle = false;
00611 
00612   while(true) {
00613     //     W.getV();
00614     //    cout <<"getuord() = " << W.getUORD() << endl;
00615     if(W.getUORD()){
00616 
00617       len0=DOWN_ONESTEP();
00618       len += len0;
00619       len1-= len0;
00620       //             printf("-- DOWN -- len =%.4f len1=%.4f\n",len,len1);
00621       //     cout << "-- DOWN -- "<<len <<" "<<len1<<endl;
00622       //     len += DOWN_ONESTEP();
00623     }else{
00624       len0=UP_ONESTEP();
00625       len += len0;
00626       len1+= len0;
00627       //            printf("--  UP  -- len =%.4f len1=%.4f\n",len,len1);
00628       //     cout << "--  UP  -- "<<len <<" "<<len1<<endl;
00629       //     len += UP_ONESTEP();
00630     }
00631     if ( EndOfCycle ) break;
00632   }
00633 #ifdef DEB
00634   if ( ! W.atOrigin() ) {  
00635     fprintf(FERR,"\nMoveHead> ### Error.\n");
00636     fprintf(FERR," ... Hasn't come back to the origin.\n");  
00637     dump();
00638     exit(0);
00639   }
00640 #endif
00641 
00642   W.remove();
00643 
00644 #ifdef MONITOR
00645   printf("\nSimulation::MoveHead> End.\n");
00646 #endif
00647   return len;
00648 
00649 }
00650 
00651 //======================================================================
00652 double Simulation::UP_ONESTEP(){
00653   double len = 0.0;
00654   double RHO = 0.0;
00655 
00656   Segment& c_S   = W.getCurrentSegment();
00657   Site&   c_Site = (* c_S.getONSITE() );
00658  
00659   Interaction** CI = c_Site.getCI();
00660   int NCI = c_Site.getNCI(); // the number of Interaction:
00661 
00662   Vertex& c_V = W.getCurrentVertex();
00663   Vertex& n_V = W.getNextVertex();
00664  
00665   double c_time;
00666   double n_time;
00667   int xinc = W.getXBEHIND();
00668   UniformInterval* UI = new UniformInterval[NCI];
00669   
00670   Ring<RegVInfo> RingRVI ;
00671   RingRVI.ROOT.V_x=&V4REF;
00672 
00673   // edit Suzuki
00674   // 逆温度に対する規格化
00675 #ifdef NORM
00676   double n_Vtime =(n_V.isTerminal())? 1.0 : n_V.time();
00677 #else
00678   double n_Vtime =(n_V.isTerminal())? LAT.BETA : n_V.time();
00679 #endif
00680   // end edit Suzuki
00681   if(c_V.isTerminal()){
00682     c_time  = 0.0;
00683     n_time  = n_Vtime ;
00684 
00685     for(int iCI = 0; iCI < NCI; iCI++){ 
00686       UI[iCI].init(CI[iCI],xinc);
00687 
00688       for(int i_body = 0; i_body <UI[iCI].nbody ; i_body++){
00689         if( &( (*CI[iCI]).site( i_body ) ) == &c_Site ) {
00690           UI[iCI].inc = 2 * i_body ;
00691           UI[iCI].n_S[i_body] = &c_S ;
00692         }else{
00693           UI[iCI].n_S[i_body] = &( (*CI[iCI]).site( i_body ).head() );
00694           Vertex& near_V = (*UI[iCI].n_S[i_body]).top();
00695           double V_time = near_V.time();
00696           if( (! near_V.isTerminal()) && ( V_time < n_time )&&( &near_V != &n_V) ){
00697             RegVInfo& RVI = TheRVIPool.pop();
00698             RVI.setRVI( &near_V, iCI, i_body, V_time);
00699             RingRVI.add_tail(RVI);
00700           }
00701         }
00702       }
00703       UI[iCI].setx();
00704       UI[iCI].setVIC();
00705 
00706       //  edit Suzuki
00707       // 逆温度に対する規格化
00708 #ifdef NORM
00709       if(UI[iCI].DefinedVIC){RHO += (*(UI[iCI].VIC)).dRHO * LAT.BETA ; }
00710 #else
00711       if(UI[iCI].DefinedVIC){RHO += (*(UI[iCI].VIC)).dRHO; }
00712 #endif
00713       //  end edit Suzuki
00714     }
00715   }else{
00716     c_time  = c_V.time();
00717     n_time  = n_Vtime-c_time;
00718 
00719     for(int iCI = 0; iCI < NCI; iCI++){
00720       UI[iCI].init(CI[iCI],xinc);
00721       for(int i_body = 0; i_body <UI[iCI].nbody ; i_body++){
00722         if( &( (*CI[iCI]).site( i_body ) ) == &c_Site ) {
00723           UI[iCI].inc = 2 * i_body ;
00724           UI[iCI].n_S[i_body] = &c_S ;
00725         }else{
00726           if( CI[iCI] ==  c_V.getONINTERACTION()  ){
00727             UI[iCI].n_S[i_body] = &( c_V.S( 2*i_body + 1 ) );
00728 
00729           }else{
00730             UI[iCI].n_S[i_body] = &( (*CI[iCI]).site( i_body ).findS( c_time ) );
00731           }
00732 
00733           Vertex& near_V = (*UI[iCI].n_S[i_body]).top();
00734           double V_time = near_V.time() - c_time;
00735           if( (! near_V.isTerminal()) && ( V_time < n_time )&&( &near_V != &n_V) ){
00736             RegVInfo& RVI = TheRVIPool.pop();
00737             RVI.setRVI( &near_V, iCI, i_body, V_time);
00738             RingRVI.add_tail(RVI);
00739           }
00740         }
00741       }
00742       UI[iCI].setx();
00743       UI[iCI].setVIC(); //どのVICかセットする
00744 
00745       //  edit Suzuki
00746       // 逆温度に対する規格化
00747 #ifdef NORM
00748       if(UI[iCI].DefinedVIC){RHO += (*(UI[iCI].VIC)).dRHO * LAT.BETA ; }
00749 #else
00750       if(UI[iCI].DefinedVIC){RHO += (*(UI[iCI].VIC)).dRHO; }
00751 #endif
00752       //  end edit Suzuki
00753     } 
00754   }
00755   double near_tau = 0.0;
00756   double try_tau  = 0.0;
00757   int out, xout;
00758 
00759   while(true){
00760 
00761     Ring<RegVInfo>::iterator it = RingRVI.sort_min();
00762   
00763     near_tau =(RingRVI.empty())? n_time - len : (*it).V_time - len ;
00764 
00765     if( RHO > EPS ){
00766       try_tau = RND.Exp()/ RHO ;
00767       if( try_tau < near_tau ){
00768         len += try_tau;
00769         double RHORND = RHO * RND.Uniform();
00770                 
00771         /*
00772           double dRHO;
00773           while(true){
00774           dRHO = (UI[i_UI].DefinedVIC)?(*UI[i_UI].VIC).dRHO : 0.0 ;
00775           if( RHORND > dRHO ){
00776           RHORND -= dRHO;
00777           ++i_UI;
00778           }else{
00779           break;
00780           }
00781           }
00782         */
00783         int i_UI=0;
00784         int s_UI=0;
00785         double iRHO = 0.0;
00786         double sRHO = 0.0;
00787         do{
00788           sRHO = iRHO;
00789           if(UI[i_UI].DefinedVIC ){
00790             //  edit Suzuki
00791             // 逆温度に対する規格化
00792 #ifdef NORM
00793             iRHO += (*UI[i_UI].VIC).dRHO * LAT.BETA ;
00794 #else
00795             iRHO += (*UI[i_UI].VIC).dRHO;
00796 #endif
00797             // end edit Suzuki
00798             s_UI = i_UI;
00799           }
00800           ++i_UI;
00801         }while( RHORND >= iRHO );
00802         RHORND -= sRHO;
00803         //  edit Suzuki
00804         // 逆温度に対する規格化
00805 #ifdef NORM
00806         ScatteringChannel& CH = (*UI[s_UI].VIC).getScatteringChannel(RHORND / LAT.BETA);
00807 
00808 #else
00809         ScatteringChannel& CH = (*UI[s_UI].VIC).getScatteringChannel(RHORND);
00810 #endif
00811         // これは必要?
00812         // end edit Suzuki
00813         
00814         double setVtime =  c_time + len ;
00815         Vertex& setV = TheVertexPool.pop();
00816         setV.BareVertex::init(UI[s_UI].I_n, setVtime, (*UI[s_UI].VP));
00817         for(int ibody = 0; ibody<UI[s_UI].nbody; ++ibody){ (* UI[s_UI].n_S[ibody]).cut(setV, ibody);}
00818         (* UI[s_UI].I_n).add_tail(setV);
00819 
00820         setV.S(UI[s_UI].inc).setX(xinc);
00821         if(! ( c_V.isTerminal() || c_V.isKink() ) ){ c_V.erase(); }
00822         W.setCurrentVertex( setV );
00823         W.setCurrentSegment( setV.S( CH.OUT ) );
00824         W.setXBEHIND( CH.XOUT );
00825         
00826         break;
00827       }
00828     }
00829     if(RingRVI.empty()){
00830       W.setCurrentVertex(n_V);
00831       int inc  = n_V.which(c_S);
00832       if ( n_V.isTerminal() ) {
00833         out = 1-inc;
00834         xout = xinc;
00835       } else{
00836         VertexInitialConfiguration& IC  = n_V.getInitialConfiguration( inc , xinc );
00837         ScatteringChannel& CH = IC.getScatteringChannel();
00838         out = CH.OUT;
00839         xout = CH.XOUT;
00840       }
00841 
00842       c_S.setX( xinc );
00843       if(! ( c_V.isTerminal() || c_V.isKink() ) ){ c_V.erase(); }
00844       if ( out == DIR::UNDEF ) {
00845         EndOfCycle = true;
00846       } else {
00847         W.setCurrentSegment( n_V.S( out ) );
00848         W.setXBEHIND( xout );
00849       }
00850       len += near_tau;
00851       break;
00852     }
00853     
00854     //UPDATE RegVI and UI
00855     Vertex* changeV = (*it).V_x;
00856     do{
00857       RegVInfo& RVImin = (*it);
00858       // +++ edit sakakura +++ //
00859 #ifdef NORM
00860       if(UI[RVImin.i_UI].DefinedVIC){RHO -= (*(UI[RVImin.i_UI].VIC)).dRHO * LAT.BETA; }
00861 #else
00862       if(UI[RVImin.i_UI].DefinedVIC){RHO -= (*(UI[RVImin.i_UI].VIC)).dRHO; }
00863 #endif
00864       // +++ edit sakakura +++ //
00865       //
00866       //
00867       // +++ edit sakakura +++ //
00868       //    if(UI[RVImin.i_UI].DefinedVIC){RHO -= (*(UI[RVImin.i_UI].VIC)).dRHO; }
00869       UI[RVImin.i_UI].n_S[RVImin.i_body] = &( (* UI[RVImin.i_UI].n_S[RVImin.i_body] ).next());//
00870       UI[RVImin.i_UI].setx(RVImin.i_body);
00871       UI[RVImin.i_UI].setVIC();
00872       // +++ edit sakakura +++ //
00873 #ifdef NORM
00874       if(UI[RVImin.i_UI].DefinedVIC){RHO += (*(UI[RVImin.i_UI].VIC)).dRHO * LAT.BETA; }
00875 #else 
00876       if(UI[RVImin.i_UI].DefinedVIC){RHO += (*(UI[RVImin.i_UI].VIC)).dRHO; }
00877 #endif 
00878       // +++ edit sakakura +++ //
00879       //kota RHOの処理で逆温度規格化しなくてよい?
00880 
00881       //add RegVI
00882       Vertex& near_V = (*UI[RVImin.i_UI].n_S[RVImin.i_body]).top();//
00883       double V_time = near_V.time() - c_time;//
00884       if( (! near_V.isTerminal()) && ( V_time < n_time )&&( &near_V != &n_V) ){
00885         RegVInfo& RVI = TheRVIPool.pop();
00886         RVI.setRVI( &near_V, RVImin.i_UI, RVImin.i_body, V_time);
00887         RingRVI.add_tail(RVI);
00888       }
00889       ++it;
00890       RVImin.erase();
00891     }while(changeV == (*it).V_x);
00892     
00893     it = RingRVI.head();
00894     len += near_tau;
00895     
00896   }
00897 
00898   while(!RingRVI.empty()){
00899     RegVInfo& DRVI = RingRVI.first();
00900     DRVI.remove();
00901     TheRVIPool.push(DRVI);
00902   };
00903   
00904   delete [] UI;
00905   return len;
00906 }
00907 //======================================================================
00908 
00909 //======================================================================
00910 double Simulation::DOWN_ONESTEP(){
00911   double len = 0.0;
00912   double RHO = 0.0;
00913 
00914   Segment& c_S   = W.getCurrentSegment();
00915   Site&   c_Site = (* c_S.getONSITE() );
00916   
00917   Interaction** CI = c_Site.getCI(); //  Interaction* CI [n_int]; which belong to curSite
00918   int NCI = c_Site.getNCI(); // the number of Interaction:
00919 
00920   Vertex& c_V = W.getCurrentVertex();
00921   Vertex& n_V = W.getNextVertex();
00922   
00923   double c_time ;
00924   double n_time ;
00925   int xinc = W.getXBEHIND();
00926   UniformInterval* UI = new UniformInterval[NCI];
00927 
00928   Ring<RegVInfo> RingRVI ;
00929   RingRVI.ROOT.V_x=&V4REF;
00930   double n_Vtime =(n_V.isTerminal())? 0.0 : n_V.time();
00931 
00932   if(c_V.isTerminal()){
00933     //  +++ add sakakura +++
00934 #ifdef NORM
00935     c_time = 1.0; //?
00936     n_time = 1.0 - n_Vtime ;
00937 #else
00938     c_time  = LAT.BETA;
00939     n_time  = LAT.BETA - n_Vtime ;
00940 #endif
00941     //  +++ add sakakura +++
00942     //    n_time  = LAT.BETA - n_Vtime ;
00943     for(int iCI = 0; iCI < NCI; iCI++){ 
00944       UI[iCI].init(CI[iCI],xinc);
00945       for(int i_body = 0; i_body <UI[iCI].nbody ; i_body++){
00946         if( &( (*CI[iCI]).site( i_body ) ) == &c_Site ) {
00947           UI[iCI].inc = 2 * i_body + 1;
00948           UI[iCI].n_S[i_body] = &c_S ;
00949         }else{
00950           UI[iCI].n_S[i_body] = &( (*CI[iCI]).site( i_body ).tail() );//
00951           Vertex& near_V = (*UI[iCI].n_S[i_body]).bottom();
00952           double V_time = c_time - near_V.time();
00953           if( (! near_V.isTerminal()) && ( V_time < n_time )&&( &near_V != &n_V) ){
00954             RegVInfo& RVI = TheRVIPool.pop();
00955             RVI.setRVI( &near_V, iCI, i_body, V_time);
00956             RingRVI.add_tail(RVI);
00957           }
00958         }
00959       }
00960       UI[iCI].setx();
00961       UI[iCI].setVIC();
00962       //  edit Suzuki
00963       // 逆温度に対する規格化
00964 #ifdef NORM
00965       if(UI[iCI].DefinedVIC){RHO += (*(UI[iCI].VIC)).dRHO * LAT.BETA ; }
00966 #else
00967       if(UI[iCI].DefinedVIC){RHO += (*(UI[iCI].VIC)).dRHO; }
00968 #endif
00969       //  end edit Suzuki
00970     }
00971     
00972   }else{
00973     c_time  = c_V.time();
00974     n_time  = c_time - n_Vtime;
00975     for(int iCI = 0; iCI < NCI; iCI++){
00976       UI[iCI].init(CI[iCI],xinc);
00977       for(int i_body = 0; i_body <UI[iCI].nbody ; i_body++){
00978         if( &( (*CI[iCI]).site( i_body ) ) == &c_Site ) {
00979           UI[iCI].inc = 2 * i_body + 1 ;//
00980           UI[iCI].n_S[i_body] = &c_S ;
00981         }else{
00982           if( CI[iCI] ==  c_V.getONINTERACTION()  ){
00983             UI[iCI].n_S[i_body] = &( c_V.S( 2*i_body ) );//
00984           }else{
00985             UI[iCI].n_S[i_body] = &( (*CI[iCI]).site( i_body ).findS( c_time ) );
00986           }
00987           Vertex& near_V = (*UI[iCI].n_S[i_body]).bottom();//
00988           double V_time = c_time - near_V.time();//
00989           if( (! near_V.isTerminal()) && ( V_time < n_time )&&( &near_V != &n_V) ){
00990             RegVInfo& RVI = TheRVIPool.pop();
00991             RVI.setRVI( &near_V, iCI, i_body, V_time);
00992             RingRVI.add_tail(RVI);
00993           }
00994         }
00995       }
00996       UI[iCI].setx();
00997       UI[iCI].setVIC();
00998       //  edit Suzuki
00999       // 逆温度に対する規格化
01000 #ifdef NORM
01001       if(UI[iCI].DefinedVIC){RHO += (*(UI[iCI].VIC)).dRHO * LAT.BETA ; }
01002 #else
01003       if(UI[iCI].DefinedVIC){RHO += (*(UI[iCI].VIC)).dRHO; }
01004 #endif 
01005       //  end edit Suzuki
01006 
01007     } 
01008   }
01009   double near_tau = 0.0;
01010   double try_tau  = 0.0;
01011   int out, xout;
01012 
01013   while(true){
01014 
01015     Ring<RegVInfo>::iterator it = RingRVI.sort_min();
01016 
01017     near_tau =(RingRVI.empty())? n_time - len : (*it).V_time - len ;
01018 
01019     if( RHO > EPS ){
01020       try_tau = RND.Exp()/ RHO ;
01021       if( try_tau < near_tau ){
01022         len += try_tau;
01023         double RHORND = RHO * RND.Uniform();
01024         /*
01025           int i_UI=0;
01026           double dRHO;
01027           while(true){
01028           dRHO = (UI[i_UI].DefinedVIC)?(*UI[i_UI].VIC).dRHO : 0.0 ;
01029           if( RHORND > dRHO ){
01030           RHORND -= dRHO;
01031           ++i_UI;
01032           }else{
01033           break;
01034           }
01035           }
01036         */
01037 
01038         int i_UI=0;
01039         int s_UI=0;
01040         double iRHO = 0.0;
01041         double sRHO = 0.0;
01042         do{
01043           sRHO = iRHO;
01044           if(UI[i_UI].DefinedVIC ){
01045             //  edit Suzuki
01046             // 逆温度に対する規格化
01047 #ifdef NORM
01048             iRHO += (*UI[i_UI].VIC).dRHO * LAT.BETA ;
01049 #else
01050             iRHO += (*UI[i_UI].VIC).dRHO;
01051 #endif
01052             // end edit Suzuki
01053             s_UI = i_UI;
01054           }
01055           ++i_UI;
01056         }while( RHORND >= iRHO );
01057         RHORND -= sRHO;
01058         //  edit Suzuki
01059         // 逆温度に対する規格化
01060 #ifdef NORM
01061         ScatteringChannel& CH = (*UI[s_UI].VIC).getScatteringChannel(RHORND / LAT.BETA);
01062 #else
01063         ScatteringChannel& CH = (*UI[s_UI].VIC).getScatteringChannel(RHORND);
01064 #endif
01065         //
01066         // end edit Suzuki
01067 
01068         double setVtime =  c_time - len ;
01069         Vertex& setV = TheVertexPool.pop();
01070         setV.BareVertex::init(UI[s_UI].I_n, setVtime, (*UI[s_UI].VP));
01071         for(int ibody = 0; ibody<UI[s_UI].nbody; ++ibody){ (* UI[s_UI].n_S[ibody]).cut(setV, ibody);}
01072         (* UI[s_UI].I_n).add_tail(setV);
01073 
01074         setV.S(UI[s_UI].inc).setX(xinc);
01075         if(! ( c_V.isTerminal() || c_V.isKink() ) ){ c_V.erase(); }
01076         W.setCurrentVertex( setV );
01077         W.setCurrentSegment( setV.S( CH.OUT ) );
01078         W.setXBEHIND( CH.XOUT );
01079         
01080         break;
01081       }
01082     }
01083     if(RingRVI.empty()){
01084       W.setCurrentVertex(n_V);
01085       int inc  = n_V.which(c_S);
01086       if ( n_V.isTerminal() ) {
01087         out = 1-inc;
01088         xout = xinc;
01089       } else {
01090 
01091         VertexInitialConfiguration& IC  = n_V.getInitialConfiguration( inc , xinc );
01092         ScatteringChannel& CH = IC.getScatteringChannel();
01093         out = CH.OUT;
01094         xout = CH.XOUT;
01095       }
01096       
01097       c_S.setX( xinc );
01098       if(! ( c_V.isTerminal() || c_V.isKink() ) ){ c_V.erase(); }
01099       if ( out == DIR::UNDEF ) {
01100         EndOfCycle = true;
01101       } else {
01102         W.setCurrentSegment( n_V.S( out ) );
01103         W.setXBEHIND( xout );
01104       }
01105       len += near_tau;
01106       break;
01107     }
01108     
01109     //UPDATE RegVI and UI
01110     Vertex* changeV = (*it).V_x;
01111     do{
01112       RegVInfo& RVImin = (*it);
01113       // +++ edit sakakura +++ //
01114 #ifdef NORM 
01115       if(UI[(*it).i_UI].DefinedVIC){RHO -= (*(UI[(*it).i_UI].VIC)).dRHO * LAT.BETA; }
01116 #else 
01117       if(UI[(*it).i_UI].DefinedVIC){RHO -= (*(UI[(*it).i_UI].VIC)).dRHO; }
01118 #endif 
01119       // +++ edit sakakura +++ //
01120       UI[(*it).i_UI].n_S[(*it).i_body] = &( (* UI[(*it).i_UI].n_S[(*it).i_body] ).prev());//
01121       UI[(*it).i_UI].setx((*it).i_body);
01122       UI[(*it).i_UI].setVIC();
01123       // +++ edit sakakura +++ //
01124 #ifdef NORM 
01125       if(UI[(*it).i_UI].DefinedVIC){RHO += (*(UI[(*it).i_UI].VIC)).dRHO * LAT.BETA; }
01126 #else 
01127       if(UI[(*it).i_UI].DefinedVIC){RHO += (*(UI[(*it).i_UI].VIC)).dRHO; }
01128 #endif 
01129       // +++ edit sakakura +++ //
01130       
01131       //add RegVI
01132       Vertex& near_V = (*UI[(*it).i_UI].n_S[(*it).i_body]).bottom();//
01133       double V_time = c_time - near_V.time();//
01134       if( (! near_V.isTerminal()) && ( V_time < n_time )&&( &near_V != &n_V) ){
01135         RegVInfo& RVI = TheRVIPool.pop();       
01136         RVI.setRVI( &near_V, (*it).i_UI, (*it).i_body, V_time);
01137         RingRVI.add_tail(RVI);
01138       }
01139       ++it;
01140       RVImin.erase();
01141     }while(changeV == (*it).V_x);
01142         
01143     it = RingRVI.head();
01144     len += near_tau;
01145         
01146   }
01147 
01148   while(!RingRVI.empty()){
01149     RegVInfo& DRVI = RingRVI.first();
01150     DRVI.remove();
01151     TheRVIPool.push(DRVI);
01152   };
01153 
01154   delete [] UI;
01155   return len;
01156 }
01157 //======================================================================
01158 
01159 inline void Simulation::dump(const char* s) {
01160 
01161   printf("\n>>> %s (iset=%d imcse=%d imcsd=%d imcs=%d icyc= %d)",
01162          s, ISET, IMCSE, IMCSD, IMCS, ICYC
01163          );
01164 
01165   LAT.dump();
01166 
01167   W.dump();
01168 
01169 }
01170 
01171 //======================================================================
01172 
01173 inline void Simulation::dump() { dump(""); }
01174 
01175 //======================================================================
01176 
01177 inline void Simulation::Check() {
01178 
01179 #ifndef DEBUG
01180   return;
01181 #endif
01182 
01183   bool ERROR = false ;
01184 
01185   for (int b=0; b<LAT.NINT; b++) {
01186     Interaction& I = LAT.I(b);
01187     if ( I.empty() ) continue;
01188     Interaction::iterator p(I);
01189     while ( ! (++p).atOrigin() ) {
01190       ERROR = ERROR || p->check();
01191     }
01192   }
01193 
01194   for (int s=0; s<LAT.NSITE; s++) {
01195     Site& S = LAT.S(s);
01196     Site::iterator p(S);
01197     while ( ! (++p).atOrigin() ) {
01198       ERROR = ERROR || p->check();
01199     }
01200   }
01201 
01202   if ( ERROR ) {
01203     dump();
01204     exit(0);
01205   }
01206 
01207 }
01208 //======================================================================
01209 // == edit sakakura ==
01210 //void Simulation::RepChange(int ix){
01211 void Simulation::RepChangeF(Algorithm& algx){
01212   //  InteractionProperty& IP4=alg[ix].getInteractionProperty(alg[ix].NSTYPE-1);
01213   InteractionProperty& IP4=algx.getInteractionProperty(algx.NSTYPE-1);
01214   InteractionProperty& IP2 = ALG.getInteractionProperty(ALG.NSTYPE-1);
01215 
01216   for (int i=0; i<ALG.NXMAX; i++){
01217     int* x = new int[2];
01218     for (int j=0; j<ALG.NXMAX; j++){
01219 
01220       x[0]=i;x[1]=j;
01221       double d= IP4.VertexDensity(x);
01222       IP2.VertexDensity(x)=d;
01223       IP2.AverageInterval(x)=1.0/d;
01224 
01225     }
01226   }
01227 
01228   IP2.EBASE=  IP4.EBASE;
01229   MSR.EBASE = IP4.EBASE * LAT.NINT;
01230   // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
01231 
01232   VertexProperty& VP2 = ALG.getVertexProperty(1);
01233   //VertexProperty& VP4 = alg[ix].getVertexProperty(1);
01234   VertexProperty& VP4 = algx.getVertexProperty(1);
01235 
01236   for (int i=0; i<VP4.NIC; i++){
01237 
01238     VP2.IC[i].NLEG =VP4.IC[i].NLEG;
01239     VP2.IC[i].INC  =VP4.IC[i].INC ;
01240     VP2.IC[i].XINC =VP4.IC[i].XINC;
01241     VP2.IC[i].NCH  =VP4.IC[i].NCH ;
01242     VP2.IC[i].dRHO =VP4.IC[i].dRHO ;
01243     VP2.IC[i].State[0] = VP4.IC[i].State[0] ; 
01244     VP2.IC[i].State[1] = VP4.IC[i].State[1] ;
01245     VP2.IC[i].State[2] = VP4.IC[i].State[2] ;
01246     VP2.IC[i].State[3] = VP4.IC[i].State[3] ;
01247 
01248     for (int j = 0; j<VP2.IC[i].NCH; j++){
01249       VP2.IC[i].CH[j].OUT =VP4.IC[i].CH[j].OUT; 
01250       VP2.IC[i].CH[j].XOUT =VP4.IC[i].CH[j].XOUT; 
01251       VP2.IC[i].CH[j].PROB =VP4.IC[i].CH[j].PROB;
01252     }
01253   }
01254 
01255   ALG.initializer(VP4.NIC);
01256 }
01257 // == edit sakakura ==
01258 
01259 void Simulation::RepChangeB(double hx){
01260 
01261   LAT.BETA=hx;
01262   for (int ix=0; ix<pow(LAT.NSITE,LAT.D);ix++){
01263     Site& SX = LAT.S(ix);
01264     SX.setBeta(1.0);
01265     }
01266 }
 全て クラス ネームスペース ファイル 関数 変数 型定義 列挙型の値 フレンド マクロ定義