ROOT logo
/**********************************************************************************************/
/* General class for solving large system of linear equations                                 */
/* Includes MINRES, FGMRES methods as well as a few precondiotiong methods                    */
/*                                                                                            */ 
/* Author: ruben.shahoyan@cern.ch                                                             */
/*                                                                                            */ 
/**********************************************************************************************/

#include "AliMinResSolve.h"
#include "AliLog.h"
#include "AliMatrixSq.h"
#include "AliMatrixSparse.h"
#include "AliSymBDMatrix.h"
#include <TStopwatch.h>
#include <float.h>
#include <TMath.h>

ClassImp(AliMinResSolve)

//______________________________________________________
AliMinResSolve::AliMinResSolve() : 
fSize(0),fPrecon(0),fMatrix(0),fRHS(0),
  fPVecY(0),fPVecR1(0),fPVecR2(0),fPVecV(0),fPVecW(0),fPVecW1(0),fPVecW2(0),
  fPvv(0),fPvz(0),fPhh(0),
  fDiagLU(0),fMatL(0),fMatU(0),fMatBD(0)
{
  // default constructor
}

//______________________________________________________
AliMinResSolve::AliMinResSolve(const AliMinResSolve& src) : 
  TObject(src),
  fSize(src.fSize),fPrecon(src.fPrecon),fMatrix(src.fMatrix),fRHS(src.fRHS),
  fPVecY(0),fPVecR1(0),fPVecR2(0),fPVecV(0),fPVecW(0),fPVecW1(0),fPVecW2(0),
  fPvv(0),fPvz(0),fPhh(0),
  fDiagLU(0),fMatL(0),fMatU(0),fMatBD(0)
{
  // copy constructor
}

//______________________________________________________
AliMinResSolve::AliMinResSolve(const AliMatrixSq *mat, const TVectorD* rhs) :
  fSize(mat->GetSize()),fPrecon(0),fMatrix((AliMatrixSq*)mat),fRHS((double*)rhs->GetMatrixArray()),
  fPVecY(0),fPVecR1(0),fPVecR2(0),fPVecV(0),fPVecW(0),fPVecW1(0),fPVecW2(0),
  fPvv(0),fPvz(0),fPhh(0),
  fDiagLU(0),fMatL(0),fMatU(0),fMatBD(0)
{
  // copy accepting equation
}

//______________________________________________________
AliMinResSolve::AliMinResSolve(const AliMatrixSq *mat, const double* rhs) :
  fSize(mat->GetSize()),fPrecon(0),fMatrix((AliMatrixSq*)mat),fRHS((double*)rhs),
  fPVecY(0),fPVecR1(0),fPVecR2(0),fPVecV(0),fPVecW(0),fPVecW1(0),fPVecW2(0),
  fPvv(0),fPvz(0),fPhh(0),
  fDiagLU(0),fMatL(0),fMatU(0),fMatBD(0)
{
  // copy accepting equation
}

//______________________________________________________
AliMinResSolve::~AliMinResSolve()
{
  // destructor
  ClearAux();
}

//______________________________________________________
AliMinResSolve& AliMinResSolve::operator=(const AliMinResSolve& src)
{
  // assignment op.
  if (this != &src) {
    fSize    = src.fSize;
    fPrecon  = src.fPrecon;
    fMatrix  = src.fMatrix;
    fRHS     = src.fRHS;
  }
  return *this;
}

//_______________________________________________________________
Int_t AliMinResSolve::BuildPrecon(Int_t prec)
{
  // preconditioner building
  //  const Double_t kTiny = 1E-12;
  fPrecon = prec;
  //
  if (fPrecon>=kPreconBD && fPrecon<kPreconILU0) { // band diagonal decomposition
    return BuildPreconBD(fPrecon - kPreconBD);     // with halfbandwidth + diagonal = fPrecon
  }
  //
  if (fPrecon>=kPreconILU0 && fPrecon<=kPreconILU10) {
    if (fMatrix->InheritsFrom("AliMatrixSparse")) return BuildPreconILUK(fPrecon-kPreconILU0);
    else                                          return BuildPreconILUKDense(fPrecon-kPreconILU0);
  }
  //
  return -1;
}

//________________________________ FGMRES METHODS ________________________________
Bool_t AliMinResSolve::SolveFGMRES(TVectorD& VecSol,Int_t precon,int itnlim,double rtol,int nkrylov)
{
  // solve by fgmres
  return SolveFGMRES(VecSol.GetMatrixArray(),precon,itnlim,rtol,nkrylov);
}

//________________________________________________________________________________
Bool_t AliMinResSolve::SolveFGMRES(double* VecSol,Int_t precon,int itnlim,double rtol,int nkrylov)
{
  // Adapted from Y.Saad fgmrs.c of ITSOL_1 package by Y.Saad: http://www-users.cs.umn.edu/~saad/software/
  /*----------------------------------------------------------------------
    |                 *** Preconditioned FGMRES ***                  
    +-----------------------------------------------------------------------
    | This is a simple version of the ARMS preconditioned FGMRES algorithm. 
    +-----------------------------------------------------------------------
    | Y. S. Dec. 2000. -- Apr. 2008  
    +-----------------------------------------------------------------------
    | VecSol  = real vector of length n containing an initial guess to the
    | precon  = precondtioner id (0 = no precon)
    | itnlim  = max n of iterations
    | rtol     = tolerance for stopping iteration
    | nkrylov = N of Krylov vectors to store
    +---------------------------------------------------------------------*/
  int l;
  double status = kTRUE;
  double t,beta,eps1=0;
  const double epsmac    = 2.22E-16;
  //
  AliInfo(Form("Solution by FGMRes: Preconditioner#%d Max.iter.: %d Tol.: %e NKrylov: %d",
	       precon,itnlim,rtol,nkrylov));
  //
  int its = 0;
  if (nkrylov<10) {
    AliInfo(Form("Changing N Krylov vectors from %d to %d",nkrylov,10));
    nkrylov = 10;
  }
  //
  if (precon>0) {
    if (precon>=kPreconsTot) {
      AliInfo(Form("Unknown preconditioner identifier %d, ignore",precon));
    }
    else {
      if  (BuildPrecon(precon)<0) {
	ClearAux();
	AliInfo("FGMRES failed to build the preconditioner");
	return kFALSE;
      }
    }
  }
  //
  if (!InitAuxFGMRES(nkrylov)) return kFALSE;
  //
  for (l=fSize;l--;) VecSol[l] = 0;
  //
  //-------------------- outer loop starts here
  TStopwatch timer; timer.Start();
  while(1) {
    //
    //-------------------- compute initial residual vector
    fMatrix->MultiplyByVec(VecSol,fPvv[0]);
    for (l=fSize;l--;) fPvv[0][l] = fRHS[l] - fPvv[0][l];    //  fPvv[0]= initial residual
    beta = 0;
    for (l=fSize;l--;) beta += fPvv[0][l]*fPvv[0][l];
    beta = TMath::Sqrt(beta);
    //
    if (beta < epsmac) break;                                      // success? 
    t = 1.0 / beta;
    //--------------------   normalize:  fPvv[0] = fPvv[0] / beta 
    for (l=fSize;l--;) fPvv[0][l] *= t;
    if (its == 0) eps1 = rtol*beta;
    //
    //    ** initialize 1-st term  of rhs of hessenberg system
    fPVecV[0] = beta;
    int i = -1;
    do {
      i++;
      its++;
      int i1 = i+1; 
      //
      //  (Right) Preconditioning Operation   z_{j} = M^{-1} v_{j}
      //
      if (precon>0) ApplyPrecon( fPvv[i], fPvz[i]);
      else          for (l=fSize;l--;)  fPvz[i][l] = fPvv[i][l];
      //
      //-------------------- matvec operation w = A z_{j} = A M^{-1} v_{j}
      fMatrix->MultiplyByVec(fPvz[i],fPvv[i1]);
      //
      // modified gram - schmidt...
      // h_{i,j} = (w,v_{i})
      // w  = w - h_{i,j} v_{i}
      //
      for (int j=0; j<=i; j++) {
	for (t=0, l=fSize;l--;) t+= fPvv[j][l]*fPvv[i1][l];
	fPhh[i][j] = t;
	for (l=fSize;l--;) fPvv[i1][l] -= t*fPvv[j][l];
      }
      // -------------------- h_{j+1,j} = ||w||_{2}
      for (t=0,l=fSize;l--;) t += fPvv[i1][l]*fPvv[i1][l]; t = TMath::Sqrt(t); 
      fPhh[i][i1] = t;
      if (t > 0) for (t=1./t, l=0; l<fSize; l++) fPvv[i1][l] *= t;  //  v_{j+1} = w / h_{j+1,j}
      //
      // done with modified gram schimdt and arnoldi step
      // now  update factorization of fPhh
      //
      // perform previous transformations  on i-th column of h
      //
      for (l=1; l<=i; l++) {
	int l1 = l-1;
	t = fPhh[i][l1];
	fPhh[i][l1] = fPVecR1[l1]*t + fPVecR2[l1]*fPhh[i][l];
	fPhh[i][l] = -fPVecR2[l1]*t + fPVecR1[l1]*fPhh[i][l];
      }
      double gam = TMath::Sqrt( fPhh[i][i]*fPhh[i][i] + fPhh[i][i1]*fPhh[i][i1]);
      //
      // if gamma is zero then any small value will do...
      // will affect only residual estimate
      if (gam < epsmac) gam = epsmac;
      //  get  next plane rotation
      fPVecR1[i] = fPhh[i][i]/gam;
      fPVecR2[i]  = fPhh[i][i1]/gam;
      fPVecV[i1]  =-fPVecR2[i]*fPVecV[i];
      fPVecV[i]  *= fPVecR1[i];
      //
      //  determine residual norm and test for convergence
      fPhh[i][i] = fPVecR1[i]*fPhh[i][i] + fPVecR2[i]*fPhh[i][i1];
      beta = TMath::Abs(fPVecV[i1]);
      //
    } while ( (i < nkrylov-1) && (beta > eps1) && (its < itnlim) );
    //
    // now compute solution. 1st, solve upper triangular system
    fPVecV[i] = fPVecV[i]/fPhh[i][i];
    for (int j=1; j<=i; j++) {
      int k=i-j;
      for (t=fPVecV[k],l=k+1; l<=i; l++) t -= fPhh[l][k]*fPVecV[l];
      fPVecV[k] = t/fPhh[k][k];
    }
    // --------------------  linear combination of v[i]'s to get sol. 
    for (int j=0; j<=i; j++) for (t=fPVecV[j],l=0; l<fSize; l++) VecSol[l] += t*fPvz[j][l];
    //
    // --------------------  restart outer loop if needed
    //    
    if (beta <= eps1) {
      timer.Stop();
      AliInfo(Form("FGMRES converged in %d iterations, CPU time: %.1f sec",its,timer.CpuTime()));
      break;   // success
    }
    //
    if (its >= itnlim) {
      timer.Stop();
      AliInfo(Form("%d iterations limit exceeded, CPU time: %.1f sec",itnlim,timer.CpuTime()));
      status = kFALSE;
      break;
    }
  }
  //
  ClearAux();
  return status;
  //
}


