ROOT logo
/**********************************************************************************************/
/* Fast symmetric matrix with dynamically expandable size.                                    */
/* Only part can be used for matrix operations. It is defined as:                             */ 
/* fNCols: rows built by constructor (GetSizeBooked)                                          */ 
/* fNRows: number of rows added dynamically (automatically added on assignment to row)        */ 
/*         GetNRowAdded                                                                       */ 
/* fNRowIndex: total size (fNCols+fNRows), GetSize                                            */ 
/* fRowLwb   : actual size to used for given operation, by default = total size, GetSizeUsed  */ 
/*                                                                                            */ 
/* Author: ruben.shahoyan@cern.ch                                                             */
/*                                                                                            */ 
/**********************************************************************************************/
#include <stdlib.h>
#include <stdio.h>
#include <iostream>
#include <float.h>
#include <string.h>
//
#include <TClass.h>
#include <TMath.h>
#include "AliSymMatrix.h"
#include "AliLog.h"
//

ClassImp(AliSymMatrix)


AliSymMatrix* AliSymMatrix::fgBuffer = 0; 
Int_t         AliSymMatrix::fgCopyCnt = 0; 
//___________________________________________________________
AliSymMatrix::AliSymMatrix() 
: fElems(0),fElemsAdd(0)
{
  // default constructor
  fSymmetric = kTRUE;
  fgCopyCnt++;
}

//___________________________________________________________
AliSymMatrix::AliSymMatrix(Int_t size)
  : AliMatrixSq(),fElems(0),fElemsAdd(0)
{
  //constructor for matrix with defined size
  fNrows = 0;
  fNrowIndex = fNcols = fRowLwb = size;
  fElems     = new Double_t[fNcols*(fNcols+1)/2];
  fSymmetric = kTRUE;
  Reset();
  fgCopyCnt++;
  //
}

//___________________________________________________________
AliSymMatrix::AliSymMatrix(const AliSymMatrix &src) 
  : AliMatrixSq(src),fElems(0),fElemsAdd(0)
{
  // copy constructor
  fNrowIndex = fNcols = src.GetSize();
  fNrows = 0;
  fRowLwb = src.GetSizeUsed();
  if (fNcols) {
    int nmainel = fNcols*(fNcols+1)/2;
    fElems     = new Double_t[nmainel];
    nmainel = src.fNcols*(src.fNcols+1)/2;
    memcpy(fElems,src.fElems,nmainel*sizeof(Double_t));
    if (src.GetSizeAdded()) { // transfer extra rows to main matrix
      Double_t *pnt = fElems + nmainel;
      int ncl = src.GetSizeBooked() + 1;
      for (int ir=0;ir<src.GetSizeAdded();ir++) {
	memcpy(pnt,src.fElemsAdd[ir],ncl*sizeof(Double_t));
	pnt += ncl;
	ncl++; 
      }
    }
  }
  else fElems = 0;
  fElemsAdd = 0;
  fgCopyCnt++;
  //
}

//___________________________________________________________
AliSymMatrix::~AliSymMatrix() 
{
  Clear();
  if (--fgCopyCnt < 1 && fgBuffer) {delete fgBuffer; fgBuffer = 0;}
}

//___________________________________________________________
AliSymMatrix&  AliSymMatrix::operator=(const AliSymMatrix& src)
{
  // assignment operator
  if (this != &src) {
    TObject::operator=(src);
    if (GetSizeBooked()!=src.GetSizeBooked() && GetSizeAdded()!=src.GetSizeAdded()) {
      // recreate the matrix
      if (fElems) delete[] fElems;
      for (int i=0;i<GetSizeAdded();i++) delete[] fElemsAdd[i]; 
      delete[] fElemsAdd;
      //
      fNrowIndex = src.GetSize(); 
      fNcols = src.GetSize();
      fNrows = 0;
      fRowLwb = src.GetSizeUsed();
      fElems     = new Double_t[GetSize()*(GetSize()+1)/2];
      int nmainel = src.GetSizeBooked()*(src.GetSizeBooked()+1);
      memcpy(fElems,src.fElems,nmainel*sizeof(Double_t));
      if (src.GetSizeAdded()) { // transfer extra rows to main matrix
	Double_t *pnt = fElems + nmainel;//*sizeof(Double_t);
	int ncl = src.GetSizeBooked() + 1;
	for (int ir=0;ir<src.GetSizeAdded();ir++) {
	  ncl += ir; 
	  memcpy(pnt,src.fElemsAdd[ir],ncl*sizeof(Double_t));
	  pnt += ncl;//*sizeof(Double_t);
	}
      }
      //
    }
    else {
      memcpy(fElems,src.fElems,GetSizeBooked()*(GetSizeBooked()+1)/2*sizeof(Double_t));
      int ncl = GetSizeBooked() + 1;
      for (int ir=0;ir<GetSizeAdded();ir++) { // dynamic rows
	ncl += ir; 
	memcpy(fElemsAdd[ir],src.fElemsAdd[ir],ncl*sizeof(Double_t));
      }
    }
  }
  //
  return *this;
}

//___________________________________________________________
AliSymMatrix& AliSymMatrix::operator+=(const AliSymMatrix& src)
{
  // add operator
  if (GetSizeUsed() != src.GetSizeUsed()) {
    AliError("Matrix sizes are different");
    return *this;
  }
  for (int i=0;i<GetSizeUsed();i++) for (int j=i;j<GetSizeUsed();j++) (*this)(j,i) += src(j,i);
  return *this;
}

