ROOT logo
// @(#)root/matrix:$Id: TMatrixDEigen.cxx 23492 2008-04-23 20:42:49Z brun $
// Authors: Fons Rademakers, Eddy Offermann  Dec 2003

/*************************************************************************
 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
 * All rights reserved.                                                  *
 *                                                                       *
 * For the licensing terms see $ROOTSYS/LICENSE.                         *
 * For the list of contributors see $ROOTSYS/README/CREDITS.             *
 *************************************************************************/

//////////////////////////////////////////////////////////////////////////
//                                                                      //
// TMatrixDEigen                                                        //
//                                                                      //
// Eigenvalues and eigenvectors of a real matrix.                       //
//                                                                      //
// If A is not symmetric, then the eigenvalue matrix D is block         //
// diagonal with the real eigenvalues in 1-by-1 blocks and any complex  //
// eigenvalues, a + i*b, in 2-by-2 blocks, [a, b; -b, a].  That is, if  //
// the complex eigenvalues look like                                    //
//                                                                      //
//     u + iv     .        .          .      .    .                     //
//       .      u - iv     .          .      .    .                     //
//       .        .      a + ib       .      .    .                     //
//       .        .        .        a - ib   .    .                     //
//       .        .        .          .      x    .                     //
//       .        .        .          .      .    y                     //
//                                                                      //
// then D looks like                                                    //
//                                                                      //
//       u        v        .          .      .    .                     //
//      -v        u        .          .      .    .                     //
//       .        .        a          b      .    .                     //
//       .        .       -b          a      .    .                     //
//       .        .        .          .      x    .                     //
//       .        .        .          .      .    y                     //
//                                                                      //
// This keeps V a real matrix in both symmetric and non-symmetric       //
// cases, and A*V = V*D.                                                //
//                                                                      //
//////////////////////////////////////////////////////////////////////////

#include "TMatrixDEigen.h"
#include "TMath.h"

ClassImp(TMatrixDEigen)

//______________________________________________________________________________
TMatrixDEigen::TMatrixDEigen(const TMatrixD &a)
{
// Constructor for eigen-problem of matrix A .

   R__ASSERT(a.IsValid());
 
   const Int_t nRows  = a.GetNrows();
   const Int_t nCols  = a.GetNcols();
   const Int_t rowLwb = a.GetRowLwb();
   const Int_t colLwb = a.GetColLwb();

   if (nRows != nCols || rowLwb != colLwb)
   {
      Error("TMatrixDEigen(TMatrixD &)","matrix should be square");
      return;
   }

   const Int_t rowUpb = rowLwb+nRows-1;
   fEigenVectors.ResizeTo(rowLwb,rowUpb,rowLwb,rowUpb);
   fEigenValuesRe.ResizeTo(rowLwb,rowUpb);
   fEigenValuesIm.ResizeTo(rowLwb,rowUpb);

   TVectorD ortho;
   Double_t work[kWorkMax];
   if (nRows > kWorkMax) ortho.ResizeTo(nRows);
   else                  ortho.Use(nRows,work);

   TMatrixD mH = a;

   // Reduce to Hessenberg form.
   MakeHessenBerg(fEigenVectors,ortho,mH);

   // Reduce Hessenberg to real Schur form.
   MakeSchurr(fEigenVectors,fEigenValuesRe,fEigenValuesIm,mH);

   // Sort eigenvalues and corresponding vectors in descending order of Re^2+Im^2
   // of the complex eigenvalues .
   Sort(fEigenVectors,fEigenValuesRe,fEigenValuesIm);
}

//______________________________________________________________________________
TMatrixDEigen::TMatrixDEigen(const TMatrixDEigen &another)
{
// Copy constructor

   *this = another;
}

//______________________________________________________________________________
void TMatrixDEigen::MakeHessenBerg(TMatrixD &v,TVectorD &ortho,TMatrixD &H)
{
// Nonsymmetric reduction to Hessenberg form.
// This is derived from the Algol procedures orthes and ortran, by Martin and Wilkinson,
// Handbook for Auto. Comp., Vol.ii-Linear Algebra, and the corresponding
// Fortran subroutines in EISPACK.

   Double_t *pV = v.GetMatrixArray();
   Double_t *pO = ortho.GetMatrixArray();
   Double_t *pH = H.GetMatrixArray();

   const Int_t n = v.GetNrows();

   const Int_t low  = 0;
   const Int_t high = n-1;

   Int_t i,j,m;
   for (m = low+1; m <= high-1; m++) {
      const Int_t off_m = m*n;

      // Scale column.

      Double_t scale = 0.0;
      for (i = m; i <= high; i++) {
         const Int_t off_i = i*n;
         scale = scale + TMath::Abs(pH[off_i+m-1]);
      }
      if (scale != 0.0) {

         // Compute Householder transformation.

         Double_t h = 0.0;
         for (i = high; i >= m; i--) {
            const Int_t off_i = i*n;
            pO[i] = pH[off_i+m-1]/scale;
            h += pO[i]*pO[i];
         }
         Double_t g = TMath::Sqrt(h);
         if (pO[m] > 0)
            g = -g;
         h = h-pO[m]*g;
         pO[m] = pO[m]-g;

         // Apply Householder similarity transformation
         // H = (I-u*u'/h)*H*(I-u*u')/h)

         for (j = m; j < n; j++) {
            Double_t f = 0.0;
            for (i = high; i >= m; i--) {
               const Int_t off_i = i*n;
               f += pO[i]*pH[off_i+j];
            }
            f = f/h;
            for (i = m; i <= high; i++) {
               const Int_t off_i = i*n;
               pH[off_i+j] -= f*pO[i];
            }
         }

         for (i = 0; i <= high; i++) {
            const Int_t off_i = i*n;
            Double_t f = 0.0;
            for (j = high; j >= m; j--)
               f += pO[j]*pH[off_i+j];
            f = f/h;
            for (j = m; j <= high; j++)
               pH[off_i+j] -= f*pO[j];
         }
         pO[m] = scale*pO[m];
         pH[off_m+m-1] = scale*g;
      }
   }

   // Accumulate transformations (Algol's ortran).

   for (i = 0; i < n; i++) {
      const Int_t off_i = i*n;
      for (j = 0; j < n; j++)
         pV[off_i+j] = (i == j ? 1.0 : 0.0);
   }

   for (m = high-1; m >= low+1; m--) {
      const Int_t off_m = m*n;
      if (pH[off_m+m-1] != 0.0) {
         for (i = m+1; i <= high; i++) {
            const Int_t off_i = i*n;
            pO[i] = pH[off_i+m-1];
         }
         for (j = m; j <= high; j++) {
            Double_t g = 0.0;
            for (i = m; i <= high; i++) {
               const Int_t off_i = i*n;
               g += pO[i]*pV[off_i+j];
            }
            // Double division avoids possible underflow
            g = (g/pO[m])/pH[off_m+m-1];
            for (i = m; i <= high; i++) {
               const Int_t off_i = i*n;
               pV[off_i+j] += g*pO[i];
            }
         }
      }
   }
}

//______________________________________________________________________________
static Double_t gCdivr, gCdivi;
static void cdiv(Double_t xr,Double_t xi,Double_t yr,Double_t yi) {
// Complex scalar division.
   Double_t r,d;
   if (TMath::Abs(yr) > TMath::Abs(yi)) {
      r = yi/yr;
      d = yr+r*yi;
      gCdivr = (xr+r*xi)/d;
      gCdivi = (xi-r*xr)/d;
   } else {
      r = yr/yi;
      d = yi+r*yr;
      gCdivr = (r*xr+xi)/d;
      gCdivi = (r*xi-xr)/d;
   }
}

