xref: /petsc/src/mat/impls/aij/seq/lusol/lusol.c (revision 38f2d2fdb3b6f522a3102c6eb796cebecf3224c0)
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 
806849ba73SBarry Smith   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
816849ba73SBarry Smith   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
826849ba73SBarry Smith   PetscErrorCode (*MatDestroy)(Mat);
834eb8e494SKris Buschelman   PetscTruth CleanUpLUSOL;
844eb8e494SKris Buschelman 
85f0c56d0fSKris Buschelman } Mat_LUSOL;
864eb8e494SKris Buschelman 
874eb8e494SKris Buschelman /*  LUSOL input/Output Parameters (Description uses C-style indexes
884eb8e494SKris Buschelman  *
894eb8e494SKris Buschelman  *  Input parameters                                        Typical value
904eb8e494SKris Buschelman  *
914eb8e494SKris Buschelman  *  luparm(0) = nout     File number for printed messages.         6
924eb8e494SKris Buschelman  *  luparm(1) = lprint   Print level.                              0
934eb8e494SKris Buschelman  *                    < 0 suppresses output.
944eb8e494SKris Buschelman  *                    = 0 gives error messages.
954eb8e494SKris Buschelman  *                    = 1 gives debug output from some of the
964eb8e494SKris Buschelman  *                        other routines in LUSOL.
974eb8e494SKris Buschelman  *                   >= 2 gives the pivot row and column and the
984eb8e494SKris Buschelman  *                        no. of rows and columns involved at
994eb8e494SKris Buschelman  *                        each elimination step in lu1fac.
1004eb8e494SKris Buschelman  *  luparm(2) = maxcol   lu1fac: maximum number of columns         5
1014eb8e494SKris Buschelman  *                        searched allowed in a Markowitz-type
1024eb8e494SKris Buschelman  *                        search for the next pivot element.
1034eb8e494SKris Buschelman  *                        For some of the factorization, the
1044eb8e494SKris Buschelman  *                        number of rows searched is
1054eb8e494SKris Buschelman  *                        maxrow = maxcol - 1.
1064eb8e494SKris Buschelman  *
1074eb8e494SKris Buschelman  *
1084eb8e494SKris Buschelman  *  Output parameters
1094eb8e494SKris Buschelman  *
1104eb8e494SKris Buschelman  *  luparm(9) = inform   Return code from last call to any LU routine.
1114eb8e494SKris Buschelman  *  luparm(10) = nsing    No. of singularities marked in the
1124eb8e494SKris Buschelman  *                        output array w(*).
1134eb8e494SKris Buschelman  *  luparm(11) = jsing    Column index of last singularity.
1144eb8e494SKris Buschelman  *  luparm(12) = minlen   Minimum recommended value for  lena.
1154eb8e494SKris Buschelman  *  luparm(13) = maxlen   ?
1164eb8e494SKris Buschelman  *  luparm(14) = nupdat   No. of updates performed by the lu8 routines.
1174eb8e494SKris Buschelman  *  luparm(15) = nrank    No. of nonempty rows of U.
1184eb8e494SKris Buschelman  *  luparm(16) = ndens1   No. of columns remaining when the density of
1194eb8e494SKris Buschelman  *                        the matrix being factorized reached dens1.
1204eb8e494SKris Buschelman  *  luparm(17) = ndens2   No. of columns remaining when the density of
1214eb8e494SKris Buschelman  *                        the matrix being factorized reached dens2.
1224eb8e494SKris Buschelman  *  luparm(18) = jumin    The column index associated with dumin.
1234eb8e494SKris Buschelman  *  luparm(19) = numl0    No. of columns in initial  L.
1244eb8e494SKris Buschelman  *  luparm(20) = lenl0    Size of initial  L  (no. of nonzeros).
1254eb8e494SKris Buschelman  *  luparm(21) = lenu0    Size of initial  U.
1264eb8e494SKris Buschelman  *  luparm(22) = lenl     Size of current  L.
1274eb8e494SKris Buschelman  *  luparm(23) = lenu     Size of current  U.
1284eb8e494SKris Buschelman  *  luparm(24) = lrow     Length of row file.
1294eb8e494SKris Buschelman  *  luparm(25) = ncp      No. of compressions of LU data structures.
1304eb8e494SKris Buschelman  *  luparm(26) = mersum   lu1fac: sum of Markowitz merit counts.
1314eb8e494SKris Buschelman  *  luparm(27) = nutri    lu1fac: triangular rows in U.
1324eb8e494SKris Buschelman  *  luparm(28) = nltri    lu1fac: triangular rows in L.
1334eb8e494SKris Buschelman  *  luparm(29) =
1344eb8e494SKris Buschelman  *
1354eb8e494SKris Buschelman  *
1364eb8e494SKris Buschelman  *  Input parameters                                        Typical value
1374eb8e494SKris Buschelman  *
1384eb8e494SKris Buschelman  *  parmlu(0) = elmax1   Max multiplier allowed in  L           10.0
1394eb8e494SKris Buschelman  *                        during factor.
1404eb8e494SKris Buschelman  *  parmlu(1) = elmax2   Max multiplier allowed in  L           10.0
1414eb8e494SKris Buschelman  *                        during updates.
1424eb8e494SKris Buschelman  *  parmlu(2) = small    Absolute tolerance for             eps**0.8
1434eb8e494SKris Buschelman  *                        treating reals as zero.     IBM double: 3.0d-13
1444eb8e494SKris Buschelman  *  parmlu(3) = utol1    Absolute tol for flagging          eps**0.66667
1454eb8e494SKris Buschelman  *                        small diagonals of U.       IBM double: 3.7d-11
1464eb8e494SKris Buschelman  *  parmlu(4) = utol2    Relative tol for flagging          eps**0.66667
1474eb8e494SKris Buschelman  *                        small diagonals of U.       IBM double: 3.7d-11
1484eb8e494SKris Buschelman  *  parmlu(5) = uspace   Factor limiting waste space in  U.      3.0
1494eb8e494SKris Buschelman  *                        In lu1fac, the row or column lists
1504eb8e494SKris Buschelman  *                        are compressed if their length
1514eb8e494SKris Buschelman  *                        exceeds uspace times the length of
1524eb8e494SKris Buschelman  *                        either file after the last compression.
1534eb8e494SKris Buschelman  *  parmlu(6) = dens1    The density at which the Markowitz      0.3
1544eb8e494SKris Buschelman  *                        strategy should search maxcol columns
1554eb8e494SKris Buschelman  *                        and no rows.
1564eb8e494SKris Buschelman  *  parmlu(7) = dens2    the density at which the Markowitz      0.6
1574eb8e494SKris Buschelman  *                        strategy should search only 1 column
1584eb8e494SKris Buschelman  *                        or (preferably) use a dense LU for
1594eb8e494SKris Buschelman  *                        all the remaining rows and columns.
1604eb8e494SKris Buschelman  *
1614eb8e494SKris Buschelman  *
1624eb8e494SKris Buschelman  *  Output parameters
1634eb8e494SKris Buschelman  *
1644eb8e494SKris Buschelman  *  parmlu(9) = amax     Maximum element in  A.
1654eb8e494SKris Buschelman  *  parmlu(10) = elmax    Maximum multiplier in current  L.
1664eb8e494SKris Buschelman  *  parmlu(11) = umax     Maximum element in current  U.
1674eb8e494SKris Buschelman  *  parmlu(12) = dumax    Maximum diagonal in  U.
1684eb8e494SKris Buschelman  *  parmlu(13) = dumin    Minimum diagonal in  U.
1694eb8e494SKris Buschelman  *  parmlu(14) =
1704eb8e494SKris Buschelman  *  parmlu(15) =
1714eb8e494SKris Buschelman  *  parmlu(16) =
1724eb8e494SKris Buschelman  *  parmlu(17) =
1734eb8e494SKris Buschelman  *  parmlu(18) =
1744eb8e494SKris Buschelman  *  parmlu(19) = resid    lu6sol: residual after solve with U or U'.
1754eb8e494SKris Buschelman  *  ...
1764eb8e494SKris Buschelman  *  parmlu(29) =
1774eb8e494SKris Buschelman  */
1784eb8e494SKris Buschelman 
1794eb8e494SKris Buschelman #define Factorization_Tolerance       1e-1
1804eb8e494SKris Buschelman #define Factorization_Pivot_Tolerance pow(2.2204460492503131E-16, 2.0 / 3.0)
1814eb8e494SKris Buschelman #define Factorization_Small_Tolerance 1e-15 /* pow(DBL_EPSILON, 0.8) */
1824eb8e494SKris Buschelman 
1832f71e704SKris Buschelman EXTERN_C_BEGIN
1842f71e704SKris Buschelman #undef __FUNCT__
1852f71e704SKris Buschelman #define __FUNCT__ "MatConvert_LUSOL_SeqAIJ"
186be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_LUSOL_SeqAIJ(Mat A,const MatType type,MatReuse reuse,Mat *newmat)
187521d7252SBarry Smith {
188dfbe8321SBarry Smith   PetscErrorCode ierr;
1892f71e704SKris Buschelman   Mat            B=*newmat;
190f0c56d0fSKris Buschelman   Mat_LUSOL      *lusol=(Mat_LUSOL *)A->spptr;
1912f71e704SKris Buschelman 
1922f71e704SKris Buschelman   PetscFunctionBegin;
193ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
1942f71e704SKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
195f0c56d0fSKris Buschelman   }
196f0c56d0fSKris Buschelman   B->ops->duplicate        = lusol->MatDuplicate;
1972f71e704SKris Buschelman   B->ops->lufactorsymbolic = lusol->MatLUFactorSymbolic;
1982f71e704SKris Buschelman   B->ops->destroy          = lusol->MatDestroy;
1992f71e704SKris Buschelman 
2002f71e704SKris Buschelman   ierr = PetscFree(lusol);CHKERRQ(ierr);
201901853e0SKris Buschelman 
202901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_lusol_C","",PETSC_NULL);CHKERRQ(ierr);
203901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_lusol_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
204901853e0SKris Buschelman 
2052f71e704SKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
2062f71e704SKris Buschelman   *newmat = B;
2072f71e704SKris Buschelman   PetscFunctionReturn(0);
2082f71e704SKris Buschelman }
2092f71e704SKris Buschelman EXTERN_C_END
2104eb8e494SKris Buschelman 
2114eb8e494SKris Buschelman #undef __FUNCT__
212f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_LUSOL"
213dfbe8321SBarry Smith PetscErrorCode MatDestroy_LUSOL(Mat A)
214dfbe8321SBarry Smith {
215dfbe8321SBarry Smith   PetscErrorCode ierr;
216f0c56d0fSKris Buschelman   Mat_LUSOL      *lusol=(Mat_LUSOL *)A->spptr;
2174eb8e494SKris Buschelman 
2184eb8e494SKris Buschelman   PetscFunctionBegin;
2194eb8e494SKris Buschelman   if (lusol->CleanUpLUSOL) {
2204eb8e494SKris Buschelman     ierr = PetscFree(lusol->ip);CHKERRQ(ierr);
2214eb8e494SKris Buschelman     ierr = PetscFree(lusol->iq);CHKERRQ(ierr);
2224eb8e494SKris Buschelman     ierr = PetscFree(lusol->lenc);CHKERRQ(ierr);
2234eb8e494SKris Buschelman     ierr = PetscFree(lusol->lenr);CHKERRQ(ierr);
2244eb8e494SKris Buschelman     ierr = PetscFree(lusol->locc);CHKERRQ(ierr);
2254eb8e494SKris Buschelman     ierr = PetscFree(lusol->locr);CHKERRQ(ierr);
2264eb8e494SKris Buschelman     ierr = PetscFree(lusol->iploc);CHKERRQ(ierr);
2274eb8e494SKris Buschelman     ierr = PetscFree(lusol->iqloc);CHKERRQ(ierr);
2284eb8e494SKris Buschelman     ierr = PetscFree(lusol->ipinv);CHKERRQ(ierr);
2294eb8e494SKris Buschelman     ierr = PetscFree(lusol->iqinv);CHKERRQ(ierr);
2304eb8e494SKris Buschelman     ierr = PetscFree(lusol->mnsw);CHKERRQ(ierr);
2314eb8e494SKris Buschelman     ierr = PetscFree(lusol->mnsv);CHKERRQ(ierr);
2324eb8e494SKris Buschelman     ierr = PetscFree(lusol->indc);CHKERRQ(ierr);
2334eb8e494SKris Buschelman   }
2344eb8e494SKris Buschelman 
235ceb03754SKris Buschelman   ierr = MatConvert_LUSOL_SeqAIJ(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);
2362f71e704SKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
2374eb8e494SKris Buschelman   PetscFunctionReturn(0);
2384eb8e494SKris Buschelman }
2394eb8e494SKris Buschelman 
2404eb8e494SKris Buschelman #undef __FUNCT__
241f0c56d0fSKris Buschelman #define __FUNCT__  "MatSolve_LUSOL"
2426849ba73SBarry Smith PetscErrorCode MatSolve_LUSOL(Mat A,Vec b,Vec x)
2436849ba73SBarry Smith {
244f0c56d0fSKris Buschelman   Mat_LUSOL *lusol=(Mat_LUSOL*)A->spptr;
2454eb8e494SKris Buschelman   double    *bb,*xx;
2464eb8e494SKris Buschelman   int       mode=5;
2476849ba73SBarry Smith   PetscErrorCode ierr;
2486849ba73SBarry Smith   int       i,m,n,nnz,status;
2494eb8e494SKris Buschelman 
2504eb8e494SKris Buschelman   PetscFunctionBegin;
2514eb8e494SKris Buschelman   ierr = VecGetArray(x, &xx);CHKERRQ(ierr);
2524eb8e494SKris Buschelman   ierr = VecGetArray(b, &bb);CHKERRQ(ierr);
2534eb8e494SKris Buschelman 
2544eb8e494SKris Buschelman   m = n = lusol->n;
2554eb8e494SKris Buschelman   nnz = lusol->nnz;
2564eb8e494SKris Buschelman 
2574eb8e494SKris Buschelman   for (i = 0; i < m; i++)
2584eb8e494SKris Buschelman     {
2594eb8e494SKris Buschelman       lusol->mnsv[i] = bb[i];
2604eb8e494SKris Buschelman     }
2614eb8e494SKris Buschelman 
2624eb8e494SKris Buschelman   LU6SOL(&mode, &m, &n, lusol->mnsv, xx, &nnz,
2634eb8e494SKris Buschelman          lusol->luparm, lusol->parmlu, lusol->data,
2644eb8e494SKris Buschelman          lusol->indc, lusol->indr, lusol->ip, lusol->iq,
2654eb8e494SKris Buschelman          lusol->lenc, lusol->lenr, lusol->locc, lusol->locr, &status);
2664eb8e494SKris Buschelman 
2674eb8e494SKris Buschelman   if (status != 0)
2684eb8e494SKris Buschelman     {
2694eb8e494SKris Buschelman       SETERRQ(PETSC_ERR_ARG_SIZ,"solve failed");
2704eb8e494SKris Buschelman     }
2714eb8e494SKris Buschelman 
2724eb8e494SKris Buschelman   ierr = VecRestoreArray(x, &xx);CHKERRQ(ierr);
2734eb8e494SKris Buschelman   ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr);
2744eb8e494SKris Buschelman   PetscFunctionReturn(0);
2754eb8e494SKris Buschelman }
2764eb8e494SKris Buschelman 
2774eb8e494SKris Buschelman #undef __FUNCT__
278f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_LUSOL"
279af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_LUSOL(Mat A,MatFactorInfo *info,Mat *F)
2806849ba73SBarry Smith {
2814eb8e494SKris Buschelman   Mat_SeqAIJ     *a;
282f0c56d0fSKris Buschelman   Mat_LUSOL      *lusol = (Mat_LUSOL*)(*F)->spptr;
2836849ba73SBarry Smith   PetscErrorCode ierr;
2844eb8e494SKris Buschelman   int            m, n, nz, nnz, status;
2856849ba73SBarry Smith   int            i, rs, re;
2864eb8e494SKris Buschelman   int            factorizations;
2874eb8e494SKris Buschelman 
2884eb8e494SKris Buschelman   PetscFunctionBegin;
2894eb8e494SKris Buschelman   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);CHKERRQ(ierr);
2904eb8e494SKris Buschelman   a = (Mat_SeqAIJ *)A->data;
2914eb8e494SKris Buschelman 
2924eb8e494SKris Buschelman   if (m != lusol->n) {
2934eb8e494SKris Buschelman     SETERRQ(PETSC_ERR_ARG_SIZ,"factorization struct inconsistent");
2944eb8e494SKris Buschelman   }
2954eb8e494SKris Buschelman 
2964eb8e494SKris Buschelman   factorizations = 0;
2974eb8e494SKris Buschelman   do
2984eb8e494SKris Buschelman     {
2994eb8e494SKris Buschelman       /*******************************************************************/
3004eb8e494SKris Buschelman       /* Check the workspace allocation.                                 */
3014eb8e494SKris Buschelman       /*******************************************************************/
3024eb8e494SKris Buschelman 
3034eb8e494SKris Buschelman       nz = a->nz;
3044eb8e494SKris Buschelman       nnz = PetscMax(lusol->nnz, (int)(lusol->elbowroom*nz));
3054eb8e494SKris Buschelman       nnz = PetscMax(nnz, 5*n);
3064eb8e494SKris Buschelman 
3074eb8e494SKris Buschelman       if (nnz < lusol->luparm[12]){
3084eb8e494SKris Buschelman         nnz = (int)(lusol->luroom * lusol->luparm[12]);
3094eb8e494SKris Buschelman       } else if ((factorizations > 0) && (lusol->luroom < 6)){
3104eb8e494SKris Buschelman         lusol->luroom += 0.1;
3114eb8e494SKris Buschelman       }
3124eb8e494SKris Buschelman 
3134eb8e494SKris Buschelman       nnz = PetscMax(nnz, (int)(lusol->luroom*(lusol->luparm[22] + lusol->luparm[23])));
3144eb8e494SKris Buschelman 
3154eb8e494SKris Buschelman       if (nnz > lusol->nnz){
3164eb8e494SKris Buschelman         ierr = PetscFree(lusol->indc);CHKERRQ(ierr);
3174eb8e494SKris Buschelman         ierr        = PetscMalloc((sizeof(double)+2*sizeof(int))*nnz,&lusol->indc);CHKERRQ(ierr);
3184eb8e494SKris Buschelman         lusol->indr = lusol->indc + nnz;
3194eb8e494SKris Buschelman         lusol->data = (double *)(lusol->indr + nnz);
3204eb8e494SKris Buschelman         lusol->nnz  = nnz;
3214eb8e494SKris Buschelman       }
3224eb8e494SKris Buschelman 
3234eb8e494SKris Buschelman       /*******************************************************************/
3244eb8e494SKris Buschelman       /* Fill in the data for the problem.      (1-based Fortran style)  */
3254eb8e494SKris Buschelman       /*******************************************************************/
3264eb8e494SKris Buschelman 
3274eb8e494SKris Buschelman       nz = 0;
3284eb8e494SKris Buschelman       for (i = 0; i < n; i++)
3294eb8e494SKris Buschelman         {
3304eb8e494SKris Buschelman           rs = a->i[i];
3314eb8e494SKris Buschelman           re = a->i[i+1];
3324eb8e494SKris Buschelman 
3334eb8e494SKris Buschelman           while (rs < re)
3344eb8e494SKris Buschelman             {
3354eb8e494SKris Buschelman               if (a->a[rs] != 0.0)
3364eb8e494SKris Buschelman                 {
3374eb8e494SKris Buschelman                   lusol->indc[nz] = i + 1;
3384eb8e494SKris Buschelman                   lusol->indr[nz] = a->j[rs] + 1;
3394eb8e494SKris Buschelman                   lusol->data[nz] = a->a[rs];
3404eb8e494SKris Buschelman                   nz++;
3414eb8e494SKris Buschelman                 }
3424eb8e494SKris Buschelman               rs++;
3434eb8e494SKris Buschelman             }
3444eb8e494SKris Buschelman         }
3454eb8e494SKris Buschelman 
3464eb8e494SKris Buschelman       /*******************************************************************/
3474eb8e494SKris Buschelman       /* Do the factorization.                                           */
3484eb8e494SKris Buschelman       /*******************************************************************/
3494eb8e494SKris Buschelman 
3504eb8e494SKris Buschelman       LU1FAC(&m, &n, &nz, &nnz,
3514eb8e494SKris Buschelman              lusol->luparm, lusol->parmlu, lusol->data,
3524eb8e494SKris Buschelman              lusol->indc, lusol->indr, lusol->ip, lusol->iq,
3534eb8e494SKris Buschelman              lusol->lenc, lusol->lenr, lusol->locc, lusol->locr,
3544eb8e494SKris Buschelman              lusol->iploc, lusol->iqloc, lusol->ipinv,
3554eb8e494SKris Buschelman              lusol->iqinv, lusol->mnsw, &status);
3564eb8e494SKris Buschelman 
3574eb8e494SKris Buschelman       switch(status)
3584eb8e494SKris Buschelman         {
3594eb8e494SKris Buschelman         case 0:		/* factored */
3604eb8e494SKris Buschelman           break;
3614eb8e494SKris Buschelman 
3624eb8e494SKris Buschelman         case 7:		/* insufficient memory */
3634eb8e494SKris Buschelman           break;
3644eb8e494SKris Buschelman 
3654eb8e494SKris Buschelman         case 1:
3664eb8e494SKris Buschelman         case -1:		/* singular */
367e005ede5SBarry Smith           SETERRQ(PETSC_ERR_LIB,"Singular matrix");
3684eb8e494SKris Buschelman 
3694eb8e494SKris Buschelman         case 3:
3704eb8e494SKris Buschelman         case 4:		/* error conditions */
371e005ede5SBarry Smith           SETERRQ(PETSC_ERR_LIB,"matrix error");
3724eb8e494SKris Buschelman 
3734eb8e494SKris Buschelman         default:		/* unknown condition */
374e005ede5SBarry Smith           SETERRQ(PETSC_ERR_LIB,"matrix unknown return code");
3754eb8e494SKris Buschelman         }
3764eb8e494SKris Buschelman 
3774eb8e494SKris Buschelman       factorizations++;
3784eb8e494SKris Buschelman     } while (status == 7);
379a8883a68SKris Buschelman   (*F)->assembled = PETSC_TRUE;
3804eb8e494SKris Buschelman   PetscFunctionReturn(0);
3814eb8e494SKris Buschelman }
3824eb8e494SKris Buschelman 
3834eb8e494SKris Buschelman #undef __FUNCT__
384f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_LUSOL"
385dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_LUSOL(Mat A, IS r, IS c,MatFactorInfo *info, Mat *F) {
3864eb8e494SKris Buschelman   /************************************************************************/
3874eb8e494SKris Buschelman   /* Input                                                                */
3884eb8e494SKris Buschelman   /*     A  - matrix to factor                                            */
3894eb8e494SKris Buschelman   /*     r  - row permutation (ignored)                                   */
3904eb8e494SKris Buschelman   /*     c  - column permutation (ignored)                                */
3914eb8e494SKris Buschelman   /*                                                                      */
3924eb8e494SKris Buschelman   /* Output                                                               */
3934eb8e494SKris Buschelman   /*     F  - matrix storing the factorization;                           */
3944eb8e494SKris Buschelman   /************************************************************************/
3954eb8e494SKris Buschelman   Mat       B;
396f0c56d0fSKris Buschelman   Mat_LUSOL *lusol;
397dfbe8321SBarry Smith   PetscErrorCode ierr;
398dfbe8321SBarry Smith   int        i, m, n, nz, nnz;
3994eb8e494SKris Buschelman 
4004eb8e494SKris Buschelman   PetscFunctionBegin;
4014eb8e494SKris Buschelman 
4024eb8e494SKris Buschelman   /************************************************************************/
4034eb8e494SKris Buschelman   /* Check the arguments.                                                 */
4044eb8e494SKris Buschelman   /************************************************************************/
4054eb8e494SKris Buschelman 
4064eb8e494SKris Buschelman   ierr = MatGetSize(A, &m, &n);CHKERRQ(ierr);
4074eb8e494SKris Buschelman   nz = ((Mat_SeqAIJ *)A->data)->nz;
4084eb8e494SKris Buschelman 
4094eb8e494SKris Buschelman   /************************************************************************/
4104eb8e494SKris Buschelman   /* Create the factorization.                                            */
4114eb8e494SKris Buschelman   /************************************************************************/
4124eb8e494SKris Buschelman 
413f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
414f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
415be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
4164eb8e494SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
4174eb8e494SKris Buschelman 
418f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatLUFactorNumeric_LUSOL;
419f0c56d0fSKris Buschelman   B->ops->solve           = MatSolve_LUSOL;
4204eb8e494SKris Buschelman   B->factor               = FACTOR_LU;
421f0c56d0fSKris Buschelman   lusol                   = (Mat_LUSOL*)(B->spptr);
4224eb8e494SKris Buschelman 
4234eb8e494SKris Buschelman   /************************************************************************/
4244eb8e494SKris Buschelman   /* Initialize parameters                                                */
4254eb8e494SKris Buschelman   /************************************************************************/
4264eb8e494SKris Buschelman 
4274eb8e494SKris Buschelman   for (i = 0; i < 30; i++)
4284eb8e494SKris Buschelman     {
4294eb8e494SKris Buschelman       lusol->luparm[i] = 0;
4304eb8e494SKris Buschelman       lusol->parmlu[i] = 0;
4314eb8e494SKris Buschelman     }
4324eb8e494SKris Buschelman 
4334eb8e494SKris Buschelman   lusol->luparm[1] = -1;
4344eb8e494SKris Buschelman   lusol->luparm[2] = 5;
4354eb8e494SKris Buschelman   lusol->luparm[7] = 1;
4364eb8e494SKris Buschelman 
4374eb8e494SKris Buschelman   lusol->parmlu[0] = 1 / Factorization_Tolerance;
4384eb8e494SKris Buschelman   lusol->parmlu[1] = 1 / Factorization_Tolerance;
4394eb8e494SKris Buschelman   lusol->parmlu[2] = Factorization_Small_Tolerance;
4404eb8e494SKris Buschelman   lusol->parmlu[3] = Factorization_Pivot_Tolerance;
4414eb8e494SKris Buschelman   lusol->parmlu[4] = Factorization_Pivot_Tolerance;
4424eb8e494SKris Buschelman   lusol->parmlu[5] = 3.0;
4434eb8e494SKris Buschelman   lusol->parmlu[6] = 0.3;
4444eb8e494SKris Buschelman   lusol->parmlu[7] = 0.6;
4454eb8e494SKris Buschelman 
4464eb8e494SKris Buschelman   /************************************************************************/
4474eb8e494SKris Buschelman   /* Allocate the workspace needed by LUSOL.                              */
4484eb8e494SKris Buschelman   /************************************************************************/
4494eb8e494SKris Buschelman 
4504eb8e494SKris Buschelman   lusol->elbowroom = PetscMax(lusol->elbowroom, info->fill);
4514eb8e494SKris Buschelman   nnz = PetscMax((int)(lusol->elbowroom*nz), 5*n);
4524eb8e494SKris Buschelman 
4534eb8e494SKris Buschelman   lusol->n = n;
4544eb8e494SKris Buschelman   lusol->nz = nz;
4554eb8e494SKris Buschelman   lusol->nnz = nnz;
4564eb8e494SKris Buschelman   lusol->luroom = 1.75;
4574eb8e494SKris Buschelman 
4584eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->ip);
4594eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->iq);
4604eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->lenc);
4614eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->lenr);
4624eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->locc);
4634eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->locr);
4644eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->iploc);
4654eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->iqloc);
4664eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->ipinv);
4674eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(int)*n,&lusol->iqinv);
4684eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(double)*n,&lusol->mnsw);
4694eb8e494SKris Buschelman   ierr = PetscMalloc(sizeof(double)*n,&lusol->mnsv);
4704eb8e494SKris Buschelman 
4714eb8e494SKris Buschelman   ierr        = PetscMalloc((sizeof(double)+2*sizeof(int))*nnz,&lusol->indc);
4724eb8e494SKris Buschelman   lusol->indr = lusol->indc + nnz;
4734eb8e494SKris Buschelman   lusol->data = (double *)(lusol->indr + nnz);
4744eb8e494SKris Buschelman   lusol->CleanUpLUSOL = PETSC_TRUE;
4754eb8e494SKris Buschelman   *F = B;
4764eb8e494SKris Buschelman   PetscFunctionReturn(0);
4774eb8e494SKris Buschelman }
4784eb8e494SKris Buschelman 
4792f71e704SKris Buschelman EXTERN_C_BEGIN
4804eb8e494SKris Buschelman #undef __FUNCT__
4812f71e704SKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_LUSOL"
482be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_LUSOL(Mat A,const MatType type,MatReuse reuse,Mat *newmat)
483521d7252SBarry Smith {
484dfbe8321SBarry Smith   PetscErrorCode ierr;
485521d7252SBarry Smith   PetscInt       m, n;
486f0c56d0fSKris Buschelman   Mat_LUSOL      *lusol;
4872f71e704SKris Buschelman   Mat            B=*newmat;
4884eb8e494SKris Buschelman 
4894eb8e494SKris Buschelman   PetscFunctionBegin;
4904eb8e494SKris Buschelman   ierr = MatGetSize(A, &m, &n);CHKERRQ(ierr);
4914eb8e494SKris Buschelman   if (m != n) {
4924eb8e494SKris Buschelman     SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square");
4934eb8e494SKris Buschelman   }
494ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
4952f71e704SKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
4962f71e704SKris Buschelman   }
4974eb8e494SKris Buschelman 
498*38f2d2fdSLisandro Dalcin   ierr                       = PetscNewLog(B,Mat_LUSOL,&lusol);CHKERRQ(ierr);
499f0c56d0fSKris Buschelman   lusol->MatDuplicate        = A->ops->duplicate;
5002f71e704SKris Buschelman   lusol->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
5012f71e704SKris Buschelman   lusol->MatDestroy          = A->ops->destroy;
5022f71e704SKris Buschelman   lusol->CleanUpLUSOL        = PETSC_FALSE;
5032f71e704SKris Buschelman 
5042f71e704SKris Buschelman   B->spptr                   = (void*)lusol;
505f0c56d0fSKris Buschelman   B->ops->duplicate          = MatDuplicate_LUSOL;
506f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic   = MatLUFactorSymbolic_LUSOL;
507f0c56d0fSKris Buschelman   B->ops->destroy            = MatDestroy_LUSOL;
5082f71e704SKris Buschelman 
509ae15b995SBarry Smith   ierr = PetscInfo(0,"Using LUSOL for LU factorization and solves.\n");CHKERRQ(ierr);
5102f71e704SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_lusol_C",
5112f71e704SKris Buschelman                                            "MatConvert_SeqAIJ_LUSOL",MatConvert_SeqAIJ_LUSOL);CHKERRQ(ierr);
5122f71e704SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_lusol_seqaij_C",
5132f71e704SKris Buschelman                                            "MatConvert_LUSOL_SeqAIJ",MatConvert_LUSOL_SeqAIJ);CHKERRQ(ierr);
5142f71e704SKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
5152f71e704SKris Buschelman   *newmat = B;
5164eb8e494SKris Buschelman   PetscFunctionReturn(0);
5174eb8e494SKris Buschelman }
5182f71e704SKris Buschelman EXTERN_C_END
5192f71e704SKris Buschelman 
520f0c56d0fSKris Buschelman #undef __FUNCT__
521f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_LUSOL"
522dfbe8321SBarry Smith PetscErrorCode MatDuplicate_LUSOL(Mat A, MatDuplicateOption op, Mat *M) {
523dfbe8321SBarry Smith   PetscErrorCode ierr;
5248f340917SKris Buschelman   Mat_LUSOL *lu=(Mat_LUSOL *)A->spptr;
525f0c56d0fSKris Buschelman   PetscFunctionBegin;
5268f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
5273f953163SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_LUSOL));CHKERRQ(ierr);
528f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
529f0c56d0fSKris Buschelman }
530f0c56d0fSKris Buschelman 
5312f71e704SKris Buschelman /*MC
532fafad747SKris Buschelman   MATLUSOL - MATLUSOL = "lusol" - A matrix type providing direct solvers (LU) for sequential matrices
5332f71e704SKris Buschelman   via the external package LUSOL.
5342f71e704SKris Buschelman 
5352f71e704SKris Buschelman   If LUSOL is installed (see the manual for
5362f71e704SKris Buschelman   instructions on how to declare the existence of external packages),
5372f71e704SKris Buschelman   a matrix type can be constructed which invokes LUSOL solvers.
5382f71e704SKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATLUSOL).
5392f71e704SKris Buschelman   This matrix type is only supported for double precision real.
5402f71e704SKris Buschelman 
5412f71e704SKris Buschelman   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
542f0c56d0fSKris Buschelman   supported for this matrix type.  MatConvert can be called for a fast inplace conversion
543f0c56d0fSKris Buschelman   to and from the MATSEQAIJ matrix type.
5442f71e704SKris Buschelman 
5452f71e704SKris Buschelman   Options Database Keys:
5460bad9183SKris Buschelman . -mat_type lusol - sets the matrix type to "lusol" during a call to MatSetFromOptions()
5472f71e704SKris Buschelman 
5482f71e704SKris Buschelman    Level: beginner
5492f71e704SKris Buschelman 
5502f71e704SKris Buschelman .seealso: PCLU
5512f71e704SKris Buschelman M*/
5524eb8e494SKris Buschelman 
5534eb8e494SKris Buschelman EXTERN_C_BEGIN
5544eb8e494SKris Buschelman #undef __FUNCT__
555f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_LUSOL"
556be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_LUSOL(Mat A)
557dfbe8321SBarry Smith {
558dfbe8321SBarry Smith   PetscErrorCode ierr;
5594eb8e494SKris Buschelman 
5604eb8e494SKris Buschelman   PetscFunctionBegin;
5614eb8e494SKris Buschelman   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
562ceb03754SKris Buschelman   ierr = MatConvert_SeqAIJ_LUSOL(A,MATLUSOL,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
5634eb8e494SKris Buschelman   PetscFunctionReturn(0);
5644eb8e494SKris Buschelman }
5654eb8e494SKris Buschelman EXTERN_C_END
566