//________________________________ MINRES METHODS ________________________________
Bool_t AliMinResSolve::SolveMinRes(TVectorD& VecSol,Int_t precon,int itnlim,double rtol)
{
  // solve by minres
  return SolveMinRes(VecSol.GetMatrixArray(), precon, itnlim, rtol);
}

//________________________________________________________________________________
Bool_t AliMinResSolve::SolveMinRes(double *VecSol,Int_t precon,int itnlim,double rtol)
{
  /*
    Adapted from author's Fortran code:
    Michael A. Saunders           na.msaunders@na-net.ornl.gov
    
    MINRES is an implementation of the algorithm described in the following reference:
    C. C. Paige and M. A. Saunders (1975),
    Solution of sparse indefinite systems of linear equations,
    SIAM J. Numer. Anal. 12(4), pp. 617-629.  

  */
  if (!fMatrix->IsSymmetric()) {
    AliInfo("MinRes cannot solve asymmetric matrices, use FGMRes instead");
    return kFALSE;
  }
  //
  ClearAux();
  const double eps    = 2.22E-16;
  double beta1;
  //
  if (precon>0) {
    if (precon>=kPreconsTot) {
      AliInfo(Form("Unknown preconditioner identifier %d, ignore",precon));
    }
    else {
      if  (BuildPrecon(precon)<0) {
	ClearAux();
	AliInfo("MinRes failed to build the preconditioner");
	return kFALSE;
      }
    }
  }
  AliInfo(Form("Solution by MinRes: Preconditioner#%d Max.iter.: %d Tol.: %e",precon,itnlim,rtol));
  //
  // ------------------------ initialization  ---------------------->>>>
  memset(VecSol,0,fSize*sizeof(double));
  int status=0, itn=0;
  double normA = 0;
  double condA = 0;
  double ynorm = 0;
  double rnorm = 0;
  double gam,gmax=1,gmin=1,gbar,oldeps,epsa,epsx,epsr,diag, delta,phi,denom,z;
  //
  if (!InitAuxMinRes()) return kFALSE;
  //
  memset(VecSol,0,fSize*sizeof(double));
  //
  // ------------ init aux -------------------------<<<<  
  //   Set up y and v for the first Lanczos vector v1.
  //   y  =  beta1 P' v1,  where  P = C**(-1). v is really P' v1.
  //
  for (int i=fSize;i--;) fPVecY[i]  = fPVecR1[i] = fRHS[i];
  //
  if ( precon>0 ) ApplyPrecon( fRHS, fPVecY);
  beta1 = 0; for (int i=fSize;i--;) beta1 += fRHS[i]*fPVecY[i]; //
  //
  if (beta1 < 0) {
    AliInfo(Form("Preconditioner is indefinite (init) (%e).",beta1));
    ClearAux();
    status = 7;
    return kFALSE;
  }
  //
  if (beta1 < eps) {
    AliInfo(Form("RHS is zero or is the nullspace of the Preconditioner: Solution is {0}"));
    ClearAux();
    return kTRUE;
  }  
  //
  beta1  = TMath::Sqrt( beta1 );       // Normalize y to get v1 later.
  //
  //      See if Msolve is symmetric. //RS: Skept
  //      See if Aprod  is symmetric. //RS: Skept
  //
  double oldb   = 0;
  double beta   = beta1;
  double dbar   = 0;
  double epsln  = 0;
  double qrnorm = beta1;
  double phibar = beta1;
  double rhs1   = beta1;
  double rhs2   = 0;
  double tnorm2 = 0;
  double ynorm2 = 0;
  double cs     = -1;
  double sn     = 0;
  for (int i=fSize;i--;) fPVecR2[i] = fPVecR1[i];
  //
  TStopwatch timer; timer.Start();
  while(status==0) { //-----------------  Main iteration loop ---------------------->>>>
    //
    itn++;
    /*-----------------------------------------------------------------
      Obtain quantities for the next Lanczos vector vk+1, k = 1, 2,...
      The general iteration is similar to the case k = 1 with v0 = 0:      
      p1      = Operator * v1  -  beta1 * v0,
      alpha1  = v1'p1,
      q2      = p2  -  alpha1 * v1,
      beta2^2 = q2'q2,
      v2      = (1/beta2) q2.      
      Again, y = betak P vk,  where  P = C**(-1).
      .... more description needed.
      -----------------------------------------------------------------*/
    //
    double s = 1./beta;                            // Normalize previous vector (in y).
    for (int i=fSize;i--;) fPVecV[i] = s*fPVecY[i];    // v = vk if P = I
    //
    fMatrix->MultiplyByVec(fPVecV,fPVecY);   //      APROD (VecV, VecY);
    //
    if (itn>=2) {
      double btrat = beta/oldb;
      for (int i=fSize;i--;) fPVecY[i] -= btrat*fPVecR1[i];
    }
    double alfa = 0; for (int i=fSize;i--;) alfa += fPVecV[i]*fPVecY[i];    //      alphak
    //
    double alf2bt = alfa/beta;
    for (int i=fSize;i--;) {
      fPVecY[i] -= alf2bt*fPVecR2[i];
      fPVecR1[i] = fPVecR2[i];
      fPVecR2[i] = fPVecY[i];
    }
    //
    if ( precon>0 ) ApplyPrecon(fPVecR2, fPVecY);
    //
    oldb  = beta;               //      oldb = betak
    beta = 0; for (int i=fSize;i--;) beta += fPVecR2[i]*fPVecY[i];    // beta = betak+1^2
    //
    if (beta<0) {
      AliInfo(Form("Preconditioner is indefinite (%e)",beta));
      status = 7;
      break;
    }
    // 
    beta   = TMath::Sqrt(beta); //            beta = betak+1
    tnorm2 += alfa*alfa + oldb*oldb + beta*beta;
    //
    if (itn == 1) {             //     Initialize a few things.
      if (beta/beta1 <= 10.0*eps) {
	status = 0; //-1   //?????  beta2 = 0 or ~ 0,  terminate later.
	AliInfo("RHS is eigenvector");
      }
      //        !tnorm2 = alfa**2
      gmax   = TMath::Abs(alfa);    //              alpha1
      gmin   = gmax;                //              alpha1
    }
    //
    /*
      Apply previous rotation Qk-1 to get
      [deltak epslnk+1] = [cs  sn][dbark    0   ]
      [gbar k dbar k+1]   [sn -cs][alfak betak+1].
    */
    //
    oldeps = epsln;
    delta  = cs * dbar  +  sn * alfa;       //  delta1 = 0         deltak
    gbar   = sn * dbar  -  cs * alfa;       //  gbar 1 = alfa1     gbar k
    epsln  =               sn * beta;       //  epsln2 = 0         epslnk+1
    dbar   =            -  cs * beta;       //  dbar 2 = beta2     dbar k+1
    //
    // Compute the next plane rotation Qk
    //
    gam    = TMath::Sqrt( gbar*gbar + beta*beta );  // gammak
    cs     = gbar / gam;                            // ck
    sn     = beta / gam;                            // sk
    phi    = cs * phibar;                           // phik
    phibar = sn * phibar;                           // phibark+1
    //
    // Update  x.
    denom = 1./gam;
    //
    for (int i=fSize;i--;) {
      fPVecW1[i] = fPVecW2[i];
      fPVecW2[i] = fPVecW[i];
      fPVecW[i]  = denom*(fPVecV[i]-oldeps*fPVecW1[i]-delta*fPVecW2[i]);
      VecSol[i] += phi*fPVecW[i];
    }
    //
    //  Go round again.
    //
    gmax = TMath::Max( gmax, gam );
    gmin = TMath::Min( gmin, gam );
    z    = rhs1 / gam;
    ynorm2 += z*z;
    rhs1   = rhs2  -  delta * z;
    rhs2   =       -  epsln * z;
    //
    //   Estimate various norms and test for convergence.
    normA  = TMath::Sqrt( tnorm2 );
    ynorm  = TMath::Sqrt( ynorm2 );
    epsa   = normA * eps;
    epsx   = normA * ynorm * eps;
    epsr   = normA * ynorm * rtol;
    diag   = gbar;
    if (diag == 0) diag = epsa;
    //
    qrnorm = phibar;
    rnorm  = qrnorm;
    /*
      Estimate  cond(A).
      In this version we look at the diagonals of  R  in the
      factorization of the lower Hessenberg matrix,  Q * H = R,
      where H is the tridiagonal matrix from Lanczos with one
      extra row, beta(k+1) e_k^T.
    */
    condA  = gmax / gmin;
    //
    // See if any of the stopping criteria are satisfied.
    // In rare cases, istop is already -1 from above (Abar = const*I).
    //
    AliDebug(2,Form("#%5d |qnrm: %+.2e Anrm:%+.2e Cnd:%+.2e Rnrm:%+.2e Ynrm:%+.2e EpsR:%+.2e EpsX:%+.2e Beta1:%+.2e",
		    itn,qrnorm, normA,condA,rnorm,ynorm,epsr,epsx,beta1));

    if (status == 0) {
      if (itn    >= itnlim    ) {status = 5; AliInfo(Form("%d iterations limit exceeded",itnlim));}
      if (condA  >= 0.1/eps   ) {status = 4; AliInfo(Form("Matrix condition nymber %e exceeds limit %e",condA,0.1/eps));}
      if (epsx   >= beta1     ) {status = 3; AliInfo(Form("Approximate convergence"));}
      if (qrnorm <= epsx      ) {status = 2; AliInfo(Form("Converged within machine precision"));}
      if (qrnorm <= epsr      ) {status = 1; AliInfo(Form("Converged"));}
    }
    //
  }  //-----------------  Main iteration loop ----------------------<<<
  //
  ClearAux();
  //
  timer.Stop();
  AliInfo(Form("Exit from MinRes: CPU time: %.2f sec\n"
	       "Status    :  %2d\n"
	       "Iterations:  %4d\n"
	       "Norm      :  %+e\n"
	       "Condition :  %+e\n"
	       "Res.Norm  :  %+e\n"
	       "Sol.Norm  :  %+e",
	       timer.CpuTime(),status,itn,normA,condA,rnorm,ynorm));
  //
  return status>=0 && status<=3;
  //
}

//______________________________________________________________
void AliMinResSolve::ApplyPrecon(const TVectorD& vecRHS, TVectorD& vecOut) const
{
  // apply precond.
  ApplyPrecon(vecRHS.GetMatrixArray(), vecOut.GetMatrixArray());
}

