1 /**********************************************************************************************/
2 /* General class for solving large system of linear equations */
3 /* Includes MINRES, FGMRES methods as well as a few precondiotiong methods */
5 /* Author: ruben.shahoyan@cern.ch */
7 /**********************************************************************************************/
9 #include "AliMinResSolve.h"
11 #include "AliMatrixSq.h"
12 #include "AliMatrixSparse.h"
13 #include "AliSymBDMatrix.h"
14 #include <TStopwatch.h>
18 ClassImp(AliMinResSolve)
20 //______________________________________________________
21 AliMinResSolve::AliMinResSolve() :
22 fSize(0),fPrecon(0),fMatrix(0),fRHS(0),
23 fPVecY(0),fPVecR1(0),fPVecR2(0),fPVecV(0),fPVecW(0),fPVecW1(0),fPVecW2(0),
24 fPvv(0),fPvz(0),fPhh(0),
25 fDiagLU(0),fMatL(0),fMatU(0),fMatBD(0)
27 // default constructor
30 //______________________________________________________
31 AliMinResSolve::AliMinResSolve(const AliMinResSolve& src) :
33 fSize(src.fSize),fPrecon(src.fPrecon),fMatrix(src.fMatrix),fRHS(src.fRHS),
34 fPVecY(0),fPVecR1(0),fPVecR2(0),fPVecV(0),fPVecW(0),fPVecW1(0),fPVecW2(0),
35 fPvv(0),fPvz(0),fPhh(0),
36 fDiagLU(0),fMatL(0),fMatU(0),fMatBD(0)
41 //______________________________________________________
42 AliMinResSolve::AliMinResSolve(const AliMatrixSq *mat, const TVectorD* rhs) :
43 fSize(mat->GetSize()),fPrecon(0),fMatrix((AliMatrixSq*)mat),fRHS((double*)rhs->GetMatrixArray()),
44 fPVecY(0),fPVecR1(0),fPVecR2(0),fPVecV(0),fPVecW(0),fPVecW1(0),fPVecW2(0),
45 fPvv(0),fPvz(0),fPhh(0),
46 fDiagLU(0),fMatL(0),fMatU(0),fMatBD(0)
48 // copy accepting equation
51 //______________________________________________________
52 AliMinResSolve::AliMinResSolve(const AliMatrixSq *mat, const double* rhs) :
53 fSize(mat->GetSize()),fPrecon(0),fMatrix((AliMatrixSq*)mat),fRHS((double*)rhs),
54 fPVecY(0),fPVecR1(0),fPVecR2(0),fPVecV(0),fPVecW(0),fPVecW1(0),fPVecW2(0),
55 fPvv(0),fPvz(0),fPhh(0),
56 fDiagLU(0),fMatL(0),fMatU(0),fMatBD(0)
58 // copy accepting equation
61 //______________________________________________________
62 AliMinResSolve::~AliMinResSolve()
68 //______________________________________________________
69 AliMinResSolve& AliMinResSolve::operator=(const AliMinResSolve& src)
74 fPrecon = src.fPrecon;
75 fMatrix = src.fMatrix;
81 //_______________________________________________________________
82 Int_t AliMinResSolve::BuildPrecon(Int_t prec)
84 // preconditioner building
85 // const Double_t kTiny = 1E-12;
88 if (fPrecon>=kPreconBD && fPrecon<kPreconILU0) { // band diagonal decomposition
89 return BuildPreconBD(fPrecon - kPreconBD); // with halfbandwidth + diagonal = fPrecon
92 if (fPrecon>=kPreconILU0 && fPrecon<=kPreconILU10) {
93 if (fMatrix->InheritsFrom("AliMatrixSparse")) return BuildPreconILUK(fPrecon-kPreconILU0);
94 else return BuildPreconILUKDense(fPrecon-kPreconILU0);
100 //________________________________ FGMRES METHODS ________________________________
101 Bool_t AliMinResSolve::SolveFGMRES(TVectorD& VecSol,Int_t precon,int itnlim,double rtol,int nkrylov)
104 return SolveFGMRES(VecSol.GetMatrixArray(),precon,itnlim,rtol,nkrylov);
107 //________________________________________________________________________________
108 Bool_t AliMinResSolve::SolveFGMRES(double* VecSol,Int_t precon,int itnlim,double rtol,int nkrylov)
110 // Adapted from Y.Saad fgmrs.c of ITSOL_1 package by Y.Saad: http://www-users.cs.umn.edu/~saad/software/
111 /*----------------------------------------------------------------------
112 | *** Preconditioned FGMRES ***
113 +-----------------------------------------------------------------------
114 | This is a simple version of the ARMS preconditioned FGMRES algorithm.
115 +-----------------------------------------------------------------------
116 | Y. S. Dec. 2000. -- Apr. 2008
117 +-----------------------------------------------------------------------
118 | VecSol = real vector of length n containing an initial guess to the
119 | precon = precondtioner id (0 = no precon)
120 | itnlim = max n of iterations
121 | rtol = tolerance for stopping iteration
122 | nkrylov = N of Krylov vectors to store
123 +---------------------------------------------------------------------*/
125 double status = kTRUE;
126 double t,beta,eps1=0;
127 const double epsmac = 2.22E-16;
129 AliInfo(Form("Solution by FGMRes: Preconditioner#%d Max.iter.: %d Tol.: %e NKrylov: %d",
130 precon,itnlim,rtol,nkrylov));
134 AliInfo(Form("Changing N Krylov vectors from %d to %d",nkrylov,10));
139 if (precon>=kPreconsTot) {
140 AliInfo(Form("Unknown preconditioner identifier %d, ignore",precon));
143 if (BuildPrecon(precon)<0) {
145 AliInfo("FGMRES failed to build the preconditioner");
151 if (!InitAuxFGMRES(nkrylov)) return kFALSE;
153 for (l=fSize;l--;) VecSol[l] = 0;
155 //-------------------- outer loop starts here
156 TStopwatch timer; timer.Start();
159 //-------------------- compute initial residual vector
160 fMatrix->MultiplyByVec(VecSol,fPvv[0]);
161 for (l=fSize;l--;) fPvv[0][l] = fRHS[l] - fPvv[0][l]; // fPvv[0]= initial residual
163 for (l=fSize;l--;) beta += fPvv[0][l]*fPvv[0][l];
164 beta = TMath::Sqrt(beta);
166 if (beta < epsmac) break; // success?
168 //-------------------- normalize: fPvv[0] = fPvv[0] / beta
169 for (l=fSize;l--;) fPvv[0][l] *= t;
170 if (its == 0) eps1 = rtol*beta;
172 // ** initialize 1-st term of rhs of hessenberg system
180 // (Right) Preconditioning Operation z_{j} = M^{-1} v_{j}
182 if (precon>0) ApplyPrecon( fPvv[i], fPvz[i]);
183 else for (l=fSize;l--;) fPvz[i][l] = fPvv[i][l];
185 //-------------------- matvec operation w = A z_{j} = A M^{-1} v_{j}
186 fMatrix->MultiplyByVec(fPvz[i],fPvv[i1]);
188 // modified gram - schmidt...
189 // h_{i,j} = (w,v_{i})
190 // w = w - h_{i,j} v_{i}
192 for (int j=0; j<=i; j++) {
193 for (t=0, l=fSize;l--;) t+= fPvv[j][l]*fPvv[i1][l];
195 for (l=fSize;l--;) fPvv[i1][l] -= t*fPvv[j][l];
197 // -------------------- h_{j+1,j} = ||w||_{2}
198 for (t=0,l=fSize;l--;) t += fPvv[i1][l]*fPvv[i1][l]; t = TMath::Sqrt(t);
200 if (t > 0) for (t=1./t, l=0; l<fSize; l++) fPvv[i1][l] *= t; // v_{j+1} = w / h_{j+1,j}
202 // done with modified gram schimdt and arnoldi step
203 // now update factorization of fPhh
205 // perform previous transformations on i-th column of h
207 for (l=1; l<=i; l++) {
210 fPhh[i][l1] = fPVecR1[l1]*t + fPVecR2[l1]*fPhh[i][l];
211 fPhh[i][l] = -fPVecR2[l1]*t + fPVecR1[l1]*fPhh[i][l];
213 double gam = TMath::Sqrt( fPhh[i][i]*fPhh[i][i] + fPhh[i][i1]*fPhh[i][i1]);
215 // if gamma is zero then any small value will do...
216 // will affect only residual estimate
217 if (gam < epsmac) gam = epsmac;
218 // get next plane rotation
219 fPVecR1[i] = fPhh[i][i]/gam;
220 fPVecR2[i] = fPhh[i][i1]/gam;
221 fPVecV[i1] =-fPVecR2[i]*fPVecV[i];
222 fPVecV[i] *= fPVecR1[i];
224 // determine residual norm and test for convergence
225 fPhh[i][i] = fPVecR1[i]*fPhh[i][i] + fPVecR2[i]*fPhh[i][i1];
226 beta = TMath::Abs(fPVecV[i1]);
228 } while ( (i < nkrylov-1) && (beta > eps1) && (its < itnlim) );
230 // now compute solution. 1st, solve upper triangular system
231 fPVecV[i] = fPVecV[i]/fPhh[i][i];
232 for (int j=1; j<=i; j++) {
234 for (t=fPVecV[k],l=k+1; l<=i; l++) t -= fPhh[l][k]*fPVecV[l];
235 fPVecV[k] = t/fPhh[k][k];
237 // -------------------- linear combination of v[i]'s to get sol.
238 for (int j=0; j<=i; j++) for (t=fPVecV[j],l=0; l<fSize; l++) VecSol[l] += t*fPvz[j][l];
240 // -------------------- restart outer loop if needed
244 AliInfo(Form("FGMRES converged in %d iterations, CPU time: %.1f sec",its,timer.CpuTime()));
250 AliInfo(Form("%d iterations limit exceeded, CPU time: %.1f sec",itnlim,timer.CpuTime()));
262 //________________________________ MINRES METHODS ________________________________
263 Bool_t AliMinResSolve::SolveMinRes(TVectorD& VecSol,Int_t precon,int itnlim,double rtol)
266 return SolveMinRes(VecSol.GetMatrixArray(), precon, itnlim, rtol);
269 //________________________________________________________________________________
270 Bool_t AliMinResSolve::SolveMinRes(double *VecSol,Int_t precon,int itnlim,double rtol)
273 Adapted from author's Fortran code:
274 Michael A. Saunders na.msaunders@na-net.ornl.gov
276 MINRES is an implementation of the algorithm described in the following reference:
277 C. C. Paige and M. A. Saunders (1975),
278 Solution of sparse indefinite systems of linear equations,
279 SIAM J. Numer. Anal. 12(4), pp. 617-629.
282 if (!fMatrix->IsSymmetric()) {
283 AliInfo("MinRes cannot solve asymmetric matrices, use FGMRes instead");
288 const double eps = 2.22E-16;
292 if (precon>=kPreconsTot) {
293 AliInfo(Form("Unknown preconditioner identifier %d, ignore",precon));
296 if (BuildPrecon(precon)<0) {
298 AliInfo("MinRes failed to build the preconditioner");
303 AliInfo(Form("Solution by MinRes: Preconditioner#%d Max.iter.: %d Tol.: %e",precon,itnlim,rtol));
305 // ------------------------ initialization ---------------------->>>>
306 memset(VecSol,0,fSize*sizeof(double));
312 double gam,gmax=1,gmin=1,gbar,oldeps,epsa,epsx,epsr,diag, delta,phi,denom,z;
314 if (!InitAuxMinRes()) return kFALSE;
316 memset(VecSol,0,fSize*sizeof(double));
318 // ------------ init aux -------------------------<<<<
319 // Set up y and v for the first Lanczos vector v1.
320 // y = beta1 P' v1, where P = C**(-1). v is really P' v1.
322 for (int i=fSize;i--;) fPVecY[i] = fPVecR1[i] = fRHS[i];
324 if ( precon>0 ) ApplyPrecon( fRHS, fPVecY);
325 beta1 = 0; for (int i=fSize;i--;) beta1 += fRHS[i]*fPVecY[i]; //
328 AliInfo(Form("Preconditioner is indefinite (init) (%e).",beta1));
335 AliInfo(Form("RHS is zero or is the nullspace of the Preconditioner: Solution is {0}"));
340 beta1 = TMath::Sqrt( beta1 ); // Normalize y to get v1 later.
342 // See if Msolve is symmetric. //RS: Skept
343 // See if Aprod is symmetric. //RS: Skept
349 double qrnorm = beta1;
350 double phibar = beta1;
357 for (int i=fSize;i--;) fPVecR2[i] = fPVecR1[i];
359 TStopwatch timer; timer.Start();
360 while(status==0) { //----------------- Main iteration loop ---------------------->>>>
363 /*-----------------------------------------------------------------
364 Obtain quantities for the next Lanczos vector vk+1, k = 1, 2,...
365 The general iteration is similar to the case k = 1 with v0 = 0:
366 p1 = Operator * v1 - beta1 * v0,
368 q2 = p2 - alpha1 * v1,
371 Again, y = betak P vk, where P = C**(-1).
372 .... more description needed.
373 -----------------------------------------------------------------*/
375 double s = 1./beta; // Normalize previous vector (in y).
376 for (int i=fSize;i--;) fPVecV[i] = s*fPVecY[i]; // v = vk if P = I
378 fMatrix->MultiplyByVec(fPVecV,fPVecY); // APROD (VecV, VecY);
381 double btrat = beta/oldb;
382 for (int i=fSize;i--;) fPVecY[i] -= btrat*fPVecR1[i];
384 double alfa = 0; for (int i=fSize;i--;) alfa += fPVecV[i]*fPVecY[i]; // alphak
386 double alf2bt = alfa/beta;
387 for (int i=fSize;i--;) {
388 fPVecY[i] -= alf2bt*fPVecR2[i];
389 fPVecR1[i] = fPVecR2[i];
390 fPVecR2[i] = fPVecY[i];
393 if ( precon>0 ) ApplyPrecon(fPVecR2, fPVecY);
395 oldb = beta; // oldb = betak
396 beta = 0; for (int i=fSize;i--;) beta += fPVecR2[i]*fPVecY[i]; // beta = betak+1^2
399 AliInfo(Form("Preconditioner is indefinite (%e)",beta));
404 beta = TMath::Sqrt(beta); // beta = betak+1
405 tnorm2 += alfa*alfa + oldb*oldb + beta*beta;
407 if (itn == 1) { // Initialize a few things.
408 if (beta/beta1 <= 10.0*eps) {
409 status = 0; //-1 //????? beta2 = 0 or ~ 0, terminate later.
410 AliInfo("RHS is eigenvector");
413 gmax = TMath::Abs(alfa); // alpha1
414 gmin = gmax; // alpha1
418 Apply previous rotation Qk-1 to get
419 [deltak epslnk+1] = [cs sn][dbark 0 ]
420 [gbar k dbar k+1] [sn -cs][alfak betak+1].
424 delta = cs * dbar + sn * alfa; // delta1 = 0 deltak
425 gbar = sn * dbar - cs * alfa; // gbar 1 = alfa1 gbar k
426 epsln = sn * beta; // epsln2 = 0 epslnk+1
427 dbar = - cs * beta; // dbar 2 = beta2 dbar k+1
429 // Compute the next plane rotation Qk
431 gam = TMath::Sqrt( gbar*gbar + beta*beta ); // gammak
432 cs = gbar / gam; // ck
433 sn = beta / gam; // sk
434 phi = cs * phibar; // phik
435 phibar = sn * phibar; // phibark+1
440 for (int i=fSize;i--;) {
441 fPVecW1[i] = fPVecW2[i];
442 fPVecW2[i] = fPVecW[i];
443 fPVecW[i] = denom*(fPVecV[i]-oldeps*fPVecW1[i]-delta*fPVecW2[i]);
444 VecSol[i] += phi*fPVecW[i];
449 gmax = TMath::Max( gmax, gam );
450 gmin = TMath::Min( gmin, gam );
453 rhs1 = rhs2 - delta * z;
456 // Estimate various norms and test for convergence.
457 normA = TMath::Sqrt( tnorm2 );
458 ynorm = TMath::Sqrt( ynorm2 );
460 epsx = normA * ynorm * eps;
461 epsr = normA * ynorm * rtol;
463 if (diag == 0) diag = epsa;
469 In this version we look at the diagonals of R in the
470 factorization of the lower Hessenberg matrix, Q * H = R,
471 where H is the tridiagonal matrix from Lanczos with one
472 extra row, beta(k+1) e_k^T.
476 // See if any of the stopping criteria are satisfied.
477 // In rare cases, istop is already -1 from above (Abar = const*I).
479 AliDebug(2,Form("#%5d |qnrm: %+.2e Anrm:%+.2e Cnd:%+.2e Rnrm:%+.2e Ynrm:%+.2e EpsR:%+.2e EpsX:%+.2e Beta1:%+.2e",
480 itn,qrnorm, normA,condA,rnorm,ynorm,epsr,epsx,beta1));
483 if (itn >= itnlim ) {status = 5; AliInfo(Form("%d iterations limit exceeded",itnlim));}
484 if (condA >= 0.1/eps ) {status = 4; AliInfo(Form("Matrix condition nymber %e exceeds limit %e",condA,0.1/eps));}
485 if (epsx >= beta1 ) {status = 3; AliInfo(Form("Approximate convergence"));}
486 if (qrnorm <= epsx ) {status = 2; AliInfo(Form("Converged within machine precision"));}
487 if (qrnorm <= epsr ) {status = 1; AliInfo(Form("Converged"));}
490 } //----------------- Main iteration loop ----------------------<<<
495 AliInfo(Form("Exit from MinRes: CPU time: %.2f sec\n"
502 timer.CpuTime(),status,itn,normA,condA,rnorm,ynorm));
504 return status>=0 && status<=3;
508 //______________________________________________________________
509 void AliMinResSolve::ApplyPrecon(const TVectorD& vecRHS, TVectorD& vecOut) const
512 ApplyPrecon(vecRHS.GetMatrixArray(), vecOut.GetMatrixArray());
515 //______________________________________________________________
516 void AliMinResSolve::ApplyPrecon(const double* vecRHS, double* vecOut) const
518 // Application of the preconditioner matrix:
519 // implicitly defines the matrix solving the M*VecOut = VecRHS
520 // const Double_t kTiny = 1E-12;
521 if (fPrecon>=kPreconBD && fPrecon<kPreconILU0) { // band diagonal decomposition
522 fMatBD->Solve(vecRHS,vecOut);
526 else if (fPrecon>=kPreconILU0 && fPrecon<=kPreconILU10) {
528 for(int i=0; i<fSize; i++) { // Block L solve
529 vecOut[i] = vecRHS[i];
530 AliVectorSparse &rowLi = *fMatL->GetRow(i);
531 int n = rowLi.GetNElems();
532 for(int j=0;j<n;j++) vecOut[i] -= vecOut[ rowLi.GetIndex(j) ] * rowLi.GetElem(j);
536 for(int i=fSize; i--; ) { // Block -- U solve
537 AliVectorSparse &rowUi = *fMatU->GetRow(i);
538 int n = rowUi.GetNElems();
539 for(int j=0;j<n;j++ ) vecOut[i] -= vecOut[ rowUi.GetIndex(j) ] * rowUi.GetElem(j);
540 vecOut[i] *= fDiagLU[i];
548 //___________________________________________________________
549 Bool_t AliMinResSolve::InitAuxMinRes()
551 // init auxiliary space for minres
552 fPVecY = new double[fSize];
553 fPVecR1 = new double[fSize];
554 fPVecR2 = new double[fSize];
555 fPVecV = new double[fSize];
556 fPVecW = new double[fSize];
557 fPVecW1 = new double[fSize];
558 fPVecW2 = new double[fSize];
560 for (int i=fSize;i--;) fPVecY[i]=fPVecR1[i]=fPVecR2[i]=fPVecV[i]=fPVecW[i]=fPVecW1[i]=fPVecW2[i]=0.0;
566 //___________________________________________________________
567 Bool_t AliMinResSolve::InitAuxFGMRES(int nkrylov)
569 // init auxiliary space for fgmres
570 fPvv = new double*[nkrylov+1];
571 fPvz = new double*[nkrylov];
572 for (int i=0; i<=nkrylov; i++) fPvv[i] = new double[fSize];
573 fPhh = new double*[nkrylov];
574 for (int i=0; i<nkrylov; i++) {
575 fPhh[i] = new double[i+2];
576 fPvz[i] = new double[fSize];
579 fPVecR1 = new double[nkrylov];
580 fPVecR2 = new double[nkrylov];
581 fPVecV = new double[nkrylov+1];
587 //___________________________________________________________
588 void AliMinResSolve::ClearAux()
591 if (fPVecY) delete[] fPVecY; fPVecY = 0;
592 if (fPVecR1) delete[] fPVecR1; fPVecR1 = 0;
593 if (fPVecR2) delete[] fPVecR2; fPVecR2 = 0;
594 if (fPVecV) delete[] fPVecV; fPVecV = 0;
595 if (fPVecW) delete[] fPVecW; fPVecW = 0;
596 if (fPVecW1) delete[] fPVecW1; fPVecW1 = 0;
597 if (fPVecW2) delete[] fPVecW2; fPVecW2 = 0;
598 if (fDiagLU) delete[] fDiagLU; fDiagLU = 0;
599 if (fMatL) delete fMatL; fMatL = 0;
600 if (fMatU) delete fMatU; fMatU = 0;
601 if (fMatBD) delete fMatBD; fMatBD = 0;
604 //___________________________________________________________
605 Int_t AliMinResSolve::BuildPreconBD(Int_t hwidth)
607 // build preconditioner
608 AliInfo(Form("Building Band-Diagonal preconditioner of half-width = %d",hwidth));
609 fMatBD = new AliSymBDMatrix( fMatrix->GetSize(), hwidth );
611 // fill the band-diagonal part of the matrix
612 if (fMatrix->InheritsFrom("AliMatrixSparse")) {
613 for (int ir=fMatrix->GetSize();ir--;) {
614 int jmin = TMath::Max(0,ir-hwidth);
615 AliVectorSparse& irow = *((AliMatrixSparse*)fMatrix)->GetRow(ir);
616 for(int j=irow.GetNElems();j--;) {
617 int jind = irow.GetIndex(j);
618 if (jind<jmin) break;
619 (*fMatBD)(ir,jind) = irow.GetElem(j);
624 for (int ir=fMatrix->GetSize();ir--;) {
625 int jmin = TMath::Max(0,ir-hwidth);
626 for(int jr=jmin;jr<=ir;jr++) (*fMatBD)(ir,jr) = fMatrix->Query(ir,jr);
630 fMatBD->DecomposeLDLT();
635 //___________________________________________________________
636 Int_t AliMinResSolve::BuildPreconILUK(Int_t lofM)
638 /*----------------------------------------------------------------------------
639 * ILUK preconditioner
640 * incomplete LU factorization with level of fill dropping
641 * Adapted from iluk.c of ITSOL_1 package by Y.Saad: http://www-users.cs.umn.edu/~saad/software/
642 *----------------------------------------------------------------------------*/
644 AliInfo(Form("Building ILU%d preconditioner",lofM));
646 TStopwatch sw; sw.Start();
647 fMatL = new AliMatrixSparse(fSize);
648 fMatU = new AliMatrixSparse(fSize);
649 fMatL->SetSymmetric(kFALSE);
650 fMatU->SetSymmetric(kFALSE);
651 fDiagLU = new Double_t[fSize];
652 AliMatrixSparse* matrix = (AliMatrixSparse*)fMatrix;
654 // symbolic factorization to calculate level of fill index arrays
655 if ( PreconILUKsymb(lofM)<0 ) {
660 Int_t *jw = new Int_t[fSize];
661 for(int j=fSize;j--;) jw[j] = -1; // set indicator array jw to -1
663 for(int i=0; i<fSize; i++ ) { // beginning of main loop
664 if ( (i%int(0.1*fSize)) == 0) {
665 AliInfo(Form("BuildPrecon: row %d of %d",i,fSize));
670 /* setup array jw[], and initial i-th row */
671 AliVectorSparse& rowLi = *fMatL->GetRow(i);
672 AliVectorSparse& rowUi = *fMatU->GetRow(i);
673 AliVectorSparse& rowM = *matrix->GetRow(i);
675 for(int j=rowLi.GetNElems(); j--;) { // initialize L part
676 int col = rowLi.GetIndex(j);
678 rowLi.GetElem(j) = 0.; // do we need this ?
681 fDiagLU[i] = 0; // initialize diagonal
683 for(int j=rowUi.GetNElems();j--;) { // initialize U part
684 int col = rowUi.GetIndex(j);
686 rowUi.GetElem(j) = 0;
688 // copy row from csmat into L,U D
689 for(int j=rowM.GetNElems(); j--;) { // L and D part
690 if (AliMatrixSq::IsZero(rowM.GetElem(j))) continue;
691 int col = rowM.GetIndex(j); // (the original matrix stores only lower triangle)
692 if( col < i ) rowLi.GetElem(jw[col]) = rowM.GetElem(j);
693 else if(col==i) fDiagLU[i] = rowM.GetElem(j);
694 else rowUi.GetElem(jw[col]) = rowM.GetElem(j);
696 if (matrix->IsSymmetric()) for (int col=i+1;col<fSize;col++) { // part of the row I on the right of diagonal is stored as
697 double vl = matrix->Query(col,i); // the lower part of the column I
698 if (AliMatrixSq::IsZero(vl)) continue;
699 rowUi.GetElem(jw[col]) = vl;
702 // eliminate previous rows
703 for(int j=0; j<rowLi.GetNElems(); j++) {
704 int jrow = rowLi.GetIndex(j);
705 // get the multiplier for row to be eliminated (jrow)
706 rowLi.GetElem(j) *= fDiagLU[jrow];
708 // combine current row and row jrow
709 AliVectorSparse& rowUj = *fMatU->GetRow(jrow);
710 for(int k=0; k<rowUj.GetNElems(); k++ ) {
711 int col = rowUj.GetIndex(k);
713 if( jpos == -1 ) continue;
714 if( col < i ) rowLi.GetElem(jpos) -= rowLi.GetElem(j) * rowUj.GetElem(k);
715 else if(col==i) fDiagLU[i] -= rowLi.GetElem(j) * rowUj.GetElem(k);
716 else rowUi.GetElem(jpos) -= rowLi.GetElem(j) * rowUj.GetElem(k);
719 // reset double-pointer to -1 ( U-part)
720 for(int j=rowLi.GetNElems(); j--;) jw[ rowLi.GetIndex(j) ] = -1;
722 for(int j=rowUi.GetNElems(); j--;) jw[ rowUi.GetIndex(j) ] = -1;
724 if( AliMatrixSq::IsZero(fDiagLU[i]) ) {
725 AliInfo(Form("Fatal error in ILIk: Zero diagonal found..."));
729 fDiagLU[i] = 1.0 / fDiagLU[i];
735 AliInfo(Form("ILU%d preconditioner OK, CPU time: %.1f sec",lofM,sw.CpuTime()));
736 AliInfo(Form("Densities: M %f L %f U %f",matrix->GetDensity(),fMatL->GetDensity(),fMatU->GetDensity()));
741 //___________________________________________________________
742 Int_t AliMinResSolve::BuildPreconILUKDense(Int_t lofM)
744 /*----------------------------------------------------------------------------
745 * ILUK preconditioner
746 * incomplete LU factorization with level of fill dropping
747 * Adapted from iluk.c of ITSOL_1 package by Y.Saad: http://www-users.cs.umn.edu/~saad/software/
748 *----------------------------------------------------------------------------*/
750 TStopwatch sw; sw.Start();
751 AliInfo(Form("Building ILU%d preconditioner for dense matrix",lofM));
753 fMatL = new AliMatrixSparse(fSize);
754 fMatU = new AliMatrixSparse(fSize);
755 fMatL->SetSymmetric(kFALSE);
756 fMatU->SetSymmetric(kFALSE);
757 fDiagLU = new Double_t[fSize];
759 // symbolic factorization to calculate level of fill index arrays
760 if ( PreconILUKsymbDense(lofM)<0 ) {
765 Int_t *jw = new Int_t[fSize];
766 for(int j=fSize;j--;) jw[j] = -1; // set indicator array jw to -1
768 for(int i=0; i<fSize; i++ ) { // beginning of main loop
769 /* setup array jw[], and initial i-th row */
770 AliVectorSparse& rowLi = *fMatL->GetRow(i);
771 AliVectorSparse& rowUi = *fMatU->GetRow(i);
773 for(int j=rowLi.GetNElems(); j--;) { // initialize L part
774 int col = rowLi.GetIndex(j);
776 rowLi.GetElem(j) = 0.; // do we need this ?
779 fDiagLU[i] = 0; // initialize diagonal
781 for(int j=rowUi.GetNElems();j--;) { // initialize U part
782 int col = rowUi.GetIndex(j);
784 rowUi.GetElem(j) = 0;
786 // copy row from csmat into L,U D
787 for(int j=fSize; j--;) { // L and D part
788 double vl = fMatrix->Query(i,j);
789 if (AliMatrixSq::IsZero(vl)) continue;
790 if( j < i ) rowLi.GetElem(jw[j]) = vl;
791 else if(j==i) fDiagLU[i] = vl;
792 else rowUi.GetElem(jw[j]) = vl;
794 // eliminate previous rows
795 for(int j=0; j<rowLi.GetNElems(); j++) {
796 int jrow = rowLi.GetIndex(j);
797 // get the multiplier for row to be eliminated (jrow)
798 rowLi.GetElem(j) *= fDiagLU[jrow];
800 // combine current row and row jrow
801 AliVectorSparse& rowUj = *fMatU->GetRow(jrow);
802 for(int k=0; k<rowUj.GetNElems(); k++ ) {
803 int col = rowUj.GetIndex(k);
805 if( jpos == -1 ) continue;
806 if( col < i ) rowLi.GetElem(jpos) -= rowLi.GetElem(j) * rowUj.GetElem(k);
807 else if(col==i) fDiagLU[i] -= rowLi.GetElem(j) * rowUj.GetElem(k);
808 else rowUi.GetElem(jpos) -= rowLi.GetElem(j) * rowUj.GetElem(k);
811 // reset double-pointer to -1 ( U-part)
812 for(int j=rowLi.GetNElems(); j--;) jw[ rowLi.GetIndex(j) ] = -1;
814 for(int j=rowUi.GetNElems(); j--;) jw[ rowUi.GetIndex(j) ] = -1;
816 if( AliMatrixSq::IsZero(fDiagLU[i])) {
817 AliInfo(Form("Fatal error in ILIk: Zero diagonal found..."));
821 fDiagLU[i] = 1.0 / fDiagLU[i];
827 AliInfo(Form("ILU%d dense preconditioner OK, CPU time: %.1f sec",lofM,sw.CpuTime()));
828 // AliInfo(Form("Densities: M %f L %f U %f",fMatrix->GetDensity(),fMatL->GetDensity(),fMatU->GetDensity()));
833 //___________________________________________________________
834 Int_t AliMinResSolve::PreconILUKsymb(Int_t lofM)
836 /*----------------------------------------------------------------------------
837 * ILUK preconditioner
838 * incomplete LU factorization with level of fill dropping
839 * Adapted from iluk.c: lofC of ITSOL_1 package by Y.Saad: http://www-users.cs.umn.edu/~saad/software/
840 *----------------------------------------------------------------------------*/
843 AliInfo("PreconILUKsymb>>");
844 AliMatrixSparse* matrix = (AliMatrixSparse*)fMatrix;
847 UChar_t **ulvl=0,*levls=0;
850 ulvl = new UChar_t*[fSize]; // stores lev-fils for U part of ILU factorization
851 levls = new UChar_t[fSize];
852 jbuf = new UShort_t[fSize];
853 iw = new Int_t[fSize];
855 for(int j=fSize; j--;) iw[j] = -1; // initialize iw
856 for(int i=0; i<fSize; i++) {
859 AliVectorSparse& row = *matrix->GetRow(i);
861 // assign lof = 0 for matrix elements
862 for(int j=0;j<row.GetNElems(); j++) {
863 int col = row.GetIndex(j);
864 if (AliMatrixSq::IsZero(row.GetElem(j))) continue; // !!!! matrix is sparse but sometimes 0 appears
865 if (col<i) { // L-part
870 else if (col>i) { // This works only for general matrix
876 if (matrix->IsSymmetric()) for (int col=i+1;col<fSize;col++) { // U-part of symmetric matrix
877 if (AliMatrixSq::IsZero(matrix->Query(col,i))) continue; // Due to the symmetry == matrix(i,col)
883 // symbolic k,i,j Gaussian elimination
885 while (++jpiv < incl) {
886 int k = jbuf[jpiv] ; // select leftmost pivot
889 for(int j=jpiv+1; j<incl; j++) if( jbuf[j]<kmin ) { kmin = jbuf[j]; jmin = j;}
891 // ------------------------------------ swap
897 int tj = levls[jpiv] ;
898 levls[jpiv] = levls[jmin];
902 // ------------------------------------ symbolic linear combinaiton of rows
903 AliVectorSparse& rowU = *fMatU->GetRow(k);
904 for(int j=0; j<rowU.GetNElems(); j++ ) {
905 int col = rowU.GetIndex(j);
906 int it = ulvl[k][j]+levls[jpiv]+1;
907 if( it > lofM ) continue;
921 else levls[ip] = TMath::Min(levls[ip], it);
924 } // end - while loop
927 for (int j=0;j<incl;j++) iw[jbuf[j]] = -1;
928 for (int j=i;j<incu;j++) iw[jbuf[j]] = -1;
931 AliVectorSparse& rowLi = *fMatL->GetRow(i);
933 if(incl>0) memcpy(rowLi.GetIndices(), jbuf, sizeof(UShort_t)*incl);
936 AliVectorSparse& rowUi = *fMatU->GetRow(i);
939 memcpy(rowUi.GetIndices(), jbuf+i, sizeof(UShort_t)*k);
940 ulvl[i] = new UChar_t[k]; // update matrix of levels
941 memcpy( ulvl[i], levls+i, k*sizeof(UChar_t) );
945 // free temp space and leave
948 for(int i=fSize; i--;) if (fMatU->GetRow(i)->GetNElems()) delete[] ulvl[i];
952 fMatL->SortIndices();
953 fMatU->SortIndices();
956 AliInfo("PreconILUKsymb<<");
961 //___________________________________________________________
962 Int_t AliMinResSolve::PreconILUKsymbDense(Int_t lofM)
964 /*----------------------------------------------------------------------------
965 * ILUK preconditioner
966 * incomplete LU factorization with level of fill dropping
967 * Adapted from iluk.c: lofC of ITSOL_1 package by Y.Saad: http://www-users.cs.umn.edu/~saad/software/
968 *----------------------------------------------------------------------------*/
970 UChar_t **ulvl=0,*levls=0;
973 ulvl = new UChar_t*[fSize]; // stores lev-fils for U part of ILU factorization
974 levls = new UChar_t[fSize];
975 jbuf = new UShort_t[fSize];
976 iw = new Int_t[fSize];
978 for(int j=fSize; j--;) iw[j] = -1; // initialize iw
979 for(int i=0; i<fSize; i++) {
983 // assign lof = 0 for matrix elements
984 for(int j=0;j<fSize; j++) {
985 if (AliMatrixSq::IsZero(fMatrix->Query(i,j))) continue;
991 else if (j>i) { // This works only for general matrix
998 // symbolic k,i,j Gaussian elimination
1000 while (++jpiv < incl) {
1001 int k = jbuf[jpiv] ; // select leftmost pivot
1004 for(int j=jpiv+1; j<incl; j++) if( jbuf[j]<kmin ) { kmin = jbuf[j]; jmin = j;}
1006 // ------------------------------------ swap
1012 int tj = levls[jpiv] ;
1013 levls[jpiv] = levls[jmin];
1017 // ------------------------------------ symbolic linear combinaiton of rows
1018 AliVectorSparse& rowU = *fMatU->GetRow(k);
1019 for(int j=0; j<rowU.GetNElems(); j++ ) {
1020 int col = rowU.GetIndex(j);
1021 int it = ulvl[k][j]+levls[jpiv]+1;
1022 if( it > lofM ) continue;
1030 else if( col > i ) {
1036 else levls[ip] = TMath::Min(levls[ip], it);
1039 } // end - while loop
1042 for (int j=0;j<incl;j++) iw[jbuf[j]] = -1;
1043 for (int j=i;j<incu;j++) iw[jbuf[j]] = -1;
1046 AliVectorSparse& rowLi = *fMatL->GetRow(i);
1048 if(incl>0) memcpy(rowLi.GetIndices(), jbuf, sizeof(UShort_t)*incl);
1051 AliVectorSparse& rowUi = *fMatU->GetRow(i);
1054 memcpy(rowUi.GetIndices(), jbuf+i, sizeof(UShort_t)*k);
1055 ulvl[i] = new UChar_t[k]; // update matrix of levels
1056 memcpy( ulvl[i], levls+i, k*sizeof(UChar_t) );
1060 // free temp space and leave
1063 for(int i=fSize; i--;) if (fMatU->GetRow(i)->GetNElems()) delete[] ulvl[i];
1067 fMatL->SortIndices();
1068 fMatU->SortIndices();