//______________________________________________________________________________
void TMatrixDEigen::MakeSchurr(TMatrixD &v,TVectorD &d,TVectorD &e,TMatrixD &H)
{
// Nonsymmetric reduction from Hessenberg to real Schur form.
// This is derived from the Algol procedure hqr2, by Martin and Wilkinson,
// Handbook for Auto. Comp., Vol.ii-Linear Algebra, and the corresponding
// Fortran subroutine in EISPACK.

   // Initialize

   const Int_t nn = v.GetNrows();
         Int_t n = nn-1;
   const Int_t low = 0;
   const Int_t high = nn-1;
   const Double_t eps = TMath::Power(2.0,-52.0);
   Double_t exshift = 0.0;
   Double_t p=0,q=0,r=0,s=0,z=0,t,w,x,y;

   Double_t *pV = v.GetMatrixArray();
   Double_t *pD = d.GetMatrixArray();
   Double_t *pE = e.GetMatrixArray();
   Double_t *pH = H.GetMatrixArray();

   // Store roots isolated by balanc and compute matrix norm

   Double_t norm = 0.0;
   Int_t i,j,k;
   for (i = 0; i < nn; i++) {
      const Int_t off_i = i*nn;
      if ((i < low) || (i > high)) {
         pD[i] = pH[off_i+i];
         pE[i] = 0.0;
      }
      for (j = TMath::Max(i-1,0); j < nn; j++)
         norm += TMath::Abs(pH[off_i+j]);
   }

   // Outer loop over eigenvalue index

   Int_t iter = 0;
   while (n >= low) {
      const Int_t off_n  = n*nn;
      const Int_t off_n1 = (n-1)*nn;

      // Look for single small sub-diagonal element

      Int_t l = n;
      while (l > low) {
         const Int_t off_l1 = (l-1)*nn;
         const Int_t off_l  = l*nn;
         s = TMath::Abs(pH[off_l1+l-1])+TMath::Abs(pH[off_l+l]);
         if (s == 0.0)
            s = norm;
         if (TMath::Abs(pH[off_l+l-1]) < eps*s)
            break;
         l--;
      }

      // Check for convergence
      // One root found

      if (l == n) {
         pH[off_n+n] = pH[off_n+n]+exshift;
         pD[n] = pH[off_n+n];
         pE[n] = 0.0;
         n--;
         iter = 0;

         // Two roots found

      } else if (l == n-1) {
         w = pH[off_n+n-1]*pH[off_n1+n];
         p = (pH[off_n1+n-1]-pH[off_n+n])/2.0;
         q = p*p+w;
         z = TMath::Sqrt(TMath::Abs(q));
         pH[off_n+n] = pH[off_n+n]+exshift;
         pH[off_n1+n-1] = pH[off_n1+n-1]+exshift;
         x = pH[off_n+n];

         // Double_t pair

         if (q >= 0) {
            if (p >= 0)
               z = p+z;
            else
               z = p-z;
            pD[n-1] = x+z;
            pD[n] = pD[n-1];
            if (z != 0.0)
               pD[n] = x-w/z;
            pE[n-1] = 0.0;
            pE[n] = 0.0;
            x = pH[off_n+n-1];
            s = TMath::Abs(x)+TMath::Abs(z);
            p = x/s;
            q = z/s;
            r = TMath::Sqrt((p*p)+(q*q));
            p = p/r;
            q = q/r;

            // Row modification

            for (j = n-1; j < nn; j++) {
               z = pH[off_n1+j];
               pH[off_n1+j] = q*z+p*pH[off_n+j];
               pH[off_n+j]  = q*pH[off_n+j]-p*z;
            }

            // Column modification

            for (i = 0; i <= n; i++) {
               const Int_t off_i = i*nn;
               z = pH[off_i+n-1];
               pH[off_i+n-1] = q*z+p*pH[off_i+n];
               pH[off_i+n]  = q*pH[off_i+n]-p*z;
            }

            // Accumulate transformations

            for (i = low; i <= high; i++) {
               const Int_t off_i = i*nn;
               z = pV[off_i+n-1];
               pV[off_i+n-1] = q*z+p*pV[off_i+n];
               pV[off_i+n]   = q*pV[off_i+n]-p*z;
            }

            // Complex pair

         } else {
            pD[n-1] = x+p;
            pD[n] = x+p;
            pE[n-1] = z;
            pE[n] = -z;
         }
         n = n-2;
         iter = 0;

         // No convergence yet

      } else {

         // Form shift

         x = pH[off_n+n];
         y = 0.0;
         w = 0.0;
         if (l < n) {
            y = pH[off_n1+n-1];
            w = pH[off_n+n-1]*pH[off_n1+n];
         }

         // Wilkinson's original ad hoc shift

         if (iter == 10) {
            exshift += x;
            for (i = low; i <= n; i++) {
               const Int_t off_i = i*nn;
               pH[off_i+i] -= x;
            }
            s = TMath::Abs(pH[off_n+n-1])+TMath::Abs(pH[off_n1+n-2]);
            x = y = 0.75*s;
            w = -0.4375*s*s;
         }

         // MATLAB's new ad hoc shift

         if (iter == 30) {
            s = (y-x)/2.0;
            s = s*s+w;
            if (s > 0) {
               s = TMath::Sqrt(s);
               if (y<x)
                  s = -s;
               s = x-w/((y-x)/2.0+s);
               for (i = low; i <= n; i++) {
                  const Int_t off_i = i*nn;
                  pH[off_i+i] -= s;
               }
               exshift += s;
               x = y = w = 0.964;
            }
         }

         if (iter++ == 50) {  // (check iteration count here.)
            Error("MakeSchurr","too many iterations");
            break;
         }

         // Look for two consecutive small sub-diagonal elements

         Int_t m = n-2;
         while (m >= l) {
            const Int_t off_m   = m*nn;
            const Int_t off_m_1 = (m-1)*nn;
            const Int_t off_m1  = (m+1)*nn;
            const Int_t off_m2  = (m+2)*nn;
            z = pH[off_m+m];
            r = x-z;
            s = y-z;
            p = (r*s-w)/pH[off_m1+m]+pH[off_m+m+1];
            q = pH[off_m1+m+1]-z-r-s;
            r = pH[off_m2+m+1];
            s = TMath::Abs(p)+TMath::Abs(q)+TMath::Abs(r);
            p = p /s;
            q = q/s;
            r = r/s;
            if (m == l)
               break;
            if (TMath::Abs(pH[off_m+m-1])*(TMath::Abs(q)+TMath::Abs(r)) <
               eps*(TMath::Abs(p)*(TMath::Abs(pH[off_m_1+m-1])+TMath::Abs(z)+
               TMath::Abs(pH[off_m1+m+1]))))
               break;
               m--;
            }

            for (i = m+2; i <= n; i++) {
               const Int_t off_i = i*nn;
               pH[off_i+i-2] = 0.0;
               if (i > m+2)
                  pH[off_i+i-3] = 0.0;
            }

            // Double QR step involving rows l:n and columns m:n

            for (k = m; k <= n-1; k++) {
               const Int_t off_k  = k*nn;
               const Int_t off_k1 = (k+1)*nn;
               const Int_t off_k2 = (k+2)*nn;
               const Int_t notlast = (k != n-1);
               if (k != m) {
                  p = pH[off_k+k-1];
                  q = pH[off_k1+k-1];
                  r = (notlast ? pH[off_k2+k-1] : 0.0);
                  x = TMath::Abs(p)+TMath::Abs(q)+TMath::Abs(r);
                  if (x != 0.0) {
                     p = p/x;
                     q = q/x;
                     r = r/x;
                  }
               }
               if (x == 0.0)
                  break;
               s = TMath::Sqrt(p*p+q*q+r*r);
               if (p < 0) {
                  s = -s;
               }
               if (s != 0) {
                  if (k != m)
                     pH[off_k+k-1] = -s*x;
                  else if (l != m)
                     pH[off_k+k-1] = -pH[off_k+k-1];
                  p = p+s;
                  x = p/s;
                  y = q/s;
                  z = r/s;
                  q = q/p;
                  r = r/p;

                  // Row modification

                  for (j = k; j < nn; j++) {
                     p = pH[off_k+j]+q*pH[off_k1+j];
                     if (notlast) {
                        p = p+r*pH[off_k2+j];
                        pH[off_k2+j] = pH[off_k2+j]-p*z;
                     }
                     pH[off_k+j]  = pH[off_k+j]-p*x;
                     pH[off_k1+j] = pH[off_k1+j]-p*y;
                  }

                  // Column modification

                  for (i = 0; i <= TMath::Min(n,k+3); i++) {
                     const Int_t off_i = i*nn;
                     p = x*pH[off_i+k]+y*pH[off_i+k+1];
                     if (notlast) {
                        p = p+z*pH[off_i+k+2];
                        pH[off_i+k+2] = pH[off_i+k+2]-p*r;
                     }
                     pH[off_i+k]   = pH[off_i+k]-p;
                     pH[off_i+k+1] = pH[off_i+k+1]-p*q;
                  }

                  // Accumulate transformations

                  for (i = low; i <= high; i++) {
                     const Int_t off_i = i*nn;
                     p = x*pV[off_i+k]+y*pV[off_i+k+1];
                     if (notlast) {
                        p = p+z*pV[off_i+k+2];
                        pV[off_i+k+2] = pV[off_i+k+2]-p*r;
                     }
                     pV[off_i+k]   = pV[off_i+k]-p;
                     pV[off_i+k+1] = pV[off_i+k+1]-p*q;
                  }
               }  // (s != 0)
            }  // k loop
         }  // check convergence
      }  // while (n >= low)

      // Backsubstitute to find vectors of upper triangular form

      if (norm == 0.0)
         return;

      for (n = nn-1; n >= 0; n--) {
         p = pD[n];
         q = pE[n];

         // Double_t vector

         const Int_t off_n = n*nn;
         if (q == 0) {
            Int_t l = n;
            pH[off_n+n] = 1.0;
            for (i = n-1; i >= 0; i--) {
               const Int_t off_i  = i*nn;
               const Int_t off_i1 = (i+1)*nn;
               w = pH[off_i+i]-p;
               r = 0.0;
               for (j = l; j <= n; j++) {
                  const Int_t off_j = j*nn;
                  r = r+pH[off_i+j]*pH[off_j+n];
               }
               if (pE[i] < 0.0) {
                  z = w;
                  s = r;
               } else {
                  l = i;
                  if (pE[i] == 0.0) {
                     if (w != 0.0)
                        pH[off_i+n] = -r/w;
                  else
                     pH[off_i+n] = -r/(eps*norm);

                  // Solve real equations

               } else {
                  x = pH[off_i+i+1];
                  y = pH[off_i1+i];
                  q = (pD[i]-p)*(pD[i]-p)+pE[i]*pE[i];
                  t = (x*s-z*r)/q;
                  pH[off_i+n] = t;
                  if (TMath::Abs(x) > TMath::Abs(z))
                     pH[i+1+n] = (-r-w*t)/x;
                  else
                     pH[i+1+n] = (-s-y*t)/z;
               }

               // Overflow control

               t = TMath::Abs(pH[off_i+n]);
               if ((eps*t)*t > 1) {
                  for (j = i; j <= n; j++) {
                     const Int_t off_j = j*nn;
                     pH[off_j+n] = pH[off_j+n]/t;
                  }
               }
            }
         }

         // Complex vector

      } else if (q < 0) {
         Int_t l = n-1;
         const Int_t off_n1 = (n-1)*nn;

         // Last vector component imaginary so matrix is triangular

         if (TMath::Abs(pH[off_n+n-1]) > TMath::Abs(pH[off_n1+n])) {
            pH[off_n1+n-1] = q/pH[off_n+n-1];
            pH[off_n1+n]   = -(pH[off_n+n]-p)/pH[off_n+n-1];
         } else {
            cdiv(0.0,-pH[off_n1+n],pH[off_n1+n-1]-p,q);
            pH[off_n1+n-1] = gCdivr;
            pH[off_n1+n]   = gCdivi;
         }
         pH[off_n+n-1] = 0.0;
         pH[off_n+n]   = 1.0;
         for (i = n-2; i >= 0; i--) {
            const Int_t off_i  = i*nn;
            const Int_t off_i1 = (i+1)*nn;
            Double_t ra = 0.0;
            Double_t sa = 0.0;
            for (j = l; j <= n; j++) {
               const Int_t off_j = j*nn;
               ra += pH[off_i+j]*pH[off_j+n-1];
               sa += pH[off_i+j]*pH[off_j+n];
            }
            w = pH[off_i+i]-p;

            if (pE[i] < 0.0) {
               z = w;
               r = ra;
               s = sa;
            } else {
               l = i;
               if (pE[i] == 0) {
                  cdiv(-ra,-sa,w,q);
                  pH[off_i+n-1] = gCdivr;
                  pH[off_i+n]   = gCdivi;
               } else {

                  // Solve complex equations

                  x = pH[off_i+i+1];
                  y = pH[off_i1+i];
                  Double_t vr = (pD[i]-p)*(pD[i]-p)+pE[i]*pE[i]-q*q;
                  Double_t vi = (pD[i]-p)*2.0*q;
                  if ((vr == 0.0) && (vi == 0.0)) {
                     vr = eps*norm*(TMath::Abs(w)+TMath::Abs(q)+
                          TMath::Abs(x)+TMath::Abs(y)+TMath::Abs(z));
                  }
                  cdiv(x*r-z*ra+q*sa,x*s-z*sa-q*ra,vr,vi);
                  pH[off_i+n-1] = gCdivr;
                  pH[off_i+n]   = gCdivi;
                  if (TMath::Abs(x) > (TMath::Abs(z)+TMath::Abs(q))) {
                     pH[off_i1+n-1] = (-ra-w*pH[off_i+n-1]+q*pH[off_i+n])/x;
                     pH[off_i1+n]   = (-sa-w*pH[off_i+n]-q*pH[off_i+n-1])/x;
                  } else {
                     cdiv(-r-y*pH[off_i+n-1],-s-y*pH[off_i+n],z,q);
                     pH[off_i1+n-1] = gCdivr;
                     pH[off_i1+n]   = gCdivi;
                  }
               }

               // Overflow control

               t = TMath::Max(TMath::Abs(pH[off_i+n-1]),TMath::Abs(pH[off_i+n]));
               if ((eps*t)*t > 1) {
                  for (j = i; j <= n; j++) {
                     const Int_t off_j = j*nn;
                     pH[off_j+n-1] = pH[off_j+n-1]/t;
                     pH[off_j+n]   = pH[off_j+n]/t;
                  }
               }
            }
         }
      }
   }

   // Vectors of isolated roots

   for (i = 0; i < nn; i++) {
      if (i < low || i > high) {
         const Int_t off_i = i*nn;
         for (j = i; j < nn; j++)
            pV[off_i+j] = pH[off_i+j];
      }
   }

   // Back transformation to get eigenvectors of original matrix

   for (j = nn-1; j >= low; j--) {
      for (i = low; i <= high; i++) {
         const Int_t off_i = i*nn;
         z = 0.0;
         for (k = low; k <= TMath::Min(j,high); k++) {
            const Int_t off_k = k*nn;
            z = z+pV[off_i+k]*pH[off_k+j];
         }
         pV[off_i+j] = z;
      }
   }

}

