xref: /petsc/src/mat/impls/aij/seq/lusol/lusol.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
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 EXTERN_C_BEGIN
234eb8e494SKris Buschelman /*
244eb8e494SKris Buschelman     Dummy symbols that the MINOS files mi25bfac.f and mi15blas.f may require
254eb8e494SKris Buschelman */
26a6dfd86eSKarl Rupp void PETSC_STDCALL M1PAGE()
27a6dfd86eSKarl Rupp {
284eb8e494SKris Buschelman   ;
294eb8e494SKris Buschelman }
30a6dfd86eSKarl Rupp void PETSC_STDCALL M5SETX()
31a6dfd86eSKarl Rupp {
324eb8e494SKris Buschelman   ;
334eb8e494SKris Buschelman }
344eb8e494SKris Buschelman 
35a6dfd86eSKarl Rupp void PETSC_STDCALL M6RDEL()
36a6dfd86eSKarl Rupp {
374eb8e494SKris Buschelman   ;
384eb8e494SKris Buschelman }
394eb8e494SKris Buschelman 
404eb8e494SKris Buschelman extern void PETSC_STDCALL LU1FAC (int *m, int *n, int *nnz, int *size, int *luparm,
414eb8e494SKris Buschelman                         double *parmlu, double *data, int *indc, int *indr,
424eb8e494SKris Buschelman                         int *rowperm, int *colperm, int *collen, int *rowlen,
434eb8e494SKris Buschelman                         int *colstart, int *rowstart, int *rploc, int *cploc,
444eb8e494SKris Buschelman                         int *rpinv, int *cpinv, double *w, int *inform);
454eb8e494SKris Buschelman 
464eb8e494SKris Buschelman extern void PETSC_STDCALL LU6SOL (int *mode, int *m, int *n, double *rhs, double *x,
474eb8e494SKris Buschelman                         int *size, int *luparm, double *parmlu, double *data,
484eb8e494SKris Buschelman                         int *indc, int *indr, int *rowperm, int *colperm,
494eb8e494SKris Buschelman                         int *collen, int *rowlen, int *colstart, int *rowstart,
504eb8e494SKris Buschelman                         int *inform);
512f71e704SKris Buschelman EXTERN_C_END
524eb8e494SKris Buschelman 
5309573ac7SBarry Smith extern PetscErrorCode MatDuplicate_LUSOL(Mat,MatDuplicateOption,Mat*);
54f0c56d0fSKris Buschelman 
55f0c56d0fSKris Buschelman typedef struct  {
564eb8e494SKris Buschelman   double *data;
574eb8e494SKris Buschelman   int    *indc;
584eb8e494SKris Buschelman   int    *indr;
594eb8e494SKris Buschelman 
604eb8e494SKris Buschelman   int    *ip;
614eb8e494SKris Buschelman   int    *iq;
624eb8e494SKris Buschelman   int    *lenc;
634eb8e494SKris Buschelman   int    *lenr;
644eb8e494SKris Buschelman   int    *locc;
654eb8e494SKris Buschelman   int    *locr;
664eb8e494SKris Buschelman   int    *iploc;
674eb8e494SKris Buschelman   int    *iqloc;
684eb8e494SKris Buschelman   int    *ipinv;
694eb8e494SKris Buschelman   int    *iqinv;
704eb8e494SKris Buschelman   double *mnsw;
714eb8e494SKris Buschelman   double *mnsv;
724eb8e494SKris Buschelman 
734eb8e494SKris Buschelman   double elbowroom;
744eb8e494SKris Buschelman   double luroom;                /* Extra space allocated when factor fails   */
754eb8e494SKris Buschelman   double parmlu[30];            /* Input/output to LUSOL                     */
764eb8e494SKris Buschelman 
774eb8e494SKris Buschelman   int n;                        /* Number of rows/columns in matrix          */
784eb8e494SKris Buschelman   int nz;                       /* Number of nonzeros                        */
794eb8e494SKris Buschelman   int nnz;                      /* Number of nonzeros allocated for factors  */
804eb8e494SKris Buschelman   int luparm[30];               /* Input/output to LUSOL                     */
814eb8e494SKris Buschelman 
82ace3abfcSBarry Smith   PetscBool CleanUpLUSOL;
834eb8e494SKris Buschelman 
84f0c56d0fSKris Buschelman } Mat_LUSOL;
854eb8e494SKris Buschelman 
864eb8e494SKris Buschelman /*  LUSOL input/Output Parameters (Description uses C-style indexes
874eb8e494SKris Buschelman  *
884eb8e494SKris Buschelman  *  Input parameters                                        Typical value
894eb8e494SKris Buschelman  *
904eb8e494SKris Buschelman  *  luparm(0) = nout     File number for printed messages.         6
914eb8e494SKris Buschelman  *  luparm(1) = lprint   Print level.                              0
924eb8e494SKris Buschelman  *                    < 0 suppresses output.
934eb8e494SKris Buschelman  *                    = 0 gives error messages.
944eb8e494SKris Buschelman  *                    = 1 gives debug output from some of the
954eb8e494SKris Buschelman  *                        other routines in LUSOL.
964eb8e494SKris Buschelman  *                   >= 2 gives the pivot row and column and the
974eb8e494SKris Buschelman  *                        no. of rows and columns involved at
984eb8e494SKris Buschelman  *                        each elimination step in lu1fac.
994eb8e494SKris Buschelman  *  luparm(2) = maxcol   lu1fac: maximum number of columns         5
1004eb8e494SKris Buschelman  *                        searched allowed in a Markowitz-type
1014eb8e494SKris Buschelman  *                        search for the next pivot element.
1024eb8e494SKris Buschelman  *                        For some of the factorization, the
1034eb8e494SKris Buschelman  *                        number of rows searched is
1044eb8e494SKris Buschelman  *                        maxrow = maxcol - 1.
1054eb8e494SKris Buschelman  *
1064eb8e494SKris Buschelman  *
1074eb8e494SKris Buschelman  *  Output parameters
1084eb8e494SKris Buschelman  *
1094eb8e494SKris Buschelman  *  luparm(9) = inform   Return code from last call to any LU routine.
1104eb8e494SKris Buschelman  *  luparm(10) = nsing    No. of singularities marked in the
1114eb8e494SKris Buschelman  *                        output array w(*).
1124eb8e494SKris Buschelman  *  luparm(11) = jsing    Column index of last singularity.
1134eb8e494SKris Buschelman  *  luparm(12) = minlen   Minimum recommended value for  lena.
1144eb8e494SKris Buschelman  *  luparm(13) = maxlen   ?
1154eb8e494SKris Buschelman  *  luparm(14) = nupdat   No. of updates performed by the lu8 routines.
1164eb8e494SKris Buschelman  *  luparm(15) = nrank    No. of nonempty rows of U.
1174eb8e494SKris Buschelman  *  luparm(16) = ndens1   No. of columns remaining when the density of
1184eb8e494SKris Buschelman  *                        the matrix being factorized reached dens1.
1194eb8e494SKris Buschelman  *  luparm(17) = ndens2   No. of columns remaining when the density of
1204eb8e494SKris Buschelman  *                        the matrix being factorized reached dens2.
1214eb8e494SKris Buschelman  *  luparm(18) = jumin    The column index associated with dumin.
1224eb8e494SKris Buschelman  *  luparm(19) = numl0    No. of columns in initial  L.
1234eb8e494SKris Buschelman  *  luparm(20) = lenl0    Size of initial  L  (no. of nonzeros).
1244eb8e494SKris Buschelman  *  luparm(21) = lenu0    Size of initial  U.
1254eb8e494SKris Buschelman  *  luparm(22) = lenl     Size of current  L.
1264eb8e494SKris Buschelman  *  luparm(23) = lenu     Size of current  U.
1274eb8e494SKris Buschelman  *  luparm(24) = lrow     Length of row file.
1284eb8e494SKris Buschelman  *  luparm(25) = ncp      No. of compressions of LU data structures.
1294eb8e494SKris Buschelman  *  luparm(26) = mersum   lu1fac: sum of Markowitz merit counts.
1304eb8e494SKris Buschelman  *  luparm(27) = nutri    lu1fac: triangular rows in U.
1314eb8e494SKris Buschelman  *  luparm(28) = nltri    lu1fac: triangular rows in L.
1324eb8e494SKris Buschelman  *  luparm(29) =
1334eb8e494SKris Buschelman  *
1344eb8e494SKris Buschelman  *
1354eb8e494SKris Buschelman  *  Input parameters                                        Typical value
1364eb8e494SKris Buschelman  *
1374eb8e494SKris Buschelman  *  parmlu(0) = elmax1   Max multiplier allowed in  L           10.0
1384eb8e494SKris Buschelman  *                        during factor.
1394eb8e494SKris Buschelman  *  parmlu(1) = elmax2   Max multiplier allowed in  L           10.0
1404eb8e494SKris Buschelman  *                        during updates.
1414eb8e494SKris Buschelman  *  parmlu(2) = small    Absolute tolerance for             eps**0.8
1424eb8e494SKris Buschelman  *                        treating reals as zero.     IBM double: 3.0d-13
1434eb8e494SKris Buschelman  *  parmlu(3) = utol1    Absolute tol for flagging          eps**0.66667
1444eb8e494SKris Buschelman  *                        small diagonals of U.       IBM double: 3.7d-11
1454eb8e494SKris Buschelman  *  parmlu(4) = utol2    Relative tol for flagging          eps**0.66667
1464eb8e494SKris Buschelman  *                        small diagonals of U.       IBM double: 3.7d-11
1474eb8e494SKris Buschelman  *  parmlu(5) = uspace   Factor limiting waste space in  U.      3.0
1484eb8e494SKris Buschelman  *                        In lu1fac, the row or column lists
1494eb8e494SKris Buschelman  *                        are compressed if their length
1504eb8e494SKris Buschelman  *                        exceeds uspace times the length of
1514eb8e494SKris Buschelman  *                        either file after the last compression.
1524eb8e494SKris Buschelman  *  parmlu(6) = dens1    The density at which the Markowitz      0.3
1534eb8e494SKris Buschelman  *                        strategy should search maxcol columns
1544eb8e494SKris Buschelman  *                        and no rows.
1554eb8e494SKris Buschelman  *  parmlu(7) = dens2    the density at which the Markowitz      0.6
1564eb8e494SKris Buschelman  *                        strategy should search only 1 column
1574eb8e494SKris Buschelman  *                        or (preferably) use a dense LU for
1584eb8e494SKris Buschelman  *                        all the remaining rows and columns.
1594eb8e494SKris Buschelman  *
1604eb8e494SKris Buschelman  *
1614eb8e494SKris Buschelman  *  Output parameters
1624eb8e494SKris Buschelman  *
1634eb8e494SKris Buschelman  *  parmlu(9) = amax     Maximum element in  A.
1644eb8e494SKris Buschelman  *  parmlu(10) = elmax    Maximum multiplier in current  L.
1654eb8e494SKris Buschelman  *  parmlu(11) = umax     Maximum element in current  U.
1664eb8e494SKris Buschelman  *  parmlu(12) = dumax    Maximum diagonal in  U.
1674eb8e494SKris Buschelman  *  parmlu(13) = dumin    Minimum diagonal in  U.
1684eb8e494SKris Buschelman  *  parmlu(14) =
1694eb8e494SKris Buschelman  *  parmlu(15) =
1704eb8e494SKris Buschelman  *  parmlu(16) =
1714eb8e494SKris Buschelman  *  parmlu(17) =
1724eb8e494SKris Buschelman  *  parmlu(18) =
1734eb8e494SKris Buschelman  *  parmlu(19) = resid    lu6sol: residual after solve with U or U'.
1744eb8e494SKris Buschelman  *  ...
1754eb8e494SKris Buschelman  *  parmlu(29) =
1764eb8e494SKris Buschelman  */
1774eb8e494SKris Buschelman 
1784eb8e494SKris Buschelman #define Factorization_Tolerance       1e-1
1794eb8e494SKris Buschelman #define Factorization_Pivot_Tolerance pow(2.2204460492503131E-16, 2.0 / 3.0)
1804eb8e494SKris Buschelman #define Factorization_Small_Tolerance 1e-15 /* pow(DBL_EPSILON, 0.8) */
1814eb8e494SKris Buschelman 
1824eb8e494SKris Buschelman #undef __FUNCT__
183f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_LUSOL"
184dfbe8321SBarry Smith PetscErrorCode MatDestroy_LUSOL(Mat A)
185dfbe8321SBarry Smith {
186dfbe8321SBarry Smith   PetscErrorCode ierr;
187f0c56d0fSKris Buschelman   Mat_LUSOL      *lusol=(Mat_LUSOL*)A->spptr;
1884eb8e494SKris Buschelman 
1894eb8e494SKris Buschelman   PetscFunctionBegin;
190bf0cc555SLisandro Dalcin   if (lusol && lusol->CleanUpLUSOL) {
1914eb8e494SKris Buschelman     ierr = PetscFree(lusol->ip);CHKERRQ(ierr);
1924eb8e494SKris Buschelman     ierr = PetscFree(lusol->iq);CHKERRQ(ierr);
1934eb8e494SKris Buschelman     ierr = PetscFree(lusol->lenc);CHKERRQ(ierr);
1944eb8e494SKris Buschelman     ierr = PetscFree(lusol->lenr);CHKERRQ(ierr);
1954eb8e494SKris Buschelman     ierr = PetscFree(lusol->locc);CHKERRQ(ierr);
1964eb8e494SKris Buschelman     ierr = PetscFree(lusol->locr);CHKERRQ(ierr);
1974eb8e494SKris Buschelman     ierr = PetscFree(lusol->iploc);CHKERRQ(ierr);
1984eb8e494SKris Buschelman     ierr = PetscFree(lusol->iqloc);CHKERRQ(ierr);
1994eb8e494SKris Buschelman     ierr = PetscFree(lusol->ipinv);CHKERRQ(ierr);
2004eb8e494SKris Buschelman     ierr = PetscFree(lusol->iqinv);CHKERRQ(ierr);
2014eb8e494SKris Buschelman     ierr = PetscFree(lusol->mnsw);CHKERRQ(ierr);
2024eb8e494SKris Buschelman     ierr = PetscFree(lusol->mnsv);CHKERRQ(ierr);
20323bdbc58SBarry Smith     ierr = PetscFree3(lusol->data,lusol->indc,lusol->indr);CHKERRQ(ierr);
2044eb8e494SKris Buschelman   }
205bf0cc555SLisandro Dalcin   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
206b24902e0SBarry Smith   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
2074eb8e494SKris Buschelman   PetscFunctionReturn(0);
2084eb8e494SKris Buschelman }
2094eb8e494SKris Buschelman 
2104eb8e494SKris Buschelman #undef __FUNCT__
211f0c56d0fSKris Buschelman #define __FUNCT__  "MatSolve_LUSOL"
2126849ba73SBarry Smith PetscErrorCode MatSolve_LUSOL(Mat A,Vec b,Vec x)
2136849ba73SBarry Smith {
214f0c56d0fSKris Buschelman   Mat_LUSOL      *lusol=(Mat_LUSOL*)A->spptr;
2154eb8e494SKris Buschelman   double         *bb,*xx;
2164eb8e494SKris Buschelman   int            mode=5;
2176849ba73SBarry Smith   PetscErrorCode ierr;
2186849ba73SBarry Smith   int            i,m,n,nnz,status;
2194eb8e494SKris Buschelman 
2204eb8e494SKris Buschelman   PetscFunctionBegin;
2214eb8e494SKris Buschelman   ierr = VecGetArray(x, &xx);CHKERRQ(ierr);
2224eb8e494SKris Buschelman   ierr = VecGetArray(b, &bb);CHKERRQ(ierr);
2234eb8e494SKris Buschelman 
2244eb8e494SKris Buschelman   m   = n = lusol->n;
2254eb8e494SKris Buschelman   nnz = lusol->nnz;
2264eb8e494SKris Buschelman 
227*2205254eSKarl Rupp   for (i = 0; i < m; i++) lusol->mnsv[i] = bb[i];
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 
23465e19b50SBarry Smith   if (status) SETERRQ1(PETSC_COMM_SELF,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"
2430481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_LUSOL(Mat F,Mat A,const MatFactorInfo *info)
2446849ba73SBarry Smith {
2454eb8e494SKris Buschelman   Mat_SeqAIJ     *a;
246719d5645SBarry 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 
256e32f2f54SBarry Smith   if (m != lusol->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"factorization struct inconsistent");
2574eb8e494SKris Buschelman 
2584eb8e494SKris Buschelman   factorizations = 0;
259*2205254eSKarl Rupp   do {
2604eb8e494SKris Buschelman     /*******************************************************************/
2614eb8e494SKris Buschelman     /* Check the workspace allocation.                                 */
2624eb8e494SKris Buschelman     /*******************************************************************/
2634eb8e494SKris Buschelman 
2644eb8e494SKris Buschelman     nz  = a->nz;
2654eb8e494SKris Buschelman     nnz = PetscMax(lusol->nnz, (int)(lusol->elbowroom*nz));
2664eb8e494SKris Buschelman     nnz = PetscMax(nnz, 5*n);
2674eb8e494SKris Buschelman 
2684eb8e494SKris Buschelman     if (nnz < lusol->luparm[12]) {
2694eb8e494SKris Buschelman       nnz = (int)(lusol->luroom * lusol->luparm[12]);
2704eb8e494SKris Buschelman     } else if ((factorizations > 0) && (lusol->luroom < 6)) {
2714eb8e494SKris Buschelman       lusol->luroom += 0.1;
2724eb8e494SKris Buschelman     }
2734eb8e494SKris Buschelman 
2744eb8e494SKris Buschelman     nnz = PetscMax(nnz, (int)(lusol->luroom*(lusol->luparm[22] + lusol->luparm[23])));
2754eb8e494SKris Buschelman 
2764eb8e494SKris Buschelman     if (nnz > lusol->nnz) {
27723bdbc58SBarry Smith       ierr       = PetscFree3(lusol->data,lusol->indc,lusol->indr);CHKERRQ(ierr);
27823bdbc58SBarry Smith       ierr       = PetscMalloc3(nnz,double,&lusol->data,nnz,PetscInt,&lusol->indc,nnz,PetscInt,&lusol->indr);CHKERRQ(ierr);
2794eb8e494SKris Buschelman       lusol->nnz = nnz;
2804eb8e494SKris Buschelman     }
2814eb8e494SKris Buschelman 
2824eb8e494SKris Buschelman     /*******************************************************************/
2834eb8e494SKris Buschelman     /* Fill in the data for the problem.      (1-based Fortran style)  */
2844eb8e494SKris Buschelman     /*******************************************************************/
2854eb8e494SKris Buschelman 
2864eb8e494SKris Buschelman     nz = 0;
287*2205254eSKarl Rupp     for (i = 0; i < n; i++) {
2884eb8e494SKris Buschelman       rs = a->i[i];
2894eb8e494SKris Buschelman       re = a->i[i+1];
2904eb8e494SKris Buschelman 
291*2205254eSKarl Rupp       while (rs < re) {
292*2205254eSKarl Rupp         if (a->a[rs] != 0.0) {
2934eb8e494SKris Buschelman           lusol->indc[nz] = i + 1;
2944eb8e494SKris Buschelman           lusol->indr[nz] = a->j[rs] + 1;
2954eb8e494SKris Buschelman           lusol->data[nz] = a->a[rs];
2964eb8e494SKris Buschelman           nz++;
2974eb8e494SKris Buschelman         }
2984eb8e494SKris Buschelman         rs++;
2994eb8e494SKris Buschelman       }
3004eb8e494SKris Buschelman     }
3014eb8e494SKris Buschelman 
3024eb8e494SKris Buschelman     /*******************************************************************/
3034eb8e494SKris Buschelman     /* Do the factorization.                                           */
3044eb8e494SKris Buschelman     /*******************************************************************/
3054eb8e494SKris Buschelman 
3064eb8e494SKris Buschelman     LU1FAC(&m, &n, &nz, &nnz,
3074eb8e494SKris Buschelman            lusol->luparm, lusol->parmlu, lusol->data,
3084eb8e494SKris Buschelman            lusol->indc, lusol->indr, lusol->ip, lusol->iq,
3094eb8e494SKris Buschelman            lusol->lenc, lusol->lenr, lusol->locc, lusol->locr,
3104eb8e494SKris Buschelman            lusol->iploc, lusol->iqloc, lusol->ipinv,
3114eb8e494SKris Buschelman            lusol->iqinv, lusol->mnsw, &status);
3124eb8e494SKris Buschelman 
313*2205254eSKarl Rupp     switch (status) {
3144eb8e494SKris Buschelman     case 0:         /* factored */
3154eb8e494SKris Buschelman       break;
3164eb8e494SKris Buschelman 
3174eb8e494SKris Buschelman     case 7:         /* insufficient memory */
3184eb8e494SKris Buschelman       break;
3194eb8e494SKris Buschelman 
3204eb8e494SKris Buschelman     case 1:
3214eb8e494SKris Buschelman     case -1:        /* singular */
322e32f2f54SBarry Smith       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Singular matrix");
3234eb8e494SKris Buschelman 
3244eb8e494SKris Buschelman     case 3:
3254eb8e494SKris Buschelman     case 4:         /* error conditions */
326e32f2f54SBarry Smith       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"matrix error");
3274eb8e494SKris Buschelman 
3284eb8e494SKris Buschelman     default:        /* unknown condition */
329e32f2f54SBarry Smith       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"matrix unknown return code");
3304eb8e494SKris Buschelman     }
3314eb8e494SKris Buschelman 
3324eb8e494SKris Buschelman     factorizations++;
3334eb8e494SKris Buschelman   } while (status == 7);
334719d5645SBarry Smith   F->ops->solve   = MatSolve_LUSOL;
335719d5645SBarry Smith   F->assembled    = PETSC_TRUE;
336719d5645SBarry Smith   F->preallocated = PETSC_TRUE;
3374eb8e494SKris Buschelman   PetscFunctionReturn(0);
3384eb8e494SKris Buschelman }
3394eb8e494SKris Buschelman 
3404eb8e494SKris Buschelman #undef __FUNCT__
341f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_LUSOL"
34235bd34faSBarry Smith PetscErrorCode MatLUFactorSymbolic_LUSOL(Mat F,Mat A, IS r, IS c,const MatFactorInfo *info)
343b24902e0SBarry Smith {
3444eb8e494SKris Buschelman   /************************************************************************/
3454eb8e494SKris Buschelman   /* Input                                                                */
3464eb8e494SKris Buschelman   /*     A  - matrix to factor                                            */
3474eb8e494SKris Buschelman   /*     r  - row permutation (ignored)                                   */
3484eb8e494SKris Buschelman   /*     c  - column permutation (ignored)                                */
3494eb8e494SKris Buschelman   /*                                                                      */
3504eb8e494SKris Buschelman   /* Output                                                               */
3514eb8e494SKris Buschelman   /*     F  - matrix storing the factorization;                           */
3524eb8e494SKris Buschelman   /************************************************************************/
353f0c56d0fSKris Buschelman   Mat_LUSOL      *lusol;
354dfbe8321SBarry Smith   PetscErrorCode ierr;
355dfbe8321SBarry Smith   int            i, m, n, nz, nnz;
3564eb8e494SKris Buschelman 
3574eb8e494SKris Buschelman   PetscFunctionBegin;
3584eb8e494SKris Buschelman   /************************************************************************/
3594eb8e494SKris Buschelman   /* Check the arguments.                                                 */
3604eb8e494SKris Buschelman   /************************************************************************/
3614eb8e494SKris Buschelman 
3624eb8e494SKris Buschelman   ierr = MatGetSize(A, &m, &n);CHKERRQ(ierr);
3634eb8e494SKris Buschelman   nz   = ((Mat_SeqAIJ*)A->data)->nz;
3644eb8e494SKris Buschelman 
3654eb8e494SKris Buschelman   /************************************************************************/
3664eb8e494SKris Buschelman   /* Create the factorization.                                            */
3674eb8e494SKris Buschelman   /************************************************************************/
3684eb8e494SKris Buschelman 
36935bd34faSBarry Smith   F->ops->lufactornumeric = MatLUFactorNumeric_LUSOL;
37035bd34faSBarry Smith   lusol                   = (Mat_LUSOL*)(F->spptr);
3714eb8e494SKris Buschelman 
3724eb8e494SKris Buschelman   /************************************************************************/
3734eb8e494SKris Buschelman   /* Initialize parameters                                                */
3744eb8e494SKris Buschelman   /************************************************************************/
3754eb8e494SKris Buschelman 
376*2205254eSKarl Rupp   for (i = 0; i < 30; i++) {
3774eb8e494SKris Buschelman     lusol->luparm[i] = 0;
3784eb8e494SKris Buschelman     lusol->parmlu[i] = 0;
3794eb8e494SKris Buschelman   }
3804eb8e494SKris Buschelman 
3814eb8e494SKris Buschelman   lusol->luparm[1] = -1;
3824eb8e494SKris Buschelman   lusol->luparm[2] = 5;
3834eb8e494SKris Buschelman   lusol->luparm[7] = 1;
3844eb8e494SKris Buschelman 
3854eb8e494SKris Buschelman   lusol->parmlu[0] = 1 / Factorization_Tolerance;
3864eb8e494SKris Buschelman   lusol->parmlu[1] = 1 / Factorization_Tolerance;
3874eb8e494SKris Buschelman   lusol->parmlu[2] = Factorization_Small_Tolerance;
3884eb8e494SKris Buschelman   lusol->parmlu[3] = Factorization_Pivot_Tolerance;
3894eb8e494SKris Buschelman   lusol->parmlu[4] = Factorization_Pivot_Tolerance;
3904eb8e494SKris Buschelman   lusol->parmlu[5] = 3.0;
3914eb8e494SKris Buschelman   lusol->parmlu[6] = 0.3;
3924eb8e494SKris Buschelman   lusol->parmlu[7] = 0.6;
3934eb8e494SKris Buschelman 
3944eb8e494SKris Buschelman   /************************************************************************/
3954eb8e494SKris Buschelman   /* Allocate the workspace needed by LUSOL.                              */
3964eb8e494SKris Buschelman   /************************************************************************/
3974eb8e494SKris Buschelman 
3984eb8e494SKris Buschelman   lusol->elbowroom = PetscMax(lusol->elbowroom, info->fill);
3994eb8e494SKris Buschelman   nnz              = PetscMax((int)(lusol->elbowroom*nz), 5*n);
4004eb8e494SKris Buschelman 
4014eb8e494SKris Buschelman   lusol->n      = n;
4024eb8e494SKris Buschelman   lusol->nz     = nz;
4034eb8e494SKris Buschelman   lusol->nnz    = nnz;
4044eb8e494SKris Buschelman   lusol->luroom = 1.75;
4054eb8e494SKris Buschelman 
4064eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->ip);
4074eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->iq);
4084eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->lenc);
4094eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->lenr);
4104eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->locc);
4114eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->locr);
4124eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->iploc);
4134eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->iqloc);
4144eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->ipinv);
4154eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->iqinv);
4164eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(double)*n,&lusol->mnsw);
4174eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(double)*n,&lusol->mnsv);
4184eb8e494SKris Buschelman 
41923bdbc58SBarry Smith   ierr = PetscMalloc3(nnz,double,&lusol->data,nnz,PetscInt,&lusol->indc,nnz,PetscInt,&lusol->indr);CHKERRQ(ierr);
420*2205254eSKarl Rupp 
4214eb8e494SKris Buschelman   lusol->CleanUpLUSOL     = PETSC_TRUE;
42235bd34faSBarry Smith   F->ops->lufactornumeric = MatLUFactorNumeric_LUSOL;
4234eb8e494SKris Buschelman   PetscFunctionReturn(0);
4244eb8e494SKris Buschelman }
4254eb8e494SKris Buschelman 
42635bd34faSBarry Smith EXTERN_C_BEGIN
42735bd34faSBarry Smith #undef __FUNCT__
42835bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_lusol"
42935bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_seqaij_lusol(Mat A,const MatSolverPackage *type)
43035bd34faSBarry Smith {
43135bd34faSBarry Smith   PetscFunctionBegin;
4322692d6eeSBarry Smith   *type = MATSOLVERLUSOL;
43335bd34faSBarry Smith   PetscFunctionReturn(0);
43435bd34faSBarry Smith }
43535bd34faSBarry Smith EXTERN_C_END
43635bd34faSBarry Smith 
4374eb8e494SKris Buschelman #undef __FUNCT__
438b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_lusol"
4395c9eb25fSBarry Smith PetscErrorCode MatGetFactor_seqaij_lusol(Mat A,MatFactorType ftype,Mat *F)
440521d7252SBarry Smith {
441b24902e0SBarry Smith   Mat            B;
442f0c56d0fSKris Buschelman   Mat_LUSOL      *lusol;
443b24902e0SBarry Smith   PetscErrorCode ierr;
44435bd34faSBarry Smith   int            m, n;
4454eb8e494SKris Buschelman 
4464eb8e494SKris Buschelman   PetscFunctionBegin;
4474eb8e494SKris Buschelman   ierr = MatGetSize(A, &m, &n);CHKERRQ(ierr);
448b24902e0SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
449b24902e0SBarry Smith   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
450b24902e0SBarry Smith   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
451b24902e0SBarry Smith   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
4524eb8e494SKris Buschelman 
45338f2d2fdSLisandro Dalcin   ierr     = PetscNewLog(B,Mat_LUSOL,&lusol);CHKERRQ(ierr);
454b24902e0SBarry Smith   B->spptr = lusol;
4552f71e704SKris Buschelman 
456f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic = MatLUFactorSymbolic_LUSOL;
457f0c56d0fSKris Buschelman   B->ops->destroy          = MatDestroy_LUSOL;
458*2205254eSKarl Rupp 
45935bd34faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_lusol",MatFactorGetSolverPackage_seqaij_lusol);CHKERRQ(ierr);
460*2205254eSKarl Rupp 
461d5f3da31SBarry Smith   B->factortype = MAT_FACTOR_LU;
462f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
463f0c56d0fSKris Buschelman }
464f0c56d0fSKris Buschelman 
4652f71e704SKris Buschelman /*MC
4662692d6eeSBarry Smith   MATSOLVERLUSOL - "lusol" - Provides direct solvers (LU) for sequential matrices
4672f71e704SKris Buschelman                          via the external package LUSOL.
4682f71e704SKris Buschelman 
4692f71e704SKris Buschelman   If LUSOL is installed (see the manual for
4702f71e704SKris Buschelman   instructions on how to declare the existence of external packages),
4712f71e704SKris Buschelman 
47241c8de11SBarry Smith   Works with MATSEQAIJ matrices
4732f71e704SKris Buschelman 
4742f71e704SKris Buschelman    Level: beginner
4752f71e704SKris Buschelman 
47641c8de11SBarry Smith .seealso: PCLU, PCFactorSetMatSolverPackage(), MatSolverPackage
47741c8de11SBarry Smith 
4782f71e704SKris Buschelman M*/
479