Trisurf Monte Carlo simulator
Samo Penic
2016-02-15 831fd28f91dd2537c99538d9664e0d62a92a791f
src/timestep.c
@@ -12,31 +12,63 @@
#include "sh.h"
#include "shcomplex.h"
#include "vesicle.h"
#include<gsl/gsl_complex.h>
#include<gsl/gsl_complex_math.h>
#include<string.h>
ts_bool run_simulation(ts_vesicle *vesicle, ts_uint mcsweeps, ts_uint inititer, ts_uint iterations, ts_uint start_iteration){
   ts_uint i, j,k;
   ts_double r0,kc1,kc2,kc3,kc4;
   ts_double l1,l2,l3,volume=0.0,area=0.0,vmsr,bfsr, vmsrt, bfsrt;
   ts_uint i, j,k,l,m;
   ts_double r0,kc1=0,kc2=0,kc3=0,kc4=0;
   ts_double l1,l2,l3,vmsr,bfsr, vmsrt, bfsrt;
   ts_ulong epochtime;
   FILE *fd1;
//    char filename[255];
   FILE *fd=fopen("statistics.csv","w");
   FILE *fd1,*fd2=NULL;
    char filename[10000];
    strcpy(filename,command_line_args.path);
    strcat(filename,"statistics.csv");
   FILE *fd=fopen(filename,"w");
   if(fd==NULL){
      fatal("Cannot open statistics.csv file for writing",1);
   }
   fprintf(fd, "Epoch OuterLoop VertexMoveSucessRate BondFlipSuccessRate Volume Area lamdba1 lambda2 lambda3 Kc(2-9) Kc(6-9) Kc(2-end) Kc(3-6)\n");
    if(vesicle->sphHarmonics!=NULL){
        strcpy(filename,command_line_args.path);
        strcat(filename,"ulm2.csv");
      fd2=fopen(filename,"w");
      if(fd2==NULL){
         fatal("Cannot open ulm2.csv file for writing",1);
      }
      fprintf(fd2, "Timestep u_00^2 u_10^2 u_11^2 u_20^2 ...\n");
   }
/* RANDOM SEED SET BY CURRENT TIME */
   epochtime=get_epoch();
   srand48(epochtime);
   centermass(vesicle);
   cell_occupation(vesicle);
   vesicle_volume(vesicle); //needed for constant volume at this moment
    vesicle_area(vesicle); //needed for constant area at this moment
   V0=vesicle->volume;
    A0=vesicle->area;
   epsvol=4.0*sqrt(2.0*M_PI)/pow(3.0,3.0/4.0)*V0/pow(vesicle->tlist->n,3.0/2.0);
    epsarea=A0/(ts_double)vesicle->tlist->n;
  //  fprintf(stderr, "DVol=%1.16f (%1.16f), V0=%1.16f\n", epsvol,0.003e-2*V0,V0);
   if(start_iteration<inititer) ts_fprintf(stdout, "Starting simulation (first %d x %d MC sweeps will not be recorded on disk)\n", inititer, mcsweeps);
   for(i=start_iteration;i<inititer+iterations;i++){
      vmsr=0.0;
      bfsr=0.0;
/*    vesicle_volume(vesicle);
    fprintf(stderr,"Volume before TS=%1.16e\n", vesicle->volume); */
      for(j=0;j<mcsweeps;j++){
         single_timestep(vesicle, &vmsrt, &bfsrt);
         vmsr+=vmsrt;
         bfsr+=bfsrt;
      }
/*
    vesicle_volume(vesicle);
    fprintf(stderr,"Volume after TS=%1.16e\n", vesicle->volume); */
      vmsr/=(ts_double)mcsweeps;
      bfsr/=(ts_double)mcsweeps;
      centermass(vesicle);
@@ -45,11 +77,11 @@
            dump_state(vesicle,i);
      if(i>=inititer){
         write_vertex_xml_file(vesicle,i-inititer);
         write_master_xml_file("test.pvd");
         write_master_xml_file(command_line_args.output_fullfilename);
         epochtime=get_epoch();         
         gyration_eigen(vesicle, &l1, &l2, &l3);
         vesicle_volume(vesicle); //calculates just volume. Area is not added to ts_vesicle yet!
         get_area_volume(vesicle, &area,&volume); //that's why I must recalculate area (and volume for no particular reason).
         vesicle_volume(vesicle); //calculates just volume.
            vesicle_area(vesicle); //calculates area.
         r0=getR0(vesicle);
            if(vesicle->sphHarmonics!=NULL){
             preparationSh(vesicle,r0);
@@ -61,8 +93,9 @@
                kc2=calculateKc(vesicle, 6,9);
                kc3=calculateKc(vesicle, 2,vesicle->sphHarmonics->l);
                kc4=calculateKc(vesicle, 3,6);
            fd1=fopen("state.dat","w");
                strcpy(filename,command_line_args.path);
                strcat(filename,"state.dat");
            fd1=fopen(filename,"w");
            fprintf(fd1,"%e %e\n",vesicle->volume, getR0(vesicle));
            for(k=0;k<vesicle->vlist->n;k++){
               fprintf(fd1,"%e %e %e %e %e\n",
@@ -74,19 +107,34 @@
               );
            }
            fclose(fd1);
         fprintf(fd2,"%u ", i);
         for(l=0;l<vesicle->sphHarmonics->l;l++){
            for(m=l;m<2*l+1;m++){
               fprintf(fd2,"%e ", gsl_complex_abs2(vesicle->sphHarmonics->ulmComplex[l][m]) );
            }
         }
            fprintf(fd2,"\n");
             fflush(fd2);
            }
         fprintf(fd, "%lu %u %e %e %1.16e %1.16e %1.16e %1.16e %1.16e %1.16e %1.16e %1.16e %1.16e\n",epochtime,i,vmsr,bfsr,volume, area,l1,l2,l3,kc1, kc2, kc3,kc4);
         fprintf(fd, "%lu %u %e %e %1.16e %1.16e %1.16e %1.16e %1.16e %1.16e %1.16e %1.16e %1.16e\n",epochtime,i,vmsr,bfsr,vesicle->volume, vesicle->area,l1,l2,l3,kc1, kc2, kc3,kc4);
          fflush(fd);   
      //   sprintf(filename,"timestep-%05d.pov",i-inititer);
      //   write_pov_file(vesicle,filename);
      }
   }
   fclose(fd);
   if(fd2!=NULL) fclose(fd2);
   return TS_SUCCESS;
}
ts_bool single_timestep(ts_vesicle *vesicle,ts_double *vmsr, ts_double *bfsr){
//    vesicle_volume(vesicle);
//    fprintf(stderr,"Volume before TS=%1.16e\n", vesicle->volume);
    ts_bool retval;
    ts_double rnvec[3];
    ts_uint i,j, b;
@@ -105,6 +153,7 @@
        //find a bond and return a pointer to a bond...
        //call single_bondflip_timestep...
        retval=single_bondflip_timestep(vesicle,vesicle->blist->bond[b],rnvec);
       //     b++; retval=TS_FAIL;
   if(retval==TS_SUCCESS) bfsrcnt++;        
    }
@@ -131,6 +180,8 @@
//   printf("Bondflip success rate in one sweep: %d/%d=%e\n", cnt,3*vesicle->blist->n,(double)cnt/(double)vesicle->blist->n/3.0);
   *vmsr=(ts_double)vmsrcnt/(ts_double)vesicle->vlist->n;
   *bfsr=(ts_double)bfsrcnt/(ts_double)vesicle->vlist->n/3.0;
//    vesicle_volume(vesicle);
//    fprintf(stderr,"Volume after TS=%1.16e\n", vesicle->volume);
    return TS_SUCCESS;
}