//______________________________________________________________________________
void TMatrixDEigen::Sort(TMatrixD &v,TVectorD &d,TVectorD &e)
{
// Sort eigenvalues and corresponding vectors in descending order of Re^2+Im^2
// of the complex eigenvalues .

   // Sort eigenvalues and corresponding vectors.
   Double_t *pV = v.GetMatrixArray();
   Double_t *pD = d.GetMatrixArray();
   Double_t *pE = e.GetMatrixArray();

   const Int_t n = v.GetNrows();

   for (Int_t i = 0; i < n-1; i++) {
      Int_t k = i;
      Double_t norm = pD[i]*pD[i]+pE[i]*pE[i];
      Int_t j;
      for (j = i+1; j < n; j++) {
         const Double_t norm_new = pD[j]*pD[j]+pE[j]*pE[j];
         if (norm_new > norm) {
            k = j;
            norm = norm_new;
         }
      }
      if (k != i) {
         Double_t tmp;
         tmp   = pD[k];
         pD[k] = pD[i];
         pD[i] = tmp;
         tmp   = pE[k];
         pE[k] = pE[i];
         pE[i] = tmp;
         for (j = 0; j < n; j++) {
            const Int_t off_j = j*n;
            tmp = pV[off_j+i];
            pV[off_j+i] = pV[off_j+k];
            pV[off_j+k] = tmp;
         }
      }
   }
}

