xref: /petsc/src/mat/impls/aij/seq/lusol/lusol.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
14eb8e494SKris Buschelman /*
24eb8e494SKris Buschelman         Provides an interface to the LUSOL package of ....
34eb8e494SKris Buschelman 
44eb8e494SKris Buschelman */
54eb8e494SKris Buschelman #include "src/mat/impls/aij/seq/aij.h"
64eb8e494SKris Buschelman 
74eb8e494SKris Buschelman #if defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
84eb8e494SKris Buschelman #define LU1FAC   lu1fac_
94eb8e494SKris Buschelman #define LU6SOL   lu6sol_
104eb8e494SKris Buschelman #define M1PAGE   m1page_
114eb8e494SKris Buschelman #define M5SETX   m5setx_
124eb8e494SKris Buschelman #define M6RDEL   m6rdel_
134eb8e494SKris Buschelman #elif !defined(PETSC_HAVE_FORTRAN_CAPS)
144eb8e494SKris Buschelman #define LU1FAC   lu1fac
154eb8e494SKris Buschelman #define LU6SOL   lu6sol
164eb8e494SKris Buschelman #define M1PAGE   m1page
174eb8e494SKris Buschelman #define M5SETX   m5setx
184eb8e494SKris Buschelman #define M6RDEL   m6rdel
194eb8e494SKris Buschelman #endif
204eb8e494SKris Buschelman 
214eb8e494SKris Buschelman EXTERN_C_BEGIN
224eb8e494SKris Buschelman /*
234eb8e494SKris Buschelman     Dummy symbols that the MINOS files mi25bfac.f and mi15blas.f may require
244eb8e494SKris Buschelman */
254eb8e494SKris Buschelman void PETSC_STDCALL M1PAGE() {
264eb8e494SKris Buschelman   ;
274eb8e494SKris Buschelman }
284eb8e494SKris Buschelman void PETSC_STDCALL M5SETX() {
294eb8e494SKris Buschelman   ;
304eb8e494SKris Buschelman }
314eb8e494SKris Buschelman 
324eb8e494SKris Buschelman void PETSC_STDCALL M6RDEL() {
334eb8e494SKris Buschelman   ;
344eb8e494SKris Buschelman }
354eb8e494SKris Buschelman 
364eb8e494SKris Buschelman extern void PETSC_STDCALL LU1FAC (int *m, int *n, int *nnz, int *size, int *luparm,
374eb8e494SKris Buschelman                         double *parmlu, double *data, int *indc, int *indr,
384eb8e494SKris Buschelman                         int *rowperm, int *colperm, int *collen, int *rowlen,
394eb8e494SKris Buschelman                         int *colstart, int *rowstart, int *rploc, int *cploc,
404eb8e494SKris Buschelman                         int *rpinv, int *cpinv, double *w, int *inform);
414eb8e494SKris Buschelman 
424eb8e494SKris Buschelman extern void PETSC_STDCALL LU6SOL (int *mode, int *m, int *n, double *rhs, double *x,
434eb8e494SKris Buschelman                         int *size, int *luparm, double *parmlu, double *data,
444eb8e494SKris Buschelman                         int *indc, int *indr, int *rowperm, int *colperm,
454eb8e494SKris Buschelman                         int *collen, int *rowlen, int *colstart, int *rowstart,
464eb8e494SKris Buschelman                         int *inform);
472f71e704SKris Buschelman EXTERN_C_END
484eb8e494SKris Buschelman 
49dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_LUSOL(Mat,MatDuplicateOption,Mat*);
50f0c56d0fSKris Buschelman 
51f0c56d0fSKris Buschelman typedef struct  {
524eb8e494SKris Buschelman   double *data;
534eb8e494SKris Buschelman   int *indc;
544eb8e494SKris Buschelman   int *indr;
554eb8e494SKris Buschelman 
564eb8e494SKris Buschelman   int *ip;
574eb8e494SKris Buschelman   int *iq;
584eb8e494SKris Buschelman   int *lenc;
594eb8e494SKris Buschelman   int *lenr;
604eb8e494SKris Buschelman   int *locc;
614eb8e494SKris Buschelman   int *locr;
624eb8e494SKris Buschelman   int *iploc;
634eb8e494SKris Buschelman   int *iqloc;
644eb8e494SKris Buschelman   int *ipinv;
654eb8e494SKris Buschelman   int *iqinv;
664eb8e494SKris Buschelman   double *mnsw;
674eb8e494SKris Buschelman   double *mnsv;
684eb8e494SKris Buschelman 
694eb8e494SKris Buschelman   double elbowroom;
704eb8e494SKris Buschelman   double luroom;		/* Extra space allocated when factor fails   */
714eb8e494SKris Buschelman   double parmlu[30];		/* Input/output to LUSOL                     */
724eb8e494SKris Buschelman 
734eb8e494SKris Buschelman   int n;			/* Number of rows/columns in matrix          */
744eb8e494SKris Buschelman   int nz;			/* Number of nonzeros                        */
754eb8e494SKris Buschelman   int nnz;			/* Number of nonzeros allocated for factors  */
764eb8e494SKris Buschelman   int luparm[30];		/* Input/output to LUSOL                     */
774eb8e494SKris Buschelman 
78*6849ba73SBarry Smith   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
79*6849ba73SBarry Smith   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
80*6849ba73SBarry Smith   PetscErrorCode (*MatDestroy)(Mat);
814eb8e494SKris Buschelman   PetscTruth CleanUpLUSOL;
824eb8e494SKris Buschelman 
83f0c56d0fSKris Buschelman } Mat_LUSOL;
844eb8e494SKris Buschelman 
854eb8e494SKris Buschelman /*  LUSOL input/Output Parameters (Description uses C-style indexes
864eb8e494SKris Buschelman  *
874eb8e494SKris Buschelman  *  Input parameters                                        Typical value
884eb8e494SKris Buschelman  *
894eb8e494SKris Buschelman  *  luparm(0) = nout     File number for printed messages.         6
904eb8e494SKris Buschelman  *  luparm(1) = lprint   Print level.                              0
914eb8e494SKris Buschelman  *                    < 0 suppresses output.
924eb8e494SKris Buschelman  *                    = 0 gives error messages.
934eb8e494SKris Buschelman  *                    = 1 gives debug output from some of the
944eb8e494SKris Buschelman  *                        other routines in LUSOL.
954eb8e494SKris Buschelman  *                   >= 2 gives the pivot row and column and the
964eb8e494SKris Buschelman  *                        no. of rows and columns involved at
974eb8e494SKris Buschelman  *                        each elimination step in lu1fac.
984eb8e494SKris Buschelman  *  luparm(2) = maxcol   lu1fac: maximum number of columns         5
994eb8e494SKris Buschelman  *                        searched allowed in a Markowitz-type
1004eb8e494SKris Buschelman  *                        search for the next pivot element.
1014eb8e494SKris Buschelman  *                        For some of the factorization, the
1024eb8e494SKris Buschelman  *                        number of rows searched is
1034eb8e494SKris Buschelman  *                        maxrow = maxcol - 1.
1044eb8e494SKris Buschelman  *
1054eb8e494SKris Buschelman  *
1064eb8e494SKris Buschelman  *  Output parameters
1074eb8e494SKris Buschelman  *
1084eb8e494SKris Buschelman  *  luparm(9) = inform   Return code from last call to any LU routine.
1094eb8e494SKris Buschelman  *  luparm(10) = nsing    No. of singularities marked in the
1104eb8e494SKris Buschelman  *                        output array w(*).
1114eb8e494SKris Buschelman  *  luparm(11) = jsing    Column index of last singularity.
1124eb8e494SKris Buschelman  *  luparm(12) = minlen   Minimum recommended value for  lena.
1134eb8e494SKris Buschelman  *  luparm(13) = maxlen   ?
1144eb8e494SKris Buschelman  *  luparm(14) = nupdat   No. of updates performed by the lu8 routines.
1154eb8e494SKris Buschelman  *  luparm(15) = nrank    No. of nonempty rows of U.
1164eb8e494SKris Buschelman  *  luparm(16) = ndens1   No. of columns remaining when the density of
1174eb8e494SKris Buschelman  *                        the matrix being factorized reached dens1.
1184eb8e494SKris Buschelman  *  luparm(17) = ndens2   No. of columns remaining when the density of
1194eb8e494SKris Buschelman  *                        the matrix being factorized reached dens2.
1204eb8e494SKris Buschelman  *  luparm(18) = jumin    The column index associated with dumin.
1214eb8e494SKris Buschelman  *  luparm(19) = numl0    No. of columns in initial  L.
1224eb8e494SKris Buschelman  *  luparm(20) = lenl0    Size of initial  L  (no. of nonzeros).
1234eb8e494SKris Buschelman  *  luparm(21) = lenu0    Size of initial  U.
1244eb8e494SKris Buschelman  *  luparm(22) = lenl     Size of current  L.
1254eb8e494SKris Buschelman  *  luparm(23) = lenu     Size of current  U.
1264eb8e494SKris Buschelman  *  luparm(24) = lrow     Length of row file.
1274eb8e494SKris Buschelman  *  luparm(25) = ncp      No. of compressions of LU data structures.
1284eb8e494SKris Buschelman  *  luparm(26) = mersum   lu1fac: sum of Markowitz merit counts.
1294eb8e494SKris Buschelman  *  luparm(27) = nutri    lu1fac: triangular rows in U.
1304eb8e494SKris Buschelman  *  luparm(28) = nltri    lu1fac: triangular rows in L.
1314eb8e494SKris Buschelman  *  luparm(29) =
1324eb8e494SKris Buschelman  *
1334eb8e494SKris Buschelman  *
1344eb8e494SKris Buschelman  *  Input parameters                                        Typical value
1354eb8e494SKris Buschelman  *
1364eb8e494SKris Buschelman  *  parmlu(0) = elmax1   Max multiplier allowed in  L           10.0
1374eb8e494SKris Buschelman  *                        during factor.
1384eb8e494SKris Buschelman  *  parmlu(1) = elmax2   Max multiplier allowed in  L           10.0
1394eb8e494SKris Buschelman  *                        during updates.
1404eb8e494SKris Buschelman  *  parmlu(2) = small    Absolute tolerance for             eps**0.8
1414eb8e494SKris Buschelman  *                        treating reals as zero.     IBM double: 3.0d-13
1424eb8e494SKris Buschelman  *  parmlu(3) = utol1    Absolute tol for flagging          eps**0.66667
1434eb8e494SKris Buschelman  *                        small diagonals of U.       IBM double: 3.7d-11
1444eb8e494SKris Buschelman  *  parmlu(4) = utol2    Relative tol for flagging          eps**0.66667
1454eb8e494SKris Buschelman  *                        small diagonals of U.       IBM double: 3.7d-11
1464eb8e494SKris Buschelman  *  parmlu(5) = uspace   Factor limiting waste space in  U.      3.0
1474eb8e494SKris Buschelman  *                        In lu1fac, the row or column lists
1484eb8e494SKris Buschelman  *                        are compressed if their length
1494eb8e494SKris Buschelman  *                        exceeds uspace times the length of
1504eb8e494SKris Buschelman  *                        either file after the last compression.
1514eb8e494SKris Buschelman  *  parmlu(6) = dens1    The density at which the Markowitz      0.3
1524eb8e494SKris Buschelman  *                        strategy should search maxcol columns
1534eb8e494SKris Buschelman  *                        and no rows.
1544eb8e494SKris Buschelman  *  parmlu(7) = dens2    the density at which the Markowitz      0.6
1554eb8e494SKris Buschelman  *                        strategy should search only 1 column
1564eb8e494SKris Buschelman  *                        or (preferably) use a dense LU for
1574eb8e494SKris Buschelman  *                        all the remaining rows and columns.
1584eb8e494SKris Buschelman  *
1594eb8e494SKris Buschelman  *
1604eb8e494SKris Buschelman  *  Output parameters
1614eb8e494SKris Buschelman  *
1624eb8e494SKris Buschelman  *  parmlu(9) = amax     Maximum element in  A.
1634eb8e494SKris Buschelman  *  parmlu(10) = elmax    Maximum multiplier in current  L.
1644eb8e494SKris Buschelman  *  parmlu(11) = umax     Maximum element in current  U.
1654eb8e494SKris Buschelman  *  parmlu(12) = dumax    Maximum diagonal in  U.
1664eb8e494SKris Buschelman  *  parmlu(13) = dumin    Minimum diagonal in  U.
1674eb8e494SKris Buschelman  *  parmlu(14) =
1684eb8e494SKris Buschelman  *  parmlu(15) =
1694eb8e494SKris Buschelman  *  parmlu(16) =
1704eb8e494SKris Buschelman  *  parmlu(17) =
1714eb8e494SKris Buschelman  *  parmlu(18) =
1724eb8e494SKris Buschelman  *  parmlu(19) = resid    lu6sol: residual after solve with U or U'.
1734eb8e494SKris Buschelman  *  ...
1744eb8e494SKris Buschelman  *  parmlu(29) =
1754eb8e494SKris Buschelman  */
1764eb8e494SKris Buschelman 
1774eb8e494SKris Buschelman #define Factorization_Tolerance       1e-1
1784eb8e494SKris Buschelman #define Factorization_Pivot_Tolerance pow(2.2204460492503131E-16, 2.0 / 3.0)
1794eb8e494SKris Buschelman #define Factorization_Small_Tolerance 1e-15 /* pow(DBL_EPSILON, 0.8) */
1804eb8e494SKris Buschelman 
1812f71e704SKris Buschelman EXTERN_C_BEGIN
1822f71e704SKris Buschelman #undef __FUNCT__
1832f71e704SKris Buschelman #define __FUNCT__ "MatConvert_LUSOL_SeqAIJ"
184dfbe8321SBarry Smith PetscErrorCode MatConvert_LUSOL_SeqAIJ(Mat A,const MatType type,Mat *newmat) {
1852f71e704SKris Buschelman   /* This routine is only called to convert an unfactored PETSc-LUSOL matrix */
1862f71e704SKris Buschelman   /* to its base PETSc type, so we will ignore 'MatType type'. */
187dfbe8321SBarry Smith   PetscErrorCode ierr;
1882f71e704SKris Buschelman   Mat       B=*newmat;
189f0c56d0fSKris Buschelman   Mat_LUSOL *lusol=(Mat_LUSOL *)A->spptr;
1902f71e704SKris Buschelman 
1912f71e704SKris Buschelman   PetscFunctionBegin;
1922f71e704SKris Buschelman   if (B != A) {
1932f71e704SKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
194f0c56d0fSKris Buschelman   }
195f0c56d0fSKris Buschelman   B->ops->duplicate        = lusol->MatDuplicate;
1962f71e704SKris Buschelman   B->ops->lufactorsymbolic = lusol->MatLUFactorSymbolic;
1972f71e704SKris Buschelman   B->ops->destroy          = lusol->MatDestroy;
1982f71e704SKris Buschelman 
1992f71e704SKris Buschelman   ierr = PetscFree(lusol);CHKERRQ(ierr);
2002f71e704SKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
2012f71e704SKris Buschelman   *newmat = B;
2022f71e704SKris Buschelman   PetscFunctionReturn(0);
2032f71e704SKris Buschelman }
2042f71e704SKris Buschelman EXTERN_C_END
2054eb8e494SKris Buschelman 
2064eb8e494SKris Buschelman #undef __FUNCT__
207f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_LUSOL"
208dfbe8321SBarry Smith PetscErrorCode MatDestroy_LUSOL(Mat A)
209dfbe8321SBarry Smith {
210dfbe8321SBarry Smith   PetscErrorCode ierr;
211f0c56d0fSKris Buschelman   Mat_LUSOL *lusol=(Mat_LUSOL *)A->spptr;
2124eb8e494SKris Buschelman 
2134eb8e494SKris Buschelman   PetscFunctionBegin;
2144eb8e494SKris Buschelman   if (lusol->CleanUpLUSOL) {
2154eb8e494SKris Buschelman     ierr = PetscFree(lusol->ip);CHKERRQ(ierr);
2164eb8e494SKris Buschelman     ierr = PetscFree(lusol->iq);CHKERRQ(ierr);
2174eb8e494SKris Buschelman     ierr = PetscFree(lusol->lenc);CHKERRQ(ierr);
2184eb8e494SKris Buschelman     ierr = PetscFree(lusol->lenr);CHKERRQ(ierr);
2194eb8e494SKris Buschelman     ierr = PetscFree(lusol->locc);CHKERRQ(ierr);
2204eb8e494SKris Buschelman     ierr = PetscFree(lusol->locr);CHKERRQ(ierr);
2214eb8e494SKris Buschelman     ierr = PetscFree(lusol->iploc);CHKERRQ(ierr);
2224eb8e494SKris Buschelman     ierr = PetscFree(lusol->iqloc);CHKERRQ(ierr);
2234eb8e494SKris Buschelman     ierr = PetscFree(lusol->ipinv);CHKERRQ(ierr);
2244eb8e494SKris Buschelman     ierr = PetscFree(lusol->iqinv);CHKERRQ(ierr);
2254eb8e494SKris Buschelman     ierr = PetscFree(lusol->mnsw);CHKERRQ(ierr);
2264eb8e494SKris Buschelman     ierr = PetscFree(lusol->mnsv);CHKERRQ(ierr);
2274eb8e494SKris Buschelman 
2284eb8e494SKris Buschelman     ierr = PetscFree(lusol->indc);CHKERRQ(ierr);
2294eb8e494SKris Buschelman   }
2304eb8e494SKris Buschelman 
2312f71e704SKris Buschelman   ierr = MatConvert_LUSOL_SeqAIJ(A,MATSEQAIJ,&A);
2322f71e704SKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
2334eb8e494SKris Buschelman   PetscFunctionReturn(0);
2344eb8e494SKris Buschelman }
2354eb8e494SKris Buschelman 
2364eb8e494SKris Buschelman #undef __FUNCT__
237f0c56d0fSKris Buschelman #define __FUNCT__  "MatSolve_LUSOL"
238*6849ba73SBarry Smith PetscErrorCode MatSolve_LUSOL(Mat A,Vec b,Vec x)
239*6849ba73SBarry Smith {
240f0c56d0fSKris Buschelman   Mat_LUSOL *lusol=(Mat_LUSOL*)A->spptr;
2414eb8e494SKris Buschelman   double    *bb,*xx;
2424eb8e494SKris Buschelman   int       mode=5;
243*6849ba73SBarry Smith   PetscErrorCode ierr;
244*6849ba73SBarry Smith   int       i,m,n,nnz,status;
2454eb8e494SKris Buschelman 
2464eb8e494SKris Buschelman   PetscFunctionBegin;
2474eb8e494SKris Buschelman   ierr = VecGetArray(x, &xx);CHKERRQ(ierr);
2484eb8e494SKris Buschelman   ierr = VecGetArray(b, &bb);CHKERRQ(ierr);
2494eb8e494SKris Buschelman 
2504eb8e494SKris Buschelman   m = n = lusol->n;
2514eb8e494SKris Buschelman   nnz = lusol->nnz;
2524eb8e494SKris Buschelman 
2534eb8e494SKris Buschelman   for (i = 0; i < m; i++)
2544eb8e494SKris Buschelman     {
2554eb8e494SKris Buschelman       lusol->mnsv[i] = bb[i];
2564eb8e494SKris Buschelman     }
2574eb8e494SKris Buschelman 
2584eb8e494SKris Buschelman   LU6SOL(&mode, &m, &n, lusol->mnsv, xx, &nnz,
2594eb8e494SKris Buschelman          lusol->luparm, lusol->parmlu, lusol->data,
2604eb8e494SKris Buschelman          lusol->indc, lusol->indr, lusol->ip, lusol->iq,
2614eb8e494SKris Buschelman          lusol->lenc, lusol->lenr, lusol->locc, lusol->locr, &status);
2624eb8e494SKris Buschelman 
2634eb8e494SKris Buschelman   if (status != 0)
2644eb8e494SKris Buschelman     {
2654eb8e494SKris Buschelman       SETERRQ(PETSC_ERR_ARG_SIZ,"solve failed");
2664eb8e494SKris Buschelman     }
2674eb8e494SKris Buschelman 
2684eb8e494SKris Buschelman   ierr = VecRestoreArray(x, &xx);CHKERRQ(ierr);
2694eb8e494SKris Buschelman   ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr);
2704eb8e494SKris Buschelman   PetscFunctionReturn(0);
2714eb8e494SKris Buschelman }
2724eb8e494SKris Buschelman 
2734eb8e494SKris Buschelman #undef __FUNCT__
274f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_LUSOL"
275*6849ba73SBarry Smith PetscErrorCode MatLUFactorNumeric_LUSOL(Mat A, Mat *F)
276*6849ba73SBarry Smith {
2774eb8e494SKris Buschelman   Mat_SeqAIJ *a;
278f0c56d0fSKris Buschelman   Mat_LUSOL  *lusol = (Mat_LUSOL*)(*F)->spptr;
279*6849ba73SBarry Smith   PetscErrorCode ierr;
2804eb8e494SKris Buschelman   int        m, n, nz, nnz, status;
281*6849ba73SBarry Smith   int        i, rs, re;
2824eb8e494SKris Buschelman   int        factorizations;
2834eb8e494SKris Buschelman 
2844eb8e494SKris Buschelman   PetscFunctionBegin;
2854eb8e494SKris Buschelman   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);CHKERRQ(ierr);
2864eb8e494SKris Buschelman   a = (Mat_SeqAIJ *)A->data;
2874eb8e494SKris Buschelman 
2884eb8e494SKris Buschelman   if (m != lusol->n) {
2894eb8e494SKris Buschelman     SETERRQ(PETSC_ERR_ARG_SIZ,"factorization struct inconsistent");
2904eb8e494SKris Buschelman   }
2914eb8e494SKris Buschelman 
2924eb8e494SKris Buschelman   factorizations = 0;
2934eb8e494SKris Buschelman   do
2944eb8e494SKris Buschelman     {
2954eb8e494SKris Buschelman       /*******************************************************************/
2964eb8e494SKris Buschelman       /* Check the workspace allocation.                                 */
2974eb8e494SKris Buschelman       /*******************************************************************/
2984eb8e494SKris Buschelman 
2994eb8e494SKris Buschelman       nz = a->nz;
3004eb8e494SKris Buschelman       nnz = PetscMax(lusol->nnz, (int)(lusol->elbowroom*nz));
3014eb8e494SKris Buschelman       nnz = PetscMax(nnz, 5*n);
3024eb8e494SKris Buschelman 
3034eb8e494SKris Buschelman       if (nnz < lusol->luparm[12]){
3044eb8e494SKris Buschelman         nnz = (int)(lusol->luroom * lusol->luparm[12]);
3054eb8e494SKris Buschelman       } else if ((factorizations > 0) && (lusol->luroom < 6)){
3064eb8e494SKris Buschelman         lusol->luroom += 0.1;
3074eb8e494SKris Buschelman       }
3084eb8e494SKris Buschelman 
3094eb8e494SKris Buschelman       nnz = PetscMax(nnz, (int)(lusol->luroom*(lusol->luparm[22] + lusol->luparm[23])));
3104eb8e494SKris Buschelman 
3114eb8e494SKris Buschelman       if (nnz > lusol->nnz){
3124eb8e494SKris Buschelman         ierr = PetscFree(lusol->indc);CHKERRQ(ierr);
3134eb8e494SKris Buschelman         ierr        = PetscMalloc((sizeof(double)+2*sizeof(int))*nnz,&lusol->indc);CHKERRQ(ierr);
3144eb8e494SKris Buschelman         lusol->indr = lusol->indc + nnz;
3154eb8e494SKris Buschelman         lusol->data = (double *)(lusol->indr + nnz);
3164eb8e494SKris Buschelman         lusol->nnz  = nnz;
3174eb8e494SKris Buschelman       }
3184eb8e494SKris Buschelman 
3194eb8e494SKris Buschelman       /*******************************************************************/
3204eb8e494SKris Buschelman       /* Fill in the data for the problem.      (1-based Fortran style)  */
3214eb8e494SKris Buschelman       /*******************************************************************/
3224eb8e494SKris Buschelman 
3234eb8e494SKris Buschelman       nz = 0;
3244eb8e494SKris Buschelman       for (i = 0; i < n; i++)
3254eb8e494SKris Buschelman         {
3264eb8e494SKris Buschelman           rs = a->i[i];
3274eb8e494SKris Buschelman           re = a->i[i+1];
3284eb8e494SKris Buschelman 
3294eb8e494SKris Buschelman           while (rs < re)
3304eb8e494SKris Buschelman             {
3314eb8e494SKris Buschelman               if (a->a[rs] != 0.0)
3324eb8e494SKris Buschelman                 {
3334eb8e494SKris Buschelman                   lusol->indc[nz] = i + 1;
3344eb8e494SKris Buschelman                   lusol->indr[nz] = a->j[rs] + 1;
3354eb8e494SKris Buschelman                   lusol->data[nz] = a->a[rs];
3364eb8e494SKris Buschelman                   nz++;
3374eb8e494SKris Buschelman                 }
3384eb8e494SKris Buschelman               rs++;
3394eb8e494SKris Buschelman             }
3404eb8e494SKris Buschelman         }
3414eb8e494SKris Buschelman 
3424eb8e494SKris Buschelman       /*******************************************************************/
3434eb8e494SKris Buschelman       /* Do the factorization.                                           */
3444eb8e494SKris Buschelman       /*******************************************************************/
3454eb8e494SKris Buschelman 
3464eb8e494SKris Buschelman       LU1FAC(&m, &n, &nz, &nnz,
3474eb8e494SKris Buschelman              lusol->luparm, lusol->parmlu, lusol->data,
3484eb8e494SKris Buschelman              lusol->indc, lusol->indr, lusol->ip, lusol->iq,
3494eb8e494SKris Buschelman              lusol->lenc, lusol->lenr, lusol->locc, lusol->locr,
3504eb8e494SKris Buschelman              lusol->iploc, lusol->iqloc, lusol->ipinv,
3514eb8e494SKris Buschelman              lusol->iqinv, lusol->mnsw, &status);
3524eb8e494SKris Buschelman 
3534eb8e494SKris Buschelman       switch(status)
3544eb8e494SKris Buschelman         {
3554eb8e494SKris Buschelman         case 0:		/* factored */
3564eb8e494SKris Buschelman           break;
3574eb8e494SKris Buschelman 
3584eb8e494SKris Buschelman         case 7:		/* insufficient memory */
3594eb8e494SKris Buschelman           break;
3604eb8e494SKris Buschelman 
3614eb8e494SKris Buschelman         case 1:
3624eb8e494SKris Buschelman         case -1:		/* singular */
3634eb8e494SKris Buschelman           SETERRQ(1,"Singular matrix");
3644eb8e494SKris Buschelman 
3654eb8e494SKris Buschelman         case 3:
3664eb8e494SKris Buschelman         case 4:		/* error conditions */
3674eb8e494SKris Buschelman           SETERRQ(1,"matrix error");
3684eb8e494SKris Buschelman 
3694eb8e494SKris Buschelman         default:		/* unknown condition */
3704eb8e494SKris Buschelman           SETERRQ(1,"matrix unknown return code");
3714eb8e494SKris Buschelman         }
3724eb8e494SKris Buschelman 
3734eb8e494SKris Buschelman       factorizations++;
3744eb8e494SKris Buschelman     } while (status == 7);
375a8883a68SKris Buschelman   (*F)->assembled = PETSC_TRUE;
3764eb8e494SKris Buschelman   PetscFunctionReturn(0);
3774eb8e494SKris Buschelman }
3784eb8e494SKris Buschelman 
3794eb8e494SKris Buschelman #undef __FUNCT__
380f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_LUSOL"
381dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_LUSOL(Mat A, IS r, IS c,MatFactorInfo *info, Mat *F) {
3824eb8e494SKris Buschelman   /************************************************************************/
3834eb8e494SKris Buschelman   /* Input                                                                */
3844eb8e494SKris Buschelman   /*     A  - matrix to factor                                            */
3854eb8e494SKris Buschelman   /*     r  - row permutation (ignored)                                   */
3864eb8e494SKris Buschelman   /*     c  - column permutation (ignored)                                */
3874eb8e494SKris Buschelman   /*                                                                      */
3884eb8e494SKris Buschelman   /* Output                                                               */
3894eb8e494SKris Buschelman   /*     F  - matrix storing the factorization;                           */
3904eb8e494SKris Buschelman   /************************************************************************/
3914eb8e494SKris Buschelman   Mat       B;
392f0c56d0fSKris Buschelman   Mat_LUSOL *lusol;
393dfbe8321SBarry Smith   PetscErrorCode ierr;
394dfbe8321SBarry Smith   int        i, m, n, nz, nnz;
3954eb8e494SKris Buschelman 
3964eb8e494SKris Buschelman   PetscFunctionBegin;
3974eb8e494SKris Buschelman 
3984eb8e494SKris Buschelman   /************************************************************************/
3994eb8e494SKris Buschelman   /* Check the arguments.                                                 */
4004eb8e494SKris Buschelman   /************************************************************************/
4014eb8e494SKris Buschelman 
4024eb8e494SKris Buschelman   ierr = MatGetSize(A, &m, &n);CHKERRQ(ierr);
4034eb8e494SKris Buschelman   nz = ((Mat_SeqAIJ *)A->data)->nz;
4044eb8e494SKris Buschelman 
4054eb8e494SKris Buschelman   /************************************************************************/
4064eb8e494SKris Buschelman   /* Create the factorization.                                            */
4074eb8e494SKris Buschelman   /************************************************************************/
4084eb8e494SKris Buschelman 
4094eb8e494SKris Buschelman   ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,m,n,&B);CHKERRQ(ierr);
410be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
4114eb8e494SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
4124eb8e494SKris Buschelman 
413f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatLUFactorNumeric_LUSOL;
414f0c56d0fSKris Buschelman   B->ops->solve           = MatSolve_LUSOL;
4154eb8e494SKris Buschelman   B->factor               = FACTOR_LU;
416f0c56d0fSKris Buschelman   lusol                   = (Mat_LUSOL*)(B->spptr);
4174eb8e494SKris Buschelman 
4184eb8e494SKris Buschelman   /************************************************************************/
4194eb8e494SKris Buschelman   /* Initialize parameters                                                */
4204eb8e494SKris Buschelman   /************************************************************************/
4214eb8e494SKris Buschelman 
4224eb8e494SKris Buschelman   for (i = 0; i < 30; i++)
4234eb8e494SKris Buschelman     {
4244eb8e494SKris Buschelman       lusol->luparm[i] = 0;
4254eb8e494SKris Buschelman       lusol->parmlu[i] = 0;
4264eb8e494SKris Buschelman     }
4274eb8e494SKris Buschelman 
4284eb8e494SKris Buschelman   lusol->luparm[1] = -1;
4294eb8e494SKris Buschelman   lusol->luparm[2] = 5;
4304eb8e494SKris Buschelman   lusol->luparm[7] = 1;
4314eb8e494SKris Buschelman 
4324eb8e494SKris Buschelman   lusol->parmlu[0] = 1 / Factorization_Tolerance;
4334eb8e494SKris Buschelman   lusol->parmlu[1] = 1 / Factorization_Tolerance;
4344eb8e494SKris Buschelman   lusol->parmlu[2] = Factorization_Small_Tolerance;
4354eb8e494SKris Buschelman   lusol->parmlu[3] = Factorization_Pivot_Tolerance;
4364eb8e494SKris Buschelman   lusol->parmlu[4] = Factorization_Pivot_Tolerance;
4374eb8e494SKris Buschelman   lusol->parmlu[5] = 3.0;
4384eb8e494SKris Buschelman   lusol->parmlu[6] = 0.3;
4394eb8e494SKris Buschelman   lusol->parmlu[7] = 0.6;
4404eb8e494SKris Buschelman 
4414eb8e494SKris Buschelman   /************************************************************************/
4424eb8e494SKris Buschelman   /* Allocate the workspace needed by LUSOL.                              */
4434eb8e494SKris Buschelman   /************************************************************************/
4444eb8e494SKris Buschelman 
4454eb8e494SKris Buschelman   lusol->elbowroom = PetscMax(lusol->elbowroom, info->fill);
4464eb8e494SKris Buschelman   nnz = PetscMax((int)(lusol->elbowroom*nz), 5*n);
4474eb8e494SKris Buschelman 
4484eb8e494SKris Buschelman   lusol->n = n;
4494eb8e494SKris Buschelman   lusol->nz = nz;
4504eb8e494SKris Buschelman   lusol->nnz = nnz;
4514eb8e494SKris Buschelman   lusol->luroom = 1.75;
4524eb8e494SKris Buschelman 
4534eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->ip);
4544eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->iq);
4554eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->lenc);
4564eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->lenr);
4574eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->locc);
4584eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->locr);
4594eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->iploc);
4604eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->iqloc);
4614eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->ipinv);
4624eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->iqinv);
4634eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(double)*n,&lusol->mnsw);
4644eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(double)*n,&lusol->mnsv);
4654eb8e494SKris Buschelman 
4664eb8e494SKris Buschelman   ierr        = PetscMalloc((sizeof(double)+2*sizeof(int))*nnz,&lusol->indc);
4674eb8e494SKris Buschelman   lusol->indr = lusol->indc + nnz;
4684eb8e494SKris Buschelman   lusol->data = (double *)(lusol->indr + nnz);
4694eb8e494SKris Buschelman   lusol->CleanUpLUSOL = PETSC_TRUE;
4704eb8e494SKris Buschelman   *F = B;
4714eb8e494SKris Buschelman   PetscFunctionReturn(0);
4724eb8e494SKris Buschelman }
4734eb8e494SKris Buschelman 
4742f71e704SKris Buschelman EXTERN_C_BEGIN
4754eb8e494SKris Buschelman #undef __FUNCT__
4762f71e704SKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_LUSOL"
477dfbe8321SBarry Smith PetscErrorCode MatConvert_SeqAIJ_LUSOL(Mat A,const MatType type,Mat *newmat) {
478dfbe8321SBarry Smith   PetscErrorCode ierr;
479dfbe8321SBarry Smith   int        m, n;
480f0c56d0fSKris Buschelman   Mat_LUSOL *lusol;
4812f71e704SKris Buschelman   Mat       B=*newmat;
4824eb8e494SKris Buschelman 
4834eb8e494SKris Buschelman   PetscFunctionBegin;
4844eb8e494SKris Buschelman   ierr = MatGetSize(A, &m, &n);CHKERRQ(ierr);
4854eb8e494SKris Buschelman   if (m != n) {
4864eb8e494SKris Buschelman     SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square");
4874eb8e494SKris Buschelman   }
4882f71e704SKris Buschelman   if (B != A) {
4892f71e704SKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
4902f71e704SKris Buschelman   }
4914eb8e494SKris Buschelman 
492f0c56d0fSKris Buschelman   ierr                       = PetscNew(Mat_LUSOL,&lusol);CHKERRQ(ierr);
493f0c56d0fSKris Buschelman   lusol->MatDuplicate        = A->ops->duplicate;
4942f71e704SKris Buschelman   lusol->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
4952f71e704SKris Buschelman   lusol->MatDestroy          = A->ops->destroy;
4962f71e704SKris Buschelman   lusol->CleanUpLUSOL        = PETSC_FALSE;
4972f71e704SKris Buschelman 
4982f71e704SKris Buschelman   B->spptr                   = (void*)lusol;
499f0c56d0fSKris Buschelman   B->ops->duplicate          = MatDuplicate_LUSOL;
500f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic   = MatLUFactorSymbolic_LUSOL;
501f0c56d0fSKris Buschelman   B->ops->destroy            = MatDestroy_LUSOL;
5022f71e704SKris Buschelman 
503f0c56d0fSKris Buschelman   PetscLogInfo(0,"Using LUSOL for LU factorization and solves.");
5042f71e704SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_lusol_C",
5052f71e704SKris Buschelman                                            "MatConvert_SeqAIJ_LUSOL",MatConvert_SeqAIJ_LUSOL);CHKERRQ(ierr);
5062f71e704SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_lusol_seqaij_C",
5072f71e704SKris Buschelman                                            "MatConvert_LUSOL_SeqAIJ",MatConvert_LUSOL_SeqAIJ);CHKERRQ(ierr);
5082f71e704SKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
5092f71e704SKris Buschelman   *newmat = B;
5104eb8e494SKris Buschelman   PetscFunctionReturn(0);
5114eb8e494SKris Buschelman }
5122f71e704SKris Buschelman EXTERN_C_END
5132f71e704SKris Buschelman 
514f0c56d0fSKris Buschelman #undef __FUNCT__
515f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_LUSOL"
516dfbe8321SBarry Smith PetscErrorCode MatDuplicate_LUSOL(Mat A, MatDuplicateOption op, Mat *M) {
517dfbe8321SBarry Smith   PetscErrorCode ierr;
5188f340917SKris Buschelman   Mat_LUSOL *lu=(Mat_LUSOL *)A->spptr;
519f0c56d0fSKris Buschelman   PetscFunctionBegin;
5208f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
5213f953163SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_LUSOL));CHKERRQ(ierr);
522f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
523f0c56d0fSKris Buschelman }
524f0c56d0fSKris Buschelman 
5252f71e704SKris Buschelman /*MC
526fafad747SKris Buschelman   MATLUSOL - MATLUSOL = "lusol" - A matrix type providing direct solvers (LU) for sequential matrices
5272f71e704SKris Buschelman   via the external package LUSOL.
5282f71e704SKris Buschelman 
5292f71e704SKris Buschelman   If LUSOL is installed (see the manual for
5302f71e704SKris Buschelman   instructions on how to declare the existence of external packages),
5312f71e704SKris Buschelman   a matrix type can be constructed which invokes LUSOL solvers.
5322f71e704SKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATLUSOL).
5332f71e704SKris Buschelman   This matrix type is only supported for double precision real.
5342f71e704SKris Buschelman 
5352f71e704SKris Buschelman   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
536f0c56d0fSKris Buschelman   supported for this matrix type.  MatConvert can be called for a fast inplace conversion
537f0c56d0fSKris Buschelman   to and from the MATSEQAIJ matrix type.
5382f71e704SKris Buschelman 
5392f71e704SKris Buschelman   Options Database Keys:
5400bad9183SKris Buschelman . -mat_type lusol - sets the matrix type to "lusol" during a call to MatSetFromOptions()
5412f71e704SKris Buschelman 
5422f71e704SKris Buschelman    Level: beginner
5432f71e704SKris Buschelman 
5442f71e704SKris Buschelman .seealso: PCLU
5452f71e704SKris Buschelman M*/
5464eb8e494SKris Buschelman 
5474eb8e494SKris Buschelman EXTERN_C_BEGIN
5484eb8e494SKris Buschelman #undef __FUNCT__
549f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_LUSOL"
550dfbe8321SBarry Smith PetscErrorCode MatCreate_LUSOL(Mat A)
551dfbe8321SBarry Smith {
552dfbe8321SBarry Smith   PetscErrorCode ierr;
5534eb8e494SKris Buschelman 
5544eb8e494SKris Buschelman   PetscFunctionBegin;
5555441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ and LUSOL types */
5565441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATLUSOL);CHKERRQ(ierr);
5574eb8e494SKris Buschelman   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
5582f71e704SKris Buschelman   ierr = MatConvert_SeqAIJ_LUSOL(A,MATLUSOL,&A);CHKERRQ(ierr);
5594eb8e494SKris Buschelman   PetscFunctionReturn(0);
5604eb8e494SKris Buschelman }
5614eb8e494SKris Buschelman EXTERN_C_END
562