#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)
{
}
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)
{
}
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)
{
}
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)
{
}
AliMinResSolve::~AliMinResSolve()
{
ClearAux();
}
AliMinResSolve& AliMinResSolve::operator=(const AliMinResSolve& src)
{
if (this != &src) {
fSize = src.fSize;
fPrecon = src.fPrecon;
fMatrix = src.fMatrix;
fRHS = src.fRHS;
}
return *this;
}
Int_t AliMinResSolve::BuildPrecon(Int_t prec)
{
fPrecon = prec;
if (fPrecon>=kPreconBD && fPrecon<kPreconILU0) {
return BuildPreconBD(fPrecon - kPreconBD);
}
if (fPrecon>=kPreconILU0 && fPrecon<=kPreconILU10) {
if (fMatrix->InheritsFrom("AliMatrixSparse")) return BuildPreconILUK(fPrecon-kPreconILU0);
else return BuildPreconILUKDense(fPrecon-kPreconILU0);
}
return -1;
}
Bool_t AliMinResSolve::SolveFGMRES(TVectorD& VecSol,Int_t precon,int itnlim,double rtol,int nkrylov)
{
return SolveFGMRES(VecSol.GetMatrixArray(),precon,itnlim,rtol,nkrylov);
}
Bool_t AliMinResSolve::SolveFGMRES(double* VecSol,Int_t precon,int itnlim,double rtol,int nkrylov)
{
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;
TStopwatch timer; timer.Start();
while(1) {
fMatrix->MultiplyByVec(VecSol,fPvv[0]);
for (l=fSize;l--;) fPvv[0][l] = fRHS[l] - fPvv[0][l];
beta = 0;
for (l=fSize;l--;) beta += fPvv[0][l]*fPvv[0][l];
beta = TMath::Sqrt(beta);
if (beta < epsmac) break;
t = 1.0 / beta;
for (l=fSize;l--;) fPvv[0][l] *= t;
if (its == 0) eps1 = rtol*beta;
fPVecV[0] = beta;
int i = -1;
do {
i++;
its++;
int i1 = i+1;
if (precon>0) ApplyPrecon( fPvv[i], fPvz[i]);
else for (l=fSize;l--;) fPvz[i][l] = fPvv[i][l];
fMatrix->MultiplyByVec(fPvz[i],fPvv[i1]);
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];
}
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;
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 (gam < epsmac) gam = epsmac;
fPVecR1[i] = fPhh[i][i]/gam;
fPVecR2[i] = fPhh[i][i1]/gam;
fPVecV[i1] =-fPVecR2[i]*fPVecV[i];
fPVecV[i] *= fPVecR1[i];
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) );
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];
}
for (int j=0; j<=i; j++) for (t=fPVecV[j],l=0; l<fSize; l++) VecSol[l] += t*fPvz[j][l];
if (beta <= eps1) {
timer.Stop();
AliInfo(Form("FGMRES converged in %d iterations, CPU time: %.1f sec",its,timer.CpuTime()));
break;
}
if (its >= itnlim) {
timer.Stop();
AliInfo(Form("%d iterations limit exceeded, CPU time: %.1f sec",itnlim,timer.CpuTime()));
status = kFALSE;
break;
}
}
ClearAux();
return status;
}
Bool_t AliMinResSolve::SolveMinRes(TVectorD& VecSol,Int_t precon,int itnlim,double rtol)
{
return SolveMinRes(VecSol.GetMatrixArray(), precon, itnlim, rtol);
}
Bool_t AliMinResSolve::SolveMinRes(double *VecSol,Int_t precon,int itnlim,double rtol)
{
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));
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));
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 );
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) {
itn++;
double s = 1./beta;
for (int i=fSize;i--;) fPVecV[i] = s*fPVecY[i];
fMatrix->MultiplyByVec(fPVecV,fPVecY);
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];
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;
beta = 0; for (int i=fSize;i--;) beta += fPVecR2[i]*fPVecY[i];
if (beta<0) {
AliInfo(Form("Preconditioner is indefinite (%e)",beta));
status = 7;
break;
}
beta = TMath::Sqrt(beta);
tnorm2 += alfa*alfa + oldb*oldb + beta*beta;
if (itn == 1) {
if (beta/beta1 <= 10.0*eps) {
status = 0;
AliInfo("RHS is eigenvector");
}
gmax = TMath::Abs(alfa);
gmin = gmax;
}
oldeps = epsln;
delta = cs * dbar + sn * alfa;
gbar = sn * dbar - cs * alfa;
epsln = sn * beta;
dbar = - cs * beta;
gam = TMath::Sqrt( gbar*gbar + beta*beta );
cs = gbar / gam;
sn = beta / gam;
phi = cs * phibar;
phibar = sn * phibar;
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];
}
gmax = TMath::Max( gmax, gam );
gmin = TMath::Min( gmin, gam );
z = rhs1 / gam;
ynorm2 += z*z;
rhs1 = rhs2 - delta * z;
rhs2 = - epsln * z;
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;
condA = gmax / gmin;
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"));}
}
}
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
{
ApplyPrecon(vecRHS.GetMatrixArray(), vecOut.GetMatrixArray());
}
void AliMinResSolve::ApplyPrecon(const double* vecRHS, double* vecOut) const
{
if (fPrecon>=kPreconBD && fPrecon<kPreconILU0) {
fMatBD->Solve(vecRHS,vecOut);
}
else if (fPrecon>=kPreconILU0 && fPrecon<=kPreconILU10) {
for(int i=0; i<fSize; i++) {
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--; ) {
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()
{
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)
{
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()
{
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)
{
AliInfo(Form("Building Band-Diagonal preconditioner of half-width = %d",hwidth));
fMatBD = new AliSymBDMatrix( fMatrix->GetSize(), hwidth );
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)
{
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;
if ( PreconILUKsymb(lofM)<0 ) {
ClearAux();
return -1;
}
Int_t *jw = new Int_t[fSize];
for(int j=fSize;j--;) jw[j] = -1;
for(int i=0; i<fSize; i++ ) {
if ( (i%int(0.1*fSize)) == 0) {
AliInfo(Form("BuildPrecon: row %d of %d",i,fSize));
sw.Stop();
sw.Print();
sw.Start(kFALSE);
}
AliVectorSparse& rowLi = *fMatL->GetRow(i);
AliVectorSparse& rowUi = *fMatU->GetRow(i);
AliVectorSparse& rowM = *matrix->GetRow(i);
for(int j=rowLi.GetNElems(); j--;) {
int col = rowLi.GetIndex(j);
jw[col] = j;
rowLi.GetElem(j) = 0.;
}
jw[i] = i;
fDiagLU[i] = 0;
for(int j=rowUi.GetNElems();j--;) {
int col = rowUi.GetIndex(j);
jw[col] = j;
rowUi.GetElem(j) = 0;
}
for(int j=rowM.GetNElems(); j--;) {
if (AliMatrixSq::IsZero(rowM.GetElem(j))) continue;
int col = rowM.GetIndex(j);
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++) {
double vl = matrix->Query(col,i);
if (AliMatrixSq::IsZero(vl)) continue;
rowUi.GetElem(jw[col]) = vl;
}
for(int j=0; j<rowLi.GetNElems(); j++) {
int jrow = rowLi.GetIndex(j);
rowLi.GetElem(j) *= fDiagLU[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);
}
}
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)
{
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];
if ( PreconILUKsymbDense(lofM)<0 ) {
ClearAux();
return -1;
}
Int_t *jw = new Int_t[fSize];
for(int j=fSize;j--;) jw[j] = -1;
for(int i=0; i<fSize; i++ ) {
AliVectorSparse& rowLi = *fMatL->GetRow(i);
AliVectorSparse& rowUi = *fMatU->GetRow(i);
for(int j=rowLi.GetNElems(); j--;) {
int col = rowLi.GetIndex(j);
jw[col] = j;
rowLi.GetElem(j) = 0.;
}
jw[i] = i;
fDiagLU[i] = 0;
for(int j=rowUi.GetNElems();j--;) {
int col = rowUi.GetIndex(j);
jw[col] = j;
rowUi.GetElem(j) = 0;
}
for(int j=fSize; j--;) {
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;
}
for(int j=0; j<rowLi.GetNElems(); j++) {
int jrow = rowLi.GetIndex(j);
rowLi.GetElem(j) *= fDiagLU[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);
}
}
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()));
return 0;
}
Int_t AliMinResSolve::PreconILUKsymb(Int_t lofM)
{
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];
levls = new UChar_t[fSize];
jbuf = new UShort_t[fSize];
iw = new Int_t[fSize];
for(int j=fSize; j--;) iw[j] = -1;
for(int i=0; i<fSize; i++) {
int incl = 0;
int incu = i;
AliVectorSparse& row = *matrix->GetRow(i);
for(int j=0;j<row.GetNElems(); j++) {
int col = row.GetIndex(j);
if (AliMatrixSq::IsZero(row.GetElem(j))) continue;
if (col<i) {
jbuf[incl] = col;
levls[incl] = 0;
iw[col] = incl++;
}
else if (col>i) {
jbuf[incu] = col;
levls[incu] = 0;
iw[col] = incu++;
}
}
if (matrix->IsSymmetric()) for (int col=i+1;col<fSize;col++) {
if (AliMatrixSq::IsZero(matrix->Query(col,i))) continue;
jbuf[incu] = col;
levls[incu] = 0;
iw[col] = incu++;
}
int jpiv = -1;
while (++jpiv < incl) {
int k = jbuf[jpiv] ;
int kmin = k;
int jmin = jpiv;
for(int j=jpiv+1; j<incl; j++) if( jbuf[j]<kmin ) { kmin = jbuf[j]; jmin = j;}
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;
}
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);
}
}
for (int j=0;j<incl;j++) iw[jbuf[j]] = -1;
for (int j=i;j<incu;j++) iw[jbuf[j]] = -1;
AliVectorSparse& rowLi = *fMatL->GetRow(i);
rowLi.ReSize(incl);
if(incl>0) memcpy(rowLi.GetIndices(), jbuf, sizeof(UShort_t)*incl);
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];
memcpy( ulvl[i], levls+i, k*sizeof(UChar_t) );
}
}
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)
{
UChar_t **ulvl=0,*levls=0;
UShort_t *jbuf=0;
Int_t *iw=0;
ulvl = new UChar_t*[fSize];
levls = new UChar_t[fSize];
jbuf = new UShort_t[fSize];
iw = new Int_t[fSize];
for(int j=fSize; j--;) iw[j] = -1;
for(int i=0; i<fSize; i++) {
int incl = 0;
int incu = i;
for(int j=0;j<fSize; j++) {
if (AliMatrixSq::IsZero(fMatrix->Query(i,j))) continue;
if (j<i) {
jbuf[incl] = j;
levls[incl] = 0;
iw[j] = incl++;
}
else if (j>i) {
jbuf[incu] = j;
levls[incu] = 0;
iw[j] = incu++;
}
}
int jpiv = -1;
while (++jpiv < incl) {
int k = jbuf[jpiv] ;
int kmin = k;
int jmin = jpiv;
for(int j=jpiv+1; j<incl; j++) if( jbuf[j]<kmin ) { kmin = jbuf[j]; jmin = j;}
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;
}
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);
}
}
for (int j=0;j<incl;j++) iw[jbuf[j]] = -1;
for (int j=i;j<incu;j++) iw[jbuf[j]] = -1;
AliVectorSparse& rowLi = *fMatL->GetRow(i);
rowLi.ReSize(incl);
if(incl>0) memcpy(rowLi.GetIndices(), jbuf, sizeof(UShort_t)*incl);
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];
memcpy( ulvl[i], levls+i, k*sizeof(UChar_t) );
}
}
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;
}