Trisurf Monte Carlo simulator
Samo Penic
2016-06-21 8db20362dcf46fe004de41ead6e52ea5b2a6f621
commit | author | age
7f6076 1 /* vim: set ts=4 sts=4 sw=4 noet : */
d7639a 2 #include<stdlib.h>
SP 3 #include<stdio.h>
aec47d 4 #include<math.h>
SP 5 //#include "io.h"
6 #include "general.h"
7 #include "timestep.h"
8 #include "vertexmove.h"
30ee9c 9 #include "bondflip.h"
d7a113 10 #include "frame.h"
SP 11 #include "io.h"
37d14a 12 #include "stats.h"
dc77e8 13 #include "sh.h"
459ff9 14 #include "shcomplex.h"
dc77e8 15 #include "vesicle.h"
5a3862 16 #include<gsl/gsl_complex.h>
M 17 #include<gsl/gsl_complex_math.h>
267db5 18 #include<string.h>
ac9826 19 #include <sys/stat.h>
SP 20
fedf2b 21
626811 22 ts_bool run_simulation(ts_vesicle *vesicle, ts_uint mcsweeps, ts_uint inititer, ts_uint iterations, ts_uint start_iteration){
5a3862 23     ts_uint i, j,k,l,m;
cfab63 24     ts_double r0,kc1=0,kc2=0,kc3=0,kc4=0;
c0ae90 25     ts_double l1,l2,l3,vmsr,bfsr, vmsrt, bfsrt;
37d14a 26     ts_ulong epochtime;
3d0247 27     FILE *fd1,*fd2=NULL,*fd3=NULL;
267db5 28      char filename[10000];
ac9826 29     //struct stat st;
SP 30     strcpy(filename,command_line_args.path);
31     strcat(filename,"statistics.csv");
32     //int result = stat(filename, &st);
33     FILE *fd;
34     if(start_iteration==0)
35         fd=fopen(filename,"w");
36     else
37         fd=fopen(filename,"a");
37d14a 38     if(fd==NULL){
SP 39         fatal("Cannot open statistics.csv file for writing",1);
40     }
ac9826 41     if(start_iteration==0)
SP 42         fprintf(fd, "Epoch OuterLoop VertexMoveSucessRate BondFlipSuccessRate Volume Area lamdba1 lambda2 lambda3 Kc(2-9) Kc(6-9) Kc(2-end) Kc(3-6)\n");
5a3862 43
M 44      if(vesicle->sphHarmonics!=NULL){
267db5 45         strcpy(filename,command_line_args.path);
SP 46         strcat(filename,"ulm2.csv"); 
ac9826 47 //    int result = stat(filename, &st);
SP 48     if(start_iteration==0)
267db5 49         fd2=fopen(filename,"w");
ac9826 50     else
SP 51         fd2=fopen(filename,"a");
5a3862 52         if(fd2==NULL){
M 53             fatal("Cannot open ulm2.csv file for writing",1);
54         }
ac9826 55         if(start_iteration==0) //file does not exist
SP 56             fprintf(fd2, "Timestep u_00^2 u_10^2 u_11^2 u_20^2 ...\n");    
5a3862 57
M 58     }
59
c60a49 60 /* RANDOM SEED SET BY CURRENT TIME */
M 61     epochtime=get_epoch();            
62     srand48(epochtime);
45c708 63 /*Nir Gov: randomly add spontaneous curvature for some vertices */
SP 64     for(i=0;i<200;i++){
65         int b=rand() % vesicle->vlist->n;
66         vesicle->vlist->vtx[b]->c=-0.1;
67     }
d7a113 68     centermass(vesicle);
SP 69     cell_occupation(vesicle);
fe5069 70     vesicle_volume(vesicle); //needed for constant volume at this moment
c0ae90 71     vesicle_area(vesicle); //needed for constant area at this moment
49981c 72     if(V0<0.000001) 
SP 73         V0=vesicle->volume; 
74     ts_fprintf(stdout,"Setting volume V0=%.17f\n",V0);
75     if(A0<0.000001)
76         A0=vesicle->area;
77         ts_fprintf(stdout,"Setting area A0=%.17f\n",A0);
a54977 78     epsvol=4.0*sqrt(2.0*M_PI)/pow(3.0,3.0/4.0)*V0/pow(vesicle->tlist->n,3.0/2.0);
c0ae90 79     epsarea=A0/(ts_double)vesicle->tlist->n;
a54977 80   //  fprintf(stderr, "DVol=%1.16f (%1.16f), V0=%1.16f\n", epsvol,0.003e-2*V0,V0);
626811 81     if(start_iteration<inititer) ts_fprintf(stdout, "Starting simulation (first %d x %d MC sweeps will not be recorded on disk)\n", inititer, mcsweeps);
SP 82     for(i=start_iteration;i<inititer+iterations;i++){
37d14a 83         vmsr=0.0;
SP 84         bfsr=0.0;
3de289 85 /*    vesicle_volume(vesicle);
SP 86     fprintf(stderr,"Volume before TS=%1.16e\n", vesicle->volume); */
d7a113 87         for(j=0;j<mcsweeps;j++){
37d14a 88             single_timestep(vesicle, &vmsrt, &bfsrt);
SP 89             vmsr+=vmsrt;
90             bfsr+=bfsrt;
d7a113 91         }
3de289 92 /*
SP 93     vesicle_volume(vesicle);
94     fprintf(stderr,"Volume after TS=%1.16e\n", vesicle->volume); */
37d14a 95         vmsr/=(ts_double)mcsweeps;
SP 96         bfsr/=(ts_double)mcsweeps;
d7a113 97         centermass(vesicle);
SP 98         cell_occupation(vesicle);
1ab449 99             dump_state(vesicle,i);
58230a 100         if(i>=inititer){
d7a113 101             write_vertex_xml_file(vesicle,i-inititer);
267db5 102             write_master_xml_file(command_line_args.output_fullfilename);
37d14a 103             epochtime=get_epoch();            
SP 104             gyration_eigen(vesicle, &l1, &l2, &l3);
c0ae90 105             vesicle_volume(vesicle); //calculates just volume. 
SP 106             vesicle_area(vesicle); //calculates area.
dc77e8 107             r0=getR0(vesicle);
632960 108             if(vesicle->sphHarmonics!=NULL){
SP 109                 preparationSh(vesicle,r0);
459ff9 110                 //calculateYlmi(vesicle);
SP 111                 calculateUlmComplex(vesicle);
112                 storeUlmComplex2(vesicle);
632960 113                 saveAvgUlm2(vesicle);
22cdfd 114                 kc1=calculateKc(vesicle, 2,9);
SP 115                 kc2=calculateKc(vesicle, 6,9);
116                 kc3=calculateKc(vesicle, 2,vesicle->sphHarmonics->l);
1665aa 117                 kc4=calculateKc(vesicle, 3,6);
267db5 118                 strcpy(filename,command_line_args.path);
SP 119                 strcat(filename,"state.dat");  
120                 fd1=fopen(filename,"w");
5bb6bb 121                 fprintf(fd1,"%e %e\n",vesicle->volume, getR0(vesicle));
M 122                 for(k=0;k<vesicle->vlist->n;k++){
123                     fprintf(fd1,"%e %e %e %e %e\n",
124                         vesicle->vlist->vtx[k]->x,
125                         vesicle->vlist->vtx[k]->y,
126                         vesicle->vlist->vtx[k]->z,
127                         vesicle->vlist->vtx[k]->solAngle,
128                         vesicle->vlist->vtx[k]->relR
129                     );
130                 }
131                 fclose(fd1);
5a3862 132         
M 133             fprintf(fd2,"%u ", i);
134             for(l=0;l<vesicle->sphHarmonics->l;l++){
135                 for(m=l;m<2*l+1;m++){
136                     fprintf(fd2,"%e ", gsl_complex_abs2(vesicle->sphHarmonics->ulmComplex[l][m]) );
137                 }
138             }
139                 fprintf(fd2,"\n");
140     
141                 fflush(fd2);    
142
632960 143             }
dc77e8 144
c0ae90 145             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);
5a3862 146
632960 147             fflush(fd);    
144784 148         //    sprintf(filename,"timestep-%05d.pov",i-inititer);
fe24d2 149         //    write_pov_file(vesicle,filename);
49981c 150         } //end if(inititer....)
SP 151         fd3=fopen(".status","w"); //write status file when everything is written to disk.
152         if(fd3==NULL){
153             fatal("Cannot open .status file for writing",1);
d7a113 154         }
49981c 155         fprintf(fd3,"%d",i);
SP 156         fclose(fd3);
157         ts_fprintf(stdout,"Done %d out of %d iterations (x %d MC sweeps).\n",i+1,inititer+iterations,mcsweeps);
d7a113 158     }
37d14a 159     fclose(fd);
5a3862 160     if(fd2!=NULL) fclose(fd2);
d7a113 161     return TS_SUCCESS;
SP 162 }
d7639a 163
37d14a 164 ts_bool single_timestep(ts_vesicle *vesicle,ts_double *vmsr, ts_double *bfsr){
3de289 165 //    vesicle_volume(vesicle);
SP 166 //    fprintf(stderr,"Volume before TS=%1.16e\n", vesicle->volume);
d7639a 167     ts_bool retval;
SP 168     ts_double rnvec[3];
fe5069 169     ts_uint i,j, b;
37d14a 170     ts_uint vmsrcnt=0;
aec47d 171     for(i=0;i<vesicle->vlist->n;i++){
d7639a 172         rnvec[0]=drand48();
SP 173         rnvec[1]=drand48();
174         rnvec[2]=drand48();
aec47d 175         retval=single_verticle_timestep(vesicle,vesicle->vlist->vtx[i],rnvec);
37d14a 176     if(retval==TS_SUCCESS) vmsrcnt++;        
d7639a 177     }
SP 178
37d14a 179     ts_int bfsrcnt=0;
fedf2b 180     for(i=0;i<3*vesicle->vlist->n;i++){
fe5069 181     b=rand() % vesicle->blist->n;
d7639a 182         //find a bond and return a pointer to a bond...
SP 183         //call single_bondflip_timestep...
fe5069 184         retval=single_bondflip_timestep(vesicle,vesicle->blist->bond[b],rnvec);
3de289 185        //     b++; retval=TS_FAIL;
37d14a 186     if(retval==TS_SUCCESS) bfsrcnt++;        
fedf2b 187     }
M 188
189     for(i=0;i<vesicle->poly_list->n;i++){
58230a 190         for(j=0;j<vesicle->poly_list->poly[i]->vlist->n;j++){
M 191             rnvec[0]=drand48();
192             rnvec[1]=drand48();
193             rnvec[2]=drand48();
194             retval=single_poly_vertex_move(vesicle,vesicle->poly_list->poly[i],vesicle->poly_list->poly[i]->vlist->vtx[j],rnvec);    
195         }
fedf2b 196     }
M 197
58230a 198
M 199     for(i=0;i<vesicle->filament_list->n;i++){
200         for(j=0;j<vesicle->filament_list->poly[i]->vlist->n;j++){
201             rnvec[0]=drand48();
202             rnvec[1]=drand48();
203             rnvec[2]=drand48();
204             retval=single_filament_vertex_move(vesicle,vesicle->filament_list->poly[i],vesicle->filament_list->poly[i]->vlist->vtx[j],rnvec);    
205         }
fedf2b 206     }
M 207  
58230a 208
fedf2b 209 //    printf("Bondflip success rate in one sweep: %d/%d=%e\n", cnt,3*vesicle->blist->n,(double)cnt/(double)vesicle->blist->n/3.0);
37d14a 210     *vmsr=(ts_double)vmsrcnt/(ts_double)vesicle->vlist->n;
SP 211     *bfsr=(ts_double)bfsrcnt/(ts_double)vesicle->vlist->n/3.0;
3de289 212 //    vesicle_volume(vesicle);
SP 213 //    fprintf(stderr,"Volume after TS=%1.16e\n", vesicle->volume);
d7639a 214     return TS_SUCCESS;
SP 215 }
216
217
218