//______________________________________________________________
void AliMinResSolve::ApplyPrecon(const double* vecRHS, double* vecOut) const
{
  // Application of the preconditioner matrix:
  // implicitly defines the matrix solving the M*VecOut = VecRHS 
  //  const Double_t kTiny = 1E-12;
  if (fPrecon>=kPreconBD && fPrecon<kPreconILU0) { // band diagonal decomposition
    fMatBD->Solve(vecRHS,vecOut);
    //    return;
  }
  //
  else if (fPrecon>=kPreconILU0 && fPrecon<=kPreconILU10) {
    //
    for(int i=0; i<fSize; i++) {     // Block L solve
      vecOut[i] = vecRHS[i];
      AliVectorSparse &rowLi = *fMatL->GetRow(i);
      int n = rowLi.GetNElems();
      for(int j=0;j<n;j++) vecOut[i] -= vecOut[ rowLi.GetIndex(j) ] * rowLi.GetElem(j);

    }
    //
    for(int i=fSize; i--; ) {    // Block -- U solve
      AliVectorSparse &rowUi = *fMatU->GetRow(i);
      int n = rowUi.GetNElems();
      for(int j=0;j<n;j++ ) vecOut[i] -= vecOut[ rowUi.GetIndex(j) ] * rowUi.GetElem(j);
      vecOut[i] *= fDiagLU[i];
    }
    //
  }
  //
}


//___________________________________________________________
Bool_t AliMinResSolve::InitAuxMinRes()
{
  // init auxiliary space for minres
  fPVecY   = new double[fSize];   
  fPVecR1  = new double[fSize];   
  fPVecR2  = new double[fSize];   
  fPVecV   = new double[fSize];   
  fPVecW   = new double[fSize];   
  fPVecW1  = new double[fSize];   
  fPVecW2  = new double[fSize];   
  //
  for (int i=fSize;i--;) fPVecY[i]=fPVecR1[i]=fPVecR2[i]=fPVecV[i]=fPVecW[i]=fPVecW1[i]=fPVecW2[i]=0.0;
  //
  return kTRUE;
}


//___________________________________________________________
Bool_t AliMinResSolve::InitAuxFGMRES(int nkrylov)
{
  // init auxiliary space for fgmres
  fPvv     = new double*[nkrylov+1];
  fPvz     = new double*[nkrylov];
  for (int i=0; i<=nkrylov; i++) fPvv[i] = new double[fSize];
  fPhh     = new double*[nkrylov];
  for (int i=0; i<nkrylov; i++) {
    fPhh[i] = new double[i+2];
    fPvz[i]  = new double[fSize];
  }
  //
  fPVecR1  = new double[nkrylov];   
  fPVecR2  = new double[nkrylov];   
  fPVecV   = new double[nkrylov+1];
  //
  return kTRUE;
}


//___________________________________________________________
void AliMinResSolve::ClearAux()
{
  // clear aux. space
  if (fPVecY)      delete[] fPVecY;  fPVecY   = 0;
  if (fPVecR1)     delete[] fPVecR1; fPVecR1  = 0; 
  if (fPVecR2)     delete[] fPVecR2; fPVecR2  = 0; 
  if (fPVecV)      delete[] fPVecV;  fPVecV   = 0;
  if (fPVecW)      delete[] fPVecW;  fPVecW   = 0;
  if (fPVecW1)     delete[] fPVecW1; fPVecW1  = 0;
  if (fPVecW2)     delete[] fPVecW2; fPVecW2  = 0;
  if (fDiagLU)   delete[] fDiagLU; fDiagLU = 0;
  if (fMatL)       delete fMatL; fMatL = 0;
  if (fMatU)       delete fMatU; fMatU = 0;
  if (fMatBD)      delete fMatBD; fMatBD = 0;
}

//___________________________________________________________
Int_t  AliMinResSolve::BuildPreconBD(Int_t hwidth)
{
  // build preconditioner
  AliInfo(Form("Building Band-Diagonal preconditioner of half-width = %d",hwidth));
  fMatBD = new AliSymBDMatrix( fMatrix->GetSize(), hwidth );
  //
  // fill the band-diagonal part of the matrix
  if (fMatrix->InheritsFrom("AliMatrixSparse")) {
    for (int ir=fMatrix->GetSize();ir--;) {
      int jmin = TMath::Max(0,ir-hwidth);
      AliVectorSparse& irow = *((AliMatrixSparse*)fMatrix)->GetRow(ir);
      for(int j=irow.GetNElems();j--;) {
	int jind = irow.GetIndex(j);
	if (jind<jmin) break;
	(*fMatBD)(ir,jind) = irow.GetElem(j);
      }
    }
  }
  else {
    for (int ir=fMatrix->GetSize();ir--;) {
      int jmin = TMath::Max(0,ir-hwidth);
      for(int jr=jmin;jr<=ir;jr++) (*fMatBD)(ir,jr) = fMatrix->Query(ir,jr);
    }
  }
  //
  fMatBD->DecomposeLDLT();
  //
  return 0;
}

//___________________________________________________________
Int_t  AliMinResSolve::BuildPreconILUK(Int_t lofM)
{
  /*----------------------------------------------------------------------------
   * ILUK preconditioner
   * incomplete LU factorization with level of fill dropping
   * Adapted from iluk.c of ITSOL_1 package by Y.Saad: http://www-users.cs.umn.edu/~saad/software/
   *----------------------------------------------------------------------------*/
  //
  AliInfo(Form("Building ILU%d preconditioner",lofM));
  //
  TStopwatch sw; sw.Start();
  fMatL = new AliMatrixSparse(fSize);
  fMatU = new AliMatrixSparse(fSize);
  fMatL->SetSymmetric(kFALSE);
  fMatU->SetSymmetric(kFALSE);
  fDiagLU = new Double_t[fSize];
  AliMatrixSparse* matrix = (AliMatrixSparse*)fMatrix;
  //
  // symbolic factorization to calculate level of fill index arrays
  if ( PreconILUKsymb(lofM)<0 ) {
    ClearAux();
    return -1;
  }
  //
  Int_t *jw = new Int_t[fSize]; 
  for(int j=fSize;j--;) jw[j] = -1;     // set indicator array jw to -1 
  //
  for(int i=0; i<fSize; i++ ) {            // beginning of main loop
    if ( (i%int(0.1*fSize)) == 0) {
      AliInfo(Form("BuildPrecon: row %d of %d",i,fSize));
      sw.Stop();
      sw.Print();
      sw.Start(kFALSE);
    }
    /* setup array jw[], and initial i-th row */
    AliVectorSparse& rowLi = *fMatL->GetRow(i);
    AliVectorSparse& rowUi = *fMatU->GetRow(i);
    AliVectorSparse& rowM  = *matrix->GetRow(i);
    //
    for(int j=rowLi.GetNElems(); j--;) {  // initialize L part
      int col = rowLi.GetIndex(j);
      jw[col] = j;
      rowLi.GetElem(j) = 0.;   // do we need this ?
    }
    jw[i] = i;
    fDiagLU[i] = 0;      // initialize diagonal
    //
    for(int j=rowUi.GetNElems();j--;)   {  // initialize U part
      int col = rowUi.GetIndex(j);
      jw[col] = j;
      rowUi.GetElem(j) = 0;
    }
    // copy row from csmat into L,U D
    for(int j=rowM.GetNElems(); j--;) {  // L and D part 
      if (AliMatrixSq::IsZero(rowM.GetElem(j))) continue;
      int col = rowM.GetIndex(j);         // (the original matrix stores only lower triangle)
      if( col < i )   rowLi.GetElem(jw[col]) = rowM.GetElem(j); 
      else if(col==i) fDiagLU[i] = rowM.GetElem(j);
      else rowUi.GetElem(jw[col]) = rowM.GetElem(j);
    }
    if (matrix->IsSymmetric()) for (int col=i+1;col<fSize;col++) {      // part of the row I on the right of diagonal is stored as 
	double vl = matrix->Query(col,i);    // the lower part of the column I
	if (AliMatrixSq::IsZero(vl)) continue;
	rowUi.GetElem(jw[col]) = vl;
      }
    //
    // eliminate previous rows
    for(int j=0; j<rowLi.GetNElems(); j++) {
      int jrow = rowLi.GetIndex(j);
      // get the multiplier for row to be eliminated (jrow)
      rowLi.GetElem(j) *= fDiagLU[jrow];
      //
      // combine current row and row jrow
      AliVectorSparse& rowUj = *fMatU->GetRow(jrow);
      for(int k=0; k<rowUj.GetNElems(); k++ ) {
	int col = rowUj.GetIndex(k);
	int jpos = jw[col];
	if( jpos == -1 ) continue;
	if( col < i )   rowLi.GetElem(jpos) -= rowLi.GetElem(j) * rowUj.GetElem(k);
	else if(col==i) fDiagLU[i] -= rowLi.GetElem(j) * rowUj.GetElem(k);
	else            rowUi.GetElem(jpos) -= rowLi.GetElem(j) * rowUj.GetElem(k);
      }
    }
    // reset double-pointer to -1 ( U-part) 
    for(int j=rowLi.GetNElems(); j--;) jw[ rowLi.GetIndex(j) ] = -1;
    jw[i] = -1;
    for(int j=rowUi.GetNElems(); j--;) jw[ rowUi.GetIndex(j) ] = -1;
    //
    if( AliMatrixSq::IsZero(fDiagLU[i]) ) {
      AliInfo(Form("Fatal error in ILIk: Zero diagonal found..."));
      delete[] jw;
      return -1;
    }
    fDiagLU[i] = 1.0 / fDiagLU[i];
  }
  //
  delete[] jw;
  //
  sw.Stop();
  AliInfo(Form("ILU%d preconditioner OK, CPU time: %.1f sec",lofM,sw.CpuTime()));
  AliInfo(Form("Densities: M %f L %f U %f",matrix->GetDensity(),fMatL->GetDensity(),fMatU->GetDensity()));
  //
  return 0;
}

