From bb207499ba467145aabb379b841094afd6d695c0 Mon Sep 17 00:00:00 2001 From: Samo Penic <samo.penic@gmail.com> Date: Sun, 15 May 2016 10:13:44 +0000 Subject: [PATCH] Fix in reading .status file. --- src/shcomplex.c | 43 ++++++++++++++++++++++++++++--------------- 1 files changed, 28 insertions(+), 15 deletions(-) diff --git a/src/shcomplex.c b/src/shcomplex.c index 9350fe5..5dc6a97 100644 --- a/src/shcomplex.c +++ b/src/shcomplex.c @@ -1,3 +1,4 @@ +/* vim: set ts=4 sts=4 sw=4 noet : */ #include<math.h> #include<stdlib.h> #include<gsl/gsl_complex.h> @@ -130,31 +131,43 @@ } -ts_bool calculateKc(ts_vesicle *vesicle){ - - ts_int i; - gsl_matrix *A=gsl_matrix_alloc(vesicle->sphHarmonics->l,2); +ts_double calculateKc(ts_vesicle *vesicle, ts_int lmin, ts_int lmax){ + ts_int min=lmin; + ts_int max=lmax; //vesicle->sphHarmonics->l-3; + ts_long i,j; + ts_double retval, bval; + gsl_matrix *A=gsl_matrix_alloc(max-min,2); gsl_vector *tau=gsl_vector_alloc(2); - gsl_vector *b=gsl_vector_alloc(vesicle->sphHarmonics->l); + gsl_vector *b=gsl_vector_alloc(max-min); gsl_vector *x=gsl_vector_alloc(2); - gsl_vector *res=gsl_vector_alloc(vesicle->sphHarmonics->l); + gsl_vector *res=gsl_vector_alloc(max-min); - //solving (A^T*A)*x=A^T*UlmSqAvg + //solving (A^T*A)*x=A^T*b //fill the data for matrix A and vector b - for(i=1;i<=vesicle->sphHarmonics->l;i++){ - gsl_matrix_set(A, i-1,0,(ts_double)((i-1)*(i+2))); - gsl_matrix_set(A, i-1,1,(ts_double)((i-1)*(i+2)*(i+1)*i)); - gsl_vector_set(b,i-1,(ts_double)vesicle->sphHarmonics->N/vesicle->sphHarmonics->sumUlm2[i-1][(i-1)*2]); + for(i=min;i<max;i++){ + gsl_matrix_set(A, i-min,0,(ts_double)((i-1)*(i+2))); + gsl_matrix_set(A, i-min,1,(ts_double)((i-1)*(i+2)*(i+1)*i)); +// fprintf(stderr,"%e %e\n", gsl_matrix_get(A,i-min,0), gsl_matrix_get(A,i-min,1)); + bval=0.0; + //average for m from 0..l (only positive m's) + for(j=0;j<=i;j++){ + bval+=vesicle->sphHarmonics->sumUlm2[i][(j+i)]; + } + bval=bval/(ts_double)vesicle->sphHarmonics->N/(ts_double)(i+1); + + gsl_vector_set(b,i-min,1.0/bval); +// fprintf(stderr,"%e\n", 1.0/gsl_vector_get(b,i-min)); } - fprintf(stderr,"b[2]=%e\n",gsl_vector_get(b,1)); +// fprintf(stderr,"b[2]=%e\n",gsl_vector_get(b,1)); gsl_linalg_QR_decomp(A,tau); gsl_linalg_QR_lssolve(A,tau,b,x,res); - fprintf(stderr,"kc=%e\n",gsl_vector_get(x,1)); - +// fprintf(stderr,"kc=%e\n",gsl_vector_get(x,1)); + retval=gsl_vector_get(x,1); gsl_matrix_free(A); gsl_vector_free(tau); gsl_vector_free(b); gsl_vector_free(x); gsl_vector_free(res); - return TS_SUCCESS; + + return retval; } -- Gitblit v1.9.3