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