00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "TMatrixDSymEigen.h"
00027 #include "TMath.h"
00028
00029 ClassImp(TMatrixDSymEigen)
00030
00031
00032 TMatrixDSymEigen::TMatrixDSymEigen(const TMatrixDSym &a)
00033 {
00034
00035
00036 R__ASSERT(a.IsValid());
00037
00038 const Int_t nRows = a.GetNrows();
00039 const Int_t rowLwb = a.GetRowLwb();
00040
00041 fEigenValues.ResizeTo(rowLwb,rowLwb+nRows-1);
00042 fEigenVectors.ResizeTo(a);
00043
00044 fEigenVectors = a;
00045
00046 TVectorD offDiag;
00047 Double_t work[kWorkMax];
00048 if (nRows > kWorkMax) offDiag.ResizeTo(nRows);
00049 else offDiag.Use(nRows,work);
00050
00051
00052 MakeTridiagonal(fEigenVectors,fEigenValues,offDiag);
00053
00054
00055 MakeEigenVectors(fEigenVectors,fEigenValues,offDiag);
00056 }
00057
00058
00059 TMatrixDSymEigen::TMatrixDSymEigen(const TMatrixDSymEigen &another)
00060 {
00061
00062
00063 *this = another;
00064 }
00065
00066
00067 void TMatrixDSymEigen::MakeTridiagonal(TMatrixD &v,TVectorD &d,TVectorD &e)
00068 {
00069
00070
00071
00072
00073 Double_t *pV = v.GetMatrixArray();
00074 Double_t *pD = d.GetMatrixArray();
00075 Double_t *pE = e.GetMatrixArray();
00076
00077 const Int_t n = v.GetNrows();
00078
00079 Int_t i,j,k;
00080 Int_t off_n1 = (n-1)*n;
00081 for (j = 0; j < n; j++)
00082 pD[j] = pV[off_n1+j];
00083
00084
00085
00086 for (i = n-1; i > 0; i--) {
00087 const Int_t off_i1 = (i-1)*n;
00088 const Int_t off_i = i*n;
00089
00090
00091
00092 Double_t scale = 0.0;
00093 Double_t h = 0.0;
00094 for (k = 0; k < i; k++)
00095 scale = scale+TMath::Abs(pD[k]);
00096 if (scale == 0.0) {
00097 pE[i] = pD[i-1];
00098 for (j = 0; j < i; j++) {
00099 const Int_t off_j = j*n;
00100 pD[j] = pV[off_i1+j];
00101 pV[off_i+j] = 0.0;
00102 pV[off_j+i] = 0.0;
00103 }
00104 } else {
00105
00106
00107
00108 for (k = 0; k < i; k++) {
00109 pD[k] /= scale;
00110 h += pD[k]*pD[k];
00111 }
00112 Double_t f = pD[i-1];
00113 Double_t g = TMath::Sqrt(h);
00114 if (f > 0)
00115 g = -g;
00116 pE[i] = scale*g;
00117 h = h-f*g;
00118 pD[i-1] = f-g;
00119 for (j = 0; j < i; j++)
00120 pE[j] = 0.0;
00121
00122
00123
00124 for (j = 0; j < i; j++) {
00125 const Int_t off_j = j*n;
00126 f = pD[j];
00127 pV[off_j+i] = f;
00128 g = pE[j]+pV[off_j+j]*f;
00129 for (k = j+1; k <= i-1; k++) {
00130 const Int_t off_k = k*n;
00131 g += pV[off_k+j]*pD[k];
00132 pE[k] += pV[off_k+j]*f;
00133 }
00134 pE[j] = g;
00135 }
00136 f = 0.0;
00137 for (j = 0; j < i; j++) {
00138 pE[j] /= h;
00139 f += pE[j]*pD[j];
00140 }
00141 Double_t hh = f/(h+h);
00142 for (j = 0; j < i; j++)
00143 pE[j] -= hh*pD[j];
00144 for (j = 0; j < i; j++) {
00145 f = pD[j];
00146 g = pE[j];
00147 for (k = j; k <= i-1; k++) {
00148 const Int_t off_k = k*n;
00149 pV[off_k+j] -= (f*pE[k]+g*pD[k]);
00150 }
00151 pD[j] = pV[off_i1+j];
00152 pV[off_i+j] = 0.0;
00153 }
00154 }
00155 pD[i] = h;
00156 }
00157
00158
00159
00160 for (i = 0; i < n-1; i++) {
00161 const Int_t off_i = i*n;
00162 pV[off_n1+i] = pV[off_i+i];
00163 pV[off_i+i] = 1.0;
00164 Double_t h = pD[i+1];
00165 if (h != 0.0) {
00166 for (k = 0; k <= i; k++) {
00167 const Int_t off_k = k*n;
00168 pD[k] = pV[off_k+i+1]/h;
00169 }
00170 for (j = 0; j <= i; j++) {
00171 Double_t g = 0.0;
00172 for (k = 0; k <= i; k++) {
00173 const Int_t off_k = k*n;
00174 g += pV[off_k+i+1]*pV[off_k+j];
00175 }
00176 for (k = 0; k <= i; k++) {
00177 const Int_t off_k = k*n;
00178 pV[off_k+j] -= g*pD[k];
00179 }
00180 }
00181 }
00182 for (k = 0; k <= i; k++) {
00183 const Int_t off_k = k*n;
00184 pV[off_k+i+1] = 0.0;
00185 }
00186 }
00187 for (j = 0; j < n; j++) {
00188 pD[j] = pV[off_n1+j];
00189 pV[off_n1+j] = 0.0;
00190 }
00191 pV[off_n1+n-1] = 1.0;
00192 pE[0] = 0.0;
00193 }
00194
00195
00196 void TMatrixDSymEigen::MakeEigenVectors(TMatrixD &v,TVectorD &d,TVectorD &e)
00197 {
00198
00199
00200
00201
00202
00203 Double_t *pV = v.GetMatrixArray();
00204 Double_t *pD = d.GetMatrixArray();
00205 Double_t *pE = e.GetMatrixArray();
00206
00207 const Int_t n = v.GetNrows();
00208
00209 Int_t i,j,k,l;
00210 for (i = 1; i < n; i++)
00211 pE[i-1] = pE[i];
00212 pE[n-1] = 0.0;
00213
00214 Double_t f = 0.0;
00215 Double_t tst1 = 0.0;
00216 Double_t eps = TMath::Power(2.0,-52.0);
00217 for (l = 0; l < n; l++) {
00218
00219
00220
00221 tst1 = TMath::Max(tst1,TMath::Abs(pD[l])+TMath::Abs(pE[l]));
00222 Int_t m = l;
00223
00224
00225 while (m < n) {
00226 if (TMath::Abs(pE[m]) <= eps*tst1) {
00227 break;
00228 }
00229 m++;
00230 }
00231
00232
00233
00234
00235 if (m > l) {
00236 Int_t iter = 0;
00237 do {
00238 if (iter++ == 30) {
00239 Error("MakeEigenVectors","too many iterations");
00240 break;
00241 }
00242
00243
00244
00245 Double_t g = pD[l];
00246 Double_t p = (pD[l+1]-g)/(2.0*pE[l]);
00247 Double_t r = TMath::Hypot(p,1.0);
00248 if (p < 0)
00249 r = -r;
00250 pD[l] = pE[l]/(p+r);
00251 pD[l+1] = pE[l]*(p+r);
00252 Double_t dl1 = pD[l+1];
00253 Double_t h = g-pD[l];
00254 for (i = l+2; i < n; i++)
00255 pD[i] -= h;
00256 f = f+h;
00257
00258
00259
00260 p = pD[m];
00261 Double_t c = 1.0;
00262 Double_t c2 = c;
00263 Double_t c3 = c;
00264 Double_t el1 = pE[l+1];
00265 Double_t s = 0.0;
00266 Double_t s2 = 0.0;
00267 for (i = m-1; i >= l; i--) {
00268 c3 = c2;
00269 c2 = c;
00270 s2 = s;
00271 g = c*pE[i];
00272 h = c*p;
00273 r = TMath::Hypot(p,pE[i]);
00274 pE[i+1] = s*r;
00275 s = pE[i]/r;
00276 c = p/r;
00277 p = c*pD[i]-s*g;
00278 pD[i+1] = h+s*(c*g+s*pD[i]);
00279
00280
00281
00282 for (k = 0; k < n; k++) {
00283 const Int_t off_k = k*n;
00284 h = pV[off_k+i+1];
00285 pV[off_k+i+1] = s*pV[off_k+i]+c*h;
00286 pV[off_k+i] = c*pV[off_k+i]-s*h;
00287 }
00288 }
00289 p = -s*s2*c3*el1*pE[l]/dl1;
00290 pE[l] = s*p;
00291 pD[l] = c*p;
00292
00293
00294
00295 } while (TMath::Abs(pE[l]) > eps*tst1);
00296 }
00297 pD[l] = pD[l]+f;
00298 pE[l] = 0.0;
00299 }
00300
00301
00302
00303 for (i = 0; i < n-1; i++) {
00304 k = i;
00305 Double_t p = pD[i];
00306 for (j = i+1; j < n; j++) {
00307 if (pD[j] > p) {
00308 k = j;
00309 p = pD[j];
00310 }
00311 }
00312 if (k != i) {
00313 pD[k] = pD[i];
00314 pD[i] = p;
00315 for (j = 0; j < n; j++) {
00316 const Int_t off_j = j*n;
00317 p = pV[off_j+i];
00318 pV[off_j+i] = pV[off_j+k];
00319 pV[off_j+k] = p;
00320 }
00321 }
00322 }
00323 }
00324
00325
00326 TMatrixDSymEigen &TMatrixDSymEigen::operator=(const TMatrixDSymEigen &source)
00327 {
00328
00329
00330 if (this != &source) {
00331 fEigenVectors.ResizeTo(source.fEigenVectors);
00332 fEigenValues.ResizeTo(source.fEigenValues);
00333 }
00334 return *this;
00335 }