Trisurf Monte Carlo simulator
Samo Penic
2021-04-19 17fe35ccc428e18dd226e07d5517c4816ef6be44
src/energy.c
@@ -9,6 +9,18 @@
#include <gsl/gsl_vector_complex.h>
#include <gsl/gsl_matrix.h>
#include <gsl/gsl_eigen.h>
int cmpfunc(const void *x, const void *y)
{
   double diff=   fabs(*(double*)x) - fabs(*(double*)y);
   if(diff<0) return 1;
   else return -1;
}
/** @brief Wrapper that calculates energy of every vertex in vesicle
 *  
 *  Function calculated energy of every vertex in vesicle. It can be used in
@@ -116,9 +128,11 @@
    ts_double We;
    ts_double Av, We_Av;
   ts_double eigenval[3];
   gsl_matrix *gsl_Sv=gsl_matrix_alloc(3,3);
   gsl_vector_complex *Sv_eigen=gsl_vector_complex_alloc(3);
   gsl_eigen_nonsymm_workspace *workspace=gsl_eigen_nonsymm_alloc(3);
   gsl_vector *Sv_eigen=gsl_vector_alloc(3);
   gsl_eigen_symm_workspace *workspace=gsl_eigen_symm_alloc(3);
   ts_double mprod[7], phi[7], he[7];
   ts_double Sv[3][3]={{0,0,0},{0,0,0},{0,0,0}};
@@ -222,7 +236,7 @@
   mprod[jj]=it->x*(k->y*edge_vector_z[jj]-edge_vector_y[jj]*k->z)-it->y*(k->x*edge_vector_z[jj]-k->z*edge_vector_x[jj])+it->z*(k->x*edge_vector_y[jj]-k->y*edge_vector_x[jj]);
   phi[jj]=copysign(acos(lm->xnorm*lp->xnorm+lm->ynorm*lp->ynorm+lm->znorm*lp->znorm-1e-8),mprod[jj])+M_PI;
   phi[jj]=copysign(acos(lm->xnorm*lp->xnorm+lm->ynorm*lp->ynorm+lm->znorm*lp->znorm-1e-15),mprod[jj])+M_PI;
//   printf("ACOS arg=%e\n", lm->xnorm*lp->xnorm+lm->ynorm*lp->ynorm+lm->znorm*lp->znorm);
   //he was multiplied with 2 before...
   he[jj]=sqrt( pow((edge_vector_x[jj]),2) + pow((edge_vector_y[jj]), 2) + pow((edge_vector_z[jj]), 2))*cos(phi[jj]/2.0);
@@ -265,21 +279,29 @@
//   printf("Se= %f, %f, %f\n    %f, %f, %f\n    %f, %f, %f\n", Se11, Se21, Se31, Se21, Se22, Se32, Se31, Se32, Se33);
//   printf("Pv= %f, %f, %f\n    %f, %f, %f\n    %f, %f, %f\n", Pv11, Pv21, Pv31, Pv21, Pv22, Pv32, Pv31, Pv32, Pv33);
   printf("Sv= %f, %f, %f\n    %f, %f, %f\n    %f, %f, %f\n", Sv[0][0], Sv[0][1], Sv[0][2], Sv[1][0], Sv[1][1], Sv[1][2], Sv[2][0], Sv[2][1], Sv[2][2]);
//   printf("Sv= %f, %f, %f\n    %f, %f, %f\n    %f, %f, %f\n", Sv[0][0], Sv[0][1], Sv[0][2], Sv[1][0], Sv[1][1], Sv[1][2], Sv[2][0], Sv[2][1], Sv[2][2]);
   gsl_eigen_nonsymm_params(0, 1, workspace);
   gsl_eigen_nonsymm(gsl_Sv, Sv_eigen, workspace);
   printf("Eigenvalues: %f, %f, %f\n",
   GSL_REAL(gsl_vector_complex_get(Sv_eigen, 0)),
   GSL_REAL(gsl_vector_complex_get(Sv_eigen, 1)),
   GSL_REAL(gsl_vector_complex_get(Sv_eigen, 2))
   );
   vtx->energy=0.0;
   gsl_eigen_symm(gsl_Sv, Sv_eigen, workspace);
//   printf("Eigenvalues: %f, %f, %f\n", gsl_vector_get(Sv_eigen, 0),gsl_vector_get(Sv_eigen, 1), gsl_vector_get(Sv_eigen, 2) );
//   printf("Eigenvalues: %f, %f, %f\n", gsl_matrix_get(evec, 0,0),gsl_matrix_get(evec, 0,1), gsl_matrix_get(evec, 0,2) );
   eigenval[0]= gsl_vector_get(Sv_eigen, 0);
   eigenval[1]= gsl_vector_get(Sv_eigen, 1);
   eigenval[2]= gsl_vector_get(Sv_eigen, 2);
   qsort(eigenval, 3, sizeof(ts_double), cmpfunc);
//   printf("Eigenvalues: %f, %f, %f\n", eigenval[0], eigenval[1], eigenval[2] );
   vtx->energy=(pow(eigenval[0]+eigenval[1],2))*Av;
   gsl_matrix_free(gsl_Sv);
   gsl_vector_complex_free(Sv_eigen);
   gsl_eigen_nonsymm_free(workspace);
   gsl_vector_free(Sv_eigen);
//   gsl_matrix_free(evec);
   gsl_eigen_symm_free(workspace);
   return TS_SUCCESS;
}