//___________________________________________________________
void AliSymMatrix::Clear(Option_t*)
{
  // clear dynamic part
  if (fElems) {delete[] fElems; fElems = 0;}
  //  
  if (fElemsAdd) {
    for (int i=0;i<GetSizeAdded();i++) delete[] fElemsAdd[i]; 
    delete[] fElemsAdd;
    fElemsAdd = 0;
  }
  fNrowIndex = fNcols = fNrows = fRowLwb = 0;
  //
}

//___________________________________________________________
Float_t AliSymMatrix::GetDensity() const
{
  // get fraction of non-zero elements
  Int_t nel = 0;
  for (int i=GetSizeUsed();i--;) for (int j=i+1;j--;) if (!IsZero(GetEl(i,j))) nel++;
  return 2.*nel/( (GetSizeUsed()+1)*GetSizeUsed() );
}

//___________________________________________________________
void AliSymMatrix::Print(Option_t* option) const
{
  // print itself
  printf("Symmetric Matrix: Size = %d (%d rows added dynamically), %d used\n",GetSize(),GetSizeAdded(),GetSizeUsed());
  TString opt = option; opt.ToLower();
  if (opt.IsNull()) return;
  opt = "%"; opt += 1+int(TMath::Log10(double(GetSize()))); opt+="d|";
  for (Int_t i=0;i<GetSizeUsed();i++) {
    printf(opt,i);
    for (Int_t j=0;j<=i;j++) printf("%+.3e|",GetEl(i,j));
    printf("\n");
  }
}

//___________________________________________________________
void AliSymMatrix::MultiplyByVec(const Double_t *vecIn,Double_t *vecOut) const
{
  // fill vecOut by matrix*vecIn
  // vector should be of the same size as the matrix
  for (int i=GetSizeUsed();i--;) {
    vecOut[i] = 0.0;
    for (int j=GetSizeUsed();j--;) vecOut[i] += vecIn[j]*GetEl(i,j);
  }
  //
}

//___________________________________________________________
AliSymMatrix* AliSymMatrix::DecomposeChol() 
{
  // Return a matrix with Choleski decomposition
  // Adopted from Numerical Recipes in C, ch.2-9, http://www.nr.com
  // consturcts Cholesky decomposition of SYMMETRIC and
  // POSITIVELY-DEFINED matrix a (a=L*Lt)
  // Only upper triangle of the matrix has to be filled.
  // In opposite to function from the book, the matrix is modified:
  // lower triangle and diagonal are refilled.
  //
  if (!fgBuffer || fgBuffer->GetSizeUsed()!=GetSizeUsed()) {
    delete fgBuffer; 
    fgBuffer = new AliSymMatrix(*this);
  }
  else (*fgBuffer) = *this;
  //
  AliSymMatrix& mchol = *fgBuffer;
  //
  for (int i=0;i<GetSizeUsed();i++) {
    Double_t *rowi = mchol.GetRow(i);
    for (int j=i;j<GetSizeUsed();j++) {
      Double_t *rowj = mchol.GetRow(j);
      double sum = rowj[i];
      for (int k=i-1;k>=0;k--) if (rowi[k]&&rowj[k]) sum -= rowi[k]*rowj[k];
      if (i == j) {
	if (sum <= 0.0) { // not positive-definite
	  AliInfo(Form("The matrix is not positive definite [%e]\n"
		       "Choleski decomposition is not possible",sum));
	  Print("l");
	  return 0;
	}
	rowi[i] = TMath::Sqrt(sum);
	//
      } else rowj[i] = sum/rowi[i];
    }
  }
  return fgBuffer;
}

//___________________________________________________________
Bool_t AliSymMatrix::InvertChol() 
{
  // Invert matrix using Choleski decomposition
  //
  AliSymMatrix* mchol = DecomposeChol();
  if (!mchol) {
    AliInfo("Failed to invert the matrix");
    return kFALSE;
  }
  //
  InvertChol(mchol);
  return kTRUE;
  //
}

//___________________________________________________________
void AliSymMatrix::InvertChol(AliSymMatrix* pmchol) 
{
  // Invert matrix using Choleski decomposition, provided the Cholseki's L matrix
  //
  Double_t sum;
  AliSymMatrix& mchol = *pmchol;
  //
  // Invert decomposed triangular L matrix (Lower triangle is filled)
  for (int i=0;i<GetSizeUsed();i++) { 
    mchol(i,i) =  1.0/mchol(i,i);
    for (int j=i+1;j<GetSizeUsed();j++) { 
      Double_t *rowj = mchol.GetRow(j);
      sum = 0.0; 
      for (int k=i;k<j;k++) if (rowj[k]) { 
	double &mki = mchol(k,i); if (mki) sum -= rowj[k]*mki;
      }
      rowj[i] = sum/rowj[j];
    } 
  }
  //
  // take product of the inverted Choleski L matrix with its transposed
  for (int i=GetSizeUsed();i--;) {
    for (int j=i+1;j--;) {
      sum = 0;
      for (int k=i;k<GetSizeUsed();k++) {
	double &mik = mchol(i,k); 
	if (mik) {
	  double &mjk = mchol(j,k);
	  if (mjk) sum += mik*mjk;
	}
      }
      (*this)(j,i) = sum;
    }
  }
  //
}


