Changeset 26892 for trunk/psLib/src/math/psMatrix.c
- Timestamp:
- Feb 10, 2010, 7:27:29 PM (16 years ago)
- File:
-
- 1 edited
-
trunk/psLib/src/math/psMatrix.c (modified) (4 diffs)
Legend:
- Unmodified
- Added
- Removed
-
trunk/psLib/src/math/psMatrix.c
r26001 r26892 1033 1033 } 1034 1034 1035 psVector *psMatrixSolveSVD(psVector *out, const psImage *matrix, const psVector *vector )1035 psVector *psMatrixSolveSVD(psVector *out, const psImage *matrix, const psVector *vector, float thresh) 1036 1036 { 1037 1037 #define psMatrixSolveSVD_EXIT {psFree(out); return NULL; } … … 1063 1063 gsl_vector_free(work); 1064 1064 1065 if (isfinite(thresh) && thresh > 0.0) { 1066 // Trim the singular values 1067 double total = 0.0; // Total of singular values 1068 for (int i = 0; i < numCols; i++) { 1069 total += gsl_vector_get(S, i); 1070 } 1071 thresh *= total; 1072 for (int i = 0; i < numCols; i++) { 1073 double value = gsl_vector_get(S, i); // Singular value 1074 if (value < thresh) { 1075 psTrace("psLib.math", 5, "Trimming singular value %d: %lg", i, value); 1076 gsl_vector_set(S, i, 0.0); 1077 #if 0 1078 for (int j = 0; j < numCols; j++) { 1079 // Being thorough; probably unnecessary 1080 gsl_matrix_set(V, j, i, 0.0); 1081 gsl_matrix_set(A, j, i, 0.0); 1082 } 1083 #endif 1084 } else { 1085 psTrace("psLib.math", 5, "Singular value %d: %lg", i, value); 1086 } 1087 } 1088 } 1089 1065 1090 // Solve system (or minimise least-squares if overconstrained): Ax = b 1066 1091 gsl_vector *b = gsl_vector_alloc(numCols); // Vector b … … 1093 1118 } 1094 1119 1095 1096 1097 1120 // This code supplied by Andy Becker (becker@astro.washington.edu) 1098 psImage *psMatrixSVD (psImage* evec, psVector* eval, const psImage* in)1121 psImage *psMatrixSVD_old(psImage* evec, psVector* eval, const psImage* in) 1099 1122 { 1100 1123 #define psMatrixSVD_EXIT {psFree(evec); psFree(eval); return NULL;} … … 1145 1168 return evec; 1146 1169 } 1170 1171 // this is basically a wrapper for the gsl function: gsl_linalg_SV_decomp() SVD decomposes 1172 // matrix A based on the following equation: A = U w V^T . This function (as usual for SVD 1173 // implementations) returns V not V^T. U and V are returned to images; w is returned to a 1174 // vector representing the diagonal of w. The input image A is not modified. U, V, and w may 1175 // be supplied as NULL or may be allocated; their lengths are set here to match the 1176 // dimensionality of A. XXX there is no error handling for the gsl functions (anywhere in 1177 // psMatrix.c) 1178 bool psMatrixSVD(psImage **U, psVector **w, psImage **V, const psImage *A) 1179 { 1180 // Error checks Missing one for eval 1181 PS_ASSERT_PTR_NON_NULL(U, false); 1182 PS_ASSERT_PTR_NON_NULL(w, false); 1183 PS_ASSERT_PTR_NON_NULL(V, false); 1184 PS_ASSERT_PTR_NON_NULL(A, false); 1185 1186 // A is provided with size Nx,Ny = numCols,numRows 1187 // U has size Nx,Ny 1188 // V has size Nx,Nx 1189 // w has size Nx 1190 1191 // Initialize data 1192 int numRows = A->numRows; 1193 int numCols = A->numCols; 1194 1195 *U = psImageRecycle(*U, numCols, numRows, A->type.type); 1196 *V = psImageRecycle(*V, numCols, numCols, A->type.type); 1197 *w = psVectorRecycle(*w, numCols, A->type.type); 1198 1199 gsl_matrix *Agsl = gsl_matrix_alloc(numRows, numCols); 1200 gsl_matrix *Vgsl = gsl_matrix_alloc(numCols, numCols); 1201 gsl_vector *Sgsl = gsl_vector_alloc(numCols); 1202 gsl_vector *work = gsl_vector_alloc(numCols); 1203 1204 // Copy psImage data into GSL matrix data 1205 matrixPStoGSL(Agsl, A); 1206 1207 // Calculate SVD decomposition 1208 gsl_linalg_SV_decomp(Agsl, Vgsl, Sgsl, work); 1209 1210 // Copy GSL matrix data to psImage data 1211 matrixGSLtoPS(*V, Vgsl); 1212 matrixGSLtoPS(*U, Agsl); // gsl_linalg_SV_decomp replaces A with U 1213 vectorGSLtoPS(*w, Sgsl); 1214 1215 // Free GSL data 1216 gsl_matrix_free(Agsl); 1217 gsl_matrix_free(Vgsl); 1218 gsl_vector_free(Sgsl); 1219 gsl_vector_free(work); 1220 1221 return true; 1222 } 1223
Note:
See TracChangeset
for help on using the changeset viewer.