//___________________________________________________________
Int_t  AliMinResSolve::BuildPreconILUKDense(Int_t lofM)
{
  /*----------------------------------------------------------------------------
   * ILUK preconditioner
   * incomplete LU factorization with level of fill dropping
   * Adapted from iluk.c of ITSOL_1 package by Y.Saad: http://www-users.cs.umn.edu/~saad/software/
   *----------------------------------------------------------------------------*/
  //
  TStopwatch sw; sw.Start();
  AliInfo(Form("Building ILU%d preconditioner for dense matrix",lofM));
  //
  fMatL = new AliMatrixSparse(fSize);
  fMatU = new AliMatrixSparse(fSize);
  fMatL->SetSymmetric(kFALSE);
  fMatU->SetSymmetric(kFALSE);
  fDiagLU = new Double_t[fSize];
  //
  // symbolic factorization to calculate level of fill index arrays
  if ( PreconILUKsymbDense(lofM)<0 ) {
    ClearAux();
    return -1;
  }
  //
  Int_t *jw = new Int_t[fSize]; 
  for(int j=fSize;j--;) jw[j] = -1;     // set indicator array jw to -1 
  //
  for(int i=0; i<fSize; i++ ) {            // beginning of main loop
    /* setup array jw[], and initial i-th row */
    AliVectorSparse& rowLi = *fMatL->GetRow(i);
    AliVectorSparse& rowUi = *fMatU->GetRow(i);
    //
    for(int j=rowLi.GetNElems(); j--;) {  // initialize L part
      int col = rowLi.GetIndex(j);
      jw[col] = j;
      rowLi.GetElem(j) = 0.;   // do we need this ?
    }
    jw[i] = i;
    fDiagLU[i] = 0;      // initialize diagonal
    //
    for(int j=rowUi.GetNElems();j--;)   {  // initialize U part
      int col = rowUi.GetIndex(j);
      jw[col] = j;
      rowUi.GetElem(j) = 0;
    }
    // copy row from csmat into L,U D
    for(int j=fSize; j--;) {  // L and D part 
      double vl = fMatrix->Query(i,j);
      if (AliMatrixSq::IsZero(vl)) continue;
      if( j < i )   rowLi.GetElem(jw[j]) = vl;
      else if(j==i) fDiagLU[i] = vl;
      else rowUi.GetElem(jw[j]) = vl;
    }
    // eliminate previous rows
    for(int j=0; j<rowLi.GetNElems(); j++) {
      int jrow = rowLi.GetIndex(j);
      // get the multiplier for row to be eliminated (jrow)
      rowLi.GetElem(j) *= fDiagLU[jrow];
      //
      // combine current row and row jrow
      AliVectorSparse& rowUj = *fMatU->GetRow(jrow);
      for(int k=0; k<rowUj.GetNElems(); k++ ) {
	int col = rowUj.GetIndex(k);
	int jpos = jw[col];
	if( jpos == -1 ) continue;
	if( col < i )   rowLi.GetElem(jpos) -= rowLi.GetElem(j) * rowUj.GetElem(k);
	else if(col==i) fDiagLU[i] -= rowLi.GetElem(j) * rowUj.GetElem(k);
	else            rowUi.GetElem(jpos) -= rowLi.GetElem(j) * rowUj.GetElem(k);
      }
    }
    // reset double-pointer to -1 ( U-part) 
    for(int j=rowLi.GetNElems(); j--;) jw[ rowLi.GetIndex(j) ] = -1;
    jw[i] = -1;
    for(int j=rowUi.GetNElems(); j--;) jw[ rowUi.GetIndex(j) ] = -1;
    //
    if( AliMatrixSq::IsZero(fDiagLU[i])) {
      AliInfo(Form("Fatal error in ILIk: Zero diagonal found..."));
      delete[] jw;
      return -1;
    }
    fDiagLU[i] = 1.0 / fDiagLU[i];
  }
  //
  delete[] jw;
  //
  sw.Stop();
  AliInfo(Form("ILU%d dense preconditioner OK, CPU time: %.1f sec",lofM,sw.CpuTime()));
  //  AliInfo(Form("Densities: M %f L %f U %f",fMatrix->GetDensity(),fMatL->GetDensity(),fMatU->GetDensity()));
  //
  return 0;
}

//___________________________________________________________
Int_t  AliMinResSolve::PreconILUKsymb(Int_t lofM)
{
  /*----------------------------------------------------------------------------
   * ILUK preconditioner
   * incomplete LU factorization with level of fill dropping
   * Adapted from iluk.c: lofC of ITSOL_1 package by Y.Saad: http://www-users.cs.umn.edu/~saad/software/
   *----------------------------------------------------------------------------*/
  //
  TStopwatch sw;
  AliInfo("PreconILUKsymb>>");
  AliMatrixSparse* matrix = (AliMatrixSparse*)fMatrix; 
  sw.Start();
  //
  UChar_t **ulvl=0,*levls=0;
  UShort_t *jbuf=0;
  Int_t    *iw=0;
  ulvl = new UChar_t*[fSize];      // stores lev-fils for U part of ILU factorization
  levls = new UChar_t[fSize];
  jbuf = new UShort_t[fSize];
  iw = new Int_t[fSize];
  //
  for(int j=fSize; j--;) iw[j] = -1;           // initialize iw 
  for(int i=0; i<fSize; i++) {
    int incl = 0;
    int incu = i; 
    AliVectorSparse& row = *matrix->GetRow(i);
    //
    // assign lof = 0 for matrix elements
    for(int j=0;j<row.GetNElems(); j++) {
      int col = row.GetIndex(j);
      if (AliMatrixSq::IsZero(row.GetElem(j))) continue;  // !!!! matrix is sparse but sometimes 0 appears 
      if (col<i) {                      // L-part
	jbuf[incl] = col;
	levls[incl] = 0;
	iw[col] = incl++;
      }
      else if (col>i) {                 // This works only for general matrix
	jbuf[incu] = col;
	levls[incu] = 0;
	iw[col] = incu++;
      }
    }
    if (matrix->IsSymmetric()) for (int col=i+1;col<fSize;col++) { // U-part of symmetric matrix
	if (AliMatrixSq::IsZero(matrix->Query(col,i))) continue;    // Due to the symmetry  == matrix(i,col)
	jbuf[incu] = col;
	levls[incu] = 0;
	iw[col] = incu++;
      }
    //
    // symbolic k,i,j Gaussian elimination
    int jpiv = -1; 
    while (++jpiv < incl) {
      int k = jbuf[jpiv] ;                        // select leftmost pivot
      int kmin = k;
      int jmin = jpiv; 
      for(int j=jpiv+1; j<incl; j++) if( jbuf[j]<kmin ) { kmin = jbuf[j]; jmin = j;}
      //
      // ------------------------------------  swap
      if(jmin!=jpiv) {
	jbuf[jpiv] = kmin; 
	jbuf[jmin] = k; 
	iw[kmin] = jpiv;
	iw[k] = jmin; 
	int tj = levls[jpiv] ;
	levls[jpiv] = levls[jmin];
	levls[jmin] = tj;
	k = kmin; 
      }
      // ------------------------------------ symbolic linear combinaiton of rows
      AliVectorSparse& rowU = *fMatU->GetRow(k);
      for(int j=0; j<rowU.GetNElems(); j++ ) {
	int col = rowU.GetIndex(j);
	int it  = ulvl[k][j]+levls[jpiv]+1; 
	if( it > lofM ) continue; 
	int ip = iw[col];
	if( ip == -1 ) {
	  if( col < i) {
	    jbuf[incl] = col;
	    levls[incl] = it;
	    iw[col] = incl++;
	  } 
	  else if( col > i ) {
	    jbuf[incu] = col;
	    levls[incu] = it;
	    iw[col] = incu++;
	  } 
	}
	else levls[ip] = TMath::Min(levls[ip], it); 
      }
      //
    } // end - while loop
    //
    // reset iw
    for (int j=0;j<incl;j++) iw[jbuf[j]] = -1;
    for (int j=i;j<incu;j++) iw[jbuf[j]] = -1;
    //
    // copy L-part
    AliVectorSparse& rowLi = *fMatL->GetRow(i);
    rowLi.ReSize(incl);
    if(incl>0) memcpy(rowLi.GetIndices(), jbuf, sizeof(UShort_t)*incl);
    // copy U-part
    int k = incu-i; 
    AliVectorSparse& rowUi = *fMatU->GetRow(i);
    rowUi.ReSize(k);
    if( k > 0 ) {
      memcpy(rowUi.GetIndices(), jbuf+i, sizeof(UShort_t)*k);
      ulvl[i] = new UChar_t[k];   // update matrix of levels 
      memcpy( ulvl[i], levls+i, k*sizeof(UChar_t) );
    }
  }
  //
  // free temp space and leave
  delete[] levls;
  delete[] jbuf;
  for(int i=fSize; i--;) if (fMatU->GetRow(i)->GetNElems()) delete[] ulvl[i]; 
  delete[] ulvl; 
  delete[] iw;
  //
  fMatL->SortIndices();
  fMatU->SortIndices();
  sw.Stop();
  sw.Print();
  AliInfo("PreconILUKsymb<<");
  return 0;
}


