Trisurf Monte Carlo simulator
Samo Penic
2014-06-15 7d62fc765adc11bbcb7db0716ac5dc19145dca80
src/cross-section.c
@@ -1,7 +1,7 @@
#include<general.h>
#include<cross-section.h>
#include<coord.h>
#include<cairo/cairo.h>
/** @brief Calculates cross-section of vesicle with plane.
 *
 * Function returns points of cross-section of vesicle with plane. Plane is described with equation $ax+by+cz+d=0$. Algorithm extracts coordinates of each vertex of a vesicle and then:
@@ -9,30 +9,43 @@
 * if a distance of point to plane (given by equation $D=\frac{ax_0+by_0+cz_0+d}{\sqrt{a^2+b^2+c^2}}$, where $x_0$, $y_0$ and $z_0$ are coordinates of a given vertex) is less than maximal allowed distance between vertices {\tt sqrt(vesicle->dmax)} than vertex is a candidate for crossection calculation.
 *
 */
ts_coord_list get_crossection_with_plane(ts_vesicle vesicle,ts_double a,ts_double b,ts_double c, ts_double d){
ts_coord_list *get_crossection_with_plane(ts_vesicle *vesicle,ts_double a,ts_double b,ts_double c, ts_double d){
   ts_uint i, j, k;
   ts_uint i, j,k,l;
   ts_double pp,Dsq; // distance from the plane squared
   ts_double ppn,Dsqn; // distance from the plane squared of a neighbor
   ts_double u; //factor to scale vector from first vector to the second to get intersection
   ts_double ppn1; // distance from the plane squared of a neighbor
   ts_vertex *vtx;
   ts_uint ntria=0; // number triangles
   ts_triangle *tria[2]; // list of triangles
   ts_coord_list *pts=init_coord_list();
   for(i=0;i<vesicle->vlist->N;i++){
   for(i=0;i<vesicle->vlist->n;i++){
      vtx=vesicle->vlist->vtx[i];
      pp=vtx->x*a+vtx->y*b+vtx->z*c+d;
      Dsq=pp*pp/(a*a+b*b+c*c);
      if(Dsq<vesicle->dmax){
         for(j=0;j<vtx->neigh_no;j++){
            ppn=vtx->neigh[j]->x*a+vtx->neigh[j]->y*b+vtx->neigh[j]->z*c+d;
            if(pp*ppn<0){ //the combination of vertices are good candidates for a crossection
               u=pp/(a*(vtx->x-vtx->neigh[j]->x)+b*(vtx->y-vtx->neigh[j]->y)+c(vtx->z-vtx->neigh[j]->z));
               add_coord(pts, vtx->x+u(vtx->neigh[j]->x - vtx->x),
                     vtx->y+u(vtx->neigh[j]->y - vtx->y),
                     vtx->z+u(vtx->neigh[j]->z - vtx->z),
                     TS_COORD_CARTESIAN);
            ppn1=vtx->neigh[j]->x*a+vtx->neigh[j]->y*b+vtx->neigh[j]->z*c+d;
            if(pp*ppn1<=0.0){ //the combination of vertices are good candidates for a crossection
               //find triangle that belongs to the two vertices
               ntria=0;
               for(k=0;k<vtx->tristar_no;k++){
                  if(vtx->tristar[k]->vertex[0]==vtx && ( vtx->tristar[k]->vertex[1]==vtx->neigh[j] || vtx->tristar[k]->vertex[2]==vtx->neigh[j]) ){
                     //triangle found.
                     tria[ntria]=vtx->tristar[k];
                     ntria++;
                  }
               }
               if(ntria==0) continue; // no need to continue
               //find the two intersections (in general) to form a intersection line
               for(l=0;l<ntria;l++){
               //we add intersection line between two points for each of the triangles found above.
                  add_crosssection_point(pts,a,b,c,d,tria[l]->vertex[0], tria[l]->vertex[1]);
                  add_crosssection_point(pts,a,b,c,d,tria[l]->vertex[0], tria[l]->vertex[2]);
                  add_crosssection_point(pts,a,b,c,d,tria[l]->vertex[1], tria[l]->vertex[2]);
               }
            }
         }
      }
@@ -41,3 +54,55 @@
}
ts_bool add_crosssection_point(ts_coord_list *pts, ts_double a, ts_double b, ts_double c, ts_double d, ts_vertex *vtx1, ts_vertex *vtx2){
   ts_double pp=vtx1->x*a+vtx1->y*b+vtx1->z*c+d;
   ts_double pp2=vtx2->x*a+vtx1->y*b+vtx2->z*c+d;
   if(pp*pp2<=0.0){
   ts_double u=pp/(a*(vtx1->x-vtx2->x)+b*(vtx1->y-vtx2->y)+c*(vtx1->z-vtx2->z));
   add_coord(pts, vtx1->x+u*(vtx2->x - vtx1->x),
      vtx1->y+u*(vtx2->y - vtx1->y),
      vtx1->z+u*(vtx2->z - vtx1->z),
      TS_COORD_CARTESIAN);
   return TS_SUCCESS;
   } else {
   return TS_FAIL;
   }
}
/** Saves calculated crossection as a png image */
ts_bool crossection_to_png(ts_coord_list *pts, char *filename){
   cairo_surface_t *surface;
   cairo_t *cr;
   ts_uint i;
   surface = cairo_image_surface_create (CAIRO_FORMAT_RGB24, 1800, 1800);
   cr = cairo_create (surface);
   cairo_rectangle(cr, 0.0, 0.0, 1800,1800);
   cairo_set_source_rgb(cr, 0.3, 0.3, 0.3);
   cairo_fill(cr);
   cairo_set_line_width (cr, 5.0/30.0);
   cairo_set_line_cap(cr, CAIRO_LINE_CAP_ROUND);
   cairo_set_line_join(cr, CAIRO_LINE_JOIN_ROUND);
   cairo_translate(cr, 900,900);
   cairo_scale (cr, 30, 30);
   cairo_set_source_rgb (cr, 1.0, 1.0, 1.0);
   for(i=0;i<pts->n;i+=2){
      cairo_move_to(cr, pts->coord[i]->e1, pts->coord[i]->e2);
      cairo_line_to(cr, pts->coord[i+1]->e1, pts->coord[i+1]->e2);
   }
   cairo_stroke(cr);
   cairo_surface_write_to_png (surface,filename);
   cairo_surface_finish (surface);
   return TS_SUCCESS;
}
ts_bool save_crossection_snapshot(ts_coord_list *pts, ts_uint timestepno){
   char filename[255];
       sprintf(filename,"timestep_%.6u.png",timestepno);
   crossection_to_png(pts,filename);
   return TS_SUCCESS;
}