1be1d678aSKris Buschelman 24eb8e494SKris Buschelman /* 34eb8e494SKris Buschelman Provides an interface to the LUSOL package of .... 44eb8e494SKris Buschelman 54eb8e494SKris Buschelman */ 6c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> 74eb8e494SKris Buschelman 84eb8e494SKris Buschelman #if defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 94eb8e494SKris Buschelman #define LU1FAC lu1fac_ 104eb8e494SKris Buschelman #define LU6SOL lu6sol_ 114eb8e494SKris Buschelman #define M1PAGE m1page_ 124eb8e494SKris Buschelman #define M5SETX m5setx_ 134eb8e494SKris Buschelman #define M6RDEL m6rdel_ 144eb8e494SKris Buschelman #elif !defined(PETSC_HAVE_FORTRAN_CAPS) 154eb8e494SKris Buschelman #define LU1FAC lu1fac 164eb8e494SKris Buschelman #define LU6SOL lu6sol 174eb8e494SKris Buschelman #define M1PAGE m1page 184eb8e494SKris Buschelman #define M5SETX m5setx 194eb8e494SKris Buschelman #define M6RDEL m6rdel 204eb8e494SKris Buschelman #endif 214eb8e494SKris Buschelman 224eb8e494SKris Buschelman /* 234eb8e494SKris Buschelman Dummy symbols that the MINOS files mi25bfac.f and mi15blas.f may require 244eb8e494SKris Buschelman */ 2519caf8f3SSatish Balay PETSC_EXTERN void M1PAGE() 26a6dfd86eSKarl Rupp { 274eb8e494SKris Buschelman ; 284eb8e494SKris Buschelman } 2919caf8f3SSatish Balay PETSC_EXTERN void M5SETX() 30a6dfd86eSKarl Rupp { 314eb8e494SKris Buschelman ; 324eb8e494SKris Buschelman } 334eb8e494SKris Buschelman 3419caf8f3SSatish Balay PETSC_EXTERN void M6RDEL() 35a6dfd86eSKarl Rupp { 364eb8e494SKris Buschelman ; 374eb8e494SKris Buschelman } 384eb8e494SKris Buschelman 3919caf8f3SSatish Balay PETSC_EXTERN void LU1FAC(int *m, int *n, int *nnz, int *size, int *luparm, 404eb8e494SKris Buschelman double *parmlu, double *data, int *indc, int *indr, 414eb8e494SKris Buschelman int *rowperm, int *colperm, int *collen, int *rowlen, 424eb8e494SKris Buschelman int *colstart, int *rowstart, int *rploc, int *cploc, 434eb8e494SKris Buschelman int *rpinv, int *cpinv, double *w, int *inform); 444eb8e494SKris Buschelman 4519caf8f3SSatish Balay PETSC_EXTERN void LU6SOL(int *mode, int *m, int *n, double *rhs, double *x, 464eb8e494SKris Buschelman int *size, int *luparm, double *parmlu, double *data, 474eb8e494SKris Buschelman int *indc, int *indr, int *rowperm, int *colperm, 484eb8e494SKris Buschelman int *collen, int *rowlen, int *colstart, int *rowstart, 494eb8e494SKris Buschelman int *inform); 504eb8e494SKris Buschelman 5109573ac7SBarry Smith extern PetscErrorCode MatDuplicate_LUSOL(Mat,MatDuplicateOption,Mat*); 52f0c56d0fSKris Buschelman 53f0c56d0fSKris Buschelman typedef struct { 544eb8e494SKris Buschelman double *data; 554eb8e494SKris Buschelman int *indc; 564eb8e494SKris Buschelman int *indr; 574eb8e494SKris Buschelman 584eb8e494SKris Buschelman int *ip; 594eb8e494SKris Buschelman int *iq; 604eb8e494SKris Buschelman int *lenc; 614eb8e494SKris Buschelman int *lenr; 624eb8e494SKris Buschelman int *locc; 634eb8e494SKris Buschelman int *locr; 644eb8e494SKris Buschelman int *iploc; 654eb8e494SKris Buschelman int *iqloc; 664eb8e494SKris Buschelman int *ipinv; 674eb8e494SKris Buschelman int *iqinv; 684eb8e494SKris Buschelman double *mnsw; 694eb8e494SKris Buschelman double *mnsv; 704eb8e494SKris Buschelman 714eb8e494SKris Buschelman double elbowroom; 724eb8e494SKris Buschelman double luroom; /* Extra space allocated when factor fails */ 734eb8e494SKris Buschelman double parmlu[30]; /* Input/output to LUSOL */ 744eb8e494SKris Buschelman 754eb8e494SKris Buschelman int n; /* Number of rows/columns in matrix */ 764eb8e494SKris Buschelman int nz; /* Number of nonzeros */ 774eb8e494SKris Buschelman int nnz; /* Number of nonzeros allocated for factors */ 784eb8e494SKris Buschelman int luparm[30]; /* Input/output to LUSOL */ 794eb8e494SKris Buschelman 80ace3abfcSBarry Smith PetscBool CleanUpLUSOL; 814eb8e494SKris Buschelman 82f0c56d0fSKris Buschelman } Mat_LUSOL; 834eb8e494SKris Buschelman 844eb8e494SKris Buschelman /* LUSOL input/Output Parameters (Description uses C-style indexes 854eb8e494SKris Buschelman * 864eb8e494SKris Buschelman * Input parameters Typical value 874eb8e494SKris Buschelman * 884eb8e494SKris Buschelman * luparm(0) = nout File number for printed messages. 6 894eb8e494SKris Buschelman * luparm(1) = lprint Print level. 0 904eb8e494SKris Buschelman * < 0 suppresses output. 914eb8e494SKris Buschelman * = 0 gives error messages. 924eb8e494SKris Buschelman * = 1 gives debug output from some of the 934eb8e494SKris Buschelman * other routines in LUSOL. 944eb8e494SKris Buschelman * >= 2 gives the pivot row and column and the 954eb8e494SKris Buschelman * no. of rows and columns involved at 964eb8e494SKris Buschelman * each elimination step in lu1fac. 974eb8e494SKris Buschelman * luparm(2) = maxcol lu1fac: maximum number of columns 5 984eb8e494SKris Buschelman * searched allowed in a Markowitz-type 994eb8e494SKris Buschelman * search for the next pivot element. 1004eb8e494SKris Buschelman * For some of the factorization, the 1014eb8e494SKris Buschelman * number of rows searched is 1024eb8e494SKris Buschelman * maxrow = maxcol - 1. 1034eb8e494SKris Buschelman * 1044eb8e494SKris Buschelman * 105*7a7aea1fSJed Brown * Output parameters: 1064eb8e494SKris Buschelman * 1074eb8e494SKris Buschelman * luparm(9) = inform Return code from last call to any LU routine. 1084eb8e494SKris Buschelman * luparm(10) = nsing No. of singularities marked in the 1094eb8e494SKris Buschelman * output array w(*). 1104eb8e494SKris Buschelman * luparm(11) = jsing Column index of last singularity. 1114eb8e494SKris Buschelman * luparm(12) = minlen Minimum recommended value for lena. 1124eb8e494SKris Buschelman * luparm(13) = maxlen ? 1134eb8e494SKris Buschelman * luparm(14) = nupdat No. of updates performed by the lu8 routines. 1144eb8e494SKris Buschelman * luparm(15) = nrank No. of nonempty rows of U. 1154eb8e494SKris Buschelman * luparm(16) = ndens1 No. of columns remaining when the density of 1164eb8e494SKris Buschelman * the matrix being factorized reached dens1. 1174eb8e494SKris Buschelman * luparm(17) = ndens2 No. of columns remaining when the density of 1184eb8e494SKris Buschelman * the matrix being factorized reached dens2. 1194eb8e494SKris Buschelman * luparm(18) = jumin The column index associated with dumin. 1204eb8e494SKris Buschelman * luparm(19) = numl0 No. of columns in initial L. 1214eb8e494SKris Buschelman * luparm(20) = lenl0 Size of initial L (no. of nonzeros). 1224eb8e494SKris Buschelman * luparm(21) = lenu0 Size of initial U. 1234eb8e494SKris Buschelman * luparm(22) = lenl Size of current L. 1244eb8e494SKris Buschelman * luparm(23) = lenu Size of current U. 1254eb8e494SKris Buschelman * luparm(24) = lrow Length of row file. 1264eb8e494SKris Buschelman * luparm(25) = ncp No. of compressions of LU data structures. 1274eb8e494SKris Buschelman * luparm(26) = mersum lu1fac: sum of Markowitz merit counts. 1284eb8e494SKris Buschelman * luparm(27) = nutri lu1fac: triangular rows in U. 1294eb8e494SKris Buschelman * luparm(28) = nltri lu1fac: triangular rows in L. 1304eb8e494SKris Buschelman * luparm(29) = 1314eb8e494SKris Buschelman * 1324eb8e494SKris Buschelman * 1334eb8e494SKris Buschelman * Input parameters Typical value 1344eb8e494SKris Buschelman * 1354eb8e494SKris Buschelman * parmlu(0) = elmax1 Max multiplier allowed in L 10.0 1364eb8e494SKris Buschelman * during factor. 1374eb8e494SKris Buschelman * parmlu(1) = elmax2 Max multiplier allowed in L 10.0 1384eb8e494SKris Buschelman * during updates. 1394eb8e494SKris Buschelman * parmlu(2) = small Absolute tolerance for eps**0.8 1404eb8e494SKris Buschelman * treating reals as zero. IBM double: 3.0d-13 1414eb8e494SKris Buschelman * parmlu(3) = utol1 Absolute tol for flagging eps**0.66667 1424eb8e494SKris Buschelman * small diagonals of U. IBM double: 3.7d-11 1434eb8e494SKris Buschelman * parmlu(4) = utol2 Relative tol for flagging eps**0.66667 1444eb8e494SKris Buschelman * small diagonals of U. IBM double: 3.7d-11 1454eb8e494SKris Buschelman * parmlu(5) = uspace Factor limiting waste space in U. 3.0 1464eb8e494SKris Buschelman * In lu1fac, the row or column lists 1474eb8e494SKris Buschelman * are compressed if their length 1484eb8e494SKris Buschelman * exceeds uspace times the length of 1494eb8e494SKris Buschelman * either file after the last compression. 1504eb8e494SKris Buschelman * parmlu(6) = dens1 The density at which the Markowitz 0.3 1514eb8e494SKris Buschelman * strategy should search maxcol columns 1524eb8e494SKris Buschelman * and no rows. 1534eb8e494SKris Buschelman * parmlu(7) = dens2 the density at which the Markowitz 0.6 1544eb8e494SKris Buschelman * strategy should search only 1 column 1554eb8e494SKris Buschelman * or (preferably) use a dense LU for 1564eb8e494SKris Buschelman * all the remaining rows and columns. 1574eb8e494SKris Buschelman * 1584eb8e494SKris Buschelman * 159*7a7aea1fSJed Brown * Output parameters: 1604eb8e494SKris Buschelman * 1614eb8e494SKris Buschelman * parmlu(9) = amax Maximum element in A. 1624eb8e494SKris Buschelman * parmlu(10) = elmax Maximum multiplier in current L. 1634eb8e494SKris Buschelman * parmlu(11) = umax Maximum element in current U. 1644eb8e494SKris Buschelman * parmlu(12) = dumax Maximum diagonal in U. 1654eb8e494SKris Buschelman * parmlu(13) = dumin Minimum diagonal in U. 1664eb8e494SKris Buschelman * parmlu(14) = 1674eb8e494SKris Buschelman * parmlu(15) = 1684eb8e494SKris Buschelman * parmlu(16) = 1694eb8e494SKris Buschelman * parmlu(17) = 1704eb8e494SKris Buschelman * parmlu(18) = 1714eb8e494SKris Buschelman * parmlu(19) = resid lu6sol: residual after solve with U or U'. 1724eb8e494SKris Buschelman * ... 1734eb8e494SKris Buschelman * parmlu(29) = 1744eb8e494SKris Buschelman */ 1754eb8e494SKris Buschelman 1764eb8e494SKris Buschelman #define Factorization_Tolerance 1e-1 1774eb8e494SKris Buschelman #define Factorization_Pivot_Tolerance pow(2.2204460492503131E-16, 2.0 / 3.0) 1784eb8e494SKris Buschelman #define Factorization_Small_Tolerance 1e-15 /* pow(DBL_EPSILON, 0.8) */ 1794eb8e494SKris Buschelman 180dfbe8321SBarry Smith PetscErrorCode MatDestroy_LUSOL(Mat A) 181dfbe8321SBarry Smith { 182dfbe8321SBarry Smith PetscErrorCode ierr; 183f0c56d0fSKris Buschelman Mat_LUSOL *lusol=(Mat_LUSOL*)A->spptr; 1844eb8e494SKris Buschelman 1854eb8e494SKris Buschelman PetscFunctionBegin; 186bf0cc555SLisandro Dalcin if (lusol && lusol->CleanUpLUSOL) { 1874eb8e494SKris Buschelman ierr = PetscFree(lusol->ip);CHKERRQ(ierr); 1884eb8e494SKris Buschelman ierr = PetscFree(lusol->iq);CHKERRQ(ierr); 1894eb8e494SKris Buschelman ierr = PetscFree(lusol->lenc);CHKERRQ(ierr); 1904eb8e494SKris Buschelman ierr = PetscFree(lusol->lenr);CHKERRQ(ierr); 1914eb8e494SKris Buschelman ierr = PetscFree(lusol->locc);CHKERRQ(ierr); 1924eb8e494SKris Buschelman ierr = PetscFree(lusol->locr);CHKERRQ(ierr); 1934eb8e494SKris Buschelman ierr = PetscFree(lusol->iploc);CHKERRQ(ierr); 1944eb8e494SKris Buschelman ierr = PetscFree(lusol->iqloc);CHKERRQ(ierr); 1954eb8e494SKris Buschelman ierr = PetscFree(lusol->ipinv);CHKERRQ(ierr); 1964eb8e494SKris Buschelman ierr = PetscFree(lusol->iqinv);CHKERRQ(ierr); 1974eb8e494SKris Buschelman ierr = PetscFree(lusol->mnsw);CHKERRQ(ierr); 1984eb8e494SKris Buschelman ierr = PetscFree(lusol->mnsv);CHKERRQ(ierr); 19923bdbc58SBarry Smith ierr = PetscFree3(lusol->data,lusol->indc,lusol->indr);CHKERRQ(ierr); 2004eb8e494SKris Buschelman } 201bf0cc555SLisandro Dalcin ierr = PetscFree(A->spptr);CHKERRQ(ierr); 202b24902e0SBarry Smith ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 2034eb8e494SKris Buschelman PetscFunctionReturn(0); 2044eb8e494SKris Buschelman } 2054eb8e494SKris Buschelman 2066849ba73SBarry Smith PetscErrorCode MatSolve_LUSOL(Mat A,Vec b,Vec x) 2076849ba73SBarry Smith { 208f0c56d0fSKris Buschelman Mat_LUSOL *lusol=(Mat_LUSOL*)A->spptr; 209d9ca1df4SBarry Smith double *xx; 210d9ca1df4SBarry Smith const double *bb; 2114eb8e494SKris Buschelman int mode=5; 2126849ba73SBarry Smith PetscErrorCode ierr; 2136849ba73SBarry Smith int i,m,n,nnz,status; 2144eb8e494SKris Buschelman 2154eb8e494SKris Buschelman PetscFunctionBegin; 2164eb8e494SKris Buschelman ierr = VecGetArray(x, &xx);CHKERRQ(ierr); 217d9ca1df4SBarry Smith ierr = VecGetArrayRead(b, &bb);CHKERRQ(ierr); 2184eb8e494SKris Buschelman 2194eb8e494SKris Buschelman m = n = lusol->n; 2204eb8e494SKris Buschelman nnz = lusol->nnz; 2214eb8e494SKris Buschelman 2222205254eSKarl Rupp for (i = 0; i < m; i++) lusol->mnsv[i] = bb[i]; 2234eb8e494SKris Buschelman 2244eb8e494SKris Buschelman LU6SOL(&mode, &m, &n, lusol->mnsv, xx, &nnz, 2254eb8e494SKris Buschelman lusol->luparm, lusol->parmlu, lusol->data, 2264eb8e494SKris Buschelman lusol->indc, lusol->indr, lusol->ip, lusol->iq, 2274eb8e494SKris Buschelman lusol->lenc, lusol->lenr, lusol->locc, lusol->locr, &status); 2284eb8e494SKris Buschelman 22965e19b50SBarry Smith if (status) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"solve failed, error code %d",status); 2304eb8e494SKris Buschelman 2314eb8e494SKris Buschelman ierr = VecRestoreArray(x, &xx);CHKERRQ(ierr); 232d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(b, &bb);CHKERRQ(ierr); 2334eb8e494SKris Buschelman PetscFunctionReturn(0); 2344eb8e494SKris Buschelman } 2354eb8e494SKris Buschelman 2360481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_LUSOL(Mat F,Mat A,const MatFactorInfo *info) 2376849ba73SBarry Smith { 2384eb8e494SKris Buschelman Mat_SeqAIJ *a; 239719d5645SBarry Smith Mat_LUSOL *lusol = (Mat_LUSOL*)F->spptr; 2406849ba73SBarry Smith PetscErrorCode ierr; 2414eb8e494SKris Buschelman int m, n, nz, nnz, status; 2426849ba73SBarry Smith int i, rs, re; 2434eb8e494SKris Buschelman int factorizations; 2444eb8e494SKris Buschelman 2454eb8e494SKris Buschelman PetscFunctionBegin; 246c3b366b1Sprj- ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 2474eb8e494SKris Buschelman a = (Mat_SeqAIJ*)A->data; 2484eb8e494SKris Buschelman 249e32f2f54SBarry Smith if (m != lusol->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"factorization struct inconsistent"); 2504eb8e494SKris Buschelman 2514eb8e494SKris Buschelman factorizations = 0; 2522205254eSKarl Rupp do { 2534eb8e494SKris Buschelman /*******************************************************************/ 2544eb8e494SKris Buschelman /* Check the workspace allocation. */ 2554eb8e494SKris Buschelman /*******************************************************************/ 2564eb8e494SKris Buschelman 2574eb8e494SKris Buschelman nz = a->nz; 2584eb8e494SKris Buschelman nnz = PetscMax(lusol->nnz, (int)(lusol->elbowroom*nz)); 2594eb8e494SKris Buschelman nnz = PetscMax(nnz, 5*n); 2604eb8e494SKris Buschelman 2614eb8e494SKris Buschelman if (nnz < lusol->luparm[12]) { 2624eb8e494SKris Buschelman nnz = (int)(lusol->luroom * lusol->luparm[12]); 2634eb8e494SKris Buschelman } else if ((factorizations > 0) && (lusol->luroom < 6)) { 2644eb8e494SKris Buschelman lusol->luroom += 0.1; 2654eb8e494SKris Buschelman } 2664eb8e494SKris Buschelman 2674eb8e494SKris Buschelman nnz = PetscMax(nnz, (int)(lusol->luroom*(lusol->luparm[22] + lusol->luparm[23]))); 2684eb8e494SKris Buschelman 2694eb8e494SKris Buschelman if (nnz > lusol->nnz) { 27023bdbc58SBarry Smith ierr = PetscFree3(lusol->data,lusol->indc,lusol->indr);CHKERRQ(ierr); 271dcca6d9dSJed Brown ierr = PetscMalloc3(nnz,&lusol->data,nnz,&lusol->indc,nnz,&lusol->indr);CHKERRQ(ierr); 2724eb8e494SKris Buschelman lusol->nnz = nnz; 2734eb8e494SKris Buschelman } 2744eb8e494SKris Buschelman 2754eb8e494SKris Buschelman /*******************************************************************/ 2764eb8e494SKris Buschelman /* Fill in the data for the problem. (1-based Fortran style) */ 2774eb8e494SKris Buschelman /*******************************************************************/ 2784eb8e494SKris Buschelman 2794eb8e494SKris Buschelman nz = 0; 2802205254eSKarl Rupp for (i = 0; i < n; i++) { 2814eb8e494SKris Buschelman rs = a->i[i]; 2824eb8e494SKris Buschelman re = a->i[i+1]; 2834eb8e494SKris Buschelman 2842205254eSKarl Rupp while (rs < re) { 2852205254eSKarl Rupp if (a->a[rs] != 0.0) { 2864eb8e494SKris Buschelman lusol->indc[nz] = i + 1; 2874eb8e494SKris Buschelman lusol->indr[nz] = a->j[rs] + 1; 2884eb8e494SKris Buschelman lusol->data[nz] = a->a[rs]; 2894eb8e494SKris Buschelman nz++; 2904eb8e494SKris Buschelman } 2914eb8e494SKris Buschelman rs++; 2924eb8e494SKris Buschelman } 2934eb8e494SKris Buschelman } 2944eb8e494SKris Buschelman 2954eb8e494SKris Buschelman /*******************************************************************/ 2964eb8e494SKris Buschelman /* Do the factorization. */ 2974eb8e494SKris Buschelman /*******************************************************************/ 2984eb8e494SKris Buschelman 2994eb8e494SKris Buschelman LU1FAC(&m, &n, &nz, &nnz, 3004eb8e494SKris Buschelman lusol->luparm, lusol->parmlu, lusol->data, 3014eb8e494SKris Buschelman lusol->indc, lusol->indr, lusol->ip, lusol->iq, 3024eb8e494SKris Buschelman lusol->lenc, lusol->lenr, lusol->locc, lusol->locr, 3034eb8e494SKris Buschelman lusol->iploc, lusol->iqloc, lusol->ipinv, 3044eb8e494SKris Buschelman lusol->iqinv, lusol->mnsw, &status); 3054eb8e494SKris Buschelman 3062205254eSKarl Rupp switch (status) { 3074eb8e494SKris Buschelman case 0: /* factored */ 3084eb8e494SKris Buschelman break; 3094eb8e494SKris Buschelman 3104eb8e494SKris Buschelman case 7: /* insufficient memory */ 3114eb8e494SKris Buschelman break; 3124eb8e494SKris Buschelman 3134eb8e494SKris Buschelman case 1: 3144eb8e494SKris Buschelman case -1: /* singular */ 315e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Singular matrix"); 3164eb8e494SKris Buschelman 3174eb8e494SKris Buschelman case 3: 3184eb8e494SKris Buschelman case 4: /* error conditions */ 319e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"matrix error"); 3204eb8e494SKris Buschelman 3214eb8e494SKris Buschelman default: /* unknown condition */ 322e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"matrix unknown return code"); 3234eb8e494SKris Buschelman } 3244eb8e494SKris Buschelman 3254eb8e494SKris Buschelman factorizations++; 3264eb8e494SKris Buschelman } while (status == 7); 327719d5645SBarry Smith F->ops->solve = MatSolve_LUSOL; 328719d5645SBarry Smith F->assembled = PETSC_TRUE; 329719d5645SBarry Smith F->preallocated = PETSC_TRUE; 3304eb8e494SKris Buschelman PetscFunctionReturn(0); 3314eb8e494SKris Buschelman } 3324eb8e494SKris Buschelman 33335bd34faSBarry Smith PetscErrorCode MatLUFactorSymbolic_LUSOL(Mat F,Mat A, IS r, IS c,const MatFactorInfo *info) 334b24902e0SBarry Smith { 3354eb8e494SKris Buschelman /************************************************************************/ 3364eb8e494SKris Buschelman /* Input */ 3374eb8e494SKris Buschelman /* A - matrix to factor */ 3384eb8e494SKris Buschelman /* r - row permutation (ignored) */ 3394eb8e494SKris Buschelman /* c - column permutation (ignored) */ 3404eb8e494SKris Buschelman /* */ 3414eb8e494SKris Buschelman /* Output */ 3424eb8e494SKris Buschelman /* F - matrix storing the factorization; */ 3434eb8e494SKris Buschelman /************************************************************************/ 344f0c56d0fSKris Buschelman Mat_LUSOL *lusol; 345dfbe8321SBarry Smith PetscErrorCode ierr; 346dfbe8321SBarry Smith int i, m, n, nz, nnz; 3474eb8e494SKris Buschelman 3484eb8e494SKris Buschelman PetscFunctionBegin; 3494eb8e494SKris Buschelman /************************************************************************/ 3504eb8e494SKris Buschelman /* Check the arguments. */ 3514eb8e494SKris Buschelman /************************************************************************/ 3524eb8e494SKris Buschelman 3534eb8e494SKris Buschelman ierr = MatGetSize(A, &m, &n);CHKERRQ(ierr); 3544eb8e494SKris Buschelman nz = ((Mat_SeqAIJ*)A->data)->nz; 3554eb8e494SKris Buschelman 3564eb8e494SKris Buschelman /************************************************************************/ 3574eb8e494SKris Buschelman /* Create the factorization. */ 3584eb8e494SKris Buschelman /************************************************************************/ 3594eb8e494SKris Buschelman 36035bd34faSBarry Smith F->ops->lufactornumeric = MatLUFactorNumeric_LUSOL; 36135bd34faSBarry Smith lusol = (Mat_LUSOL*)(F->spptr); 3624eb8e494SKris Buschelman 3634eb8e494SKris Buschelman /************************************************************************/ 3644eb8e494SKris Buschelman /* Initialize parameters */ 3654eb8e494SKris Buschelman /************************************************************************/ 3664eb8e494SKris Buschelman 3672205254eSKarl Rupp for (i = 0; i < 30; i++) { 3684eb8e494SKris Buschelman lusol->luparm[i] = 0; 3694eb8e494SKris Buschelman lusol->parmlu[i] = 0; 3704eb8e494SKris Buschelman } 3714eb8e494SKris Buschelman 3724eb8e494SKris Buschelman lusol->luparm[1] = -1; 3734eb8e494SKris Buschelman lusol->luparm[2] = 5; 3744eb8e494SKris Buschelman lusol->luparm[7] = 1; 3754eb8e494SKris Buschelman 3764eb8e494SKris Buschelman lusol->parmlu[0] = 1 / Factorization_Tolerance; 3774eb8e494SKris Buschelman lusol->parmlu[1] = 1 / Factorization_Tolerance; 3784eb8e494SKris Buschelman lusol->parmlu[2] = Factorization_Small_Tolerance; 3794eb8e494SKris Buschelman lusol->parmlu[3] = Factorization_Pivot_Tolerance; 3804eb8e494SKris Buschelman lusol->parmlu[4] = Factorization_Pivot_Tolerance; 3814eb8e494SKris Buschelman lusol->parmlu[5] = 3.0; 3824eb8e494SKris Buschelman lusol->parmlu[6] = 0.3; 3834eb8e494SKris Buschelman lusol->parmlu[7] = 0.6; 3844eb8e494SKris Buschelman 3854eb8e494SKris Buschelman /************************************************************************/ 3864eb8e494SKris Buschelman /* Allocate the workspace needed by LUSOL. */ 3874eb8e494SKris Buschelman /************************************************************************/ 3884eb8e494SKris Buschelman 3894eb8e494SKris Buschelman lusol->elbowroom = PetscMax(lusol->elbowroom, info->fill); 3904eb8e494SKris Buschelman nnz = PetscMax((int)(lusol->elbowroom*nz), 5*n); 3914eb8e494SKris Buschelman 3924eb8e494SKris Buschelman lusol->n = n; 3934eb8e494SKris Buschelman lusol->nz = nz; 3944eb8e494SKris Buschelman lusol->nnz = nnz; 3954eb8e494SKris Buschelman lusol->luroom = 1.75; 3964eb8e494SKris Buschelman 3974eb8e494SKris Buschelman ierr = PetscMalloc(sizeof(int)*n,&lusol->ip); 3984eb8e494SKris Buschelman ierr = PetscMalloc(sizeof(int)*n,&lusol->iq); 3994eb8e494SKris Buschelman ierr = PetscMalloc(sizeof(int)*n,&lusol->lenc); 4004eb8e494SKris Buschelman ierr = PetscMalloc(sizeof(int)*n,&lusol->lenr); 4014eb8e494SKris Buschelman ierr = PetscMalloc(sizeof(int)*n,&lusol->locc); 4024eb8e494SKris Buschelman ierr = PetscMalloc(sizeof(int)*n,&lusol->locr); 4034eb8e494SKris Buschelman ierr = PetscMalloc(sizeof(int)*n,&lusol->iploc); 4044eb8e494SKris Buschelman ierr = PetscMalloc(sizeof(int)*n,&lusol->iqloc); 4054eb8e494SKris Buschelman ierr = PetscMalloc(sizeof(int)*n,&lusol->ipinv); 4064eb8e494SKris Buschelman ierr = PetscMalloc(sizeof(int)*n,&lusol->iqinv); 4074eb8e494SKris Buschelman ierr = PetscMalloc(sizeof(double)*n,&lusol->mnsw); 4084eb8e494SKris Buschelman ierr = PetscMalloc(sizeof(double)*n,&lusol->mnsv); 4094eb8e494SKris Buschelman 410dcca6d9dSJed Brown ierr = PetscMalloc3(nnz,&lusol->data,nnz,&lusol->indc,nnz,&lusol->indr);CHKERRQ(ierr); 4112205254eSKarl Rupp 4124eb8e494SKris Buschelman lusol->CleanUpLUSOL = PETSC_TRUE; 41335bd34faSBarry Smith F->ops->lufactornumeric = MatLUFactorNumeric_LUSOL; 4144eb8e494SKris Buschelman PetscFunctionReturn(0); 4154eb8e494SKris Buschelman } 4164eb8e494SKris Buschelman 417ea799195SBarry Smith PetscErrorCode MatFactorGetSolverType_seqaij_lusol(Mat A,MatSolverType *type) 41835bd34faSBarry Smith { 41935bd34faSBarry Smith PetscFunctionBegin; 4202692d6eeSBarry Smith *type = MATSOLVERLUSOL; 42135bd34faSBarry Smith PetscFunctionReturn(0); 42235bd34faSBarry Smith } 42335bd34faSBarry Smith 4248cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_lusol(Mat A,MatFactorType ftype,Mat *F) 425521d7252SBarry Smith { 426b24902e0SBarry Smith Mat B; 427f0c56d0fSKris Buschelman Mat_LUSOL *lusol; 428b24902e0SBarry Smith PetscErrorCode ierr; 42935bd34faSBarry Smith int m, n; 4304eb8e494SKris Buschelman 4314eb8e494SKris Buschelman PetscFunctionBegin; 4324eb8e494SKris Buschelman ierr = MatGetSize(A, &m, &n);CHKERRQ(ierr); 433ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 434b24902e0SBarry Smith ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr); 435b24902e0SBarry Smith ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 4360298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr); 4374eb8e494SKris Buschelman 438b00a9115SJed Brown ierr = PetscNewLog(B,&lusol);CHKERRQ(ierr); 439b24902e0SBarry Smith B->spptr = lusol; 4402f71e704SKris Buschelman 441f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = MatLUFactorSymbolic_LUSOL; 442f0c56d0fSKris Buschelman B->ops->destroy = MatDestroy_LUSOL; 4432205254eSKarl Rupp 4443ca39a21SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_seqaij_lusol);CHKERRQ(ierr); 4452205254eSKarl Rupp 446d5f3da31SBarry Smith B->factortype = MAT_FACTOR_LU; 44700c67f3bSHong Zhang ierr = PetscFree(B->solvertype);CHKERRQ(ierr); 44800c67f3bSHong Zhang ierr = PetscStrallocpy(MATSOLVERLUSOL,&B->solvertype);CHKERRQ(ierr); 44900c67f3bSHong Zhang 450f0c56d0fSKris Buschelman PetscFunctionReturn(0); 451f0c56d0fSKris Buschelman } 452f0c56d0fSKris Buschelman 4533ca39a21SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_Lusol(void) 45442c9c57cSBarry Smith { 45542c9c57cSBarry Smith PetscErrorCode ierr; 45642c9c57cSBarry Smith 45742c9c57cSBarry Smith PetscFunctionBegin; 4583ca39a21SBarry Smith ierr = MatSolverTypeRegister(MATSOLVERLUSOL,MATSEQAIJ, MAT_FACTOR_LU,MatGetFactor_seqaij_lusol);CHKERRQ(ierr); 45942c9c57cSBarry Smith PetscFunctionReturn(0); 46042c9c57cSBarry Smith } 46142c9c57cSBarry Smith 4622f71e704SKris Buschelman /*MC 4632692d6eeSBarry Smith MATSOLVERLUSOL - "lusol" - Provides direct solvers (LU) for sequential matrices 4642f71e704SKris Buschelman via the external package LUSOL. 4652f71e704SKris Buschelman 4662f71e704SKris Buschelman If LUSOL is installed (see the manual for 4672f71e704SKris Buschelman instructions on how to declare the existence of external packages), 4682f71e704SKris Buschelman 46941c8de11SBarry Smith Works with MATSEQAIJ matrices 4702f71e704SKris Buschelman 4712f71e704SKris Buschelman Level: beginner 4722f71e704SKris Buschelman 4733ca39a21SBarry Smith .seealso: PCLU, PCFactorSetMatSolverType(), MatSolverType 47441c8de11SBarry Smith 4752f71e704SKris Buschelman M*/ 476