Trisurf Monte Carlo simulator
mihaf
2014-03-24 ea1cce4c80ddb68c0b46a6873ef3b64cc4a958a1
commit | author | age
7958e9 1 #include<stdlib.h>
SP 2 #include<math.h>
3 #include<stdio.h>
4 #include "general.h"
5 #include "vertex.h"
6 #include "bond.h"
7 #include "vesicle.h"
8 #include "vertex.h"
9 #include "triangle.h"
10 #include "initial_distribution.h"
f74313 11 #include "energy.h"
1ab449 12 #include "poly.h"
8a6614 13 #include "io.h"
7958e9 14
SP 15 ts_vesicle *initial_distribution_dipyramid(ts_uint nshell, ts_uint ncmax1, ts_uint ncmax2, ts_uint ncmax3, ts_double stepsize){
1ab449 16     ts_fprintf(stdout,"Starting initial_distribution on vesicle with %u shells!...\n",nshell);
7958e9 17     ts_bool retval;
1ab449 18     ts_uint no_vertices=5*nshell*nshell+2;    
SP 19     ts_vesicle *vesicle=init_vesicle(no_vertices,ncmax1,ncmax2,ncmax3,stepsize);
20     vesicle->nshell=nshell;
21     //retval = vtx_set_global_values(vesicle);
22     retval = pentagonal_dipyramid_vertex_distribution(vesicle->vlist);
23     retval = init_vertex_neighbours(vesicle->vlist);
24     vesicle->vlist = init_sort_neighbours(vesicle->blist,vesicle->vlist);
b01cc1 25    // retval = init_vesicle_bonds(vesicle); // bonds are created in sort_neigh
1ab449 26     retval = init_triangles(vesicle);
SP 27     retval = init_triangle_neighbours(vesicle);
28     retval = init_common_vertex_triangle_neighbours(vesicle);
29     retval = init_normal_vectors(vesicle->tlist);
30     retval = mean_curvature_and_energy(vesicle);
31     ts_fprintf(stdout,"initial_distribution finished!\n");
41a035 32     if(retval);
7958e9 33     return vesicle;
SP 34
35
36
1ab449 37
SP 38 ts_vesicle *create_vesicle_from_tape(ts_tape *tape){
39     ts_vesicle *vesicle;
bcf455 40     ts_vertex *vtx;
M 41
1ab449 42     vesicle=initial_distribution_dipyramid(tape->nshell,tape->ncxmax,tape->ncymax,tape->nczmax,tape->stepsize);
58230a 43     // Nucleus:
M 44     vesicle->R_nucleus=tape->R_nucleus;
bcf455 45
58230a 46     //Initialize grafted polymers (brush):
624f81 47     vesicle->poly_list=init_poly_list(tape->npoly,tape->nmono, vesicle->vlist, vesicle);
1ab449 48     vesicle->spring_constant=tape->kspring;
SP 49     poly_assign_spring_const(vesicle);
bcf455 50
58230a 51     //Initialize filaments (polymers inside the vesicle):
M 52     vesicle->filament_list=init_poly_list(tape->nfil,tape->nfono, NULL, vesicle);
bcf455 53     poly_assign_filament_xi(vesicle,tape);
58230a 54
bcf455 55     ts_uint i,j;
M 56     for(i=0;i<vesicle->filament_list->n;i++){
57         for(j=0;j<vesicle->filament_list->poly[i]->blist->n;j++){
58             bond_vector(vesicle->filament_list->poly[i]->blist->bond[j]);
59             vesicle->filament_list->poly[i]->blist->bond[j]->bond_length = sqrt(vtx_distance_sq(vesicle->filament_list->poly[i]->blist->bond[j]->vtx1,vesicle->filament_list->poly[i]->blist->bond[j]->vtx2));
60         }
58230a 61     }
bcf455 62
M 63     for(i=0;i<vesicle->filament_list->n;i++){
64         for(j=0;j<vesicle->filament_list->poly[i]->vlist->n;j++){
65             vtx = vesicle->filament_list->poly[i]->vlist->vtx[j];
66             if(vtx->bond_no == 2){
67             vtx->energy = -(vtx->bond[0]->x*vtx->bond[1]->x + vtx->bond[0]->y*vtx->bond[1]->y + vtx->bond[0]->z*vtx->bond[1]->z)/vtx->bond[0]->bond_length/vtx->bond[1]->bond_length;
68             }
69         }
58230a 70     }
bcf455 71
ea1cce 72     for(i=0;i<vesicle->filament_list->n;i++){
M 73         vertex_list_assign_id(vesicle->filament_list->poly[i]->vlist,TS_ID_FILAMENT);
74     }
bcf455 75
58230a 76 //    vesicle->spring_constant=tape->kspring;
M 77 //    poly_assign_spring_const(vesicle);
78
1ab449 79     
SP 80     vesicle->nshell=tape->nshell;
81     vesicle->dmax=tape->dmax*tape->dmax; /* dmax^2 in the vesicle dmax variable */
ea1cce 82     vesicle->clist->dmin_interspecies = tape->dmin_interspecies*tape->dmin_interspecies;
1ab449 83     vesicle->bending_rigidity=tape->xk0;
SP 84     vtx_set_global_values(vesicle); /* make xk0 default value for every vertex */ 
85     ts_fprintf(stdout, "Tape setting: xk0=%e\n",tape->xk0);
86     vesicle->stepsize=tape->stepsize;
87     vesicle->clist->ncmax[0]=tape->ncxmax;
88     vesicle->clist->ncmax[1]=tape->ncymax;
89     vesicle->clist->ncmax[2]=tape->nczmax;
90     vesicle->clist->max_occupancy=8; /* hard coded max occupancy? */
91
92     vesicle->pressure= tape->pressure;
93     vesicle->pswitch=tape->pswitch;
94
95     return vesicle;
96
97 }
98
99
100
101
102
7958e9 103 ts_bool pentagonal_dipyramid_vertex_distribution(ts_vertex_list *vlist){
SP 104     /* Some often used relations */
105     const ts_double s1= sin(2.0*M_PI/5.0);
106     const ts_double s2= sin(4.0*M_PI/5.0);
107     const ts_double c1= cos(2.0*M_PI/5.0);
108     const ts_double c2= cos(4.0*M_PI/5.0);
109
110     /* Calculates projection lenght of an edge bond to pentagram plane */
111     const ts_double xl0=A0/(2.0*sin(M_PI/5.0));
112 #ifdef TS_DOUBLE_DOUBLE
113     const ts_double z0=sqrt(pow(A0,2)-pow(xl0,2));
114 #endif
115 #ifdef TS_DOUBLE_FLOAT
116     const ts_double z0=sqrtf(powf(A0,2)-powf(xl0,2));
117 #endif
118 #ifdef TS_DOUBLE_LONGDOUBLE
119     const ts_double z0=sqrtl(powl(A0,2)-powl(xl0,2));
120 #endif
121 //    const z0=sqrt(A0*A0 -xl0*xl0); /* I could use pow function but if pow is used make a check on the float type. If float then powf, if long double use powl */
122
123 /*placeholder for the pointer to vertex datastructure list... DIRTY: actual pointer points towards invalid address, one position before actual beginning of the list... This is to solve the difference between 1 based indexing in original program in fortran and 0 based indexing in C. All algorithms remain unchanged because of this!*/
124     ts_vertex **vtx=vlist->vtx -1 ; 
125
126
127     ts_uint nshell=(ts_uint)( sqrt((ts_double)(vlist->n-2)/5));
128 //    printf("nshell=%u\n",nshell);
129     ts_uint i,n0; // some for loop prereq
130     ts_int j,k;
131     ts_double dx,dy; // end loop prereq
132
133     /* topmost vertex */
8f6a69 134     vtx[1]->x=0.0;
SP 135     vtx[1]->y=0.0;
136     vtx[1]->z=z0*(ts_double)nshell;
7958e9 137     
SP 138     /* starting from to in circular order on pentagrams */    
139     for(i=1;i<=nshell;i++){
140         n0=2+5*i*(i-1)/2; //-1 would be for the reason that C index starts from 0 
8f6a69 141         vtx[n0]->x=0.0;
SP 142         vtx[n0]->y=(ts_double)i*xl0;
143         vtx[n0+i]->x=vtx[n0]->y*s1;
144         vtx[n0+i]->y=vtx[n0]->y*c1;
145         vtx[n0+2*i]->x=vtx[n0]->y*s2;
146         vtx[n0+2*i]->y=vtx[n0]->y*c2;
147         vtx[n0+3*i]->x=-vtx[n0+2*i]->x;
148         vtx[n0+3*i]->y=vtx[n0+2*i]->y;
149         vtx[n0+4*i]->x=-vtx[n0+i]->x;
150         vtx[n0+4*i]->y=vtx[n0+i]->y;
7958e9 151     }
SP 152
153     /* vertexes on the faces of the dipyramid */
154     for(i=1;i<=nshell;i++){
155         n0=2+5*i*(i-1)/2; // -1 would be because of C!
156         for(j=1;j<=i-1;j++){
8f6a69 157             dx=(vtx[n0]->x-vtx[n0+4*i]->x)/(ts_double)i;
SP 158             dy=(vtx[n0]->y-vtx[n0+4*i]->y)/(ts_double)i;
159             vtx[n0+4*i+j]->x=(ts_double)j*dx+vtx[n0+4*i]->x;
160             vtx[n0+4*i+j]->y=(ts_double)j*dy+vtx[n0+4*i]->y;
7958e9 161         }
SP 162         for(k=0;k<=3;k++){ // I would be worried about zero starting of for
8f6a69 163             dx=(vtx[n0+(k+1)*i]->x - vtx[n0+k*i]->x)/(ts_double) i;
SP 164             dy=(vtx[n0+(k+1)*i]->y - vtx[n0+k*i]->y)/(ts_double) i;
7958e9 165             for(j=1; j<=i-1;j++){
8f6a69 166                 vtx[n0+k*i+j]->x= (ts_double)j*dx+vtx[n0+k*i]->x;
SP 167                 vtx[n0+k*i+j]->y= (ts_double)j*dy+vtx[n0+k*i]->y;
7958e9 168             } 
SP 169         } 
170     }
171
172     for(i=1;i<=nshell;i++){
173         n0= 2+ 5*i*(i-1)/2;
174         for(j=0;j<=5*i-1;j++){
8f6a69 175         vtx[n0+j]->z= z0*(ts_double)(nshell-i);   // I would be worried about zero starting of for
7958e9 176         }
SP 177     }
178
179 /* for botom part of dipyramide we calculate the positions of vertices */
180     for(i=2+5*nshell*(nshell+1)/2;i<=vlist->n;i++){
8f6a69 181         vtx[i]->x=vtx[vlist->n - i +1]->x;
SP 182         vtx[i]->y=vtx[vlist->n - i +1]->y;
183         vtx[i]->z=-vtx[vlist->n - i +1]->z;
7958e9 184     }
SP 185
186     for(i=1;i<=vlist->n;i++){
187         for(j=1;j<=vlist->n;j++){
188             if(i!=j && vtx_distance_sq(vtx[i],vtx[j])<0.001){
189                 printf("Vertices %u and %u are the same!\n",i,j);
190             }
191         }
192     }
193     return TS_SUCCESS;
194 }
195
196
197
198 ts_bool init_vertex_neighbours(ts_vertex_list *vlist){
199     ts_vertex **vtx=vlist->vtx -1; // take a look at dipyramid function for comment.
200     const ts_double eps=0.001; //TODO: find out if you can use EPS from math.h
201     ts_uint i,j;
202     ts_double dist2; // Square of distance of neighbours
203     /*this is not required if we zero all data in vertex structure at initialization */
204     /*if we force zeroing at initialization this for loop can safely be deleted */
205     //for(i=1;i<=vlist->n;i++){
206     //    vtx[i].neigh_no=0;
207     //}
208     for(i=1;i<=vlist->n;i++){
209         for(j=1;j<=vlist->n;j++){
210             dist2=vtx_distance_sq(vtx[i],vtx[j]);
211             if( (dist2>eps) && (dist2<(A0*A0+eps))){ 
212     //if it is close enough, but not too much close (solves problem of comparing when i==j)
213                 vtx_add_neighbour(vtx[i],vtx[j]);
214             }
215         }
216     //        printf ("vertex %u ima %u sosedov!\n",i,vtx[i]->data->neigh_no);
217     }
218
219     return TS_SUCCESS;
220 }
221
b01cc1 222 // TODO: with new datastructure can be rewritten. Partially it is done, but it is complicated.
SP 223 ts_vertex_list *init_sort_neighbours(ts_bond_list *blist,ts_vertex_list *vlist){
7958e9 224     ts_vertex **vtx=vlist->vtx -1; // take a look at dipyramid function for comment.
SP 225     ts_uint i,l,j,jj,jjj,k=0;   
226     ts_double eps=0.001; // Take a look if EPS from math.h can be used
227
228 /*lets initialize memory for temporary vertex_list. Should we write a function instead */
b01cc1 229     ts_vertex_list *tvlist=vertex_list_copy(vlist);
7958e9 230     ts_vertex **tvtx=tvlist->vtx -1;  /* again to compensate for 0-indexing */
SP 231
232     ts_double dist2; // Square of distance of neighbours
233     ts_double direct; // Something, dont know what, but could be normal of some kind
234     for(i=1;i<=vlist->n;i++){
235         k++; // WHY i IS NOT GOOD??
8f6a69 236            vtx_add_cneighbour(blist,tvtx[k], tvtx[vtx[i]->neigh[0]->idx+1]); //always add 1st
7958e9 237            jjj=1;
SP 238            jj=1;
8f6a69 239            for(l=2;l<=vtx[i]->neigh_no;l++){
SP 240                for(j=2;j<=vtx[i]->neigh_no;j++){
241                    dist2=vtx_distance_sq(vtx[i]->neigh[j-1],vtx[i]->neigh[jj-1]);
242                    direct=vtx_direct(vtx[i],vtx[i]->neigh[j-1],vtx[i]->neigh[jj-1]);
243 // TODO: check if fabs can be used with all floating point types!!
7958e9 244                    if( (fabs(dist2-A0*A0)<=eps) && (direct>0.0) && (j!=jjj) ){
8f6a69 245                        vtx_add_cneighbour(blist,tvtx[k],tvtx[vtx[i]->neigh[j-1]->idx+1]);
7958e9 246                        jjj=jj;
SP 247                        jj=j;
248                        break;
249                    }
250                }
251            }    
252     }
b01cc1 253 /* We use the temporary vertex for our main vertices and we abandon main
SP 254  * vertices, because their neighbours are not correctly ordered */
255    // tvtx=vlist->vtx;
256    // vlist->vtx=tvtx;
257    // tvlist->vtx=vtx;
258     vtx_list_free(vlist);
259 /* Let's make a check if the number of bonds is correct */
260     if((blist->n)!=3*(tvlist->n-2)){
261         ts_fprintf(stderr,"Number of bonds is %u should be %u!\n", blist->n, 3*(tvlist->n-2));
262         fatal("Number of bonds is not 3*(no_vertex-2).",4);
7958e9 263     }
SP 264
b01cc1 265     return tvlist;
7958e9 266 }
SP 267
268
269 ts_bool init_vesicle_bonds(ts_vesicle *vesicle){
270     ts_vertex_list *vlist=vesicle->vlist;
271     ts_bond_list *blist=vesicle->blist;
272     ts_vertex **vtx=vesicle->vlist->vtx - 1; // Because of 0 indexing
273 /* lets make correct clockwise ordering of in nearest neighbour list */
274     ts_uint i,j,k;
275     for(i=1;i<=vlist->n;i++){
276         for(j=i+1;j<=vlist->n;j++){
8f6a69 277             for(k=0;k<vtx[i]->neigh_no;k++){ // has changed 0 to < instead of 1 and <=
SP 278                 if(vtx[i]->neigh[k]==vtx[j]){  //if addresses matches it is the same
7958e9 279                     bond_add(blist,vtx[i],vtx[j]);
SP 280                     break;
281                 }
282             }
283         }
284     } 
285 /* Let's make a check if the number of bonds is correct */
286     if((blist->n)!=3*(vlist->n-2)){
287         ts_fprintf(stderr,"Number of bonds is %u should be %u!\n", blist->n, 3*(vlist->n-2));
288         fatal("Number of bonds is not 3*(no_vertex-2).",4);
289     }
290     return TS_SUCCESS;
291 }
292
293
294
295 ts_bool init_triangles(ts_vesicle *vesicle){
296     ts_uint i,j,jj,k;
297     ts_vertex **vtx=vesicle->vlist->vtx -1; // difference between 0 indexing and 1 indexing
298     ts_triangle_list *tlist=vesicle->tlist;
299     ts_double dist, direct;
300     ts_double eps=0.001; // can we use EPS from math.h?
301     k=0;
302     for(i=1;i<=vesicle->vlist->n;i++){
8f6a69 303         for(j=1;j<=vtx[i]->neigh_no;j++){
SP 304             for(jj=1;jj<=vtx[i]->neigh_no;jj++){
7958e9 305         //        ts_fprintf(stderr,"%u: (%u,%u) neigh_no=%u ",i,j,jj,vtx[i].neigh_no);
SP 306         //      ts_fprintf(stderr,"%e, %e",vtx[i].neigh[j-1]->x,vtx[i].neigh[jj-1]->x);
8f6a69 307                 dist=vtx_distance_sq(vtx[i]->neigh[j-1],vtx[i]->neigh[jj-1]);
SP 308                 direct=vtx_direct(vtx[i],vtx[i]->neigh[j-1],vtx[i]->neigh[jj-1]);                
309 // TODO: same as above                
310                 if(fabs(dist-A0*A0)<=eps && direct < 0.0 && vtx[i]->neigh[j-1]->idx+1 > i && vtx[i]->neigh[jj-1]->idx+1 >i){
311                     triangle_add(tlist,vtx[i],vtx[i]->neigh[j-1],vtx[i]->neigh[jj-1]);
7958e9 312                 }    
SP 313             }    
314         }
315     }
316 /* We check if all triangles have 3 vertices and if the number of triangles
317  * matches the theoretical value.
318  */
319     for(i=0;i<tlist->n;i++){
320         k=0;
321         for(j=0;j<3;j++){
41a035 322             if(tlist->tria[i]->vertex[j]!=NULL)
7958e9 323             k++;
SP 324         }
325             if(k!=3){
8f6a69 326                 fatal("Some triangles have less than 3 vertices..",4);
7958e9 327             }   
SP 328     } 
329     if(tlist->n!=2*(vesicle->vlist->n -2)){
330         ts_fprintf(stderr,"The number of triangles is %u but should be %u!\n",tlist->n,2*(vesicle->vlist->n -2));
331         fatal("The number of triangles doesn't match 2*(no_vertex -2).",4);
332     }
333     return TS_SUCCESS;
334 }
335
336
337
338 ts_bool init_triangle_neighbours(ts_vesicle *vesicle){
339     ts_uint i,j,nobo;
340     ts_vertex *i1,*i2,*i3,*j1,*j2,*j3;
341 //    ts_vertex **vtx=vesicle->vlist->vtx -1; // difference between 0 indexing and 1 indexing
342     ts_triangle_list *tlist=vesicle->tlist;
343     ts_triangle **tria=tlist->tria -1;
344     nobo=0;
345     for(i=1;i<=tlist->n;i++){
41a035 346         i1=tria[i]->vertex[0]; 
SP 347         i2=tria[i]->vertex[1]; 
348         i3=tria[i]->vertex[2]; 
7958e9 349         for(j=1;j<=tlist->n;j++){
SP 350             if(j==i) continue;
41a035 351             j1=tria[j]->vertex[0]; 
SP 352             j2=tria[j]->vertex[1]; 
353             j3=tria[j]->vertex[2]; 
7958e9 354             if((i1==j1 && i3==j2) || (i1==j2 && i3==j3) || (i1==j3 && i3==j1)){
SP 355                     triangle_add_neighbour(tria[i],tria[j]);
356                     nobo++;
357             }
358         }
359     }
360     for(i=1;i<=tlist->n;i++){
41a035 361         i1=tria[i]->vertex[0]; 
SP 362         i2=tria[i]->vertex[1]; 
363         i3=tria[i]->vertex[2]; 
7958e9 364         for(j=1;j<=tlist->n;j++){
SP 365             if(j==i) continue;
41a035 366             j1=tria[j]->vertex[0]; 
SP 367             j2=tria[j]->vertex[1]; 
368             j3=tria[j]->vertex[2]; 
7958e9 369             if((i1==j1 && i2==j3) || (i1==j3 && i2==j2) || (i1==j2 && i2==j1)){
SP 370                 triangle_add_neighbour(tria[i],tria[j]);
371                 nobo++;
372             }
373         }
374     }
375     for(i=1;i<=tlist->n;i++){
41a035 376         i1=tria[i]->vertex[0]; 
SP 377         i2=tria[i]->vertex[1]; 
378         i3=tria[i]->vertex[2]; 
7958e9 379         for(j=1;j<=tlist->n;j++){
SP 380             if(j==i) continue;
41a035 381             j1=tria[j]->vertex[0]; 
SP 382             j2=tria[j]->vertex[1]; 
383             j3=tria[j]->vertex[2]; 
7958e9 384             if((i2==j1 && i3==j3) || (i2==j3 && i3==j2) || (i2==j2 && i3==j1)){
SP 385                 triangle_add_neighbour(tria[i],tria[j]);
386                 nobo++;
387             }
388         }
389     }
390     if(nobo != vesicle->blist->n*2) {
391             ts_fprintf(stderr,"Number of triangles= %u, number of bonds= %u\n",nobo/2, vesicle->blist->n);
392             fatal("Number of triangle neighbour pairs differs from double the number of bonds!",4);
393     }
394     return TS_SUCCESS;
395 }
396
397
398 ts_bool init_common_vertex_triangle_neighbours(ts_vesicle *vesicle){
399     ts_uint i,j,jp,k;
400     ts_vertex *k1,*k2,*k3,*k4,*k5;
401     ts_vertex **vtx=vesicle->vlist->vtx -1; // difference between 0 indexing and 1 indexing
402     ts_triangle_list *tlist=vesicle->tlist;
403     ts_triangle **tria=tlist->tria -1;
404
405     for(i=1;i<=vesicle->vlist->n;i++){
8f6a69 406         for(j=1;j<=vtx[i]->neigh_no;j++){
SP 407             k1=vtx[i]->neigh[j-1];
7958e9 408             jp=j+1;
8f6a69 409             if(j == vtx[i]->neigh_no) jp=1;
SP 410             k2=vtx[i]->neigh[jp-1];
7958e9 411             for(k=1;k<=tlist->n;k++){        // VERY NON-OPTIMAL!!! too many loops (vlist.n * vtx.neigh * tlist.n )!
41a035 412                 k3=tria[k]->vertex[0];
SP 413                 k4=tria[k]->vertex[1];
414                 k5=tria[k]->vertex[2];
7958e9 415 //                ts_fprintf(stderr,"%u %u: k=(%u %u %u)\n",k1,k2,k3,k4,k5);
SP 416                 if((vtx[i]==k3 && k1==k4 && k2==k5) ||
417                 (vtx[i]==k4 && k1==k5 && k2==k3) ||
418                 (vtx[i]==k5 && k1==k3 && k2==k4)){
b01cc1 419
SP 420 //TODO: probably something wrong with neighbour distribution.
421 //                if(vtx[i]==k3 || vtx[i]==k4 || vtx[i]==k5){
dac2e5 422     //                    if(i==6) ts_fprintf(stdout, "Vtx[%u] > Added to tristar!\n",i);
7958e9 423                     vertex_add_tristar(vtx[i],tria[k]);
SP 424                 }
425             }
426         }
427 /*        ts_fprintf(stderr,"TRISTAR for %u (%u):",i-1,vtx[i].tristar_no);
428         for(j=0;j<vtx[i].tristar_no;j++){
429             ts_fprintf(stderr," %u,",vtx[i].tristar[j]->idx);
430         }
431         ts_fprintf(stderr,"\n"); */
432     }
433     return TS_SUCCESS;
434 }
435
436
437 ts_bool init_normal_vectors(ts_triangle_list *tlist){
438     /* Normals point INSIDE vesicle */
439     ts_uint k;
440     ts_triangle **tria=tlist->tria -1; //for 0 indexing
441     for(k=1;k<=tlist->n;k++){
442         triangle_normal_vector(tria[k]);    
443     }
444     return TS_SUCCESS;
445 }