//___________________________________________________________
Int_t  AliMinResSolve::PreconILUKsymbDense(Int_t lofM)
{
  /*----------------------------------------------------------------------------
   * ILUK preconditioner
   * incomplete LU factorization with level of fill dropping
   * Adapted from iluk.c: lofC of ITSOL_1 package by Y.Saad: http://www-users.cs.umn.edu/~saad/software/
   *----------------------------------------------------------------------------*/
  //
  UChar_t **ulvl=0,*levls=0;
  UShort_t *jbuf=0;
  Int_t    *iw=0;
  ulvl = new UChar_t*[fSize];      // stores lev-fils for U part of ILU factorization
  levls = new UChar_t[fSize];
  jbuf = new UShort_t[fSize];
  iw = new Int_t[fSize];
  //
  for(int j=fSize; j--;) iw[j] = -1;           // initialize iw 
  for(int i=0; i<fSize; i++) {
    int incl = 0;
    int incu = i; 
    //
    // assign lof = 0 for matrix elements
    for(int j=0;j<fSize; j++) {
      if (AliMatrixSq::IsZero(fMatrix->Query(i,j))) continue;
      if (j<i) {                      // L-part
	jbuf[incl] = j;
	levls[incl] = 0;
	iw[j] = incl++;
      }
      else if (j>i) {                 // This works only for general matrix
	jbuf[incu] = j;
	levls[incu] = 0;
	iw[j] = incu++;
      }
    }
    //
    // symbolic k,i,j Gaussian elimination
    int jpiv = -1; 
    while (++jpiv < incl) {
      int k = jbuf[jpiv] ;                        // select leftmost pivot
      int kmin = k;
      int jmin = jpiv; 
      for(int j=jpiv+1; j<incl; j++) if( jbuf[j]<kmin ) { kmin = jbuf[j]; jmin = j;}
      //
      // ------------------------------------  swap
      if(jmin!=jpiv) {
	jbuf[jpiv] = kmin; 
	jbuf[jmin] = k; 
	iw[kmin] = jpiv;
	iw[k] = jmin; 
	int tj = levls[jpiv] ;
	levls[jpiv] = levls[jmin];
	levls[jmin] = tj;
	k = kmin; 
      }
      // ------------------------------------ symbolic linear combinaiton of rows
      AliVectorSparse& rowU = *fMatU->GetRow(k);
      for(int j=0; j<rowU.GetNElems(); j++ ) {
	int col = rowU.GetIndex(j);
	int it  = ulvl[k][j]+levls[jpiv]+1; 
	if( it > lofM ) continue; 
	int ip = iw[col];
	if( ip == -1 ) {
	  if( col < i) {
	    jbuf[incl] = col;
	    levls[incl] = it;
	    iw[col] = incl++;
	  } 
	  else if( col > i ) {
	    jbuf[incu] = col;
	    levls[incu] = it;
	    iw[col] = incu++;
	  } 
	}
	else levls[ip] = TMath::Min(levls[ip], it); 
      }
      //
    } // end - while loop
    //
    // reset iw
    for (int j=0;j<incl;j++) iw[jbuf[j]] = -1;
    for (int j=i;j<incu;j++) iw[jbuf[j]] = -1;
    //
    // copy L-part
    AliVectorSparse& rowLi = *fMatL->GetRow(i);
    rowLi.ReSize(incl);
    if(incl>0) memcpy(rowLi.GetIndices(), jbuf, sizeof(UShort_t)*incl);
    // copy U-part
    int k = incu-i; 
    AliVectorSparse& rowUi = *fMatU->GetRow(i);
    rowUi.ReSize(k);
    if( k > 0 ) {
      memcpy(rowUi.GetIndices(), jbuf+i, sizeof(UShort_t)*k);
      ulvl[i] = new UChar_t[k];   // update matrix of levels 
      memcpy( ulvl[i], levls+i, k*sizeof(UChar_t) );
    }
  }
  //
  // free temp space and leave
  delete[] levls;
  delete[] jbuf;
  for(int i=fSize; i--;) if (fMatU->GetRow(i)->GetNElems()) delete[] ulvl[i]; 
  delete[] ulvl; 
  delete[] iw;
  //
  fMatL->SortIndices();
  fMatU->SortIndices();
  return 0;
}
 AliMinResSolve.cxx:1
 AliMinResSolve.cxx:2
 AliMinResSolve.cxx:3
 AliMinResSolve.cxx:4
 AliMinResSolve.cxx:5
 AliMinResSolve.cxx:6
 AliMinResSolve.cxx:7
 AliMinResSolve.cxx:8
 AliMinResSolve.cxx:9
 AliMinResSolve.cxx:10
 AliMinResSolve.cxx:11
 AliMinResSolve.cxx:12
 AliMinResSolve.cxx:13
 AliMinResSolve.cxx:14
 AliMinResSolve.cxx:15
 AliMinResSolve.cxx:16
 AliMinResSolve.cxx:17
 AliMinResSolve.cxx:18
 AliMinResSolve.cxx:19
 AliMinResSolve.cxx:20
 AliMinResSolve.cxx:21
 AliMinResSolve.cxx:22
 AliMinResSolve.cxx:23
 AliMinResSolve.cxx:24
 AliMinResSolve.cxx:25
 AliMinResSolve.cxx:26
 AliMinResSolve.cxx:27
 AliMinResSolve.cxx:28
 AliMinResSolve.cxx:29
 AliMinResSolve.cxx:30
 AliMinResSolve.cxx:31
 AliMinResSolve.cxx:32
 AliMinResSolve.cxx:33
 AliMinResSolve.cxx:34
 AliMinResSolve.cxx:35
 AliMinResSolve.cxx:36
 AliMinResSolve.cxx:37
 AliMinResSolve.cxx:38
 AliMinResSolve.cxx:39
 AliMinResSolve.cxx:40
 AliMinResSolve.cxx:41
 AliMinResSolve.cxx:42
 AliMinResSolve.cxx:43
 AliMinResSolve.cxx:44
 AliMinResSolve.cxx:45
 AliMinResSolve.cxx:46
 AliMinResSolve.cxx:47
 AliMinResSolve.cxx:48
 AliMinResSolve.cxx:49
 AliMinResSolve.cxx:50
 AliMinResSolve.cxx:51
 AliMinResSolve.cxx:52
 AliMinResSolve.cxx:53
 AliMinResSolve.cxx:54
 AliMinResSolve.cxx:55
 AliMinResSolve.cxx:56
 AliMinResSolve.cxx:57
 AliMinResSolve.cxx:58
 AliMinResSolve.cxx:59
 AliMinResSolve.cxx:60
 AliMinResSolve.cxx:61
 AliMinResSolve.cxx:62
 AliMinResSolve.cxx:63
 AliMinResSolve.cxx:64
 AliMinResSolve.cxx:65
 AliMinResSolve.cxx:66
 AliMinResSolve.cxx:67
 AliMinResSolve.cxx:68
 AliMinResSolve.cxx:69
 AliMinResSolve.cxx:70
 AliMinResSolve.cxx:71
 AliMinResSolve.cxx:72
 AliMinResSolve.cxx:73
 AliMinResSolve.cxx:74
 AliMinResSolve.cxx:75
 AliMinResSolve.cxx:76
 AliMinResSolve.cxx:77
 AliMinResSolve.cxx:78
 AliMinResSolve.cxx:79
 AliMinResSolve.cxx:80
 AliMinResSolve.cxx:81
 AliMinResSolve.cxx:82
 AliMinResSolve.cxx:83
 AliMinResSolve.cxx:84
 AliMinResSolve.cxx:85
 AliMinResSolve.cxx:86
 AliMinResSolve.cxx:87
 AliMinResSolve.cxx:88
 AliMinResSolve.cxx:89
 AliMinResSolve.cxx:90
 AliMinResSolve.cxx:91
 AliMinResSolve.cxx:92
 AliMinResSolve.cxx:93
 AliMinResSolve.cxx:94
 AliMinResSolve.cxx:95
 AliMinResSolve.cxx:96
 AliMinResSolve.cxx:97
 AliMinResSolve.cxx:98
 AliMinResSolve.cxx:99
 AliMinResSolve.cxx:100
 AliMinResSolve.cxx:101
 AliMinResSolve.cxx:102
 AliMinResSolve.cxx:103
 AliMinResSolve.cxx:104
 AliMinResSolve.cxx:105
 AliMinResSolve.cxx:106
 AliMinResSolve.cxx:107
 AliMinResSolve.cxx:108
 AliMinResSolve.cxx:109
 AliMinResSolve.cxx:110
 AliMinResSolve.cxx:111
 AliMinResSolve.cxx:112
 AliMinResSolve.cxx:113
 AliMinResSolve.cxx:114
 AliMinResSolve.cxx:115
 AliMinResSolve.cxx:116
 AliMinResSolve.cxx:117
 AliMinResSolve.cxx:118
 AliMinResSolve.cxx:119
 AliMinResSolve.cxx:120
 AliMinResSolve.cxx:121
 AliMinResSolve.cxx:122
 AliMinResSolve.cxx:123
 AliMinResSolve.cxx:124
 AliMinResSolve.cxx:125
 AliMinResSolve.cxx:126
 AliMinResSolve.cxx:127
 AliMinResSolve.cxx:128
 AliMinResSolve.cxx:129
 AliMinResSolve.cxx:130
 AliMinResSolve.cxx:131
 AliMinResSolve.cxx:132
 AliMinResSolve.cxx:133
 AliMinResSolve.cxx:134
 AliMinResSolve.cxx:135
 AliMinResSolve.cxx:136
 AliMinResSolve.cxx:137
 AliMinResSolve.cxx:138
 AliMinResSolve.cxx:139
 AliMinResSolve.cxx:140
 AliMinResSolve.cxx:141
 AliMinResSolve.cxx:142
 AliMinResSolve.cxx:143
 AliMinResSolve.cxx:144
 AliMinResSolve.cxx:145
 AliMinResSolve.cxx:146
 AliMinResSolve.cxx:147
 AliMinResSolve.cxx:148
 AliMinResSolve.cxx:149
 AliMinResSolve.cxx:150
 AliMinResSolve.cxx:151
 AliMinResSolve.cxx:152
 AliMinResSolve.cxx:153
 AliMinResSolve.cxx:154
 AliMinResSolve.cxx:155
 AliMinResSolve.cxx:156
 AliMinResSolve.cxx:157
 AliMinResSolve.cxx:158
 AliMinResSolve.cxx:159
 AliMinResSolve.cxx:160
 AliMinResSolve.cxx:161
 AliMinResSolve.cxx:162
 AliMinResSolve.cxx:163
 AliMinResSolve.cxx:164
 AliMinResSolve.cxx:165
 AliMinResSolve.cxx:166
 AliMinResSolve.cxx:167
 AliMinResSolve.cxx:168
 AliMinResSolve.cxx:169
 AliMinResSolve.cxx:170
 AliMinResSolve.cxx:171
 AliMinResSolve.cxx:172
 AliMinResSolve.cxx:173
 AliMinResSolve.cxx:174
 AliMinResSolve.cxx:175
 AliMinResSolve.cxx:176
 AliMinResSolve.cxx:177
 AliMinResSolve.cxx:178
 AliMinResSolve.cxx:179
 AliMinResSolve.cxx:180
 AliMinResSolve.cxx:181
 AliMinResSolve.cxx:182
 AliMinResSolve.cxx:183
 AliMinResSolve.cxx:184
 AliMinResSolve.cxx:185
 AliMinResSolve.cxx:186
 AliMinResSolve.cxx:187
 AliMinResSolve.cxx:188
 AliMinResSolve.cxx:189
 AliMinResSolve.cxx:190
 AliMinResSolve.cxx:191
 AliMinResSolve.cxx:192
 AliMinResSolve.cxx:193
 AliMinResSolve.cxx:194
 AliMinResSolve.cxx:195
 AliMinResSolve.cxx:196
 AliMinResSolve.cxx:197
 AliMinResSolve.cxx:198
 AliMinResSolve.cxx:199
 AliMinResSolve.cxx:200
 AliMinResSolve.cxx:201
 AliMinResSolve.cxx:202
 AliMinResSolve.cxx:203
 AliMinResSolve.cxx:204
 AliMinResSolve.cxx:205
 AliMinResSolve.cxx:206
 AliMinResSolve.cxx:207
 AliMinResSolve.cxx:208
 AliMinResSolve.cxx:209
 AliMinResSolve.cxx:210
 AliMinResSolve.cxx:211
 AliMinResSolve.cxx:212
 AliMinResSolve.cxx:213
 AliMinResSolve.cxx:214
 AliMinResSolve.cxx:215
 AliMinResSolve.cxx:216
 AliMinResSolve.cxx:217
 AliMinResSolve.cxx:218
 AliMinResSolve.cxx:219
 AliMinResSolve.cxx:220
 AliMinResSolve.cxx:221
 AliMinResSolve.cxx:222
 AliMinResSolve.cxx:223
 AliMinResSolve.cxx:224
 AliMinResSolve.cxx:225
 AliMinResSolve.cxx:226
 AliMinResSolve.cxx:227
 AliMinResSolve.cxx:228
 AliMinResSolve.cxx:229
 AliMinResSolve.cxx:230
 AliMinResSolve.cxx:231
 AliMinResSolve.cxx:232
 AliMinResSolve.cxx:233
 AliMinResSolve.cxx:234
 AliMinResSolve.cxx:235
 AliMinResSolve.cxx:236
 AliMinResSolve.cxx:237
 AliMinResSolve.cxx:238
 AliMinResSolve.cxx:239
 AliMinResSolve.cxx:240
 AliMinResSolve.cxx:241
 AliMinResSolve.cxx:242
 AliMinResSolve.cxx:243
 AliMinResSolve.cxx:244
 AliMinResSolve.cxx:245
 AliMinResSolve.cxx:246
 AliMinResSolve.cxx:247
 AliMinResSolve.cxx:248
 AliMinResSolve.cxx:249
 AliMinResSolve.cxx:250
 AliMinResSolve.cxx:251
 AliMinResSolve.cxx:252
 AliMinResSolve.cxx:253
 AliMinResSolve.cxx:254
 AliMinResSolve.cxx:255
 AliMinResSolve.cxx:256
 AliMinResSolve.cxx:257
 AliMinResSolve.cxx:258
 AliMinResSolve.cxx:259
 AliMinResSolve.cxx:260
 AliMinResSolve.cxx:261
 AliMinResSolve.cxx:262
 AliMinResSolve.cxx:263
 AliMinResSolve.cxx:264
 AliMinResSolve.cxx:265
 AliMinResSolve.cxx:266
 AliMinResSolve.cxx:267
 AliMinResSolve.cxx:268
 AliMinResSolve.cxx:269
 AliMinResSolve.cxx:270
 AliMinResSolve.cxx:271
 AliMinResSolve.cxx:272
 AliMinResSolve.cxx:273
 AliMinResSolve.cxx:274
 AliMinResSolve.cxx:275
 AliMinResSolve.cxx:276
 AliMinResSolve.cxx:277
 AliMinResSolve.cxx:278
 AliMinResSolve.cxx:279
 AliMinResSolve.cxx:280
 AliMinResSolve.cxx:281
 AliMinResSolve.cxx:282
 AliMinResSolve.cxx:283
 AliMinResSolve.cxx:284
 AliMinResSolve.cxx:285
 AliMinResSolve.cxx:286
 AliMinResSolve.cxx:287
 AliMinResSolve.cxx:288
 AliMinResSolve.cxx:289
 AliMinResSolve.cxx:290
 AliMinResSolve.cxx:291
 AliMinResSolve.cxx:292
 AliMinResSolve.cxx:293
 AliMinResSolve.cxx:294
 AliMinResSolve.cxx:295
 AliMinResSolve.cxx:296
 AliMinResSolve.cxx:297
 AliMinResSolve.cxx:298
 AliMinResSolve.cxx:299
 AliMinResSolve.cxx:300
 AliMinResSolve.cxx:301
 AliMinResSolve.cxx:302
 AliMinResSolve.cxx:303
 AliMinResSolve.cxx:304
 AliMinResSolve.cxx:305
 AliMinResSolve.cxx:306
 AliMinResSolve.cxx:307
 AliMinResSolve.cxx:308
 AliMinResSolve.cxx:309
 AliMinResSolve.cxx:310
 AliMinResSolve.cxx:311
 AliMinResSolve.cxx:312
 AliMinResSolve.cxx:313
 AliMinResSolve.cxx:314
 AliMinResSolve.cxx:315
 AliMinResSolve.cxx:316
 AliMinResSolve.cxx:317
 AliMinResSolve.cxx:318
 AliMinResSolve.cxx:319
 AliMinResSolve.cxx:320
 AliMinResSolve.cxx:321
 AliMinResSolve.cxx:322
 AliMinResSolve.cxx:323
 AliMinResSolve.cxx:324
 AliMinResSolve.cxx:325
 AliMinResSolve.cxx:326
 AliMinResSolve.cxx:327
 AliMinResSolve.cxx:328
 AliMinResSolve.cxx:329
 AliMinResSolve.cxx:330
 AliMinResSolve.cxx:331
 AliMinResSolve.cxx:332
 AliMinResSolve.cxx:333
 AliMinResSolve.cxx:334
 AliMinResSolve.cxx:335
 AliMinResSolve.cxx:336
 AliMinResSolve.cxx:337
 AliMinResSolve.cxx:338
 AliMinResSolve.cxx:339
 AliMinResSolve.cxx:340
 AliMinResSolve.cxx:341
 AliMinResSolve.cxx:342
 AliMinResSolve.cxx:343
 AliMinResSolve.cxx:344
 AliMinResSolve.cxx:345
 AliMinResSolve.cxx:346
 AliMinResSolve.cxx:347
 AliMinResSolve.cxx:348
 AliMinResSolve.cxx:349
 AliMinResSolve.cxx:350
 AliMinResSolve.cxx:351
 AliMinResSolve.cxx:352
 AliMinResSolve.cxx:353
 AliMinResSolve.cxx:354
 AliMinResSolve.cxx:355
 AliMinResSolve.cxx:356
 AliMinResSolve.cxx:357
 AliMinResSolve.cxx:358
 AliMinResSolve.cxx:359
 AliMinResSolve.cxx:360
 AliMinResSolve.cxx:361
 AliMinResSolve.cxx:362
 AliMinResSolve.cxx:363
 AliMinResSolve.cxx:364
 AliMinResSolve.cxx:365
 AliMinResSolve.cxx:366
 AliMinResSolve.cxx:367
 AliMinResSolve.cxx:368
 AliMinResSolve.cxx:369
 AliMinResSolve.cxx:370
 AliMinResSolve.cxx:371
 AliMinResSolve.cxx:372
 AliMinResSolve.cxx:373
 AliMinResSolve.cxx:374
 AliMinResSolve.cxx:375
 AliMinResSolve.cxx:376
 AliMinResSolve.cxx:377
 AliMinResSolve.cxx:378
 AliMinResSolve.cxx:379
 AliMinResSolve.cxx:380
 AliMinResSolve.cxx:381
 AliMinResSolve.cxx:382
 AliMinResSolve.cxx:383
 AliMinResSolve.cxx:384
 AliMinResSolve.cxx:385
 AliMinResSolve.cxx:386
 AliMinResSolve.cxx:387
 AliMinResSolve.cxx:388
 AliMinResSolve.cxx:389
 AliMinResSolve.cxx:390
 AliMinResSolve.cxx:391
 AliMinResSolve.cxx:392
 AliMinResSolve.cxx:393
 AliMinResSolve.cxx:394
 AliMinResSolve.cxx:395
 AliMinResSolve.cxx:396
 AliMinResSolve.cxx:397
 AliMinResSolve.cxx:398
 AliMinResSolve.cxx:399
 AliMinResSolve.cxx:400
 AliMinResSolve.cxx:401
 AliMinResSolve.cxx:402
 AliMinResSolve.cxx:403
 AliMinResSolve.cxx:404
 AliMinResSolve.cxx:405
 AliMinResSolve.cxx:406
 AliMinResSolve.cxx:407
 AliMinResSolve.cxx:408
 AliMinResSolve.cxx:409
 AliMinResSolve.cxx:410
 AliMinResSolve.cxx:411
 AliMinResSolve.cxx:412
 AliMinResSolve.cxx:413
 AliMinResSolve.cxx:414
 AliMinResSolve.cxx:415
 AliMinResSolve.cxx:416
 AliMinResSolve.cxx:417
 AliMinResSolve.cxx:418
 AliMinResSolve.cxx:419
 AliMinResSolve.cxx:420
 AliMinResSolve.cxx:421
 AliMinResSolve.cxx:422
 AliMinResSolve.cxx:423
 AliMinResSolve.cxx:424
 AliMinResSolve.cxx:425
 AliMinResSolve.cxx:426
 AliMinResSolve.cxx:427
 AliMinResSolve.cxx:428
 AliMinResSolve.cxx:429
 AliMinResSolve.cxx:430
 AliMinResSolve.cxx:431
 AliMinResSolve.cxx:432
 AliMinResSolve.cxx:433
 AliMinResSolve.cxx:434
 AliMinResSolve.cxx:435
 AliMinResSolve.cxx:436
 AliMinResSolve.cxx:437
 AliMinResSolve.cxx:438
 AliMinResSolve.cxx:439
 AliMinResSolve.cxx:440
 AliMinResSolve.cxx:441
 AliMinResSolve.cxx:442
 AliMinResSolve.cxx:443
 AliMinResSolve.cxx:444
 AliMinResSolve.cxx:445
 AliMinResSolve.cxx:446
 AliMinResSolve.cxx:447
 AliMinResSolve.cxx:448
 AliMinResSolve.cxx:449
 AliMinResSolve.cxx:450
 AliMinResSolve.cxx:451
 AliMinResSolve.cxx:452
 AliMinResSolve.cxx:453
 AliMinResSolve.cxx:454
 AliMinResSolve.cxx:455
 AliMinResSolve.cxx:456
 AliMinResSolve.cxx:457
 AliMinResSolve.cxx:458
 AliMinResSolve.cxx:459
 AliMinResSolve.cxx:460
 AliMinResSolve.cxx:461
 AliMinResSolve.cxx:462
 AliMinResSolve.cxx:463
 AliMinResSolve.cxx:464
 AliMinResSolve.cxx:465
 AliMinResSolve.cxx:466
 AliMinResSolve.cxx:467
 AliMinResSolve.cxx:468
 AliMinResSolve.cxx:469
 AliMinResSolve.cxx:470
 AliMinResSolve.cxx:471
 AliMinResSolve.cxx:472
 AliMinResSolve.cxx:473
 AliMinResSolve.cxx:474
 AliMinResSolve.cxx:475
 AliMinResSolve.cxx:476
 AliMinResSolve.cxx:477
 AliMinResSolve.cxx:478
 AliMinResSolve.cxx:479
 AliMinResSolve.cxx:480
 AliMinResSolve.cxx:481
 AliMinResSolve.cxx:482
 AliMinResSolve.cxx:483
 AliMinResSolve.cxx:484
 AliMinResSolve.cxx:485
 AliMinResSolve.cxx:486
 AliMinResSolve.cxx:487
 AliMinResSolve.cxx:488
 AliMinResSolve.cxx:489
 AliMinResSolve.cxx:490
 AliMinResSolve.cxx:491
 AliMinResSolve.cxx:492
 AliMinResSolve.cxx:493
 AliMinResSolve.cxx:494
 AliMinResSolve.cxx:495
 AliMinResSolve.cxx:496
 AliMinResSolve.cxx:497
 AliMinResSolve.cxx:498
 AliMinResSolve.cxx:499
 AliMinResSolve.cxx:500
 AliMinResSolve.cxx:501
 AliMinResSolve.cxx:502
 AliMinResSolve.cxx:503
 AliMinResSolve.cxx:504
 AliMinResSolve.cxx:505
 AliMinResSolve.cxx:506
 AliMinResSolve.cxx:507
 AliMinResSolve.cxx:508
 AliMinResSolve.cxx:509
 AliMinResSolve.cxx:510
 AliMinResSolve.cxx:511
 AliMinResSolve.cxx:512
 AliMinResSolve.cxx:513
 AliMinResSolve.cxx:514
 AliMinResSolve.cxx:515
 AliMinResSolve.cxx:516
 AliMinResSolve.cxx:517
 AliMinResSolve.cxx:518
 AliMinResSolve.cxx:519
 AliMinResSolve.cxx:520
 AliMinResSolve.cxx:521
 AliMinResSolve.cxx:522
 AliMinResSolve.cxx:523
 AliMinResSolve.cxx:524
 AliMinResSolve.cxx:525
 AliMinResSolve.cxx:526
 AliMinResSolve.cxx:527
 AliMinResSolve.cxx:528
 AliMinResSolve.cxx:529
 AliMinResSolve.cxx:530
 AliMinResSolve.cxx:531
 AliMinResSolve.cxx:532
 AliMinResSolve.cxx:533
 AliMinResSolve.cxx:534
 AliMinResSolve.cxx:535
 AliMinResSolve.cxx:536
 AliMinResSolve.cxx:537
 AliMinResSolve.cxx:538
 AliMinResSolve.cxx:539
 AliMinResSolve.cxx:540
 AliMinResSolve.cxx:541
 AliMinResSolve.cxx:542
 AliMinResSolve.cxx:543
 AliMinResSolve.cxx:544
 AliMinResSolve.cxx:545
 AliMinResSolve.cxx:546
 AliMinResSolve.cxx:547
 AliMinResSolve.cxx:548
 AliMinResSolve.cxx:549
 AliMinResSolve.cxx:550
 AliMinResSolve.cxx:551
 AliMinResSolve.cxx:552
 AliMinResSolve.cxx:553
 AliMinResSolve.cxx:554
 AliMinResSolve.cxx:555
 AliMinResSolve.cxx:556
 AliMinResSolve.cxx:557
 AliMinResSolve.cxx:558
 AliMinResSolve.cxx:559
 AliMinResSolve.cxx:560
 AliMinResSolve.cxx:561
 AliMinResSolve.cxx:562
 AliMinResSolve.cxx:563
 AliMinResSolve.cxx:564
 AliMinResSolve.cxx:565
 AliMinResSolve.cxx:566
 AliMinResSolve.cxx:567
 AliMinResSolve.cxx:568
 AliMinResSolve.cxx:569
 AliMinResSolve.cxx:570
 AliMinResSolve.cxx:571
 AliMinResSolve.cxx:572
 AliMinResSolve.cxx:573
 AliMinResSolve.cxx:574
 AliMinResSolve.cxx:575
 AliMinResSolve.cxx:576
 AliMinResSolve.cxx:577
 AliMinResSolve.cxx:578
 AliMinResSolve.cxx:579
 AliMinResSolve.cxx:580
 AliMinResSolve.cxx:581
 AliMinResSolve.cxx:582
 AliMinResSolve.cxx:583
 AliMinResSolve.cxx:584
 AliMinResSolve.cxx:585
 AliMinResSolve.cxx:586
 AliMinResSolve.cxx:587
 AliMinResSolve.cxx:588
 AliMinResSolve.cxx:589
 AliMinResSolve.cxx:590
 AliMinResSolve.cxx:591
 AliMinResSolve.cxx:592
 AliMinResSolve.cxx:593
 AliMinResSolve.cxx:594
 AliMinResSolve.cxx:595
 AliMinResSolve.cxx:596
 AliMinResSolve.cxx:597
 AliMinResSolve.cxx:598
 AliMinResSolve.cxx:599
 AliMinResSolve.cxx:600
 AliMinResSolve.cxx:601
 AliMinResSolve.cxx:602
 AliMinResSolve.cxx:603
 AliMinResSolve.cxx:604
 AliMinResSolve.cxx:605
 AliMinResSolve.cxx:606
 AliMinResSolve.cxx:607
 AliMinResSolve.cxx:608
 AliMinResSolve.cxx:609
 AliMinResSolve.cxx:610
 AliMinResSolve.cxx:611
 AliMinResSolve.cxx:612
 AliMinResSolve.cxx:613
 AliMinResSolve.cxx:614
 AliMinResSolve.cxx:615
 AliMinResSolve.cxx:616
 AliMinResSolve.cxx:617
 AliMinResSolve.cxx:618
 AliMinResSolve.cxx:619
 AliMinResSolve.cxx:620
 AliMinResSolve.cxx:621
 AliMinResSolve.cxx:622
 AliMinResSolve.cxx:623
 AliMinResSolve.cxx:624
 AliMinResSolve.cxx:625
 AliMinResSolve.cxx:626
 AliMinResSolve.cxx:627
 AliMinResSolve.cxx:628
 AliMinResSolve.cxx:629
 AliMinResSolve.cxx:630
 AliMinResSolve.cxx:631
 AliMinResSolve.cxx:632
 AliMinResSolve.cxx:633
 AliMinResSolve.cxx:634
 AliMinResSolve.cxx:635
 AliMinResSolve.cxx:636
 AliMinResSolve.cxx:637
 AliMinResSolve.cxx:638
 AliMinResSolve.cxx:639
 AliMinResSolve.cxx:640
 AliMinResSolve.cxx:641
 AliMinResSolve.cxx:642
 AliMinResSolve.cxx:643
 AliMinResSolve.cxx:644
 AliMinResSolve.cxx:645
 AliMinResSolve.cxx:646
 AliMinResSolve.cxx:647
 AliMinResSolve.cxx:648
 AliMinResSolve.cxx:649
 AliMinResSolve.cxx:650
 AliMinResSolve.cxx:651
 AliMinResSolve.cxx:652
 AliMinResSolve.cxx:653
 AliMinResSolve.cxx:654
 AliMinResSolve.cxx:655
 AliMinResSolve.cxx:656
 AliMinResSolve.cxx:657
 AliMinResSolve.cxx:658
 AliMinResSolve.cxx:659
 AliMinResSolve.cxx:660
 AliMinResSolve.cxx:661
 AliMinResSolve.cxx:662
 AliMinResSolve.cxx:663
 AliMinResSolve.cxx:664
 AliMinResSolve.cxx:665
 AliMinResSolve.cxx:666
 AliMinResSolve.cxx:667
 AliMinResSolve.cxx:668
 AliMinResSolve.cxx:669
 AliMinResSolve.cxx:670
 AliMinResSolve.cxx:671
 AliMinResSolve.cxx:672
 AliMinResSolve.cxx:673
 AliMinResSolve.cxx:674
 AliMinResSolve.cxx:675
 AliMinResSolve.cxx:676
 AliMinResSolve.cxx:677
 AliMinResSolve.cxx:678
 AliMinResSolve.cxx:679
 AliMinResSolve.cxx:680
 AliMinResSolve.cxx:681
 AliMinResSolve.cxx:682
 AliMinResSolve.cxx:683
 AliMinResSolve.cxx:684
 AliMinResSolve.cxx:685
 AliMinResSolve.cxx:686
 AliMinResSolve.cxx:687
 AliMinResSolve.cxx:688
 AliMinResSolve.cxx:689
 AliMinResSolve.cxx:690
 AliMinResSolve.cxx:691
 AliMinResSolve.cxx:692
 AliMinResSolve.cxx:693
 AliMinResSolve.cxx:694
 AliMinResSolve.cxx:695
 AliMinResSolve.cxx:696
 AliMinResSolve.cxx:697
 AliMinResSolve.cxx:698
 AliMinResSolve.cxx:699
 AliMinResSolve.cxx:700
 AliMinResSolve.cxx:701
 AliMinResSolve.cxx:702
 AliMinResSolve.cxx:703
 AliMinResSolve.cxx:704
 AliMinResSolve.cxx:705
 AliMinResSolve.cxx:706
 AliMinResSolve.cxx:707
 AliMinResSolve.cxx:708
 AliMinResSolve.cxx:709
 AliMinResSolve.cxx:710
 AliMinResSolve.cxx:711
 AliMinResSolve.cxx:712
 AliMinResSolve.cxx:713
 AliMinResSolve.cxx:714
 AliMinResSolve.cxx:715
 AliMinResSolve.cxx:716
 AliMinResSolve.cxx:717
 AliMinResSolve.cxx:718
 AliMinResSolve.cxx:719
 AliMinResSolve.cxx:720
 AliMinResSolve.cxx:721
 AliMinResSolve.cxx:722
 AliMinResSolve.cxx:723
 AliMinResSolve.cxx:724
 AliMinResSolve.cxx:725
 AliMinResSolve.cxx:726
 AliMinResSolve.cxx:727
 AliMinResSolve.cxx:728
 AliMinResSolve.cxx:729
 AliMinResSolve.cxx:730
 AliMinResSolve.cxx:731
 AliMinResSolve.cxx:732
 AliMinResSolve.cxx:733
 AliMinResSolve.cxx:734
 AliMinResSolve.cxx:735
 AliMinResSolve.cxx:736
 AliMinResSolve.cxx:737
 AliMinResSolve.cxx:738
 AliMinResSolve.cxx:739
 AliMinResSolve.cxx:740
 AliMinResSolve.cxx:741
 AliMinResSolve.cxx:742
 AliMinResSolve.cxx:743
 AliMinResSolve.cxx:744
 AliMinResSolve.cxx:745
 AliMinResSolve.cxx:746
 AliMinResSolve.cxx:747
 AliMinResSolve.cxx:748
 AliMinResSolve.cxx:749
 AliMinResSolve.cxx:750
 AliMinResSolve.cxx:751
 AliMinResSolve.cxx:752
 AliMinResSolve.cxx:753
 AliMinResSolve.cxx:754
 AliMinResSolve.cxx:755
 AliMinResSolve.cxx:756
 AliMinResSolve.cxx:757
 AliMinResSolve.cxx:758
 AliMinResSolve.cxx:759
 AliMinResSolve.cxx:760
 AliMinResSolve.cxx:761
 AliMinResSolve.cxx:762
 AliMinResSolve.cxx:763
 AliMinResSolve.cxx:764
 AliMinResSolve.cxx:765
 AliMinResSolve.cxx:766
 AliMinResSolve.cxx:767
 AliMinResSolve.cxx:768
 AliMinResSolve.cxx:769
 AliMinResSolve.cxx:770
 AliMinResSolve.cxx:771
 AliMinResSolve.cxx:772
 AliMinResSolve.cxx:773
 AliMinResSolve.cxx:774
 AliMinResSolve.cxx:775
 AliMinResSolve.cxx:776
 AliMinResSolve.cxx:777
 AliMinResSolve.cxx:778
 AliMinResSolve.cxx:779
 AliMinResSolve.cxx:780
 AliMinResSolve.cxx:781
 AliMinResSolve.cxx:782
 AliMinResSolve.cxx:783
 AliMinResSolve.cxx:784
 AliMinResSolve.cxx:785
 AliMinResSolve.cxx:786
 AliMinResSolve.cxx:787
 AliMinResSolve.cxx:788
 AliMinResSolve.cxx:789
 AliMinResSolve.cxx:790
 AliMinResSolve.cxx:791
 AliMinResSolve.cxx:792
 AliMinResSolve.cxx:793
 AliMinResSolve.cxx:794
 AliMinResSolve.cxx:795
 AliMinResSolve.cxx:796
 AliMinResSolve.cxx:797
 AliMinResSolve.cxx:798
 AliMinResSolve.cxx:799
 AliMinResSolve.cxx:800
 AliMinResSolve.cxx:801
 AliMinResSolve.cxx:802
 AliMinResSolve.cxx:803
 AliMinResSolve.cxx:804
 AliMinResSolve.cxx:805
 AliMinResSolve.cxx:806
 AliMinResSolve.cxx:807
 AliMinResSolve.cxx:808
 AliMinResSolve.cxx:809
 AliMinResSolve.cxx:810
 AliMinResSolve.cxx:811
 AliMinResSolve.cxx:812
 AliMinResSolve.cxx:813
 AliMinResSolve.cxx:814
 AliMinResSolve.cxx:815
 AliMinResSolve.cxx:816
 AliMinResSolve.cxx:817
 AliMinResSolve.cxx:818
 AliMinResSolve.cxx:819
 AliMinResSolve.cxx:820
 AliMinResSolve.cxx:821
 AliMinResSolve.cxx:822
 AliMinResSolve.cxx:823
 AliMinResSolve.cxx:824
 AliMinResSolve.cxx:825
 AliMinResSolve.cxx:826
 AliMinResSolve.cxx:827
 AliMinResSolve.cxx:828
 AliMinResSolve.cxx:829
 AliMinResSolve.cxx:830
 AliMinResSolve.cxx:831
 AliMinResSolve.cxx:832
 AliMinResSolve.cxx:833
 AliMinResSolve.cxx:834
 AliMinResSolve.cxx:835
 AliMinResSolve.cxx:836
 AliMinResSolve.cxx:837
 AliMinResSolve.cxx:838
 AliMinResSolve.cxx:839
 AliMinResSolve.cxx:840
 AliMinResSolve.cxx:841
 AliMinResSolve.cxx:842
 AliMinResSolve.cxx:843
 AliMinResSolve.cxx:844
 AliMinResSolve.cxx:845
 AliMinResSolve.cxx:846
 AliMinResSolve.cxx:847
 AliMinResSolve.cxx:848
 AliMinResSolve.cxx:849
 AliMinResSolve.cxx:850
 AliMinResSolve.cxx:851
 AliMinResSolve.cxx:852
 AliMinResSolve.cxx:853
 AliMinResSolve.cxx:854
 AliMinResSolve.cxx:855
 AliMinResSolve.cxx:856
 AliMinResSolve.cxx:857
 AliMinResSolve.cxx:858
 AliMinResSolve.cxx:859
 AliMinResSolve.cxx:860
 AliMinResSolve.cxx:861
 AliMinResSolve.cxx:862
 AliMinResSolve.cxx:863
 AliMinResSolve.cxx:864
 AliMinResSolve.cxx:865
 AliMinResSolve.cxx:866
 AliMinResSolve.cxx:867
 AliMinResSolve.cxx:868
 AliMinResSolve.cxx:869
 AliMinResSolve.cxx:870
 AliMinResSolve.cxx:871
 AliMinResSolve.cxx:872
 AliMinResSolve.cxx:873
 AliMinResSolve.cxx:874
 AliMinResSolve.cxx:875
 AliMinResSolve.cxx:876
 AliMinResSolve.cxx:877
 AliMinResSolve.cxx:878
 AliMinResSolve.cxx:879
 AliMinResSolve.cxx:880
 AliMinResSolve.cxx:881
 AliMinResSolve.cxx:882
 AliMinResSolve.cxx:883
 AliMinResSolve.cxx:884
 AliMinResSolve.cxx:885
 AliMinResSolve.cxx:886
 AliMinResSolve.cxx:887
 AliMinResSolve.cxx:888
 AliMinResSolve.cxx:889
 AliMinResSolve.cxx:890
 AliMinResSolve.cxx:891
 AliMinResSolve.cxx:892
 AliMinResSolve.cxx:893
 AliMinResSolve.cxx:894
 AliMinResSolve.cxx:895
 AliMinResSolve.cxx:896
 AliMinResSolve.cxx:897
 AliMinResSolve.cxx:898
 AliMinResSolve.cxx:899
 AliMinResSolve.cxx:900
 AliMinResSolve.cxx:901
 AliMinResSolve.cxx:902
 AliMinResSolve.cxx:903
 AliMinResSolve.cxx:904
 AliMinResSolve.cxx:905
 AliMinResSolve.cxx:906
 AliMinResSolve.cxx:907
 AliMinResSolve.cxx:908
 AliMinResSolve.cxx:909
 AliMinResSolve.cxx:910
 AliMinResSolve.cxx:911
 AliMinResSolve.cxx:912
 AliMinResSolve.cxx:913
 AliMinResSolve.cxx:914
 AliMinResSolve.cxx:915
 AliMinResSolve.cxx:916
 AliMinResSolve.cxx:917
 AliMinResSolve.cxx:918
 AliMinResSolve.cxx:919
 AliMinResSolve.cxx:920
 AliMinResSolve.cxx:921
 AliMinResSolve.cxx:922
 AliMinResSolve.cxx:923
 AliMinResSolve.cxx:924
 AliMinResSolve.cxx:925
 AliMinResSolve.cxx:926
 AliMinResSolve.cxx:927
 AliMinResSolve.cxx:928
 AliMinResSolve.cxx:929
 AliMinResSolve.cxx:930
 AliMinResSolve.cxx:931
 AliMinResSolve.cxx:932
 AliMinResSolve.cxx:933
 AliMinResSolve.cxx:934
 AliMinResSolve.cxx:935
 AliMinResSolve.cxx:936
 AliMinResSolve.cxx:937
 AliMinResSolve.cxx:938
 AliMinResSolve.cxx:939
 AliMinResSolve.cxx:940
 AliMinResSolve.cxx:941
 AliMinResSolve.cxx:942
 AliMinResSolve.cxx:943
 AliMinResSolve.cxx:944
 AliMinResSolve.cxx:945
 AliMinResSolve.cxx:946
 AliMinResSolve.cxx:947
 AliMinResSolve.cxx:948
 AliMinResSolve.cxx:949
 AliMinResSolve.cxx:950
 AliMinResSolve.cxx:951
 AliMinResSolve.cxx:952
 AliMinResSolve.cxx:953
 AliMinResSolve.cxx:954
 AliMinResSolve.cxx:955
 AliMinResSolve.cxx:956
 AliMinResSolve.cxx:957
 AliMinResSolve.cxx:958
 AliMinResSolve.cxx:959
 AliMinResSolve.cxx:960
 AliMinResSolve.cxx:961
 AliMinResSolve.cxx:962
 AliMinResSolve.cxx:963
 AliMinResSolve.cxx:964
 AliMinResSolve.cxx:965
 AliMinResSolve.cxx:966
 AliMinResSolve.cxx:967
 AliMinResSolve.cxx:968
 AliMinResSolve.cxx:969
 AliMinResSolve.cxx:970
 AliMinResSolve.cxx:971
 AliMinResSolve.cxx:972
 AliMinResSolve.cxx:973
 AliMinResSolve.cxx:974
 AliMinResSolve.cxx:975
 AliMinResSolve.cxx:976
 AliMinResSolve.cxx:977
 AliMinResSolve.cxx:978
 AliMinResSolve.cxx:979
 AliMinResSolve.cxx:980
 AliMinResSolve.cxx:981
 AliMinResSolve.cxx:982
 AliMinResSolve.cxx:983
 AliMinResSolve.cxx:984
 AliMinResSolve.cxx:985
 AliMinResSolve.cxx:986
 AliMinResSolve.cxx:987
 AliMinResSolve.cxx:988
 AliMinResSolve.cxx:989
 AliMinResSolve.cxx:990
 AliMinResSolve.cxx:991
 AliMinResSolve.cxx:992
 AliMinResSolve.cxx:993
 AliMinResSolve.cxx:994
 AliMinResSolve.cxx:995
 AliMinResSolve.cxx:996
 AliMinResSolve.cxx:997
 AliMinResSolve.cxx:998
 AliMinResSolve.cxx:999
 AliMinResSolve.cxx:1000
 AliMinResSolve.cxx:1001
 AliMinResSolve.cxx:1002
 AliMinResSolve.cxx:1003
 AliMinResSolve.cxx:1004
 AliMinResSolve.cxx:1005
 AliMinResSolve.cxx:1006
 AliMinResSolve.cxx:1007
 AliMinResSolve.cxx:1008
 AliMinResSolve.cxx:1009
 AliMinResSolve.cxx:1010
 AliMinResSolve.cxx:1011
 AliMinResSolve.cxx:1012
 AliMinResSolve.cxx:1013
 AliMinResSolve.cxx:1014
 AliMinResSolve.cxx:1015
 AliMinResSolve.cxx:1016
 AliMinResSolve.cxx:1017
 AliMinResSolve.cxx:1018
 AliMinResSolve.cxx:1019
 AliMinResSolve.cxx:1020
 AliMinResSolve.cxx:1021
 AliMinResSolve.cxx:1022
 AliMinResSolve.cxx:1023
 AliMinResSolve.cxx:1024
 AliMinResSolve.cxx:1025
 AliMinResSolve.cxx:1026
 AliMinResSolve.cxx:1027
 AliMinResSolve.cxx:1028
 AliMinResSolve.cxx:1029
 AliMinResSolve.cxx:1030
 AliMinResSolve.cxx:1031
 AliMinResSolve.cxx:1032
 AliMinResSolve.cxx:1033
 AliMinResSolve.cxx:1034
 AliMinResSolve.cxx:1035
 AliMinResSolve.cxx:1036
 AliMinResSolve.cxx:1037
 AliMinResSolve.cxx:1038
 AliMinResSolve.cxx:1039
 AliMinResSolve.cxx:1040
 AliMinResSolve.cxx:1041
 AliMinResSolve.cxx:1042
 AliMinResSolve.cxx:1043
 AliMinResSolve.cxx:1044
 AliMinResSolve.cxx:1045
 AliMinResSolve.cxx:1046
 AliMinResSolve.cxx:1047
 AliMinResSolve.cxx:1048
 AliMinResSolve.cxx:1049
 AliMinResSolve.cxx:1050
 AliMinResSolve.cxx:1051
 AliMinResSolve.cxx:1052
 AliMinResSolve.cxx:1053
 AliMinResSolve.cxx:1054
 AliMinResSolve.cxx:1055
 AliMinResSolve.cxx:1056
 AliMinResSolve.cxx:1057
 AliMinResSolve.cxx:1058
 AliMinResSolve.cxx:1059
 AliMinResSolve.cxx:1060
 AliMinResSolve.cxx:1061
 AliMinResSolve.cxx:1062
 AliMinResSolve.cxx:1063
 AliMinResSolve.cxx:1064
 AliMinResSolve.cxx:1065
 AliMinResSolve.cxx:1066
 AliMinResSolve.cxx:1067
 AliMinResSolve.cxx:1068
 AliMinResSolve.cxx:1069
 AliMinResSolve.cxx:1070