//______________________________________________________________________________
TMatrixDEigen &TMatrixDEigen::operator=(const TMatrixDEigen &source)
{
// Assignment operator

   if (this != &source) {
      fEigenVectors.ResizeTo(source.fEigenVectors);
      fEigenValuesRe.ResizeTo(source.fEigenValuesRe);
      fEigenValuesIm.ResizeTo(source.fEigenValuesIm);
   }
   return *this;
}

//______________________________________________________________________________
const TMatrixD TMatrixDEigen::GetEigenValues() const
{
// Computes the block diagonal eigenvalue matrix.
// If the original matrix A is not symmetric, then the eigenvalue
// matrix D is block diagonal with the real eigenvalues in 1-by-1
// blocks and any complex eigenvalues,
//    a + i*b, in 2-by-2 blocks, [a, b; -b, a].
//  That is, if the complex eigenvalues look like
//
//     u + iv     .        .          .      .    .
//       .      u - iv     .          .      .    .
//       .        .      a + ib       .      .    .
//       .        .        .        a - ib   .    .
//       .        .        .          .      x    .
//       .        .        .          .      .    y
//
// then D looks like
//
//     u        v        .          .      .    .
//    -v        u        .          .      .    .
//     .        .        a          b      .    .
//     .        .       -b          a      .    .
//     .        .        .          .      x    .
//     .        .        .          .      .    y
//
// This keeps V a real matrix in both symmetric and non-symmetric
// cases, and A*V = V*D.
//
// Indexing:
//  If matrix A has the index/shape (rowLwb,rowUpb,rowLwb,rowUpb)
//  each eigen-vector must have the shape (rowLwb,rowUpb) .
//  For convinience, the column index of the eigen-vector matrix
//  also runs from rowLwb to rowUpb so that the returned matrix
//  has also index/shape (rowLwb,rowUpb,rowLwb,rowUpb) .
//
   const Int_t nrows  = fEigenVectors.GetNrows();
   const Int_t rowLwb = fEigenVectors.GetRowLwb();
   const Int_t rowUpb = rowLwb+nrows-1;

   TMatrixD mD(rowLwb,rowUpb,rowLwb,rowUpb);

   Double_t *pD = mD.GetMatrixArray();
   const Double_t * const pd = fEigenValuesRe.GetMatrixArray();
   const Double_t * const pe = fEigenValuesIm.GetMatrixArray();

   for (Int_t i = 0; i < nrows; i++) {
      const Int_t off_i = i*nrows;
      for (Int_t j = 0; j < nrows; j++)
         pD[off_i+j] = 0.0;
      pD[off_i+i] = pd[i];
      if (pe[i] > 0) {
         pD[off_i+i+1] = pe[i];
      } else if (pe[i] < 0) {
         pD[off_i+i-1] = pe[i];
      }
   }

   return mD;
}
 TMatrixDEigen.cxx:1
 TMatrixDEigen.cxx:2
 TMatrixDEigen.cxx:3
 TMatrixDEigen.cxx:4
 TMatrixDEigen.cxx:5
 TMatrixDEigen.cxx:6
 TMatrixDEigen.cxx:7
 TMatrixDEigen.cxx:8
 TMatrixDEigen.cxx:9
 TMatrixDEigen.cxx:10
 TMatrixDEigen.cxx:11
 TMatrixDEigen.cxx:12
 TMatrixDEigen.cxx:13
 TMatrixDEigen.cxx:14
 TMatrixDEigen.cxx:15
 TMatrixDEigen.cxx:16
 TMatrixDEigen.cxx:17
 TMatrixDEigen.cxx:18
 TMatrixDEigen.cxx:19
 TMatrixDEigen.cxx:20
 TMatrixDEigen.cxx:21
 TMatrixDEigen.cxx:22
 TMatrixDEigen.cxx:23
 TMatrixDEigen.cxx:24
 TMatrixDEigen.cxx:25
 TMatrixDEigen.cxx:26
 TMatrixDEigen.cxx:27
 TMatrixDEigen.cxx:28
 TMatrixDEigen.cxx:29
 TMatrixDEigen.cxx:30
 TMatrixDEigen.cxx:31
 TMatrixDEigen.cxx:32
 TMatrixDEigen.cxx:33
 TMatrixDEigen.cxx:34
 TMatrixDEigen.cxx:35
 TMatrixDEigen.cxx:36
 TMatrixDEigen.cxx:37
 TMatrixDEigen.cxx:38
 TMatrixDEigen.cxx:39
 TMatrixDEigen.cxx:40
 TMatrixDEigen.cxx:41
 TMatrixDEigen.cxx:42
 TMatrixDEigen.cxx:43
 TMatrixDEigen.cxx:44
 TMatrixDEigen.cxx:45
 TMatrixDEigen.cxx:46
 TMatrixDEigen.cxx:47
 TMatrixDEigen.cxx:48
 TMatrixDEigen.cxx:49
 TMatrixDEigen.cxx:50
 TMatrixDEigen.cxx:51
 TMatrixDEigen.cxx:52
 TMatrixDEigen.cxx:53
 TMatrixDEigen.cxx:54
 TMatrixDEigen.cxx:55
 TMatrixDEigen.cxx:56
 TMatrixDEigen.cxx:57
 TMatrixDEigen.cxx:58
 TMatrixDEigen.cxx:59
 TMatrixDEigen.cxx:60
 TMatrixDEigen.cxx:61
 TMatrixDEigen.cxx:62
 TMatrixDEigen.cxx:63
 TMatrixDEigen.cxx:64
 TMatrixDEigen.cxx:65
 TMatrixDEigen.cxx:66
 TMatrixDEigen.cxx:67
 TMatrixDEigen.cxx:68
 TMatrixDEigen.cxx:69
 TMatrixDEigen.cxx:70
 TMatrixDEigen.cxx:71
 TMatrixDEigen.cxx:72
 TMatrixDEigen.cxx:73
 TMatrixDEigen.cxx:74
 TMatrixDEigen.cxx:75
 TMatrixDEigen.cxx:76
 TMatrixDEigen.cxx:77
 TMatrixDEigen.cxx:78
 TMatrixDEigen.cxx:79
 TMatrixDEigen.cxx:80
 TMatrixDEigen.cxx:81
 TMatrixDEigen.cxx:82
 TMatrixDEigen.cxx:83
 TMatrixDEigen.cxx:84
 TMatrixDEigen.cxx:85
 TMatrixDEigen.cxx:86
 TMatrixDEigen.cxx:87
 TMatrixDEigen.cxx:88
 TMatrixDEigen.cxx:89
 TMatrixDEigen.cxx:90
 TMatrixDEigen.cxx:91
 TMatrixDEigen.cxx:92
 TMatrixDEigen.cxx:93
 TMatrixDEigen.cxx:94
 TMatrixDEigen.cxx:95
 TMatrixDEigen.cxx:96
 TMatrixDEigen.cxx:97
 TMatrixDEigen.cxx:98
 TMatrixDEigen.cxx:99
 TMatrixDEigen.cxx:100
 TMatrixDEigen.cxx:101
 TMatrixDEigen.cxx:102
 TMatrixDEigen.cxx:103
 TMatrixDEigen.cxx:104
 TMatrixDEigen.cxx:105
 TMatrixDEigen.cxx:106
 TMatrixDEigen.cxx:107
 TMatrixDEigen.cxx:108
 TMatrixDEigen.cxx:109
 TMatrixDEigen.cxx:110
 TMatrixDEigen.cxx:111
 TMatrixDEigen.cxx:112
 TMatrixDEigen.cxx:113
 TMatrixDEigen.cxx:114
 TMatrixDEigen.cxx:115
 TMatrixDEigen.cxx:116
 TMatrixDEigen.cxx:117
 TMatrixDEigen.cxx:118
 TMatrixDEigen.cxx:119
 TMatrixDEigen.cxx:120
 TMatrixDEigen.cxx:121
 TMatrixDEigen.cxx:122
 TMatrixDEigen.cxx:123
 TMatrixDEigen.cxx:124
 TMatrixDEigen.cxx:125
 TMatrixDEigen.cxx:126
 TMatrixDEigen.cxx:127
 TMatrixDEigen.cxx:128
 TMatrixDEigen.cxx:129
 TMatrixDEigen.cxx:130
 TMatrixDEigen.cxx:131
 TMatrixDEigen.cxx:132
 TMatrixDEigen.cxx:133
 TMatrixDEigen.cxx:134
 TMatrixDEigen.cxx:135
 TMatrixDEigen.cxx:136
 TMatrixDEigen.cxx:137
 TMatrixDEigen.cxx:138
 TMatrixDEigen.cxx:139
 TMatrixDEigen.cxx:140
 TMatrixDEigen.cxx:141
 TMatrixDEigen.cxx:142
 TMatrixDEigen.cxx:143
 TMatrixDEigen.cxx:144
 TMatrixDEigen.cxx:145
 TMatrixDEigen.cxx:146
 TMatrixDEigen.cxx:147
 TMatrixDEigen.cxx:148
 TMatrixDEigen.cxx:149
 TMatrixDEigen.cxx:150
 TMatrixDEigen.cxx:151
 TMatrixDEigen.cxx:152
 TMatrixDEigen.cxx:153
 TMatrixDEigen.cxx:154
 TMatrixDEigen.cxx:155
 TMatrixDEigen.cxx:156
 TMatrixDEigen.cxx:157
 TMatrixDEigen.cxx:158
 TMatrixDEigen.cxx:159
 TMatrixDEigen.cxx:160
 TMatrixDEigen.cxx:161
 TMatrixDEigen.cxx:162
 TMatrixDEigen.cxx:163
 TMatrixDEigen.cxx:164
 TMatrixDEigen.cxx:165
 TMatrixDEigen.cxx:166
 TMatrixDEigen.cxx:167
 TMatrixDEigen.cxx:168
 TMatrixDEigen.cxx:169
 TMatrixDEigen.cxx:170
 TMatrixDEigen.cxx:171
 TMatrixDEigen.cxx:172
 TMatrixDEigen.cxx:173
 TMatrixDEigen.cxx:174
 TMatrixDEigen.cxx:175
 TMatrixDEigen.cxx:176
 TMatrixDEigen.cxx:177
 TMatrixDEigen.cxx:178
 TMatrixDEigen.cxx:179
 TMatrixDEigen.cxx:180
 TMatrixDEigen.cxx:181
 TMatrixDEigen.cxx:182
 TMatrixDEigen.cxx:183
 TMatrixDEigen.cxx:184
 TMatrixDEigen.cxx:185
 TMatrixDEigen.cxx:186
 TMatrixDEigen.cxx:187
 TMatrixDEigen.cxx:188
 TMatrixDEigen.cxx:189
 TMatrixDEigen.cxx:190
 TMatrixDEigen.cxx:191
 TMatrixDEigen.cxx:192
 TMatrixDEigen.cxx:193
 TMatrixDEigen.cxx:194
 TMatrixDEigen.cxx:195
 TMatrixDEigen.cxx:196
 TMatrixDEigen.cxx:197
 TMatrixDEigen.cxx:198
 TMatrixDEigen.cxx:199
 TMatrixDEigen.cxx:200
 TMatrixDEigen.cxx:201
 TMatrixDEigen.cxx:202
 TMatrixDEigen.cxx:203
 TMatrixDEigen.cxx:204
 TMatrixDEigen.cxx:205
 TMatrixDEigen.cxx:206
 TMatrixDEigen.cxx:207
 TMatrixDEigen.cxx:208
 TMatrixDEigen.cxx:209
 TMatrixDEigen.cxx:210
 TMatrixDEigen.cxx:211
 TMatrixDEigen.cxx:212
 TMatrixDEigen.cxx:213
 TMatrixDEigen.cxx:214
 TMatrixDEigen.cxx:215
 TMatrixDEigen.cxx:216
 TMatrixDEigen.cxx:217
 TMatrixDEigen.cxx:218
 TMatrixDEigen.cxx:219
 TMatrixDEigen.cxx:220
 TMatrixDEigen.cxx:221
 TMatrixDEigen.cxx:222
 TMatrixDEigen.cxx:223
 TMatrixDEigen.cxx:224
 TMatrixDEigen.cxx:225
 TMatrixDEigen.cxx:226
 TMatrixDEigen.cxx:227
 TMatrixDEigen.cxx:228
 TMatrixDEigen.cxx:229
 TMatrixDEigen.cxx:230
 TMatrixDEigen.cxx:231
 TMatrixDEigen.cxx:232
 TMatrixDEigen.cxx:233
 TMatrixDEigen.cxx:234
 TMatrixDEigen.cxx:235
 TMatrixDEigen.cxx:236
 TMatrixDEigen.cxx:237
 TMatrixDEigen.cxx:238
 TMatrixDEigen.cxx:239
 TMatrixDEigen.cxx:240
 TMatrixDEigen.cxx:241
 TMatrixDEigen.cxx:242
 TMatrixDEigen.cxx:243
 TMatrixDEigen.cxx:244
 TMatrixDEigen.cxx:245
 TMatrixDEigen.cxx:246
 TMatrixDEigen.cxx:247
 TMatrixDEigen.cxx:248
 TMatrixDEigen.cxx:249
 TMatrixDEigen.cxx:250
 TMatrixDEigen.cxx:251
 TMatrixDEigen.cxx:252
 TMatrixDEigen.cxx:253
 TMatrixDEigen.cxx:254
 TMatrixDEigen.cxx:255
 TMatrixDEigen.cxx:256
 TMatrixDEigen.cxx:257
 TMatrixDEigen.cxx:258
 TMatrixDEigen.cxx:259
 TMatrixDEigen.cxx:260
 TMatrixDEigen.cxx:261
 TMatrixDEigen.cxx:262
 TMatrixDEigen.cxx:263
 TMatrixDEigen.cxx:264
 TMatrixDEigen.cxx:265
 TMatrixDEigen.cxx:266
 TMatrixDEigen.cxx:267
 TMatrixDEigen.cxx:268
 TMatrixDEigen.cxx:269
 TMatrixDEigen.cxx:270
 TMatrixDEigen.cxx:271
 TMatrixDEigen.cxx:272
 TMatrixDEigen.cxx:273
 TMatrixDEigen.cxx:274
 TMatrixDEigen.cxx:275
 TMatrixDEigen.cxx:276
 TMatrixDEigen.cxx:277
 TMatrixDEigen.cxx:278
 TMatrixDEigen.cxx:279
 TMatrixDEigen.cxx:280
 TMatrixDEigen.cxx:281
 TMatrixDEigen.cxx:282
 TMatrixDEigen.cxx:283
 TMatrixDEigen.cxx:284
 TMatrixDEigen.cxx:285
 TMatrixDEigen.cxx:286
 TMatrixDEigen.cxx:287
 TMatrixDEigen.cxx:288
 TMatrixDEigen.cxx:289
 TMatrixDEigen.cxx:290
 TMatrixDEigen.cxx:291
 TMatrixDEigen.cxx:292
 TMatrixDEigen.cxx:293
 TMatrixDEigen.cxx:294
 TMatrixDEigen.cxx:295
 TMatrixDEigen.cxx:296
 TMatrixDEigen.cxx:297
 TMatrixDEigen.cxx:298
 TMatrixDEigen.cxx:299
 TMatrixDEigen.cxx:300
 TMatrixDEigen.cxx:301
 TMatrixDEigen.cxx:302
 TMatrixDEigen.cxx:303
 TMatrixDEigen.cxx:304
 TMatrixDEigen.cxx:305
 TMatrixDEigen.cxx:306
 TMatrixDEigen.cxx:307
 TMatrixDEigen.cxx:308
 TMatrixDEigen.cxx:309
 TMatrixDEigen.cxx:310
 TMatrixDEigen.cxx:311
 TMatrixDEigen.cxx:312
 TMatrixDEigen.cxx:313
 TMatrixDEigen.cxx:314
 TMatrixDEigen.cxx:315
 TMatrixDEigen.cxx:316
 TMatrixDEigen.cxx:317
 TMatrixDEigen.cxx:318
 TMatrixDEigen.cxx:319
 TMatrixDEigen.cxx:320
 TMatrixDEigen.cxx:321
 TMatrixDEigen.cxx:322
 TMatrixDEigen.cxx:323
 TMatrixDEigen.cxx:324
 TMatrixDEigen.cxx:325
 TMatrixDEigen.cxx:326
 TMatrixDEigen.cxx:327
 TMatrixDEigen.cxx:328
 TMatrixDEigen.cxx:329
 TMatrixDEigen.cxx:330
 TMatrixDEigen.cxx:331
 TMatrixDEigen.cxx:332
 TMatrixDEigen.cxx:333
 TMatrixDEigen.cxx:334
 TMatrixDEigen.cxx:335
 TMatrixDEigen.cxx:336
 TMatrixDEigen.cxx:337
 TMatrixDEigen.cxx:338
 TMatrixDEigen.cxx:339
 TMatrixDEigen.cxx:340
 TMatrixDEigen.cxx:341
 TMatrixDEigen.cxx:342
 TMatrixDEigen.cxx:343
 TMatrixDEigen.cxx:344
 TMatrixDEigen.cxx:345
 TMatrixDEigen.cxx:346
 TMatrixDEigen.cxx:347
 TMatrixDEigen.cxx:348
 TMatrixDEigen.cxx:349
 TMatrixDEigen.cxx:350
 TMatrixDEigen.cxx:351
 TMatrixDEigen.cxx:352
 TMatrixDEigen.cxx:353
 TMatrixDEigen.cxx:354
 TMatrixDEigen.cxx:355
 TMatrixDEigen.cxx:356
 TMatrixDEigen.cxx:357
 TMatrixDEigen.cxx:358
 TMatrixDEigen.cxx:359
 TMatrixDEigen.cxx:360
 TMatrixDEigen.cxx:361
 TMatrixDEigen.cxx:362
 TMatrixDEigen.cxx:363
 TMatrixDEigen.cxx:364
 TMatrixDEigen.cxx:365
 TMatrixDEigen.cxx:366
 TMatrixDEigen.cxx:367
 TMatrixDEigen.cxx:368
 TMatrixDEigen.cxx:369
 TMatrixDEigen.cxx:370
 TMatrixDEigen.cxx:371
 TMatrixDEigen.cxx:372
 TMatrixDEigen.cxx:373
 TMatrixDEigen.cxx:374
 TMatrixDEigen.cxx:375
 TMatrixDEigen.cxx:376
 TMatrixDEigen.cxx:377
 TMatrixDEigen.cxx:378
 TMatrixDEigen.cxx:379
 TMatrixDEigen.cxx:380
 TMatrixDEigen.cxx:381
 TMatrixDEigen.cxx:382
 TMatrixDEigen.cxx:383
 TMatrixDEigen.cxx:384
 TMatrixDEigen.cxx:385
 TMatrixDEigen.cxx:386
 TMatrixDEigen.cxx:387
 TMatrixDEigen.cxx:388
 TMatrixDEigen.cxx:389
 TMatrixDEigen.cxx:390
 TMatrixDEigen.cxx:391
 TMatrixDEigen.cxx:392
 TMatrixDEigen.cxx:393
 TMatrixDEigen.cxx:394
 TMatrixDEigen.cxx:395
 TMatrixDEigen.cxx:396
 TMatrixDEigen.cxx:397
 TMatrixDEigen.cxx:398
 TMatrixDEigen.cxx:399
 TMatrixDEigen.cxx:400
 TMatrixDEigen.cxx:401
 TMatrixDEigen.cxx:402
 TMatrixDEigen.cxx:403
 TMatrixDEigen.cxx:404
 TMatrixDEigen.cxx:405
 TMatrixDEigen.cxx:406
 TMatrixDEigen.cxx:407
 TMatrixDEigen.cxx:408
 TMatrixDEigen.cxx:409
 TMatrixDEigen.cxx:410
 TMatrixDEigen.cxx:411
 TMatrixDEigen.cxx:412
 TMatrixDEigen.cxx:413
 TMatrixDEigen.cxx:414
 TMatrixDEigen.cxx:415
 TMatrixDEigen.cxx:416
 TMatrixDEigen.cxx:417
 TMatrixDEigen.cxx:418
 TMatrixDEigen.cxx:419
 TMatrixDEigen.cxx:420
 TMatrixDEigen.cxx:421
 TMatrixDEigen.cxx:422
 TMatrixDEigen.cxx:423
 TMatrixDEigen.cxx:424
 TMatrixDEigen.cxx:425
 TMatrixDEigen.cxx:426
 TMatrixDEigen.cxx:427
 TMatrixDEigen.cxx:428
 TMatrixDEigen.cxx:429
 TMatrixDEigen.cxx:430
 TMatrixDEigen.cxx:431
 TMatrixDEigen.cxx:432
 TMatrixDEigen.cxx:433
 TMatrixDEigen.cxx:434
 TMatrixDEigen.cxx:435
 TMatrixDEigen.cxx:436
 TMatrixDEigen.cxx:437
 TMatrixDEigen.cxx:438
 TMatrixDEigen.cxx:439
 TMatrixDEigen.cxx:440
 TMatrixDEigen.cxx:441
 TMatrixDEigen.cxx:442
 TMatrixDEigen.cxx:443
 TMatrixDEigen.cxx:444
 TMatrixDEigen.cxx:445
 TMatrixDEigen.cxx:446
 TMatrixDEigen.cxx:447
 TMatrixDEigen.cxx:448
 TMatrixDEigen.cxx:449
 TMatrixDEigen.cxx:450
 TMatrixDEigen.cxx:451
 TMatrixDEigen.cxx:452
 TMatrixDEigen.cxx:453
 TMatrixDEigen.cxx:454
 TMatrixDEigen.cxx:455
 TMatrixDEigen.cxx:456
 TMatrixDEigen.cxx:457
 TMatrixDEigen.cxx:458
 TMatrixDEigen.cxx:459
 TMatrixDEigen.cxx:460
 TMatrixDEigen.cxx:461
 TMatrixDEigen.cxx:462
 TMatrixDEigen.cxx:463
 TMatrixDEigen.cxx:464
 TMatrixDEigen.cxx:465
 TMatrixDEigen.cxx:466
 TMatrixDEigen.cxx:467
 TMatrixDEigen.cxx:468
 TMatrixDEigen.cxx:469
 TMatrixDEigen.cxx:470
 TMatrixDEigen.cxx:471
 TMatrixDEigen.cxx:472
 TMatrixDEigen.cxx:473
 TMatrixDEigen.cxx:474
 TMatrixDEigen.cxx:475
 TMatrixDEigen.cxx:476
 TMatrixDEigen.cxx:477
 TMatrixDEigen.cxx:478
 TMatrixDEigen.cxx:479
 TMatrixDEigen.cxx:480
 TMatrixDEigen.cxx:481
 TMatrixDEigen.cxx:482
 TMatrixDEigen.cxx:483
 TMatrixDEigen.cxx:484
 TMatrixDEigen.cxx:485
 TMatrixDEigen.cxx:486
 TMatrixDEigen.cxx:487
 TMatrixDEigen.cxx:488
 TMatrixDEigen.cxx:489
 TMatrixDEigen.cxx:490
 TMatrixDEigen.cxx:491
 TMatrixDEigen.cxx:492
 TMatrixDEigen.cxx:493
 TMatrixDEigen.cxx:494
 TMatrixDEigen.cxx:495
 TMatrixDEigen.cxx:496
 TMatrixDEigen.cxx:497
 TMatrixDEigen.cxx:498
 TMatrixDEigen.cxx:499
 TMatrixDEigen.cxx:500
 TMatrixDEigen.cxx:501
 TMatrixDEigen.cxx:502
 TMatrixDEigen.cxx:503
 TMatrixDEigen.cxx:504
 TMatrixDEigen.cxx:505
 TMatrixDEigen.cxx:506
 TMatrixDEigen.cxx:507
 TMatrixDEigen.cxx:508
 TMatrixDEigen.cxx:509
 TMatrixDEigen.cxx:510
 TMatrixDEigen.cxx:511
 TMatrixDEigen.cxx:512
 TMatrixDEigen.cxx:513
 TMatrixDEigen.cxx:514
 TMatrixDEigen.cxx:515
 TMatrixDEigen.cxx:516
 TMatrixDEigen.cxx:517
 TMatrixDEigen.cxx:518
 TMatrixDEigen.cxx:519
 TMatrixDEigen.cxx:520
 TMatrixDEigen.cxx:521
 TMatrixDEigen.cxx:522
 TMatrixDEigen.cxx:523
 TMatrixDEigen.cxx:524
 TMatrixDEigen.cxx:525
 TMatrixDEigen.cxx:526
 TMatrixDEigen.cxx:527
 TMatrixDEigen.cxx:528
 TMatrixDEigen.cxx:529
 TMatrixDEigen.cxx:530
 TMatrixDEigen.cxx:531
 TMatrixDEigen.cxx:532
 TMatrixDEigen.cxx:533
 TMatrixDEigen.cxx:534
 TMatrixDEigen.cxx:535
 TMatrixDEigen.cxx:536
 TMatrixDEigen.cxx:537
 TMatrixDEigen.cxx:538
 TMatrixDEigen.cxx:539
 TMatrixDEigen.cxx:540
 TMatrixDEigen.cxx:541
 TMatrixDEigen.cxx:542
 TMatrixDEigen.cxx:543
 TMatrixDEigen.cxx:544
 TMatrixDEigen.cxx:545
 TMatrixDEigen.cxx:546
 TMatrixDEigen.cxx:547
 TMatrixDEigen.cxx:548
 TMatrixDEigen.cxx:549
 TMatrixDEigen.cxx:550
 TMatrixDEigen.cxx:551
 TMatrixDEigen.cxx:552
 TMatrixDEigen.cxx:553
 TMatrixDEigen.cxx:554
 TMatrixDEigen.cxx:555
 TMatrixDEigen.cxx:556
 TMatrixDEigen.cxx:557
 TMatrixDEigen.cxx:558
 TMatrixDEigen.cxx:559
 TMatrixDEigen.cxx:560
 TMatrixDEigen.cxx:561
 TMatrixDEigen.cxx:562
 TMatrixDEigen.cxx:563
 TMatrixDEigen.cxx:564
 TMatrixDEigen.cxx:565
 TMatrixDEigen.cxx:566
 TMatrixDEigen.cxx:567
 TMatrixDEigen.cxx:568
 TMatrixDEigen.cxx:569
 TMatrixDEigen.cxx:570
 TMatrixDEigen.cxx:571
 TMatrixDEigen.cxx:572
 TMatrixDEigen.cxx:573
 TMatrixDEigen.cxx:574
 TMatrixDEigen.cxx:575
 TMatrixDEigen.cxx:576
 TMatrixDEigen.cxx:577
 TMatrixDEigen.cxx:578
 TMatrixDEigen.cxx:579
 TMatrixDEigen.cxx:580
 TMatrixDEigen.cxx:581
 TMatrixDEigen.cxx:582
 TMatrixDEigen.cxx:583
 TMatrixDEigen.cxx:584
 TMatrixDEigen.cxx:585
 TMatrixDEigen.cxx:586
 TMatrixDEigen.cxx:587
 TMatrixDEigen.cxx:588
 TMatrixDEigen.cxx:589
 TMatrixDEigen.cxx:590
 TMatrixDEigen.cxx:591
 TMatrixDEigen.cxx:592
 TMatrixDEigen.cxx:593
 TMatrixDEigen.cxx:594
 TMatrixDEigen.cxx:595
 TMatrixDEigen.cxx:596
 TMatrixDEigen.cxx:597
 TMatrixDEigen.cxx:598
 TMatrixDEigen.cxx:599
 TMatrixDEigen.cxx:600
 TMatrixDEigen.cxx:601
 TMatrixDEigen.cxx:602
 TMatrixDEigen.cxx:603
 TMatrixDEigen.cxx:604
 TMatrixDEigen.cxx:605
 TMatrixDEigen.cxx:606
 TMatrixDEigen.cxx:607
 TMatrixDEigen.cxx:608
 TMatrixDEigen.cxx:609
 TMatrixDEigen.cxx:610
 TMatrixDEigen.cxx:611
 TMatrixDEigen.cxx:612
 TMatrixDEigen.cxx:613
 TMatrixDEigen.cxx:614
 TMatrixDEigen.cxx:615
 TMatrixDEigen.cxx:616
 TMatrixDEigen.cxx:617
 TMatrixDEigen.cxx:618
 TMatrixDEigen.cxx:619
 TMatrixDEigen.cxx:620
 TMatrixDEigen.cxx:621
 TMatrixDEigen.cxx:622
 TMatrixDEigen.cxx:623
 TMatrixDEigen.cxx:624
 TMatrixDEigen.cxx:625
 TMatrixDEigen.cxx:626
 TMatrixDEigen.cxx:627
 TMatrixDEigen.cxx:628
 TMatrixDEigen.cxx:629
 TMatrixDEigen.cxx:630
 TMatrixDEigen.cxx:631
 TMatrixDEigen.cxx:632
 TMatrixDEigen.cxx:633
 TMatrixDEigen.cxx:634
 TMatrixDEigen.cxx:635
 TMatrixDEigen.cxx:636
 TMatrixDEigen.cxx:637
 TMatrixDEigen.cxx:638
 TMatrixDEigen.cxx:639
 TMatrixDEigen.cxx:640
 TMatrixDEigen.cxx:641
 TMatrixDEigen.cxx:642
 TMatrixDEigen.cxx:643
 TMatrixDEigen.cxx:644
 TMatrixDEigen.cxx:645
 TMatrixDEigen.cxx:646
 TMatrixDEigen.cxx:647
 TMatrixDEigen.cxx:648
 TMatrixDEigen.cxx:649
 TMatrixDEigen.cxx:650
 TMatrixDEigen.cxx:651
 TMatrixDEigen.cxx:652
 TMatrixDEigen.cxx:653
 TMatrixDEigen.cxx:654
 TMatrixDEigen.cxx:655
 TMatrixDEigen.cxx:656
 TMatrixDEigen.cxx:657
 TMatrixDEigen.cxx:658
 TMatrixDEigen.cxx:659
 TMatrixDEigen.cxx:660
 TMatrixDEigen.cxx:661
 TMatrixDEigen.cxx:662
 TMatrixDEigen.cxx:663
 TMatrixDEigen.cxx:664
 TMatrixDEigen.cxx:665
 TMatrixDEigen.cxx:666
 TMatrixDEigen.cxx:667
 TMatrixDEigen.cxx:668
 TMatrixDEigen.cxx:669
 TMatrixDEigen.cxx:670
 TMatrixDEigen.cxx:671
 TMatrixDEigen.cxx:672
 TMatrixDEigen.cxx:673
 TMatrixDEigen.cxx:674
 TMatrixDEigen.cxx:675
 TMatrixDEigen.cxx:676
 TMatrixDEigen.cxx:677
 TMatrixDEigen.cxx:678
 TMatrixDEigen.cxx:679
 TMatrixDEigen.cxx:680
 TMatrixDEigen.cxx:681
 TMatrixDEigen.cxx:682
 TMatrixDEigen.cxx:683
 TMatrixDEigen.cxx:684
 TMatrixDEigen.cxx:685
 TMatrixDEigen.cxx:686
 TMatrixDEigen.cxx:687
 TMatrixDEigen.cxx:688
 TMatrixDEigen.cxx:689
 TMatrixDEigen.cxx:690
 TMatrixDEigen.cxx:691
 TMatrixDEigen.cxx:692
 TMatrixDEigen.cxx:693
 TMatrixDEigen.cxx:694
 TMatrixDEigen.cxx:695
 TMatrixDEigen.cxx:696
 TMatrixDEigen.cxx:697
 TMatrixDEigen.cxx:698
 TMatrixDEigen.cxx:699
 TMatrixDEigen.cxx:700
 TMatrixDEigen.cxx:701
 TMatrixDEigen.cxx:702
 TMatrixDEigen.cxx:703
 TMatrixDEigen.cxx:704
 TMatrixDEigen.cxx:705
 TMatrixDEigen.cxx:706
 TMatrixDEigen.cxx:707
 TMatrixDEigen.cxx:708
 TMatrixDEigen.cxx:709
 TMatrixDEigen.cxx:710
 TMatrixDEigen.cxx:711
 TMatrixDEigen.cxx:712
 TMatrixDEigen.cxx:713
 TMatrixDEigen.cxx:714
 TMatrixDEigen.cxx:715
 TMatrixDEigen.cxx:716
 TMatrixDEigen.cxx:717
 TMatrixDEigen.cxx:718
 TMatrixDEigen.cxx:719
 TMatrixDEigen.cxx:720
 TMatrixDEigen.cxx:721
 TMatrixDEigen.cxx:722
 TMatrixDEigen.cxx:723
 TMatrixDEigen.cxx:724
 TMatrixDEigen.cxx:725
 TMatrixDEigen.cxx:726
 TMatrixDEigen.cxx:727
 TMatrixDEigen.cxx:728
 TMatrixDEigen.cxx:729
 TMatrixDEigen.cxx:730
 TMatrixDEigen.cxx:731
 TMatrixDEigen.cxx:732
 TMatrixDEigen.cxx:733
 TMatrixDEigen.cxx:734
 TMatrixDEigen.cxx:735
 TMatrixDEigen.cxx:736
 TMatrixDEigen.cxx:737
 TMatrixDEigen.cxx:738
 TMatrixDEigen.cxx:739
 TMatrixDEigen.cxx:740
 TMatrixDEigen.cxx:741
 TMatrixDEigen.cxx:742
 TMatrixDEigen.cxx:743
 TMatrixDEigen.cxx:744
 TMatrixDEigen.cxx:745
 TMatrixDEigen.cxx:746
 TMatrixDEigen.cxx:747
 TMatrixDEigen.cxx:748
 TMatrixDEigen.cxx:749
 TMatrixDEigen.cxx:750
 TMatrixDEigen.cxx:751
 TMatrixDEigen.cxx:752
 TMatrixDEigen.cxx:753
 TMatrixDEigen.cxx:754
 TMatrixDEigen.cxx:755
 TMatrixDEigen.cxx:756
 TMatrixDEigen.cxx:757
 TMatrixDEigen.cxx:758
 TMatrixDEigen.cxx:759
 TMatrixDEigen.cxx:760
 TMatrixDEigen.cxx:761
 TMatrixDEigen.cxx:762
 TMatrixDEigen.cxx:763
 TMatrixDEigen.cxx:764
 TMatrixDEigen.cxx:765
 TMatrixDEigen.cxx:766
 TMatrixDEigen.cxx:767
 TMatrixDEigen.cxx:768
 TMatrixDEigen.cxx:769
 TMatrixDEigen.cxx:770
 TMatrixDEigen.cxx:771
 TMatrixDEigen.cxx:772
 TMatrixDEigen.cxx:773
 TMatrixDEigen.cxx:774
 TMatrixDEigen.cxx:775
 TMatrixDEigen.cxx:776
 TMatrixDEigen.cxx:777
 TMatrixDEigen.cxx:778
 TMatrixDEigen.cxx:779
 TMatrixDEigen.cxx:780
 TMatrixDEigen.cxx:781
 TMatrixDEigen.cxx:782
 TMatrixDEigen.cxx:783
 TMatrixDEigen.cxx:784
 TMatrixDEigen.cxx:785
 TMatrixDEigen.cxx:786
 TMatrixDEigen.cxx:787
 TMatrixDEigen.cxx:788
 TMatrixDEigen.cxx:789
 TMatrixDEigen.cxx:790
 TMatrixDEigen.cxx:791
 TMatrixDEigen.cxx:792
 TMatrixDEigen.cxx:793
 TMatrixDEigen.cxx:794
 TMatrixDEigen.cxx:795
 TMatrixDEigen.cxx:796
 TMatrixDEigen.cxx:797
 TMatrixDEigen.cxx:798
 TMatrixDEigen.cxx:799
 TMatrixDEigen.cxx:800
 TMatrixDEigen.cxx:801
 TMatrixDEigen.cxx:802
 TMatrixDEigen.cxx:803