00001 #ifndef __GSL_MATRIX_CC__
00002 #define __GSL_MATRIX_CC__
00003 #include "gslmatrix.h"
00004
00005
00006 CGSLMatrix_Real::CGSLMatrix_Real(int dimset){
00007 dim=dimset;
00008 w=NULL;
00009 eval=NULL;
00010 evec=NULL;
00011 g=NULL;
00012 p=NULL;
00013 U=NULL;
00014 m=NULL;
00015 v=NULL;
00016 }
00017
00018 CGSLMatrix_Real::~CGSLMatrix_Real(){
00019 if(w!=NULL) gsl_eigen_symmv_free(w);
00020 if(eval!=NULL) gsl_vector_free(eval);
00021 if(evec!=NULL) gsl_matrix_free(evec);
00022 if(g!=NULL) gsl_vector_free(g);
00023 if(p!=NULL) gsl_permutation_free(p);
00024 if(m!=NULL) gsl_matrix_free(m);
00025 if(v!=NULL) gsl_vector_free(v);
00026 if(U!=NULL){
00027 for(int i=0;i<dim;i++) delete [] U[i];
00028 delete [] U;
00029 }
00030 }
00031
00032 void CGSLMatrix_Real::SolveLinearEqs(double *y,double **A,double *x){
00033 int i,j,s;
00034 if(v==NULL) v=gsl_vector_alloc(dim);
00035
00036 if(m==NULL) m=gsl_matrix_alloc(dim,dim);
00037
00038 for(i=0;i<dim;i++){
00039 gsl_vector_set(v,i,y[i]);
00040 for(j=0;j<dim;j++) gsl_matrix_set(m,i,j,A[i][j]);
00041 }
00042
00043 if(g==NULL) g = gsl_vector_alloc (dim);
00044 if(p==NULL) p = gsl_permutation_alloc (dim);
00045
00046 gsl_linalg_LU_decomp (m,p,&s);
00047 gsl_linalg_LU_solve (m, p,v,g);
00048
00049 for(i=0;i<dim;i++) x[i]=gsl_vector_get(g,i);
00050 }
00051
00052 void CGSLMatrix_Real::EigenFind(double **A,double **eigenvec,
00053 double *eigenval){
00054 int i,j;
00055 if(m==NULL) m=gsl_matrix_alloc(dim,dim);
00056 for(i=0;i<dim;i++){
00057 for(j=0;j<dim;j++) gsl_matrix_set(m,i,j,A[i][j]);
00058 }
00059
00060 if(eval==NULL) eval=gsl_vector_alloc(dim);
00061 if(evec==NULL) evec=gsl_matrix_alloc(dim,dim);
00062 if(w==NULL) w=gsl_eigen_symmv_alloc(dim);
00063
00064 gsl_eigen_symmv(m,eval,evec,w);
00065
00066
00067 for(i=0;i<dim;i++){
00068 eigenval[i]=gsl_vector_get(eval,i);
00069 for(j=0;j<dim;j++) eigenvec[i][j]=gsl_matrix_get(evec,i,j);
00070 }
00071
00072 }
00073
00074 void CGSLMatrix_Real::Invert(double **A,double **Ainv){
00075
00076 int i,j;
00077 if(m==NULL) m=gsl_matrix_alloc(dim,dim);
00078 for(i=0;i<dim;i++){
00079 for(j=0;j<dim;j++) gsl_matrix_set(m,i,j,A[i][j]);
00080 }
00081
00082 if(eval==NULL) eval=gsl_vector_alloc(dim);
00083 if(evec==NULL) evec=gsl_matrix_alloc(dim,dim);
00084 if(w==NULL) w=gsl_eigen_symmv_alloc(dim);
00085
00086 gsl_eigen_symmv(m,eval,evec,w);
00087
00088
00089 if(U==NULL){
00090 U=new double*[dim];
00091 for(i=0;i<dim;i++) U[i]=new double[dim];
00092 }
00093
00094 for(i=0;i<dim;i++){
00095 for(j=0;j<dim;j++) {
00096 U[i][j]=gsl_matrix_get(evec,j,i);
00097 Ainv[i][j]=0.0;
00098 }
00099 }
00100 int k;
00101 for(i=0;i<dim;i++){
00102 for(j=0;j<dim;j++){
00103 for(k=0;k<dim;k++) Ainv[i][j]+=U[k][i]*U[k][j]/gsl_vector_get(eval,k);
00104 }
00105 }
00106 }
00107
00108
00109
00110
00111 CGSLMatrix_Complex::CGSLMatrix_Complex(int dimset){
00112 dim=dimset;
00113 eval=NULL;
00114 evec=NULL;
00115 w=NULL;
00116 g=NULL;
00117 p=NULL;
00118 U=NULL;
00119 m=NULL;
00120 v=NULL;
00121 }
00122
00123 CGSLMatrix_Complex::~CGSLMatrix_Complex(){
00124 if(eval!=NULL) gsl_vector_free(eval);
00125 if(evec!=NULL) gsl_matrix_complex_free(evec);
00126 if(w!=NULL) gsl_eigen_hermv_free(w);
00127 if(g!=NULL) gsl_vector_complex_free(g);
00128 if(p!=NULL) gsl_permutation_free(p);
00129 if(m!=NULL) gsl_matrix_complex_free(m);
00130 if(v!=NULL) gsl_vector_complex_free(v);
00131 if(U!=NULL){
00132 for(int i=0;i<dim;i++) delete [] U[i];
00133 delete [] U;
00134 }
00135 }
00136
00137 void CGSLMatrix_Complex::EigenFind(complex<double> **A,
00138 complex<double> **eigenvec,double *eigenval){
00139 complex<double> ci(0.0,1.0);
00140 gsl_complex z;
00141
00142 if(m==NULL) m=gsl_matrix_complex_alloc(dim,dim);
00143 int i,j;
00144 for(i=0;i<dim;i++){
00145 for(j=0;j<dim;j++){
00146 GSL_SET_COMPLEX(&z,real(A[i][j]),imag(A[i][j]));
00147 gsl_matrix_complex_set(m,i,j,z);
00148 }
00149 }
00150
00151 if(eval==NULL) eval=gsl_vector_alloc(dim);
00152 if(evec==NULL) evec=gsl_matrix_complex_alloc(dim,dim);
00153 if(w==NULL) w=gsl_eigen_hermv_alloc(dim);
00154
00155 gsl_eigen_hermv(m,eval,evec,w);
00156
00157
00158
00159 for(i=0;i<dim;i++){
00160 eigenval[i]=gsl_vector_get(eval,i);
00161 for(j=0;j<dim;j++){
00162 z=gsl_matrix_complex_get(evec,i,j);
00163 eigenvec[i][j]=GSL_REAL(z)+ci*GSL_IMAG(z);
00164 }
00165 }
00166
00167 }
00168
00169 void CGSLMatrix_Complex::Invert(complex<double> **A,complex<double> **Ainv){
00170 complex<double> ci(0.0,1.0);
00171 gsl_complex z;
00172
00173 if(m==NULL) m=gsl_matrix_complex_alloc(dim,dim);
00174 int i,j;
00175 for(i=0;i<dim;i++){
00176 for(j=0;j<dim;j++){
00177 GSL_SET_COMPLEX(&z,real(A[i][j]),imag(A[i][j]));
00178 gsl_matrix_complex_set(m,i,j,z);
00179 }
00180 }
00181
00182 if(eval==NULL) eval=gsl_vector_alloc(dim);
00183 if(evec==NULL) evec=gsl_matrix_complex_alloc(dim,dim);
00184 if(w==NULL) w=gsl_eigen_hermv_alloc(dim);
00185
00186 gsl_eigen_hermv(m,eval,evec,w);
00187
00188
00189 if(U==NULL){
00190 U=new complex<double> *[dim];
00191 for(i=0;i<dim;i++) U[i]=new complex<double>[dim];
00192 }
00193
00194 for(i=0;i<dim;i++){
00195 for(j=0;j<dim;j++) {
00196 z=gsl_matrix_complex_get(evec,j,i);
00197 U[i][j]=GSL_REAL(z)+ci*GSL_IMAG(z);
00198 Ainv[i][j]=0.0;
00199 }
00200 }
00201 int k;
00202 for(i=0;i<dim;i++){
00203 for(j=0;j<dim;j++){
00204 for(k=0;k<dim;k++)
00205 Ainv[i][j]+=U[k][i]*conj(U[k][j])/gsl_vector_get(eval,k);
00206 }
00207 }
00208
00209 }
00210
00211 void CGSLMatrix_Complex::SolveLinearEqs(complex<double> *y,
00212 complex<double> **A,complex<double> *x){
00213 complex<double> ci(0.0,1.0);
00214 int i,j,s;
00215 gsl_complex z;
00216
00217 if(m==NULL) m=gsl_matrix_complex_alloc(dim,dim);
00218 if(v==NULL) v=gsl_vector_complex_alloc(dim);
00219 for(i=0;i<dim;i++){
00220 GSL_SET_COMPLEX(&z,real(y[i]),imag(y[i]));
00221 gsl_vector_complex_set(v,i,z);
00222 for(j=0;j<dim;j++){
00223 GSL_SET_COMPLEX(&z,real(A[i][j]),imag(A[i][j]));
00224 gsl_matrix_complex_set(m,i,j,z);
00225 }
00226 }
00227
00228 if(g==NULL) g = gsl_vector_complex_alloc (dim);
00229 if(p==NULL) p = gsl_permutation_alloc (dim);
00230
00231
00232 gsl_linalg_complex_LU_decomp (m, p, &s);
00233 gsl_linalg_complex_LU_solve (m, p, v, g);
00234
00235 for(i=0;i<dim;i++){
00236 z=gsl_vector_complex_get(g,i);
00237 x[i]=GSL_REAL(z)+ci*GSL_IMAG(z);
00238 }
00239 }
00240
00241 #endif