xref: /petsc/src/mat/impls/aij/seq/lusol/lusol.c (revision 719d5645761d844e4357b7ee00a3296dffe0b787)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
34eb8e494SKris Buschelman /*
44eb8e494SKris Buschelman         Provides an interface to the LUSOL package of ....
54eb8e494SKris Buschelman 
64eb8e494SKris Buschelman */
74eb8e494SKris Buschelman #include "src/mat/impls/aij/seq/aij.h"
84eb8e494SKris Buschelman 
94eb8e494SKris Buschelman #if defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
104eb8e494SKris Buschelman #define LU1FAC   lu1fac_
114eb8e494SKris Buschelman #define LU6SOL   lu6sol_
124eb8e494SKris Buschelman #define M1PAGE   m1page_
134eb8e494SKris Buschelman #define M5SETX   m5setx_
144eb8e494SKris Buschelman #define M6RDEL   m6rdel_
154eb8e494SKris Buschelman #elif !defined(PETSC_HAVE_FORTRAN_CAPS)
164eb8e494SKris Buschelman #define LU1FAC   lu1fac
174eb8e494SKris Buschelman #define LU6SOL   lu6sol
184eb8e494SKris Buschelman #define M1PAGE   m1page
194eb8e494SKris Buschelman #define M5SETX   m5setx
204eb8e494SKris Buschelman #define M6RDEL   m6rdel
214eb8e494SKris Buschelman #endif
224eb8e494SKris Buschelman 
234eb8e494SKris Buschelman EXTERN_C_BEGIN
244eb8e494SKris Buschelman /*
254eb8e494SKris Buschelman     Dummy symbols that the MINOS files mi25bfac.f and mi15blas.f may require
264eb8e494SKris Buschelman */
274eb8e494SKris Buschelman void PETSC_STDCALL M1PAGE() {
284eb8e494SKris Buschelman   ;
294eb8e494SKris Buschelman }
304eb8e494SKris Buschelman void PETSC_STDCALL M5SETX() {
314eb8e494SKris Buschelman   ;
324eb8e494SKris Buschelman }
334eb8e494SKris Buschelman 
344eb8e494SKris Buschelman void PETSC_STDCALL M6RDEL() {
354eb8e494SKris Buschelman   ;
364eb8e494SKris Buschelman }
374eb8e494SKris Buschelman 
384eb8e494SKris Buschelman extern void PETSC_STDCALL LU1FAC (int *m, int *n, int *nnz, int *size, int *luparm,
394eb8e494SKris Buschelman                         double *parmlu, double *data, int *indc, int *indr,
404eb8e494SKris Buschelman                         int *rowperm, int *colperm, int *collen, int *rowlen,
414eb8e494SKris Buschelman                         int *colstart, int *rowstart, int *rploc, int *cploc,
424eb8e494SKris Buschelman                         int *rpinv, int *cpinv, double *w, int *inform);
434eb8e494SKris Buschelman 
444eb8e494SKris Buschelman extern void PETSC_STDCALL LU6SOL (int *mode, int *m, int *n, double *rhs, double *x,
454eb8e494SKris Buschelman                         int *size, int *luparm, double *parmlu, double *data,
464eb8e494SKris Buschelman                         int *indc, int *indr, int *rowperm, int *colperm,
474eb8e494SKris Buschelman                         int *collen, int *rowlen, int *colstart, int *rowstart,
484eb8e494SKris Buschelman                         int *inform);
492f71e704SKris Buschelman EXTERN_C_END
504eb8e494SKris Buschelman 
51dfbe8321SBarry 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 
804eb8e494SKris Buschelman   PetscTruth 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  *
1054eb8e494SKris Buschelman  *  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  *
1594eb8e494SKris Buschelman  *  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 
1804eb8e494SKris Buschelman #undef __FUNCT__
181f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_LUSOL"
182dfbe8321SBarry Smith PetscErrorCode MatDestroy_LUSOL(Mat A)
183dfbe8321SBarry Smith {
184dfbe8321SBarry Smith   PetscErrorCode ierr;
185f0c56d0fSKris Buschelman   Mat_LUSOL      *lusol=(Mat_LUSOL *)A->spptr;
1864eb8e494SKris Buschelman 
1874eb8e494SKris Buschelman   PetscFunctionBegin;
1884eb8e494SKris Buschelman   if (lusol->CleanUpLUSOL) {
1894eb8e494SKris Buschelman     ierr = PetscFree(lusol->ip);CHKERRQ(ierr);
1904eb8e494SKris Buschelman     ierr = PetscFree(lusol->iq);CHKERRQ(ierr);
1914eb8e494SKris Buschelman     ierr = PetscFree(lusol->lenc);CHKERRQ(ierr);
1924eb8e494SKris Buschelman     ierr = PetscFree(lusol->lenr);CHKERRQ(ierr);
1934eb8e494SKris Buschelman     ierr = PetscFree(lusol->locc);CHKERRQ(ierr);
1944eb8e494SKris Buschelman     ierr = PetscFree(lusol->locr);CHKERRQ(ierr);
1954eb8e494SKris Buschelman     ierr = PetscFree(lusol->iploc);CHKERRQ(ierr);
1964eb8e494SKris Buschelman     ierr = PetscFree(lusol->iqloc);CHKERRQ(ierr);
1974eb8e494SKris Buschelman     ierr = PetscFree(lusol->ipinv);CHKERRQ(ierr);
1984eb8e494SKris Buschelman     ierr = PetscFree(lusol->iqinv);CHKERRQ(ierr);
1994eb8e494SKris Buschelman     ierr = PetscFree(lusol->mnsw);CHKERRQ(ierr);
2004eb8e494SKris Buschelman     ierr = PetscFree(lusol->mnsv);CHKERRQ(ierr);
2014eb8e494SKris Buschelman     ierr = PetscFree(lusol->indc);CHKERRQ(ierr);
2024eb8e494SKris Buschelman   }
203b24902e0SBarry Smith   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
2044eb8e494SKris Buschelman   PetscFunctionReturn(0);
2054eb8e494SKris Buschelman }
2064eb8e494SKris Buschelman 
2074eb8e494SKris Buschelman #undef __FUNCT__
208f0c56d0fSKris Buschelman #define __FUNCT__  "MatSolve_LUSOL"
2096849ba73SBarry Smith PetscErrorCode MatSolve_LUSOL(Mat A,Vec b,Vec x)
2106849ba73SBarry Smith {
211f0c56d0fSKris Buschelman   Mat_LUSOL      *lusol=(Mat_LUSOL*)A->spptr;
2124eb8e494SKris Buschelman   double         *bb,*xx;
2134eb8e494SKris Buschelman   int            mode=5;
2146849ba73SBarry Smith   PetscErrorCode ierr;
2156849ba73SBarry Smith   int            i,m,n,nnz,status;
2164eb8e494SKris Buschelman 
2174eb8e494SKris Buschelman   PetscFunctionBegin;
2184eb8e494SKris Buschelman   ierr = VecGetArray(x, &xx);CHKERRQ(ierr);
2194eb8e494SKris Buschelman   ierr = VecGetArray(b, &bb);CHKERRQ(ierr);
2204eb8e494SKris Buschelman 
2214eb8e494SKris Buschelman   m = n = lusol->n;
2224eb8e494SKris Buschelman   nnz = lusol->nnz;
2234eb8e494SKris Buschelman 
2244eb8e494SKris Buschelman   for (i = 0; i < m; i++)
2254eb8e494SKris Buschelman     {
2264eb8e494SKris Buschelman       lusol->mnsv[i] = bb[i];
2274eb8e494SKris Buschelman     }
2284eb8e494SKris Buschelman 
2294eb8e494SKris Buschelman   LU6SOL(&mode, &m, &n, lusol->mnsv, xx, &nnz,
2304eb8e494SKris Buschelman          lusol->luparm, lusol->parmlu, lusol->data,
2314eb8e494SKris Buschelman          lusol->indc, lusol->indr, lusol->ip, lusol->iq,
2324eb8e494SKris Buschelman          lusol->lenc, lusol->lenr, lusol->locc, lusol->locr, &status);
2334eb8e494SKris Buschelman 
234b24902e0SBarry Smith   if (status != 0) SETERRQ1(PETSC_ERR_ARG_SIZ,"solve failed, error code %d",status);
2354eb8e494SKris Buschelman 
2364eb8e494SKris Buschelman   ierr = VecRestoreArray(x, &xx);CHKERRQ(ierr);
2374eb8e494SKris Buschelman   ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr);
2384eb8e494SKris Buschelman   PetscFunctionReturn(0);
2394eb8e494SKris Buschelman }
2404eb8e494SKris Buschelman 
2414eb8e494SKris Buschelman #undef __FUNCT__
242f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_LUSOL"
243*719d5645SBarry Smith PetscErrorCode MatLUFactorNumeric_LUSOL(Mat F,Mat A,MatFactorInfo *info)
2446849ba73SBarry Smith {
2454eb8e494SKris Buschelman   Mat_SeqAIJ     *a;
246*719d5645SBarry Smith   Mat_LUSOL      *lusol = (Mat_LUSOL*)F->spptr;
2476849ba73SBarry Smith   PetscErrorCode ierr;
2484eb8e494SKris Buschelman   int            m, n, nz, nnz, status;
2496849ba73SBarry Smith   int            i, rs, re;
2504eb8e494SKris Buschelman   int            factorizations;
2514eb8e494SKris Buschelman 
2524eb8e494SKris Buschelman   PetscFunctionBegin;
2534eb8e494SKris Buschelman   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);CHKERRQ(ierr);
2544eb8e494SKris Buschelman   a = (Mat_SeqAIJ *)A->data;
2554eb8e494SKris Buschelman 
256b24902e0SBarry Smith   if (m != lusol->n) SETERRQ(PETSC_ERR_ARG_SIZ,"factorization struct inconsistent");
2574eb8e494SKris Buschelman 
2584eb8e494SKris Buschelman   factorizations = 0;
2594eb8e494SKris Buschelman   do
2604eb8e494SKris Buschelman     {
2614eb8e494SKris Buschelman       /*******************************************************************/
2624eb8e494SKris Buschelman       /* Check the workspace allocation.                                 */
2634eb8e494SKris Buschelman       /*******************************************************************/
2644eb8e494SKris Buschelman 
2654eb8e494SKris Buschelman       nz = a->nz;
2664eb8e494SKris Buschelman       nnz = PetscMax(lusol->nnz, (int)(lusol->elbowroom*nz));
2674eb8e494SKris Buschelman       nnz = PetscMax(nnz, 5*n);
2684eb8e494SKris Buschelman 
2694eb8e494SKris Buschelman       if (nnz < lusol->luparm[12]){
2704eb8e494SKris Buschelman         nnz = (int)(lusol->luroom * lusol->luparm[12]);
2714eb8e494SKris Buschelman       } else if ((factorizations > 0) && (lusol->luroom < 6)){
2724eb8e494SKris Buschelman         lusol->luroom += 0.1;
2734eb8e494SKris Buschelman       }
2744eb8e494SKris Buschelman 
2754eb8e494SKris Buschelman       nnz = PetscMax(nnz, (int)(lusol->luroom*(lusol->luparm[22] + lusol->luparm[23])));
2764eb8e494SKris Buschelman 
2774eb8e494SKris Buschelman       if (nnz > lusol->nnz){
2784eb8e494SKris Buschelman         ierr = PetscFree(lusol->indc);CHKERRQ(ierr);
2794eb8e494SKris Buschelman         ierr        = PetscMalloc((sizeof(double)+2*sizeof(int))*nnz,&lusol->indc);CHKERRQ(ierr);
2804eb8e494SKris Buschelman         lusol->indr = lusol->indc + nnz;
2814eb8e494SKris Buschelman         lusol->data = (double *)(lusol->indr + nnz);
2824eb8e494SKris Buschelman         lusol->nnz  = nnz;
2834eb8e494SKris Buschelman       }
2844eb8e494SKris Buschelman 
2854eb8e494SKris Buschelman       /*******************************************************************/
2864eb8e494SKris Buschelman       /* Fill in the data for the problem.      (1-based Fortran style)  */
2874eb8e494SKris Buschelman       /*******************************************************************/
2884eb8e494SKris Buschelman 
2894eb8e494SKris Buschelman       nz = 0;
2904eb8e494SKris Buschelman       for (i = 0; i < n; i++)
2914eb8e494SKris Buschelman         {
2924eb8e494SKris Buschelman           rs = a->i[i];
2934eb8e494SKris Buschelman           re = a->i[i+1];
2944eb8e494SKris Buschelman 
2954eb8e494SKris Buschelman           while (rs < re)
2964eb8e494SKris Buschelman             {
2974eb8e494SKris Buschelman               if (a->a[rs] != 0.0)
2984eb8e494SKris Buschelman                 {
2994eb8e494SKris Buschelman                   lusol->indc[nz] = i + 1;
3004eb8e494SKris Buschelman                   lusol->indr[nz] = a->j[rs] + 1;
3014eb8e494SKris Buschelman                   lusol->data[nz] = a->a[rs];
3024eb8e494SKris Buschelman                   nz++;
3034eb8e494SKris Buschelman                 }
3044eb8e494SKris Buschelman               rs++;
3054eb8e494SKris Buschelman             }
3064eb8e494SKris Buschelman         }
3074eb8e494SKris Buschelman 
3084eb8e494SKris Buschelman       /*******************************************************************/
3094eb8e494SKris Buschelman       /* Do the factorization.                                           */
3104eb8e494SKris Buschelman       /*******************************************************************/
3114eb8e494SKris Buschelman 
3124eb8e494SKris Buschelman       LU1FAC(&m, &n, &nz, &nnz,
3134eb8e494SKris Buschelman              lusol->luparm, lusol->parmlu, lusol->data,
3144eb8e494SKris Buschelman              lusol->indc, lusol->indr, lusol->ip, lusol->iq,
3154eb8e494SKris Buschelman              lusol->lenc, lusol->lenr, lusol->locc, lusol->locr,
3164eb8e494SKris Buschelman              lusol->iploc, lusol->iqloc, lusol->ipinv,
3174eb8e494SKris Buschelman              lusol->iqinv, lusol->mnsw, &status);
3184eb8e494SKris Buschelman 
3194eb8e494SKris Buschelman       switch(status)
3204eb8e494SKris Buschelman         {
3214eb8e494SKris Buschelman         case 0:		/* factored */
3224eb8e494SKris Buschelman           break;
3234eb8e494SKris Buschelman 
3244eb8e494SKris Buschelman         case 7:		/* insufficient memory */
3254eb8e494SKris Buschelman           break;
3264eb8e494SKris Buschelman 
3274eb8e494SKris Buschelman         case 1:
3284eb8e494SKris Buschelman         case -1:		/* singular */
329e005ede5SBarry Smith           SETERRQ(PETSC_ERR_LIB,"Singular matrix");
3304eb8e494SKris Buschelman 
3314eb8e494SKris Buschelman         case 3:
3324eb8e494SKris Buschelman         case 4:		/* error conditions */
333e005ede5SBarry Smith           SETERRQ(PETSC_ERR_LIB,"matrix error");
3344eb8e494SKris Buschelman 
3354eb8e494SKris Buschelman         default:		/* unknown condition */
336e005ede5SBarry Smith           SETERRQ(PETSC_ERR_LIB,"matrix unknown return code");
3374eb8e494SKris Buschelman         }
3384eb8e494SKris Buschelman 
3394eb8e494SKris Buschelman       factorizations++;
3404eb8e494SKris Buschelman     } while (status == 7);
341*719d5645SBarry Smith   F->ops->solve   = MatSolve_LUSOL;
342*719d5645SBarry Smith   F->assembled    = PETSC_TRUE;
343*719d5645SBarry Smith   F->preallocated = PETSC_TRUE;
3444eb8e494SKris Buschelman   PetscFunctionReturn(0);
3454eb8e494SKris Buschelman }
3464eb8e494SKris Buschelman 
3474eb8e494SKris Buschelman #undef __FUNCT__
348f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_LUSOL"
349b24902e0SBarry Smith PetscErrorCode MatLUFactorSymbolic_LUSOL(Mat A, IS r, IS c,MatFactorInfo *info, Mat *F)
350b24902e0SBarry Smith {
3514eb8e494SKris Buschelman   /************************************************************************/
3524eb8e494SKris Buschelman   /* Input                                                                */
3534eb8e494SKris Buschelman   /*     A  - matrix to factor                                            */
3544eb8e494SKris Buschelman   /*     r  - row permutation (ignored)                                   */
3554eb8e494SKris Buschelman   /*     c  - column permutation (ignored)                                */
3564eb8e494SKris Buschelman   /*                                                                      */
3574eb8e494SKris Buschelman   /* Output                                                               */
3584eb8e494SKris Buschelman   /*     F  - matrix storing the factorization;                           */
3594eb8e494SKris Buschelman   /************************************************************************/
3604eb8e494SKris Buschelman   Mat            B;
361f0c56d0fSKris Buschelman   Mat_LUSOL      *lusol;
362dfbe8321SBarry Smith   PetscErrorCode ierr;
363dfbe8321SBarry Smith   int            i, m, n, nz, nnz;
3644eb8e494SKris Buschelman 
3654eb8e494SKris Buschelman   PetscFunctionBegin;
3664eb8e494SKris Buschelman 
3674eb8e494SKris Buschelman   /************************************************************************/
3684eb8e494SKris Buschelman   /* Check the arguments.                                                 */
3694eb8e494SKris Buschelman   /************************************************************************/
3704eb8e494SKris Buschelman 
3714eb8e494SKris Buschelman   ierr = MatGetSize(A, &m, &n);CHKERRQ(ierr);
3724eb8e494SKris Buschelman   nz = ((Mat_SeqAIJ *)A->data)->nz;
3734eb8e494SKris Buschelman 
3744eb8e494SKris Buschelman   /************************************************************************/
3754eb8e494SKris Buschelman   /* Create the factorization.                                            */
3764eb8e494SKris Buschelman   /************************************************************************/
3774eb8e494SKris Buschelman 
378f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatLUFactorNumeric_LUSOL;
379f0c56d0fSKris Buschelman   lusol                   = (Mat_LUSOL*)(B->spptr);
3804eb8e494SKris Buschelman 
3814eb8e494SKris Buschelman   /************************************************************************/
3824eb8e494SKris Buschelman   /* Initialize parameters                                                */
3834eb8e494SKris Buschelman   /************************************************************************/
3844eb8e494SKris Buschelman 
3854eb8e494SKris Buschelman   for (i = 0; i < 30; i++)
3864eb8e494SKris Buschelman     {
3874eb8e494SKris Buschelman       lusol->luparm[i] = 0;
3884eb8e494SKris Buschelman       lusol->parmlu[i] = 0;
3894eb8e494SKris Buschelman     }
3904eb8e494SKris Buschelman 
3914eb8e494SKris Buschelman   lusol->luparm[1] = -1;
3924eb8e494SKris Buschelman   lusol->luparm[2] = 5;
3934eb8e494SKris Buschelman   lusol->luparm[7] = 1;
3944eb8e494SKris Buschelman 
3954eb8e494SKris Buschelman   lusol->parmlu[0] = 1 / Factorization_Tolerance;
3964eb8e494SKris Buschelman   lusol->parmlu[1] = 1 / Factorization_Tolerance;
3974eb8e494SKris Buschelman   lusol->parmlu[2] = Factorization_Small_Tolerance;
3984eb8e494SKris Buschelman   lusol->parmlu[3] = Factorization_Pivot_Tolerance;
3994eb8e494SKris Buschelman   lusol->parmlu[4] = Factorization_Pivot_Tolerance;
4004eb8e494SKris Buschelman   lusol->parmlu[5] = 3.0;
4014eb8e494SKris Buschelman   lusol->parmlu[6] = 0.3;
4024eb8e494SKris Buschelman   lusol->parmlu[7] = 0.6;
4034eb8e494SKris Buschelman 
4044eb8e494SKris Buschelman   /************************************************************************/
4054eb8e494SKris Buschelman   /* Allocate the workspace needed by LUSOL.                              */
4064eb8e494SKris Buschelman   /************************************************************************/
4074eb8e494SKris Buschelman 
4084eb8e494SKris Buschelman   lusol->elbowroom = PetscMax(lusol->elbowroom, info->fill);
4094eb8e494SKris Buschelman   nnz = PetscMax((int)(lusol->elbowroom*nz), 5*n);
4104eb8e494SKris Buschelman 
4114eb8e494SKris Buschelman   lusol->n = n;
4124eb8e494SKris Buschelman   lusol->nz = nz;
4134eb8e494SKris Buschelman   lusol->nnz = nnz;
4144eb8e494SKris Buschelman   lusol->luroom = 1.75;
4154eb8e494SKris Buschelman 
4164eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->ip);
4174eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->iq);
4184eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->lenc);
4194eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->lenr);
4204eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->locc);
4214eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->locr);
4224eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->iploc);
4234eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->iqloc);
4244eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->ipinv);
4254eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->iqinv);
4264eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(double)*n,&lusol->mnsw);
4274eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(double)*n,&lusol->mnsv);
4284eb8e494SKris Buschelman 
4294eb8e494SKris Buschelman   ierr        = PetscMalloc((sizeof(double)+2*sizeof(int))*nnz,&lusol->indc);
4304eb8e494SKris Buschelman   lusol->indr = lusol->indc + nnz;
4314eb8e494SKris Buschelman   lusol->data = (double *)(lusol->indr + nnz);
4324eb8e494SKris Buschelman   lusol->CleanUpLUSOL = PETSC_TRUE;
433db4efbfdSBarry Smith   B->ops->lufactornumeric  = MatLUFactorNumeric_LUSOL;
434db4efbfdSBarry Smith   B->ops->solve            = MatSolve_LUSOL;
4354eb8e494SKris Buschelman   *F = B;
4364eb8e494SKris Buschelman   PetscFunctionReturn(0);
4374eb8e494SKris Buschelman }
4384eb8e494SKris Buschelman 
4394eb8e494SKris Buschelman #undef __FUNCT__
440b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_lusol"
4415c9eb25fSBarry Smith PetscErrorCode MatGetFactor_seqaij_lusol(Mat A,MatFactorType ftype,Mat *F)
442521d7252SBarry Smith {
443b24902e0SBarry Smith   Mat            B;
444f0c56d0fSKris Buschelman   Mat_LUSOL      *lusol;
445b24902e0SBarry Smith   PetscErrorCode ierr;
446b24902e0SBarry Smith   int            i, m, n, nz, nnz;
4474eb8e494SKris Buschelman 
4484eb8e494SKris Buschelman   PetscFunctionBegin;
4494eb8e494SKris Buschelman   ierr = MatGetSize(A, &m, &n);CHKERRQ(ierr);
450b24902e0SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
451b24902e0SBarry Smith   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
452b24902e0SBarry Smith   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
453b24902e0SBarry Smith   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
4544eb8e494SKris Buschelman 
45538f2d2fdSLisandro Dalcin   ierr = PetscNewLog(B,Mat_LUSOL,&lusol);CHKERRQ(ierr);
456b24902e0SBarry Smith   B->spptr                 = lusol;
4572f71e704SKris Buschelman 
458f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic = MatLUFactorSymbolic_LUSOL;
459f0c56d0fSKris Buschelman   B->ops->destroy          = MatDestroy_LUSOL;
4605c9eb25fSBarry Smith   B->factor                = MAT_FACTOR_LU;
461f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
462f0c56d0fSKris Buschelman }
463f0c56d0fSKris Buschelman 
4642f71e704SKris Buschelman /*MC
465fafad747SKris Buschelman   MATLUSOL - MATLUSOL = "lusol" - A matrix type providing direct solvers (LU) for sequential matrices
4662f71e704SKris Buschelman   via the external package LUSOL.
4672f71e704SKris Buschelman 
4682f71e704SKris Buschelman   If LUSOL is installed (see the manual for
4692f71e704SKris Buschelman   instructions on how to declare the existence of external packages),
4702f71e704SKris Buschelman   a matrix type can be constructed which invokes LUSOL solvers.
4712f71e704SKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATLUSOL).
4722f71e704SKris Buschelman   This matrix type is only supported for double precision real.
4732f71e704SKris Buschelman 
4742f71e704SKris Buschelman   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
475f0c56d0fSKris Buschelman   supported for this matrix type.  MatConvert can be called for a fast inplace conversion
476f0c56d0fSKris Buschelman   to and from the MATSEQAIJ matrix type.
4772f71e704SKris Buschelman 
4782f71e704SKris Buschelman   Options Database Keys:
4790bad9183SKris Buschelman . -mat_type lusol - sets the matrix type to "lusol" during a call to MatSetFromOptions()
4802f71e704SKris Buschelman 
4812f71e704SKris Buschelman    Level: beginner
4822f71e704SKris Buschelman 
4832f71e704SKris Buschelman .seealso: PCLU
4842f71e704SKris Buschelman M*/
485