00001
00002
00003
00004
00005
00006
00007
00008
00016 #ifndef __F_MATRIX_H__
00017 #define __F_MATRIX_H__
00018
00019 #include <dims_common.h>
00020 #include <dims_clapack.h>
00021 #if defined(__ia64) || defined(__x68_64)
00022 #include <cmath>
00023 #else
00024 #include <math.h>
00025 #endif
00026
00027 #include <iostream>
00028 using namespace std;
00029
00030 #include <assert.h>
00031
00032 class fVec;
00033 class fMat;
00034
00035 fMat inv_svd(const fMat& mat, int lwork = -1);
00036 fMat p_inv_svd(const fMat& mat, int lwork = -1);
00037 fMat sr_inv_svd(const fMat& mat, fVec& w_err, fVec& w_norm, double k, int lwork = -1);
00038 fMat lineq_svd(const fMat& A, const fMat& b, int lwork = -1);
00039
00040
00046 class fMat
00047 {
00048 public:
00050
00053 fMat() {
00054 p_data = 0;
00055 n_col = 0;
00056 n_row = 0;
00057 m_info = 0;
00058 temp = 0.0;
00059 }
00061
00067 fMat(int m, int n, double* ini = 0) {
00068 p_data = 0;
00069 n_col = 0;
00070 n_row = 0;
00071 m_info = 0;
00072 temp = 0.0;
00073 resize(m, n);
00074 if(ini)
00075 {
00076 int k = m * n;
00077 for(int i=0; i<k; i++) p_data[i] = ini[i];
00078 }
00079 }
00081
00087 fMat(int m, int n, double ini) {
00088 p_data = 0;
00089 n_col = 0;
00090 n_row = 0;
00091 m_info = 0;
00092 temp = 0.0;
00093 resize(m, n);
00094 int k = m * n;
00095 for(int i=0; i<k; i++) p_data[i] = ini;
00096 }
00098
00102 fMat(const fMat& ini) {
00103 p_data = 0;
00104 n_col = 0;
00105 n_row = 0;
00106 m_info = ini.m_info;
00107 temp = ini.temp;
00108 resize(ini.row(), ini.col());
00109 int k = ini.row() * ini.col();
00110 for(int i=0; i<k; i++) p_data[i] = ini.data()[i];
00111 }
00113
00117 fMat(const fVec& ini);
00118
00120 ~fMat() {
00121 if(p_data) delete[] p_data;
00122 }
00123
00125
00133 int resize(int i, int j) {
00134 if(n_row == i && n_col == j) return 0;
00135 if(p_data) delete[] p_data;
00136 p_data = 0;
00137 n_row = i;
00138 n_col = j;
00139 if(n_row > 0 && n_col > 0)
00140 {
00141 int nn = n_row * n_col;
00142 p_data = new double[nn];
00143 if(!p_data)
00144 {
00145 n_row = 0;
00146 n_col = 0;
00147 m_info = 1;
00148 return -1;
00149 }
00150 }
00151 return 0;
00152 }
00154 int col() const {
00155 return n_col;
00156 }
00158 int row() const {
00159 return n_row;
00160 }
00162 int info() {
00163 return m_info;
00164 }
00166 void setinfo(int _info) {
00167 m_info = _info;
00168 }
00169
00171 double& operator() (int i, int j) {
00172 assert(i>=0 && i<n_row && j>=0 && j<n_col);
00173 return *(p_data + i + j*n_row);
00174 }
00175
00177 double operator() (int i, int j) const {
00178 assert(i>=0 && i<n_row && j>=0 && j<n_col);
00179 return *(p_data + i + j*n_row);
00180 }
00181
00183 fMat operator = (const fMat& mat);
00185 void operator = (double d);
00187 void set(double* _d) {
00188 int n_elements = n_row * n_col;
00189 dims_copy(_d, p_data, n_elements);
00190 }
00191
00193 double* data() const {
00194 return p_data;
00195 }
00197 operator double *() {
00198 return p_data;
00199 }
00200
00202 friend ostream& operator << (ostream& ost, const fMat& mat) {
00203 int i, j;
00204 ost << "[" << mat.row() << ", " << mat.col() << "]" << endl;
00205 for(i=0; i<mat.row(); i++)
00206 {
00207 for(j=0; j<mat.col(); j++)
00208 {
00209 ost << *(mat.data() + i + j*mat.row());
00210 if(j != mat.col()-1) ost << "\t";
00211 }
00212 if(i == mat.row()-1) ost << flush;
00213 else ost << endl;
00214 }
00215 return ost;
00216 }
00217
00219
00225 void set_submat(int row_start, int col_start, const fMat& subM);
00226
00228
00234 void get_submat(int row_start, int col_start, const fMat& allM);
00235
00237 friend fMat tran(const fMat& mat);
00238
00240 void tran(const fMat& mat);
00241
00243 void tran();
00244
00246 void identity();
00247
00249 void zero() {
00250 int ndata = n_row * n_col;
00251 for(int i=0; i<ndata; i++) p_data[i] = 0.0;
00252 }
00253
00264
00265 friend fMat inv(const fMat& mat);
00267 friend fMat p_inv(const fMat& mat);
00269 friend fMat sr_inv(const fMat& mat, fVec& w_err, fVec& w_norm, double k);
00271 friend fMat lineq(const fMat& A, const fMat& b);
00273 friend fMat lineq_sr(const fMat& A, fVec& w_err, fVec& w_norm, double k, const fMat& b);
00280
00281 int inv(const fMat&);
00283 int p_inv(const fMat&);
00285 int sr_inv(const fMat& mat, fVec& w_err, fVec& w_norm, double k);
00287 int lineq(const fMat& A, const fMat& b);
00289 int lineq_sr(const fMat& A, fVec& w_err, fVec& w_norm, double k, const fMat& b);
00292
00303
00304 friend fMat inv_svd(const fMat& mat, int lwork);
00306 friend fMat p_inv_svd(const fMat& mat, int lwork);
00308 friend fMat sr_inv_svd(const fMat& mat, fVec& w_err, fVec& w_norm, double k, int lwork);
00310 friend fMat lineq_svd(const fMat& A, const fMat& b, int lwork);
00317
00318 int inv_svd(const fMat&, int lwork = -1);
00320 int p_inv_svd(const fMat&, int lwork = -1);
00322 int sr_inv_svd(const fMat& mat, fVec& w_err, fVec& w_norm, double k, int lwork = -1);
00324 int lineq_svd(const fMat& A, const fMat& b, int lwork = -1);
00327
00329 int inv_posv(const fMat&);
00331 int lineq_posv(const fMat& A, const fMat& b);
00333 int inv_porfs(const fMat&);
00335 int lineq_porfs(const fMat& A, const fMat& b);
00336
00342 friend fMat operator - (const fMat& mat);
00343 void operator += (const fMat& mat);
00344 void operator -= (const fMat& mat);
00345 void operator *= (double d);
00346 void operator /= (double d);
00347
00348 friend fMat operator + (const fMat& mat1, const fMat& mat2);
00349 friend fMat operator - (const fMat& mat1, const fMat& mat2);
00350 friend fMat operator * (const fMat& mat1, const fMat& mat2);
00351 friend fVec operator * (const fMat& mat, const fVec& vec);
00352 friend fMat operator * (double d, const fMat& mat);
00353 friend fMat operator * (const fMat& mat, double d);
00354 friend fMat operator / (const fMat& mat, double d);
00362 void set(const fMat& mat);
00363 void neg(const fMat& mat);
00364 void add(const fMat& mat1, const fMat& mat2);
00365 void add(const fMat& mat);
00366 void sub(const fMat& mat1, const fMat& mat2);
00367 void mul(const fMat& mat1, const fMat& mat2);
00368 void mul(double d, const fMat& mat);
00369 void mul(const fMat& mat, double d);
00370 void mul_tran(const fMat& mat1, int trans_first);
00371 void div(const fMat& mat, double d);
00374
00375 void dmul(const fMat& mat1, const fMat& mat2);
00377 void ddiv(const fMat& mat1, const fMat& mat2);
00378
00380
00386 int svd(fMat& U, fVec& Sigma, fMat& VT);
00387
00389
00393 void symmetric(char t = 'U');
00394
00396 double det(void);
00397
00399 friend fMat eigs(const fMat& mat, double *w);
00401 friend fMat eigs2(const fMat& mat, double *w);
00402
00404
00409 int eig(fVec& wr, fVec& wi);
00410
00412
00418 int eig(fVec& wr, fVec& wi, fMat& v);
00419
00421
00428 int eig(fVec& wr, fVec& wi, fMat& vr, fMat& vi);
00429
00431
00438 friend int inv_enlarge(const fMat& m12, const fMat& m21, const fMat& m22, const fMat& P, fMat& X, fMat& y, fMat& z, fMat& w);
00439
00441
00449 friend int inv_shrink(const fMat& P, const fMat& q, const fMat& r, const fMat& s, fMat& X);
00450
00452
00463 friend int inv_row_replaced(const fMat& P, const fMat& q, const fMat& m2d, fMat& X, fMat& y);
00464
00466
00477 friend int inv_col_replaced(const fMat& P, const fMat& q, const fMat& m2d, fMat& X, fMat& y);
00478
00479 protected:
00480 double* p_data;
00481 double temp;
00482 int m_info;
00483 int n_col;
00484 int n_row;
00485 };
00486
00491 class fVec
00492 : public fMat
00493 {
00494 public:
00495 fVec() : fMat() {
00496 }
00497 fVec(int m, double* ini = 0) : fMat(m, 1, ini) {
00498 }
00499 fVec(int m, double ini) : fMat(m, 1, ini) {
00500 }
00501 fVec(const fVec& ini) : fMat(ini) {
00502 }
00503 ~fVec() {
00504 }
00505
00507 int size() const {
00508 return n_row;
00509 }
00511 void resize(int i) {
00512 if(n_row == i) return;
00513 if(p_data) delete[] p_data;
00514 p_data = 0;
00515 n_row = i;
00516 n_col = 1;
00517 if(n_row > 0)
00518 {
00519 p_data = new double[n_row];
00520 }
00521 }
00522
00524 void get_subvec(int start, const fVec& allV);
00526 void set_subvec(int start, const fVec& subV);
00527
00529 void operator = (double d);
00530 fVec operator = (const fVec& vec);
00531
00533 void set(double* _d) {
00534 dims_copy(_d, p_data, n_row);
00535 }
00536
00538 double& operator () (int i) {
00539 assert(i >= 0 && i < n_row);
00540 return *(p_data + i);
00541 }
00543 double operator () (int i) const {
00544 assert(i >= 0 && i < n_row);
00545 return *(p_data + i);
00546 }
00547
00549 friend double length(const fVec& v) {
00550 return sqrt(v*v);
00551 }
00553 double length() const {
00554 return sqrt((*this) * (*this));
00555 }
00557 friend fVec unit(const fVec& v) {
00558 fVec ret(v.n_col);
00559 double len = v.length();
00560 ret = v / len;
00561 return ret;
00562 }
00564 void unit() {
00565 double len = length();
00566 (*this) /= len;
00567 }
00569 double sum() {
00570 double ret = 0.0;
00571 for(int i=0; i<n_row; i++) ret += p_data[i];
00572 return ret;
00573 }
00575 void zero() {
00576 for(int i=0; i<n_row; i++) p_data[i] = 0;
00577 }
00578
00580 int max_index();
00582 double max_value();
00584 int min_index();
00586 double min_value();
00587
00589 void abs();
00590
00592 int lineq(const fMat& A, const fVec& b);
00593 int lineq_svd(const fMat& A, const fVec& b, int lwork = -1);
00594 int lineq_posv(const fMat& A, const fVec& b);
00595 int lineq_porfs(const fMat& A, const fVec& b);
00596 int lineq_sr(const fMat& A, fVec& w_err, fVec& w_norm, double k, const fVec& b);
00597
00603 friend fVec operator - (const fVec& vec);
00604 void operator += (const fVec& vec);
00605 void operator -= (const fVec& vec);
00606 void operator *= (double d);
00607 void operator /= (double d);
00608
00609 friend fVec operator + (const fVec& vec1, const fVec& vec2);
00610 friend fVec operator - (const fVec& vec1, const fVec& vec2);
00611 friend fVec operator / (const fVec& vec, double d);
00612 friend fVec operator * (double d, const fVec& vec);
00613 friend fVec operator * (const fVec& vec, double d);
00614 friend double operator * (const fVec& vec1, const fVec& vec2);
00622 void set(const fVec& vec);
00623 void neg(const fVec& vec);
00624 void add(const fVec& vec1, const fVec& vec2);
00625 void add(const fVec& vec);
00626 void sub(const fVec& vec1, const fVec& vec2);
00627 void div(const fVec& vec, double d);
00628 void mul(const fVec& vec, double d);
00629 void mul(double d, const fVec& vec);
00631 void mul(const fMat& mat, const fVec& vec);
00633 void mul(const fVec& vec, const fMat& mat);
00635 protected:
00636 };
00637
00638 #endif