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>
20 ClassImp(AliMinResSolve)
22 //______________________________________________________
23 AliMinResSolve::AliMinResSolve() :
24 fSize(0),fPrecon(0),fMatrix(0),fRHS(0),
25 fPVecY(0),fPVecR1(0),fPVecR2(0),fPVecV(0),fPVecW(0),fPVecW1(0),fPVecW2(0),
26 fPvv(0),fPvz(0),fPhh(0),
27 fDiagLU(0),fMatL(0),fMatU(0),fMatBD(0)
29 // default constructor
32 //______________________________________________________
33 AliMinResSolve::AliMinResSolve(const AliMinResSolve& src) :
35 fSize(src.fSize),fPrecon(src.fPrecon),fMatrix(src.fMatrix),fRHS(src.fRHS),
36 fPVecY(0),fPVecR1(0),fPVecR2(0),fPVecV(0),fPVecW(0),fPVecW1(0),fPVecW2(0),
37 fPvv(0),fPvz(0),fPhh(0),
38 fDiagLU(0),fMatL(0),fMatU(0),fMatBD(0)
43 //______________________________________________________
44 AliMinResSolve::AliMinResSolve(const AliMatrixSq *mat, const TVectorD* rhs) :
45 fSize(mat->GetSize()),fPrecon(0),fMatrix((AliMatrixSq*)mat),fRHS((double*)rhs->GetMatrixArray()),
46 fPVecY(0),fPVecR1(0),fPVecR2(0),fPVecV(0),fPVecW(0),fPVecW1(0),fPVecW2(0),
47 fPvv(0),fPvz(0),fPhh(0),
48 fDiagLU(0),fMatL(0),fMatU(0),fMatBD(0)
50 // copy accepting equation
53 //______________________________________________________
54 AliMinResSolve::AliMinResSolve(const AliMatrixSq *mat, const double* rhs) :
55 fSize(mat->GetSize()),fPrecon(0),fMatrix((AliMatrixSq*)mat),fRHS((double*)rhs),
56 fPVecY(0),fPVecR1(0),fPVecR2(0),fPVecV(0),fPVecW(0),fPVecW1(0),fPVecW2(0),
57 fPvv(0),fPvz(0),fPhh(0),
58 fDiagLU(0),fMatL(0),fMatU(0),fMatBD(0)
60 // copy accepting equation
63 //______________________________________________________
64 AliMinResSolve::~AliMinResSolve()
70 //______________________________________________________
71 AliMinResSolve& AliMinResSolve::operator=(const AliMinResSolve& src)
76 fPrecon = src.fPrecon;
77 fMatrix = src.fMatrix;
83 //_______________________________________________________________
84 Int_t AliMinResSolve::BuildPrecon(Int_t prec)
86 // preconditioner building
87 // const Double_t kTiny = 1E-12;
90 if (fPrecon>=kPreconBD && fPrecon<kPreconILU0) { // band diagonal decomposition
91 return BuildPreconBD(fPrecon - kPreconBD); // with halfbandwidth + diagonal = fPrecon
94 if (fPrecon>=kPreconILU0 && fPrecon<=kPreconILU10) {
95 if (fMatrix->InheritsFrom("AliMatrixSparse")) return BuildPreconILUK(fPrecon-kPreconILU0);
96 else return BuildPreconILUKDense(fPrecon-kPreconILU0);
102 //________________________________ FGMRES METHODS ________________________________
103 Bool_t AliMinResSolve::SolveFGMRES(TVectorD& VecSol,Int_t precon,int itnlim,double rtol,int nkrylov)
106 return SolveFGMRES(VecSol.GetMatrixArray(),precon,itnlim,rtol,nkrylov);
109 //________________________________________________________________________________
110 Bool_t AliMinResSolve::SolveFGMRES(double* VecSol,Int_t precon,int itnlim,double rtol,int nkrylov)
112 // Adapted from Y.Saad fgmrs.c of ITSOL_1 package by Y.Saad: http://www-users.cs.umn.edu/~saad/software/
113 /*----------------------------------------------------------------------
114 | *** Preconditioned FGMRES ***
115 +-----------------------------------------------------------------------
116 | This is a simple version of the ARMS preconditioned FGMRES algorithm.
117 +-----------------------------------------------------------------------
118 | Y. S. Dec. 2000. -- Apr. 2008
119 +-----------------------------------------------------------------------
120 | VecSol = real vector of length n containing an initial guess to the
121 | precon = precondtioner id (0 = no precon)
122 | itnlim = max n of iterations
123 | rtol = tolerance for stopping iteration
124 | nkrylov = N of Krylov vectors to store
125 +---------------------------------------------------------------------*/
127 double status = kTRUE;
128 double t,beta,eps1=0;
129 const double epsmac = 2.22E-16;
131 AliInfo(Form("Solution by FGMRes: Preconditioner#%d Max.iter.: %d Tol.: %e NKrylov: %d",
132 precon,itnlim,rtol,nkrylov));
136 AliInfo(Form("Changing N Krylov vectors from %d to %d",nkrylov,10));
141 if (precon>=kPreconsTot) {
142 AliInfo(Form("Unknown preconditioner identifier %d, ignore",precon));
145 if (BuildPrecon(precon)<0) {
147 AliInfo("FGMRES failed to build the preconditioner");
153 if (!InitAuxFGMRES(nkrylov)) return kFALSE;
155 for (l=fSize;l--;) VecSol[l] = 0;
157 //-------------------- outer loop starts here
158 TStopwatch timer; timer.Start();
161 //-------------------- compute initial residual vector
162 fMatrix->MultiplyByVec(VecSol,fPvv[0]);
163 for (l=fSize;l--;) fPvv[0][l] = fRHS[l] - fPvv[0][l]; // fPvv[0]= initial residual
165 for (l=fSize;l--;) beta += fPvv[0][l]*fPvv[0][l];
166 beta = TMath::Sqrt(beta);
168 if (beta < epsmac) break; // success?
170 //-------------------- normalize: fPvv[0] = fPvv[0] / beta
171 for (l=fSize;l--;) fPvv[0][l] *= t;
172 if (its == 0) eps1 = rtol*beta;
174 // ** initialize 1-st term of rhs of hessenberg system
182 // (Right) Preconditioning Operation z_{j} = M^{-1} v_{j}
184 if (precon>0) ApplyPrecon( fPvv[i], fPvz[i]);
185 else for (l=fSize;l--;) fPvz[i][l] = fPvv[i][l];
187 //-------------------- matvec operation w = A z_{j} = A M^{-1} v_{j}
188 fMatrix->MultiplyByVec(fPvz[i],fPvv[i1]);
190 // modified gram - schmidt...
191 // h_{i,j} = (w,v_{i})
192 // w = w - h_{i,j} v_{i}
194 for (int j=0; j<=i; j++) {
195 for (t=0, l=fSize;l--;) t+= fPvv[j][l]*fPvv[i1][l];
197 for (l=fSize;l--;) fPvv[i1][l] -= t*fPvv[j][l];
199 // -------------------- h_{j+1,j} = ||w||_{2}
200 for (t=0,l=fSize;l--;) t += fPvv[i1][l]*fPvv[i1][l]; t = TMath::Sqrt(t);
202 if (t > 0) for (t=1./t, l=0; l<fSize; l++) fPvv[i1][l] *= t; // v_{j+1} = w / h_{j+1,j}
204 // done with modified gram schimdt and arnoldi step
205 // now update factorization of fPhh
207 // perform previous transformations on i-th column of h
209 for (l=1; l<=i; l++) {
212 fPhh[i][l1] = fPVecR1[l1]*t + fPVecR2[l1]*fPhh[i][l];
213 fPhh[i][l] = -fPVecR2[l1]*t + fPVecR1[l1]*fPhh[i][l];
215 double gam = TMath::Sqrt( fPhh[i][i]*fPhh[i][i] + fPhh[i][i1]*fPhh[i][i1]);
217 // if gamma is zero then any small value will do...
218 // will affect only residual estimate
219 if (gam < epsmac) gam = epsmac;
220 // get next plane rotation
221 fPVecR1[i] = fPhh[i][i]/gam;
222 fPVecR2[i] = fPhh[i][i1]/gam;
223 fPVecV[i1] =-fPVecR2[i]*fPVecV[i];
224 fPVecV[i] *= fPVecR1[i];
226 // determine residual norm and test for convergence
227 fPhh[i][i] = fPVecR1[i]*fPhh[i][i] + fPVecR2[i]*fPhh[i][i1];
228 beta = TMath::Abs(fPVecV[i1]);
230 } while ( (i < nkrylov-1) && (beta > eps1) && (its < itnlim) );
232 // now compute solution. 1st, solve upper triangular system
233 fPVecV[i] = fPVecV[i]/fPhh[i][i];
234 for (int j=1; j<=i; j++) {
236 for (t=fPVecV[k],l=k+1; l<=i; l++) t -= fPhh[l][k]*fPVecV[l];
237 fPVecV[k] = t/fPhh[k][k];
239 // -------------------- linear combination of v[i]'s to get sol.
240 for (int j=0; j<=i; j++) for (t=fPVecV[j],l=0; l<fSize; l++) VecSol[l] += t*fPvz[j][l];
242 // -------------------- restart outer loop if needed
246 AliInfo(Form("FGMRES converged in %d iterations, CPU time: %.1f sec",its,timer.CpuTime()));
252 AliInfo(Form("%d iterations limit exceeded, CPU time: %.1f sec",itnlim,timer.CpuTime()));
264 //________________________________ MINRES METHODS ________________________________
265 Bool_t AliMinResSolve::SolveMinRes(TVectorD& VecSol,Int_t precon,int itnlim,double rtol)
268 return SolveMinRes(VecSol.GetMatrixArray(), precon, itnlim, rtol);
271 //________________________________________________________________________________
272 Bool_t AliMinResSolve::SolveMinRes(double *VecSol,Int_t precon,int itnlim,double rtol)
275 Adapted from author's Fortran code:
276 Michael A. Saunders na.msaunders@na-net.ornl.gov
278 MINRES is an implementation of the algorithm described in the following reference:
279 C. C. Paige and M. A. Saunders (1975),
280 Solution of sparse indefinite systems of linear equations,
281 SIAM J. Numer. Anal. 12(4), pp. 617-629.
284 if (!fMatrix->IsSymmetric()) {
285 AliInfo("MinRes cannot solve asymmetric matrices, use FGMRes instead");
290 const double eps = 2.22E-16;
294 if (precon>=kPreconsTot) {
295 AliInfo(Form("Unknown preconditioner identifier %d, ignore",precon));
298 if (BuildPrecon(precon)<0) {
300 AliInfo("MinRes failed to build the preconditioner");
305 AliInfo(Form("Solution by MinRes: Preconditioner#%d Max.iter.: %d Tol.: %e",precon,itnlim,rtol));
307 // ------------------------ initialization ---------------------->>>>
308 memset(VecSol,0,fSize*sizeof(double));
314 double gam,gmax=1,gmin=1,gbar,oldeps,epsa,epsx,epsr,diag, delta,phi,denom,z;
316 if (!InitAuxMinRes()) return kFALSE;
318 memset(VecSol,0,fSize*sizeof(double));
320 // ------------ init aux -------------------------<<<<
321 // Set up y and v for the first Lanczos vector v1.
322 // y = beta1 P' v1, where P = C**(-1). v is really P' v1.
324 for (int i=fSize;i--;) fPVecY[i] = fPVecR1[i] = fRHS[i];
326 if ( precon>0 ) ApplyPrecon( fRHS, fPVecY);
327 beta1 = 0; for (int i=fSize;i--;) beta1 += fRHS[i]*fPVecY[i]; //
330 AliInfo(Form("Preconditioner is indefinite (init) (%e).",beta1));
337 AliInfo(Form("RHS is zero or is the nullspace of the Preconditioner: Solution is {0}"));
342 beta1 = TMath::Sqrt( beta1 ); // Normalize y to get v1 later.
344 // See if Msolve is symmetric. //RS: Skept
345 // See if Aprod is symmetric. //RS: Skept
351 double qrnorm = beta1;
352 double phibar = beta1;
359 for (int i=fSize;i--;) fPVecR2[i] = fPVecR1[i];
361 TStopwatch timer; timer.Start();
362 while(status==0) { //----------------- Main iteration loop ---------------------->>>>
365 /*-----------------------------------------------------------------
366 Obtain quantities for the next Lanczos vector vk+1, k = 1, 2,...
367 The general iteration is similar to the case k = 1 with v0 = 0:
368 p1 = Operator * v1 - beta1 * v0,
370 q2 = p2 - alpha1 * v1,
373 Again, y = betak P vk, where P = C**(-1).
374 .... more description needed.
375 -----------------------------------------------------------------*/
377 double s = 1./beta; // Normalize previous vector (in y).
378 for (int i=fSize;i--;) fPVecV[i] = s*fPVecY[i]; // v = vk if P = I
380 fMatrix->MultiplyByVec(fPVecV,fPVecY); // APROD (VecV, VecY);
383 double btrat = beta/oldb;
384 for (int i=fSize;i--;) fPVecY[i] -= btrat*fPVecR1[i];
386 double alfa = 0; for (int i=fSize;i--;) alfa += fPVecV[i]*fPVecY[i]; // alphak
388 double alf2bt = alfa/beta;
389 for (int i=fSize;i--;) {
390 fPVecY[i] -= alf2bt*fPVecR2[i];
391 fPVecR1[i] = fPVecR2[i];
392 fPVecR2[i] = fPVecY[i];
395 if ( precon>0 ) ApplyPrecon(fPVecR2, fPVecY);
397 oldb = beta; // oldb = betak
398 beta = 0; for (int i=fSize;i--;) beta += fPVecR2[i]*fPVecY[i]; // beta = betak+1^2
401 AliInfo(Form("Preconditioner is indefinite (%e)",beta));
406 beta = TMath::Sqrt(beta); // beta = betak+1
407 tnorm2 += alfa*alfa + oldb*oldb + beta*beta;
409 if (itn == 1) { // Initialize a few things.
410 if (beta/beta1 <= 10.0*eps) {
411 status = 0; //-1 //????? beta2 = 0 or ~ 0, terminate later.
412 AliInfo("RHS is eigenvector");
415 gmax = TMath::Abs(alfa); // alpha1
416 gmin = gmax; // alpha1
420 Apply previous rotation Qk-1 to get
421 [deltak epslnk+1] = [cs sn][dbark 0 ]
422 [gbar k dbar k+1] [sn -cs][alfak betak+1].
426 delta = cs * dbar + sn * alfa; // delta1 = 0 deltak
427 gbar = sn * dbar - cs * alfa; // gbar 1 = alfa1 gbar k
428 epsln = sn * beta; // epsln2 = 0 epslnk+1
429 dbar = - cs * beta; // dbar 2 = beta2 dbar k+1
431 // Compute the next plane rotation Qk
433 gam = TMath::Sqrt( gbar*gbar + beta*beta ); // gammak
434 cs = gbar / gam; // ck
435 sn = beta / gam; // sk
436 phi = cs * phibar; // phik
437 phibar = sn * phibar; // phibark+1
442 for (int i=fSize;i--;) {
443 fPVecW1[i] = fPVecW2[i];
444 fPVecW2[i] = fPVecW[i];
445 fPVecW[i] = denom*(fPVecV[i]-oldeps*fPVecW1[i]-delta*fPVecW2[i]);
446 VecSol[i] += phi*fPVecW[i];
451 gmax = TMath::Max( gmax, gam );
452 gmin = TMath::Min( gmin, gam );
455 rhs1 = rhs2 - delta * z;
458 // Estimate various norms and test for convergence.
459 normA = TMath::Sqrt( tnorm2 );
460 ynorm = TMath::Sqrt( ynorm2 );
462 epsx = normA * ynorm * eps;
463 epsr = normA * ynorm * rtol;
465 if (diag == 0) diag = epsa;
471 In this version we look at the diagonals of R in the
472 factorization of the lower Hessenberg matrix, Q * H = R,
473 where H is the tridiagonal matrix from Lanczos with one
474 extra row, beta(k+1) e_k^T.
478 // See if any of the stopping criteria are satisfied.
479 // In rare cases, istop is already -1 from above (Abar = const*I).
481 AliDebug(2,Form("#%5d |qnrm: %+.2e Anrm:%+.2e Cnd:%+.2e Rnrm:%+.2e Ynrm:%+.2e EpsR:%+.2e EpsX:%+.2e Beta1:%+.2e",
482 itn,qrnorm, normA,condA,rnorm,ynorm,epsr,epsx,beta1));
485 if (itn >= itnlim ) {status = 5; AliInfo(Form("%d iterations limit exceeded",itnlim));}
486 if (condA >= 0.1/eps ) {status = 4; AliInfo(Form("Matrix condition nymber %e exceeds limit %e",condA,0.1/eps));}
487 if (epsx >= beta1 ) {status = 3; AliInfo(Form("Approximate convergence"));}
488 if (qrnorm <= epsx ) {status = 2; AliInfo(Form("Converged within machine precision"));}
489 if (qrnorm <= epsr ) {status = 1; AliInfo(Form("Converged"));}
492 } //----------------- Main iteration loop ----------------------<<<
497 AliInfo(Form("Exit from MinRes: CPU time: %.2f sec\n"
504 timer.CpuTime(),status,itn,normA,condA,rnorm,ynorm));
506 return status>=0 && status<=3;
510 //______________________________________________________________
511 void AliMinResSolve::ApplyPrecon(const TVectorD& vecRHS, TVectorD& vecOut) const
514 ApplyPrecon(vecRHS.GetMatrixArray(), vecOut.GetMatrixArray());
517 //______________________________________________________________
518 void AliMinResSolve::ApplyPrecon(const double* vecRHS, double* vecOut) const
520 // Application of the preconditioner matrix:
521 // implicitly defines the matrix solving the M*VecOut = VecRHS
522 // const Double_t kTiny = 1E-12;
523 if (fPrecon>=kPreconBD && fPrecon<kPreconILU0) { // band diagonal decomposition
524 fMatBD->Solve(vecRHS,vecOut);
528 else if (fPrecon>=kPreconILU0 && fPrecon<=kPreconILU10) {
530 for(int i=0; i<fSize; i++) { // Block L solve
531 vecOut[i] = vecRHS[i];
532 AliVectorSparse &rowLi = *fMatL->GetRow(i);
533 int n = rowLi.GetNElems();
534 for(int j=0;j<n;j++) vecOut[i] -= vecOut[ rowLi.GetIndex(j) ] * rowLi.GetElem(j);
538 for(int i=fSize; i--; ) { // Block -- U solve
539 AliVectorSparse &rowUi = *fMatU->GetRow(i);
540 int n = rowUi.GetNElems();
541 for(int j=0;j<n;j++ ) vecOut[i] -= vecOut[ rowUi.GetIndex(j) ] * rowUi.GetElem(j);
542 vecOut[i] *= fDiagLU[i];
550 //___________________________________________________________
551 Bool_t AliMinResSolve::InitAuxMinRes()
553 // init auxiliary space for minres
555 fPVecY = new double[fSize];
556 fPVecR1 = new double[fSize];
557 fPVecR2 = new double[fSize];
558 fPVecV = new double[fSize];
559 fPVecW = new double[fSize];
560 fPVecW1 = new double[fSize];
561 fPVecW2 = new double[fSize];
564 AliInfo(Form("Failed to allocate the memory for auxialiary arrays for %d equations",fSize));
569 for (int i=fSize;i--;) fPVecY[i]=fPVecR1[i]=fPVecR2[i]=fPVecV[i]=fPVecW[i]=fPVecW1[i]=fPVecW2[i]=0.0;
575 //___________________________________________________________
576 Bool_t AliMinResSolve::InitAuxFGMRES(int nkrylov)
578 // init auxiliary space for fgmres
580 fPvv = new double*[nkrylov+1];
581 fPvz = new double*[nkrylov];
582 for (int i=0; i<=nkrylov; i++) fPvv[i] = new double[fSize];
583 fPhh = new double*[nkrylov];
584 for (int i=0; i<nkrylov; i++) {
585 fPhh[i] = new double[i+2];
586 fPvz[i] = new double[fSize];
589 fPVecR1 = new double[nkrylov];
590 fPVecR2 = new double[nkrylov];
591 fPVecV = new double[nkrylov+1];
594 AliInfo(Form("Failed to allocate the memory for auxialiary arrays for %d equations",fSize));
603 //___________________________________________________________
604 void AliMinResSolve::ClearAux()
607 if (fPVecY) delete[] fPVecY; fPVecY = 0;
608 if (fPVecR1) delete[] fPVecR1; fPVecR1 = 0;
609 if (fPVecR2) delete[] fPVecR2; fPVecR2 = 0;
610 if (fPVecV) delete[] fPVecV; fPVecV = 0;
611 if (fPVecW) delete[] fPVecW; fPVecW = 0;
612 if (fPVecW1) delete[] fPVecW1; fPVecW1 = 0;
613 if (fPVecW2) delete[] fPVecW2; fPVecW2 = 0;
614 if (fDiagLU) delete[] fDiagLU; fDiagLU = 0;
615 if (fMatL) delete fMatL; fMatL = 0;
616 if (fMatU) delete fMatU; fMatU = 0;
617 if (fMatBD) delete fMatBD; fMatBD = 0;
620 //___________________________________________________________
621 Int_t AliMinResSolve::BuildPreconBD(Int_t hwidth)
623 // build preconditioner
624 AliInfo(Form("Building Band-Diagonal preconditioner of half-width = %d",hwidth));
625 fMatBD = new AliSymBDMatrix( fMatrix->GetSize(), hwidth );
627 // fill the band-diagonal part of the matrix
628 if (fMatrix->InheritsFrom("AliMatrixSparse")) {
629 for (int ir=fMatrix->GetSize();ir--;) {
630 int jmin = TMath::Max(0,ir-hwidth);
631 AliVectorSparse& irow = *((AliMatrixSparse*)fMatrix)->GetRow(ir);
632 for(int j=irow.GetNElems();j--;) {
633 int jind = irow.GetIndex(j);
634 if (jind<jmin) break;
635 (*fMatBD)(ir,jind) = irow.GetElem(j);
640 for (int ir=fMatrix->GetSize();ir--;) {
641 int jmin = TMath::Max(0,ir-hwidth);
642 for(int jr=jmin;jr<=ir;jr++) (*fMatBD)(ir,jr) = fMatrix->Query(ir,jr);
646 fMatBD->DecomposeLDLT();
651 //___________________________________________________________
652 Int_t AliMinResSolve::BuildPreconILUK(Int_t lofM)
654 /*----------------------------------------------------------------------------
655 * ILUK preconditioner
656 * incomplete LU factorization with level of fill dropping
657 * Adapted from iluk.c of ITSOL_1 package by Y.Saad: http://www-users.cs.umn.edu/~saad/software/
658 *----------------------------------------------------------------------------*/
660 AliInfo(Form("Building ILU%d preconditioner",lofM));
662 TStopwatch sw; sw.Start();
663 fMatL = new AliMatrixSparse(fSize);
664 fMatU = new AliMatrixSparse(fSize);
665 fMatL->SetSymmetric(kFALSE);
666 fMatU->SetSymmetric(kFALSE);
667 fDiagLU = new Double_t[fSize];
668 AliMatrixSparse* matrix = (AliMatrixSparse*)fMatrix;
670 // symbolic factorization to calculate level of fill index arrays
671 if ( PreconILUKsymb(lofM)<0 ) {
676 Int_t *jw = new Int_t[fSize];
677 for(int j=fSize;j--;) jw[j] = -1; // set indicator array jw to -1
679 for(int i=0; i<fSize; i++ ) { // beginning of main loop
680 if ( (i%int(0.1*fSize)) == 0) {
681 AliInfo(Form("BuildPrecon: row %d of %d",i,fSize));
686 /* setup array jw[], and initial i-th row */
687 AliVectorSparse& rowLi = *fMatL->GetRow(i);
688 AliVectorSparse& rowUi = *fMatU->GetRow(i);
689 AliVectorSparse& rowM = *matrix->GetRow(i);
691 for(int j=rowLi.GetNElems(); j--;) { // initialize L part
692 int col = rowLi.GetIndex(j);
694 rowLi.GetElem(j) = 0.; // do we need this ?
697 fDiagLU[i] = 0; // initialize diagonal
699 for(int j=rowUi.GetNElems();j--;) { // initialize U part
700 int col = rowUi.GetIndex(j);
702 rowUi.GetElem(j) = 0;
704 // copy row from csmat into L,U D
705 for(int j=rowM.GetNElems(); j--;) { // L and D part
706 if (AliMatrixSq::IsZero(rowM.GetElem(j))) continue;
707 int col = rowM.GetIndex(j); // (the original matrix stores only lower triangle)
708 if( col < i ) rowLi.GetElem(jw[col]) = rowM.GetElem(j);
709 else if(col==i) fDiagLU[i] = rowM.GetElem(j);
710 else rowUi.GetElem(jw[col]) = rowM.GetElem(j);
712 if (matrix->IsSymmetric()) for (int col=i+1;col<fSize;col++) { // part of the row I on the right of diagonal is stored as
713 double vl = matrix->Query(col,i); // the lower part of the column I
714 if (AliMatrixSq::IsZero(vl)) continue;
715 rowUi.GetElem(jw[col]) = vl;
718 // eliminate previous rows
719 for(int j=0; j<rowLi.GetNElems(); j++) {
720 int jrow = rowLi.GetIndex(j);
721 // get the multiplier for row to be eliminated (jrow)
722 rowLi.GetElem(j) *= fDiagLU[jrow];
724 // combine current row and row jrow
725 AliVectorSparse& rowUj = *fMatU->GetRow(jrow);
726 for(int k=0; k<rowUj.GetNElems(); k++ ) {
727 int col = rowUj.GetIndex(k);
729 if( jpos == -1 ) continue;
730 if( col < i ) rowLi.GetElem(jpos) -= rowLi.GetElem(j) * rowUj.GetElem(k);
731 else if(col==i) fDiagLU[i] -= rowLi.GetElem(j) * rowUj.GetElem(k);
732 else rowUi.GetElem(jpos) -= rowLi.GetElem(j) * rowUj.GetElem(k);
735 // reset double-pointer to -1 ( U-part)
736 for(int j=rowLi.GetNElems(); j--;) jw[ rowLi.GetIndex(j) ] = -1;
738 for(int j=rowUi.GetNElems(); j--;) jw[ rowUi.GetIndex(j) ] = -1;
740 if( AliMatrixSq::IsZero(fDiagLU[i]) ) {
741 AliInfo(Form("Fatal error in ILIk: Zero diagonal found..."));
745 fDiagLU[i] = 1.0 / fDiagLU[i];
751 AliInfo(Form("ILU%d preconditioner OK, CPU time: %.1f sec",lofM,sw.CpuTime()));
752 AliInfo(Form("Densities: M %f L %f U %f",matrix->GetDensity(),fMatL->GetDensity(),fMatU->GetDensity()));
757 //___________________________________________________________
758 Int_t AliMinResSolve::BuildPreconILUKDense(Int_t lofM)
760 /*----------------------------------------------------------------------------
761 * ILUK preconditioner
762 * incomplete LU factorization with level of fill dropping
763 * Adapted from iluk.c of ITSOL_1 package by Y.Saad: http://www-users.cs.umn.edu/~saad/software/
764 *----------------------------------------------------------------------------*/
766 TStopwatch sw; sw.Start();
767 AliInfo(Form("Building ILU%d preconditioner for dense matrix",lofM));
769 fMatL = new AliMatrixSparse(fSize);
770 fMatU = new AliMatrixSparse(fSize);
771 fMatL->SetSymmetric(kFALSE);
772 fMatU->SetSymmetric(kFALSE);
773 fDiagLU = new Double_t[fSize];
775 // symbolic factorization to calculate level of fill index arrays
776 if ( PreconILUKsymbDense(lofM)<0 ) {
781 Int_t *jw = new Int_t[fSize];
782 for(int j=fSize;j--;) jw[j] = -1; // set indicator array jw to -1
784 for(int i=0; i<fSize; i++ ) { // beginning of main loop
785 /* setup array jw[], and initial i-th row */
786 AliVectorSparse& rowLi = *fMatL->GetRow(i);
787 AliVectorSparse& rowUi = *fMatU->GetRow(i);
789 for(int j=rowLi.GetNElems(); j--;) { // initialize L part
790 int col = rowLi.GetIndex(j);
792 rowLi.GetElem(j) = 0.; // do we need this ?
795 fDiagLU[i] = 0; // initialize diagonal
797 for(int j=rowUi.GetNElems();j--;) { // initialize U part
798 int col = rowUi.GetIndex(j);
800 rowUi.GetElem(j) = 0;
802 // copy row from csmat into L,U D
803 for(int j=fSize; j--;) { // L and D part
804 double vl = fMatrix->Query(i,j);
805 if (AliMatrixSq::IsZero(vl)) continue;
806 if( j < i ) rowLi.GetElem(jw[j]) = vl;
807 else if(j==i) fDiagLU[i] = vl;
808 else rowUi.GetElem(jw[j]) = vl;
810 // eliminate previous rows
811 for(int j=0; j<rowLi.GetNElems(); j++) {
812 int jrow = rowLi.GetIndex(j);
813 // get the multiplier for row to be eliminated (jrow)
814 rowLi.GetElem(j) *= fDiagLU[jrow];
816 // combine current row and row jrow
817 AliVectorSparse& rowUj = *fMatU->GetRow(jrow);
818 for(int k=0; k<rowUj.GetNElems(); k++ ) {
819 int col = rowUj.GetIndex(k);
821 if( jpos == -1 ) continue;
822 if( col < i ) rowLi.GetElem(jpos) -= rowLi.GetElem(j) * rowUj.GetElem(k);
823 else if(col==i) fDiagLU[i] -= rowLi.GetElem(j) * rowUj.GetElem(k);
824 else rowUi.GetElem(jpos) -= rowLi.GetElem(j) * rowUj.GetElem(k);
827 // reset double-pointer to -1 ( U-part)
828 for(int j=rowLi.GetNElems(); j--;) jw[ rowLi.GetIndex(j) ] = -1;
830 for(int j=rowUi.GetNElems(); j--;) jw[ rowUi.GetIndex(j) ] = -1;
832 if( AliMatrixSq::IsZero(fDiagLU[i])) {
833 AliInfo(Form("Fatal error in ILIk: Zero diagonal found..."));
837 fDiagLU[i] = 1.0 / fDiagLU[i];
843 AliInfo(Form("ILU%d dense preconditioner OK, CPU time: %.1f sec",lofM,sw.CpuTime()));
844 // AliInfo(Form("Densities: M %f L %f U %f",fMatrix->GetDensity(),fMatL->GetDensity(),fMatU->GetDensity()));
849 //___________________________________________________________
850 Int_t AliMinResSolve::PreconILUKsymb(Int_t lofM)
852 /*----------------------------------------------------------------------------
853 * ILUK preconditioner
854 * incomplete LU factorization with level of fill dropping
855 * Adapted from iluk.c: lofC of ITSOL_1 package by Y.Saad: http://www-users.cs.umn.edu/~saad/software/
856 *----------------------------------------------------------------------------*/
859 AliInfo("PreconILUKsymb>>");
860 AliMatrixSparse* matrix = (AliMatrixSparse*)fMatrix;
863 UChar_t **ulvl=0,*levls=0;
867 ulvl = new UChar_t*[fSize]; // stores lev-fils for U part of ILU factorization
868 levls = new UChar_t[fSize];
869 jbuf = new UShort_t[fSize];
870 iw = new Int_t[fSize];
874 AliInfo("Failed to allocate the memory in AliMinResSolve::PreconILUKsymb");
875 if (ulvl) delete[] ulvl;
876 if (levls) delete[] levls;
877 if (jbuf) delete[] jbuf;
883 for(int j=fSize; j--;) iw[j] = -1; // initialize iw
884 for(int i=0; i<fSize; i++) {
887 AliVectorSparse& row = *matrix->GetRow(i);
889 // assign lof = 0 for matrix elements
890 for(int j=0;j<row.GetNElems(); j++) {
891 int col = row.GetIndex(j);
892 if (AliMatrixSq::IsZero(row.GetElem(j))) continue; // !!!! matrix is sparse but sometimes 0 appears
893 if (col<i) { // L-part
898 else if (col>i) { // This works only for general matrix
904 if (matrix->IsSymmetric()) for (int col=i+1;col<fSize;col++) { // U-part of symmetric matrix
905 if (AliMatrixSq::IsZero(matrix->Query(col,i))) continue; // Due to the symmetry == matrix(i,col)
911 // symbolic k,i,j Gaussian elimination
913 while (++jpiv < incl) {
914 int k = jbuf[jpiv] ; // select leftmost pivot
917 for(int j=jpiv+1; j<incl; j++) if( jbuf[j]<kmin ) { kmin = jbuf[j]; jmin = j;}
919 // ------------------------------------ swap
925 int tj = levls[jpiv] ;
926 levls[jpiv] = levls[jmin];
930 // ------------------------------------ symbolic linear combinaiton of rows
931 AliVectorSparse& rowU = *fMatU->GetRow(k);
932 for(int j=0; j<rowU.GetNElems(); j++ ) {
933 int col = rowU.GetIndex(j);
934 int it = ulvl[k][j]+levls[jpiv]+1;
935 if( it > lofM ) continue;
949 else levls[ip] = TMath::Min(levls[ip], it);
952 } // end - while loop
955 for (int j=0;j<incl;j++) iw[jbuf[j]] = -1;
956 for (int j=i;j<incu;j++) iw[jbuf[j]] = -1;
959 AliVectorSparse& rowLi = *fMatL->GetRow(i);
961 if(incl>0) memcpy(rowLi.GetIndices(), jbuf, sizeof(UShort_t)*incl);
964 AliVectorSparse& rowUi = *fMatU->GetRow(i);
967 memcpy(rowUi.GetIndices(), jbuf+i, sizeof(UShort_t)*k);
968 ulvl[i] = new UChar_t[k]; // update matrix of levels
969 memcpy( ulvl[i], levls+i, k*sizeof(UChar_t) );
973 // free temp space and leave
976 for(int i=fSize; i--;) if (fMatU->GetRow(i)->GetNElems()) delete[] ulvl[i];
979 fMatL->SortIndices();
980 fMatU->SortIndices();
983 AliInfo("PreconILUKsymb<<");
988 //___________________________________________________________
989 Int_t AliMinResSolve::PreconILUKsymbDense(Int_t lofM)
991 /*----------------------------------------------------------------------------
992 * ILUK preconditioner
993 * incomplete LU factorization with level of fill dropping
994 * Adapted from iluk.c: lofC of ITSOL_1 package by Y.Saad: http://www-users.cs.umn.edu/~saad/software/
995 *----------------------------------------------------------------------------*/
997 UChar_t **ulvl=0,*levls=0;
1001 ulvl = new UChar_t*[fSize]; // stores lev-fils for U part of ILU factorization
1002 levls = new UChar_t[fSize];
1003 jbuf = new UShort_t[fSize];
1004 iw = new Int_t[fSize];
1008 AliInfo("Failed to allocate the memory in AliMinResSolve::PreconILUKsymb");
1009 if (ulvl) delete[] ulvl;
1010 if (levls) delete[] levls;
1011 if (jbuf) delete[] jbuf;
1012 if (iw) delete[] iw;
1017 for(int j=fSize; j--;) iw[j] = -1; // initialize iw
1018 for(int i=0; i<fSize; i++) {
1022 // assign lof = 0 for matrix elements
1023 for(int j=0;j<fSize; j++) {
1024 if (AliMatrixSq::IsZero(fMatrix->Query(i,j))) continue;
1025 if (j<i) { // L-part
1030 else if (j>i) { // This works only for general matrix
1037 // symbolic k,i,j Gaussian elimination
1039 while (++jpiv < incl) {
1040 int k = jbuf[jpiv] ; // select leftmost pivot
1043 for(int j=jpiv+1; j<incl; j++) if( jbuf[j]<kmin ) { kmin = jbuf[j]; jmin = j;}
1045 // ------------------------------------ swap
1051 int tj = levls[jpiv] ;
1052 levls[jpiv] = levls[jmin];
1056 // ------------------------------------ symbolic linear combinaiton of rows
1057 AliVectorSparse& rowU = *fMatU->GetRow(k);
1058 for(int j=0; j<rowU.GetNElems(); j++ ) {
1059 int col = rowU.GetIndex(j);
1060 int it = ulvl[k][j]+levls[jpiv]+1;
1061 if( it > lofM ) continue;
1069 else if( col > i ) {
1075 else levls[ip] = TMath::Min(levls[ip], it);
1078 } // end - while loop
1081 for (int j=0;j<incl;j++) iw[jbuf[j]] = -1;
1082 for (int j=i;j<incu;j++) iw[jbuf[j]] = -1;
1085 AliVectorSparse& rowLi = *fMatL->GetRow(i);
1087 if(incl>0) memcpy(rowLi.GetIndices(), jbuf, sizeof(UShort_t)*incl);
1090 AliVectorSparse& rowUi = *fMatU->GetRow(i);
1093 memcpy(rowUi.GetIndices(), jbuf+i, sizeof(UShort_t)*k);
1094 ulvl[i] = new UChar_t[k]; // update matrix of levels
1095 memcpy( ulvl[i], levls+i, k*sizeof(UChar_t) );
1099 // free temp space and leave
1102 for(int i=fSize; i--;) if (fMatU->GetRow(i)->GetNElems()) delete[] ulvl[i];
1105 fMatL->SortIndices();
1106 fMatU->SortIndices();