//___________________________________________________________
Bool_t AliSymMatrix::SolveChol(Double_t *b, Bool_t invert) 
{
  // Adopted from Numerical Recipes in C, ch.2-9, http://www.nr.com
  // Solves the set of n linear equations A x = b, 
  // where a is a positive-definite symmetric matrix. 
  // a[1..n][1..n] is the output of the routine CholDecomposw. 
  // Only the lower triangle of a is accessed. b[1..n] is input as the 
  // right-hand side vector. The solution vector is returned in b[1..n].
  //
  Int_t i,k;
  Double_t sum;
  //
  AliSymMatrix *pmchol = DecomposeChol();
  if (!pmchol) {
    AliInfo("SolveChol failed");
    //    Print("l");
    return kFALSE;
  }
  AliSymMatrix& mchol = *pmchol;
  //
  for (i=0;i<GetSizeUsed();i++) {
    Double_t *rowi = mchol.GetRow(i);
    for (sum=b[i],k=i-1;k>=0;k--) if (rowi[k]&&b[k]) sum -= rowi[k]*b[k];
    b[i]=sum/rowi[i];
  }
  //
  for (i=GetSizeUsed()-1;i>=0;i--) {
    for (sum=b[i],k=i+1;k<GetSizeUsed();k++) if (b[k]) {
      double &mki=mchol(k,i); if (mki) sum -= mki*b[k];
    }
    b[i]=sum/mchol(i,i);
  }
  //
  if (invert) InvertChol(pmchol);
  return kTRUE;
  //
}

//___________________________________________________________
Bool_t AliSymMatrix::SolveChol(TVectorD &b, Bool_t invert) 
{
  return SolveChol((Double_t*)b.GetMatrixArray(),invert);
}


//___________________________________________________________
Bool_t AliSymMatrix::SolveChol(Double_t *brhs, Double_t *bsol,Bool_t invert) 
{
  memcpy(bsol,brhs,GetSizeUsed()*sizeof(Double_t));
  return SolveChol(bsol,invert);  
}

//___________________________________________________________
Bool_t AliSymMatrix::SolveChol(const TVectorD &brhs, TVectorD &bsol,Bool_t invert) 
{
  bsol = brhs;
  return SolveChol(bsol,invert);
}

//___________________________________________________________
void AliSymMatrix::AddRows(int nrows)
{
  // add empty rows
  if (nrows<1) return;
  Double_t **pnew = new Double_t*[nrows+fNrows];
  for (int ir=0;ir<fNrows;ir++) pnew[ir] = fElemsAdd[ir]; // copy old extra rows
  for (int ir=0;ir<nrows;ir++) {
    int ncl = GetSize()+1;
    pnew[fNrows] = new Double_t[ncl];
    memset(pnew[fNrows],0,ncl*sizeof(Double_t));
    fNrows++;
    fNrowIndex++;
    fRowLwb++;
  }
  delete[] fElemsAdd;
  fElemsAdd = pnew;
  //
}

//___________________________________________________________
void AliSymMatrix::Reset()
{
  // if additional rows exist, regularize it
  if (fElemsAdd) {
    delete[] fElems;
    for (int i=0;i<fNrows;i++) delete[] fElemsAdd[i]; 
    delete[] fElemsAdd; fElemsAdd = 0;
    fNcols = fRowLwb = fNrowIndex;
    fElems = new Double_t[GetSize()*(GetSize()+1)/2];
    fNrows = 0;
  }
  if (fElems) memset(fElems,0,GetSize()*(GetSize()+1)/2*sizeof(Double_t));
  //
}

//___________________________________________________________
/*
void AliSymMatrix::AddToRow(Int_t r, Double_t *valc,Int_t *indc,Int_t n)
{
  //   for (int i=n;i--;) {
  //     (*this)(indc[i],r) += valc[i];
  //   }
  //   return;

  double *row;
  if (r>=fNrowIndex) {
    AddRows(r-fNrowIndex+1); 
    row = &((fElemsAdd[r-fNcols])[0]);
  }
  else row = &fElems[GetIndex(r,0)];
  //
  int nadd = 0;
  for (int i=n;i--;) {
    if (indc[i]>r) continue;
    row[indc[i]] += valc[i];
    nadd++;
  }
  if (nadd == n) return;
  //
  // add to col>row
  for (int i=n;i--;) {
    if (indc[i]>r) (*this)(indc[i],r) += valc[i];
  }
  //
}
*/

//___________________________________________________________
Double_t* AliSymMatrix::GetRow(Int_t r)
{
  // get pointer on the row
  if (r>=GetSize()) {
    int nn = GetSize();
    AddRows(r-GetSize()+1); 
    AliDebug(2,Form("create %d of %d\n",r, nn));
    return &((fElemsAdd[r-GetSizeBooked()])[0]);
  }
  else return &fElems[GetIndex(r,0)];
}


