Trisurf Monte Carlo simulator
Samo Penic
2016-04-07 47d80d88e986fc039866c3a627ffe1140538d4b2
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
/* vim: set ts=4 sts=4 sw=4 noet : */
#include<stdlib.h>
#include<stdio.h>
#include<math.h>
//#include "io.h"
#include "general.h"
#include "timestep.h"
#include "vertexmove.h"
#include "bondflip.h"
#include "frame.h"
#include "io.h"
#include "stats.h"
#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,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,*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);
        cell_occupation(vesicle);
        ts_fprintf(stdout,"Done %d out of %d iterations (x %d MC sweeps).\n",i+1,inititer+iterations,mcsweeps);
            dump_state(vesicle,i);
        if(i>=inititer){
            write_vertex_xml_file(vesicle,i-inititer);
            write_master_xml_file(command_line_args.output_fullfilename);
            epochtime=get_epoch();            
            gyration_eigen(vesicle, &l1, &l2, &l3);
            vesicle_volume(vesicle); //calculates just volume. 
            vesicle_area(vesicle); //calculates area.
            r0=getR0(vesicle);
            if(vesicle->sphHarmonics!=NULL){
                preparationSh(vesicle,r0);
                //calculateYlmi(vesicle);
                calculateUlmComplex(vesicle);
                storeUlmComplex2(vesicle);
                saveAvgUlm2(vesicle);
                kc1=calculateKc(vesicle, 2,9);
                kc2=calculateKc(vesicle, 6,9);
                kc3=calculateKc(vesicle, 2,vesicle->sphHarmonics->l);
                kc4=calculateKc(vesicle, 3,6);
                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",
                        vesicle->vlist->vtx[k]->x,
                        vesicle->vlist->vtx[k]->y,
                        vesicle->vlist->vtx[k]->z,
                        vesicle->vlist->vtx[k]->solAngle,
                        vesicle->vlist->vtx[k]->relR
                    );
                }
                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,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;
    ts_uint vmsrcnt=0;
    for(i=0;i<vesicle->vlist->n;i++){
        rnvec[0]=drand48();
        rnvec[1]=drand48();
        rnvec[2]=drand48();
        retval=single_verticle_timestep(vesicle,vesicle->vlist->vtx[i],rnvec);
    if(retval==TS_SUCCESS) vmsrcnt++;        
    }
 
    ts_int bfsrcnt=0;
    for(i=0;i<3*vesicle->vlist->n;i++){
    b=rand() % vesicle->blist->n;
        //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++;        
    }
 
    for(i=0;i<vesicle->poly_list->n;i++){
        for(j=0;j<vesicle->poly_list->poly[i]->vlist->n;j++){
            rnvec[0]=drand48();
            rnvec[1]=drand48();
            rnvec[2]=drand48();
            retval=single_poly_vertex_move(vesicle,vesicle->poly_list->poly[i],vesicle->poly_list->poly[i]->vlist->vtx[j],rnvec);    
        }
    }
 
 
    for(i=0;i<vesicle->filament_list->n;i++){
        for(j=0;j<vesicle->filament_list->poly[i]->vlist->n;j++){
            rnvec[0]=drand48();
            rnvec[1]=drand48();
            rnvec[2]=drand48();
            retval=single_filament_vertex_move(vesicle,vesicle->filament_list->poly[i],vesicle->filament_list->poly[i]->vlist->vtx[j],rnvec);    
        }
    }
 
 
//    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;
}