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