//___________________________________________________________
int AliSymMatrix::SolveSpmInv(double *vecB, Bool_t stabilize)
{
  //   Solution a la MP1: gaussian eliminations
  ///  Obtain solution of a system of linear equations with symmetric matrix 
  ///  and the inverse (using 'singular-value friendly' GAUSS pivot)
  //

  Int_t nRank = 0;
  int iPivot;
  double vPivot = 0.;
  double eps = 1e-14;
  int nGlo = GetSizeUsed();
  bool   *bUnUsed = new bool[nGlo];
  double *rowMax,*colMax=0;
  rowMax  = new double[nGlo];
  //
  if (stabilize) {
    colMax   = new double[nGlo];
    for (Int_t i=nGlo; i--;) rowMax[i] = colMax[i] = 0.0;
    for (Int_t i=nGlo; i--;) for (Int_t j=i+1;j--;) { 
	double vl = TMath::Abs(Query(i,j));
	if (IsZero(vl)) continue;
	if (vl > rowMax[i]) rowMax[i] = vl; // Max elemt of row i
	if (vl > colMax[j]) colMax[j] = vl; // Max elemt of column j
	if (i==j) continue;
	if (vl > rowMax[j]) rowMax[j] = vl; // Max elemt of row j
	if (vl > colMax[i]) colMax[i] = vl; // Max elemt of column i
      }
    //
    for (Int_t i=nGlo; i--;) {
      if (!IsZero(rowMax[i])) rowMax[i] = 1./rowMax[i]; // Max elemt of row i
      if (!IsZero(colMax[i])) colMax[i] = 1./colMax[i]; // Max elemt of column i
    }
    //
  }
  //
  for (Int_t i=nGlo; i--;) bUnUsed[i] = true;
  //  
  if (!fgBuffer || fgBuffer->GetSizeUsed()!=GetSizeUsed()) {
    delete fgBuffer; 
    fgBuffer = new AliSymMatrix(*this);
  }
  else (*fgBuffer) = *this;
  //
  if (stabilize) for (int i=0;i<nGlo; i++) { // Small loop for matrix equilibration (gives a better conditioning) 
      for (int j=0;j<=i; j++) {
	double vl = Query(i,j);
	if (!IsZero(vl)) SetEl(i,j, TMath::Sqrt(rowMax[i])*vl*TMath::Sqrt(colMax[j]) ); // Equilibrate the V matrix
      }
      for (int j=i+1;j<nGlo;j++) {
	double vl = Query(j,i);
	if (!IsZero(vl)) fgBuffer->SetEl(j,i,TMath::Sqrt(rowMax[i])*vl*TMath::Sqrt(colMax[j]) ); // Equilibrate the V matrix
      }
    }
  //
  for (Int_t j=nGlo; j--;) fgBuffer->DiagElem(j) = TMath::Abs(QueryDiag(j)); // save diagonal elem absolute values 
  //
  for (Int_t i=0; i<nGlo; i++) {
    vPivot = 0.0;
    iPivot = -1;
    //
    for (Int_t j=0; j<nGlo; j++) { // First look for the pivot, ie max unused diagonal element       
      double vl;
      if (bUnUsed[j] && (TMath::Abs(vl=QueryDiag(j))>TMath::Max(TMath::Abs(vPivot),eps*fgBuffer->QueryDiag(j)))) {    
	vPivot = vl;
	iPivot = j;
      }
    }
    //
    if (iPivot >= 0) {   // pivot found          
      nRank++;
      bUnUsed[iPivot] = false; // This value is used
      vPivot = 1.0/vPivot;
      DiagElem(iPivot) = -vPivot; // Replace pivot by its inverse
      //
      for (Int_t j=0; j<nGlo; j++) {      
	for (Int_t jj=0; jj<nGlo; jj++) {  
	  if (j != iPivot && jj != iPivot) {// Other elements (!!! do them first as you use old matV[k][j]'s !!!)	  
	    double &r = j>=jj ? (*this)(j,jj) : (*fgBuffer)(jj,j);
	    r -= vPivot* ( j>iPivot  ? Query(j,iPivot)  : fgBuffer->Query(iPivot,j) )
	      *          ( iPivot>jj ? Query(iPivot,jj) : fgBuffer->Query(jj,iPivot));
	  }
	}
      }
      //
      for (Int_t j=0; j<nGlo; j++) if (j != iPivot) { // Pivot row or column elements 
	  (*this)(j,iPivot)     *= vPivot;
	  (*fgBuffer)(iPivot,j) *= vPivot;
	}
      //
    }
    else {  // No more pivot value (clear those elements)
      for (Int_t j=0; j<nGlo; j++) {
	if (bUnUsed[j]) {
	  vecB[j] = 0.0;
	  for (Int_t k=0; k<nGlo; k++) {
	    (*this)(j,k) = 0.;
	    if (j!=k) (*fgBuffer)(j,k) = 0;
	  }
	}
      }
      break;  // No more pivots anyway, stop here
    }
  }
  //
  if (stabilize) for (Int_t i=0; i<nGlo; i++) for (Int_t j=0; j<nGlo; j++) {
	double vl = TMath::Sqrt(colMax[i])*TMath::Sqrt(rowMax[j]); // Correct matrix V
	if (i>=j) (*this)(i,j) *= vl;
	else      (*fgBuffer)(j,i) *= vl;
      }
  //
  for (Int_t j=0; j<nGlo; j++) {
    rowMax[j] = 0.0;
    for (Int_t jj=0; jj<nGlo; jj++) { // Reverse matrix elements
      double vl;
      if (j>=jj) vl = (*this)(j,jj)     = -Query(j,jj);
      else       vl = (*fgBuffer)(j,jj) = -fgBuffer->Query(j,jj);
      rowMax[j] += vl*vecB[jj];
    }		
  }
  
  for (Int_t j=0; j<nGlo; j++) {
    vecB[j] = rowMax[j]; // The final result
  }
  //
  delete [] bUnUsed;
  delete [] rowMax;
  if (stabilize) delete [] colMax;

  return nRank;
}


 AliSymMatrix.cxx:1
 AliSymMatrix.cxx:2
 AliSymMatrix.cxx:3
 AliSymMatrix.cxx:4
 AliSymMatrix.cxx:5
 AliSymMatrix.cxx:6
 AliSymMatrix.cxx:7
 AliSymMatrix.cxx:8
 AliSymMatrix.cxx:9
 AliSymMatrix.cxx:10
 AliSymMatrix.cxx:11
 AliSymMatrix.cxx:12
 AliSymMatrix.cxx:13
 AliSymMatrix.cxx:14
 AliSymMatrix.cxx:15
 AliSymMatrix.cxx:16
 AliSymMatrix.cxx:17
 AliSymMatrix.cxx:18
 AliSymMatrix.cxx:19
 AliSymMatrix.cxx:20
 AliSymMatrix.cxx:21
 AliSymMatrix.cxx:22
 AliSymMatrix.cxx:23
 AliSymMatrix.cxx:24
 AliSymMatrix.cxx:25
 AliSymMatrix.cxx:26
 AliSymMatrix.cxx:27
 AliSymMatrix.cxx:28
 AliSymMatrix.cxx:29
 AliSymMatrix.cxx:30
 AliSymMatrix.cxx:31
 AliSymMatrix.cxx:32
 AliSymMatrix.cxx:33
 AliSymMatrix.cxx:34
 AliSymMatrix.cxx:35
 AliSymMatrix.cxx:36
 AliSymMatrix.cxx:37
 AliSymMatrix.cxx:38
 AliSymMatrix.cxx:39
 AliSymMatrix.cxx:40
 AliSymMatrix.cxx:41
 AliSymMatrix.cxx:42
 AliSymMatrix.cxx:43
 AliSymMatrix.cxx:44
 AliSymMatrix.cxx:45
 AliSymMatrix.cxx:46
 AliSymMatrix.cxx:47
 AliSymMatrix.cxx:48
 AliSymMatrix.cxx:49
 AliSymMatrix.cxx:50
 AliSymMatrix.cxx:51
 AliSymMatrix.cxx:52
 AliSymMatrix.cxx:53
 AliSymMatrix.cxx:54
 AliSymMatrix.cxx:55
 AliSymMatrix.cxx:56
 AliSymMatrix.cxx:57
 AliSymMatrix.cxx:58
 AliSymMatrix.cxx:59
 AliSymMatrix.cxx:60
 AliSymMatrix.cxx:61
 AliSymMatrix.cxx:62
 AliSymMatrix.cxx:63
 AliSymMatrix.cxx:64
 AliSymMatrix.cxx:65
 AliSymMatrix.cxx:66
 AliSymMatrix.cxx:67
 AliSymMatrix.cxx:68
 AliSymMatrix.cxx:69
 AliSymMatrix.cxx:70
 AliSymMatrix.cxx:71
 AliSymMatrix.cxx:72
 AliSymMatrix.cxx:73
 AliSymMatrix.cxx:74
 AliSymMatrix.cxx:75
 AliSymMatrix.cxx:76
 AliSymMatrix.cxx:77
 AliSymMatrix.cxx:78
 AliSymMatrix.cxx:79
 AliSymMatrix.cxx:80
 AliSymMatrix.cxx:81
 AliSymMatrix.cxx:82
 AliSymMatrix.cxx:83
 AliSymMatrix.cxx:84
 AliSymMatrix.cxx:85
 AliSymMatrix.cxx:86
 AliSymMatrix.cxx:87
 AliSymMatrix.cxx:88
 AliSymMatrix.cxx:89
 AliSymMatrix.cxx:90
 AliSymMatrix.cxx:91
 AliSymMatrix.cxx:92
 AliSymMatrix.cxx:93
 AliSymMatrix.cxx:94
 AliSymMatrix.cxx:95
 AliSymMatrix.cxx:96
 AliSymMatrix.cxx:97
 AliSymMatrix.cxx:98
 AliSymMatrix.cxx:99
 AliSymMatrix.cxx:100
 AliSymMatrix.cxx:101
 AliSymMatrix.cxx:102
 AliSymMatrix.cxx:103
 AliSymMatrix.cxx:104
 AliSymMatrix.cxx:105
 AliSymMatrix.cxx:106
 AliSymMatrix.cxx:107
 AliSymMatrix.cxx:108
 AliSymMatrix.cxx:109
 AliSymMatrix.cxx:110
 AliSymMatrix.cxx:111
 AliSymMatrix.cxx:112
 AliSymMatrix.cxx:113
 AliSymMatrix.cxx:114
 AliSymMatrix.cxx:115
 AliSymMatrix.cxx:116
 AliSymMatrix.cxx:117
 AliSymMatrix.cxx:118
 AliSymMatrix.cxx:119
 AliSymMatrix.cxx:120
 AliSymMatrix.cxx:121
 AliSymMatrix.cxx:122
 AliSymMatrix.cxx:123
 AliSymMatrix.cxx:124
 AliSymMatrix.cxx:125
 AliSymMatrix.cxx:126
 AliSymMatrix.cxx:127
 AliSymMatrix.cxx:128
 AliSymMatrix.cxx:129
 AliSymMatrix.cxx:130
 AliSymMatrix.cxx:131
 AliSymMatrix.cxx:132
 AliSymMatrix.cxx:133
 AliSymMatrix.cxx:134
 AliSymMatrix.cxx:135
 AliSymMatrix.cxx:136
 AliSymMatrix.cxx:137
 AliSymMatrix.cxx:138
 AliSymMatrix.cxx:139
 AliSymMatrix.cxx:140
 AliSymMatrix.cxx:141
 AliSymMatrix.cxx:142
 AliSymMatrix.cxx:143
 AliSymMatrix.cxx:144
 AliSymMatrix.cxx:145
 AliSymMatrix.cxx:146
 AliSymMatrix.cxx:147
 AliSymMatrix.cxx:148
 AliSymMatrix.cxx:149
 AliSymMatrix.cxx:150
 AliSymMatrix.cxx:151
 AliSymMatrix.cxx:152
 AliSymMatrix.cxx:153
 AliSymMatrix.cxx:154
 AliSymMatrix.cxx:155
 AliSymMatrix.cxx:156
 AliSymMatrix.cxx:157
 AliSymMatrix.cxx:158
 AliSymMatrix.cxx:159
 AliSymMatrix.cxx:160
 AliSymMatrix.cxx:161
 AliSymMatrix.cxx:162
 AliSymMatrix.cxx:163
 AliSymMatrix.cxx:164
 AliSymMatrix.cxx:165
 AliSymMatrix.cxx:166
 AliSymMatrix.cxx:167
 AliSymMatrix.cxx:168
 AliSymMatrix.cxx:169
 AliSymMatrix.cxx:170
 AliSymMatrix.cxx:171
 AliSymMatrix.cxx:172
 AliSymMatrix.cxx:173
 AliSymMatrix.cxx:174
 AliSymMatrix.cxx:175
 AliSymMatrix.cxx:176
 AliSymMatrix.cxx:177
 AliSymMatrix.cxx:178
 AliSymMatrix.cxx:179
 AliSymMatrix.cxx:180
 AliSymMatrix.cxx:181
 AliSymMatrix.cxx:182
 AliSymMatrix.cxx:183
 AliSymMatrix.cxx:184
 AliSymMatrix.cxx:185
 AliSymMatrix.cxx:186
 AliSymMatrix.cxx:187
 AliSymMatrix.cxx:188
 AliSymMatrix.cxx:189
 AliSymMatrix.cxx:190
 AliSymMatrix.cxx:191
 AliSymMatrix.cxx:192
 AliSymMatrix.cxx:193
 AliSymMatrix.cxx:194
 AliSymMatrix.cxx:195
 AliSymMatrix.cxx:196
 AliSymMatrix.cxx:197
 AliSymMatrix.cxx:198
 AliSymMatrix.cxx:199
 AliSymMatrix.cxx:200
 AliSymMatrix.cxx:201
 AliSymMatrix.cxx:202
 AliSymMatrix.cxx:203
 AliSymMatrix.cxx:204
 AliSymMatrix.cxx:205
 AliSymMatrix.cxx:206
 AliSymMatrix.cxx:207
 AliSymMatrix.cxx:208
 AliSymMatrix.cxx:209
 AliSymMatrix.cxx:210
 AliSymMatrix.cxx:211
 AliSymMatrix.cxx:212
 AliSymMatrix.cxx:213
 AliSymMatrix.cxx:214
 AliSymMatrix.cxx:215
 AliSymMatrix.cxx:216
 AliSymMatrix.cxx:217
 AliSymMatrix.cxx:218
 AliSymMatrix.cxx:219
 AliSymMatrix.cxx:220
 AliSymMatrix.cxx:221
 AliSymMatrix.cxx:222
 AliSymMatrix.cxx:223
 AliSymMatrix.cxx:224
 AliSymMatrix.cxx:225
 AliSymMatrix.cxx:226
 AliSymMatrix.cxx:227
 AliSymMatrix.cxx:228
 AliSymMatrix.cxx:229
 AliSymMatrix.cxx:230
 AliSymMatrix.cxx:231
 AliSymMatrix.cxx:232
 AliSymMatrix.cxx:233
 AliSymMatrix.cxx:234
 AliSymMatrix.cxx:235
 AliSymMatrix.cxx:236
 AliSymMatrix.cxx:237
 AliSymMatrix.cxx:238
 AliSymMatrix.cxx:239
 AliSymMatrix.cxx:240
 AliSymMatrix.cxx:241
 AliSymMatrix.cxx:242
 AliSymMatrix.cxx:243
 AliSymMatrix.cxx:244
 AliSymMatrix.cxx:245
 AliSymMatrix.cxx:246
 AliSymMatrix.cxx:247
 AliSymMatrix.cxx:248
 AliSymMatrix.cxx:249
 AliSymMatrix.cxx:250
 AliSymMatrix.cxx:251
 AliSymMatrix.cxx:252
 AliSymMatrix.cxx:253
 AliSymMatrix.cxx:254
 AliSymMatrix.cxx:255
 AliSymMatrix.cxx:256
 AliSymMatrix.cxx:257
 AliSymMatrix.cxx:258
 AliSymMatrix.cxx:259
 AliSymMatrix.cxx:260
 AliSymMatrix.cxx:261
 AliSymMatrix.cxx:262
 AliSymMatrix.cxx:263
 AliSymMatrix.cxx:264
 AliSymMatrix.cxx:265
 AliSymMatrix.cxx:266
 AliSymMatrix.cxx:267
 AliSymMatrix.cxx:268
 AliSymMatrix.cxx:269
 AliSymMatrix.cxx:270
 AliSymMatrix.cxx:271
 AliSymMatrix.cxx:272
 AliSymMatrix.cxx:273
 AliSymMatrix.cxx:274
 AliSymMatrix.cxx:275
 AliSymMatrix.cxx:276
 AliSymMatrix.cxx:277
 AliSymMatrix.cxx:278
 AliSymMatrix.cxx:279
 AliSymMatrix.cxx:280
 AliSymMatrix.cxx:281
 AliSymMatrix.cxx:282
 AliSymMatrix.cxx:283
 AliSymMatrix.cxx:284
 AliSymMatrix.cxx:285
 AliSymMatrix.cxx:286
 AliSymMatrix.cxx:287
 AliSymMatrix.cxx:288
 AliSymMatrix.cxx:289
 AliSymMatrix.cxx:290
 AliSymMatrix.cxx:291
 AliSymMatrix.cxx:292
 AliSymMatrix.cxx:293
 AliSymMatrix.cxx:294
 AliSymMatrix.cxx:295
 AliSymMatrix.cxx:296
 AliSymMatrix.cxx:297
 AliSymMatrix.cxx:298
 AliSymMatrix.cxx:299
 AliSymMatrix.cxx:300
 AliSymMatrix.cxx:301
 AliSymMatrix.cxx:302
 AliSymMatrix.cxx:303
 AliSymMatrix.cxx:304
 AliSymMatrix.cxx:305
 AliSymMatrix.cxx:306
 AliSymMatrix.cxx:307
 AliSymMatrix.cxx:308
 AliSymMatrix.cxx:309
 AliSymMatrix.cxx:310
 AliSymMatrix.cxx:311
 AliSymMatrix.cxx:312
 AliSymMatrix.cxx:313
 AliSymMatrix.cxx:314
 AliSymMatrix.cxx:315
 AliSymMatrix.cxx:316
 AliSymMatrix.cxx:317
 AliSymMatrix.cxx:318
 AliSymMatrix.cxx:319
 AliSymMatrix.cxx:320
 AliSymMatrix.cxx:321
 AliSymMatrix.cxx:322
 AliSymMatrix.cxx:323
 AliSymMatrix.cxx:324
 AliSymMatrix.cxx:325
 AliSymMatrix.cxx:326
 AliSymMatrix.cxx:327
 AliSymMatrix.cxx:328
 AliSymMatrix.cxx:329
 AliSymMatrix.cxx:330
 AliSymMatrix.cxx:331
 AliSymMatrix.cxx:332
 AliSymMatrix.cxx:333
 AliSymMatrix.cxx:334
 AliSymMatrix.cxx:335
 AliSymMatrix.cxx:336
 AliSymMatrix.cxx:337
 AliSymMatrix.cxx:338
 AliSymMatrix.cxx:339
 AliSymMatrix.cxx:340
 AliSymMatrix.cxx:341
 AliSymMatrix.cxx:342
 AliSymMatrix.cxx:343
 AliSymMatrix.cxx:344
 AliSymMatrix.cxx:345
 AliSymMatrix.cxx:346
 AliSymMatrix.cxx:347
 AliSymMatrix.cxx:348
 AliSymMatrix.cxx:349
 AliSymMatrix.cxx:350
 AliSymMatrix.cxx:351
 AliSymMatrix.cxx:352
 AliSymMatrix.cxx:353
 AliSymMatrix.cxx:354
 AliSymMatrix.cxx:355
 AliSymMatrix.cxx:356
 AliSymMatrix.cxx:357
 AliSymMatrix.cxx:358
 AliSymMatrix.cxx:359
 AliSymMatrix.cxx:360
 AliSymMatrix.cxx:361
 AliSymMatrix.cxx:362
 AliSymMatrix.cxx:363
 AliSymMatrix.cxx:364
 AliSymMatrix.cxx:365
 AliSymMatrix.cxx:366
 AliSymMatrix.cxx:367
 AliSymMatrix.cxx:368
 AliSymMatrix.cxx:369
 AliSymMatrix.cxx:370
 AliSymMatrix.cxx:371
 AliSymMatrix.cxx:372
 AliSymMatrix.cxx:373
 AliSymMatrix.cxx:374
 AliSymMatrix.cxx:375
 AliSymMatrix.cxx:376
 AliSymMatrix.cxx:377
 AliSymMatrix.cxx:378
 AliSymMatrix.cxx:379
 AliSymMatrix.cxx:380
 AliSymMatrix.cxx:381
 AliSymMatrix.cxx:382
 AliSymMatrix.cxx:383
 AliSymMatrix.cxx:384
 AliSymMatrix.cxx:385
 AliSymMatrix.cxx:386
 AliSymMatrix.cxx:387
 AliSymMatrix.cxx:388
 AliSymMatrix.cxx:389
 AliSymMatrix.cxx:390
 AliSymMatrix.cxx:391
 AliSymMatrix.cxx:392
 AliSymMatrix.cxx:393
 AliSymMatrix.cxx:394
 AliSymMatrix.cxx:395
 AliSymMatrix.cxx:396
 AliSymMatrix.cxx:397
 AliSymMatrix.cxx:398
 AliSymMatrix.cxx:399
 AliSymMatrix.cxx:400
 AliSymMatrix.cxx:401
 AliSymMatrix.cxx:402
 AliSymMatrix.cxx:403
 AliSymMatrix.cxx:404
 AliSymMatrix.cxx:405
 AliSymMatrix.cxx:406
 AliSymMatrix.cxx:407
 AliSymMatrix.cxx:408
 AliSymMatrix.cxx:409
 AliSymMatrix.cxx:410
 AliSymMatrix.cxx:411
 AliSymMatrix.cxx:412
 AliSymMatrix.cxx:413
 AliSymMatrix.cxx:414
 AliSymMatrix.cxx:415
 AliSymMatrix.cxx:416
 AliSymMatrix.cxx:417
 AliSymMatrix.cxx:418
 AliSymMatrix.cxx:419
 AliSymMatrix.cxx:420
 AliSymMatrix.cxx:421
 AliSymMatrix.cxx:422
 AliSymMatrix.cxx:423
 AliSymMatrix.cxx:424
 AliSymMatrix.cxx:425
 AliSymMatrix.cxx:426
 AliSymMatrix.cxx:427
 AliSymMatrix.cxx:428
 AliSymMatrix.cxx:429
 AliSymMatrix.cxx:430
 AliSymMatrix.cxx:431
 AliSymMatrix.cxx:432
 AliSymMatrix.cxx:433
 AliSymMatrix.cxx:434
 AliSymMatrix.cxx:435
 AliSymMatrix.cxx:436
 AliSymMatrix.cxx:437
 AliSymMatrix.cxx:438
 AliSymMatrix.cxx:439
 AliSymMatrix.cxx:440
 AliSymMatrix.cxx:441
 AliSymMatrix.cxx:442
 AliSymMatrix.cxx:443
 AliSymMatrix.cxx:444
 AliSymMatrix.cxx:445
 AliSymMatrix.cxx:446
 AliSymMatrix.cxx:447
 AliSymMatrix.cxx:448
 AliSymMatrix.cxx:449
 AliSymMatrix.cxx:450
 AliSymMatrix.cxx:451
 AliSymMatrix.cxx:452
 AliSymMatrix.cxx:453
 AliSymMatrix.cxx:454
 AliSymMatrix.cxx:455
 AliSymMatrix.cxx:456
 AliSymMatrix.cxx:457
 AliSymMatrix.cxx:458
 AliSymMatrix.cxx:459
 AliSymMatrix.cxx:460
 AliSymMatrix.cxx:461
 AliSymMatrix.cxx:462
 AliSymMatrix.cxx:463
 AliSymMatrix.cxx:464
 AliSymMatrix.cxx:465
 AliSymMatrix.cxx:466
 AliSymMatrix.cxx:467
 AliSymMatrix.cxx:468
 AliSymMatrix.cxx:469
 AliSymMatrix.cxx:470
 AliSymMatrix.cxx:471
 AliSymMatrix.cxx:472
 AliSymMatrix.cxx:473
 AliSymMatrix.cxx:474
 AliSymMatrix.cxx:475
 AliSymMatrix.cxx:476
 AliSymMatrix.cxx:477
 AliSymMatrix.cxx:478
 AliSymMatrix.cxx:479
 AliSymMatrix.cxx:480
 AliSymMatrix.cxx:481
 AliSymMatrix.cxx:482
 AliSymMatrix.cxx:483
 AliSymMatrix.cxx:484
 AliSymMatrix.cxx:485
 AliSymMatrix.cxx:486
 AliSymMatrix.cxx:487
 AliSymMatrix.cxx:488
 AliSymMatrix.cxx:489
 AliSymMatrix.cxx:490
 AliSymMatrix.cxx:491
 AliSymMatrix.cxx:492
 AliSymMatrix.cxx:493
 AliSymMatrix.cxx:494
 AliSymMatrix.cxx:495
 AliSymMatrix.cxx:496
 AliSymMatrix.cxx:497
 AliSymMatrix.cxx:498
 AliSymMatrix.cxx:499
 AliSymMatrix.cxx:500
 AliSymMatrix.cxx:501
 AliSymMatrix.cxx:502
 AliSymMatrix.cxx:503
 AliSymMatrix.cxx:504
 AliSymMatrix.cxx:505
 AliSymMatrix.cxx:506
 AliSymMatrix.cxx:507
 AliSymMatrix.cxx:508
 AliSymMatrix.cxx:509
 AliSymMatrix.cxx:510
 AliSymMatrix.cxx:511
 AliSymMatrix.cxx:512
 AliSymMatrix.cxx:513
 AliSymMatrix.cxx:514
 AliSymMatrix.cxx:515
 AliSymMatrix.cxx:516
 AliSymMatrix.cxx:517
 AliSymMatrix.cxx:518
 AliSymMatrix.cxx:519
 AliSymMatrix.cxx:520
 AliSymMatrix.cxx:521
 AliSymMatrix.cxx:522
 AliSymMatrix.cxx:523
 AliSymMatrix.cxx:524
 AliSymMatrix.cxx:525
 AliSymMatrix.cxx:526
 AliSymMatrix.cxx:527
 AliSymMatrix.cxx:528
 AliSymMatrix.cxx:529
 AliSymMatrix.cxx:530
 AliSymMatrix.cxx:531
 AliSymMatrix.cxx:532
 AliSymMatrix.cxx:533
 AliSymMatrix.cxx:534
 AliSymMatrix.cxx:535
 AliSymMatrix.cxx:536
 AliSymMatrix.cxx:537
 AliSymMatrix.cxx:538
 AliSymMatrix.cxx:539
 AliSymMatrix.cxx:540
 AliSymMatrix.cxx:541
 AliSymMatrix.cxx:542
 AliSymMatrix.cxx:543
 AliSymMatrix.cxx:544
 AliSymMatrix.cxx:545
 AliSymMatrix.cxx:546
 AliSymMatrix.cxx:547
 AliSymMatrix.cxx:548
 AliSymMatrix.cxx:549
 AliSymMatrix.cxx:550
 AliSymMatrix.cxx:551
 AliSymMatrix.cxx:552
 AliSymMatrix.cxx:553
 AliSymMatrix.cxx:554
 AliSymMatrix.cxx:555
 AliSymMatrix.cxx:556
 AliSymMatrix.cxx:557
 AliSymMatrix.cxx:558
 AliSymMatrix.cxx:559
 AliSymMatrix.cxx:560
 AliSymMatrix.cxx:561
 AliSymMatrix.cxx:562
 AliSymMatrix.cxx:563
 AliSymMatrix.cxx:564
 AliSymMatrix.cxx:565