xref: /petsc/src/mat/impls/aij/seq/lusol/lusol.c (revision 3e7ff0edd3573be01c8c0fa32db97c3db8fa5c8d)
1 
2 /*
3         Provides an interface to the LUSOL package of ....
4 
5 */
6 #include <../src/mat/impls/aij/seq/aij.h>
7 
8 #if defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
9 #define LU1FAC   lu1fac_
10 #define LU6SOL   lu6sol_
11 #define M1PAGE   m1page_
12 #define M5SETX   m5setx_
13 #define M6RDEL   m6rdel_
14 #elif !defined(PETSC_HAVE_FORTRAN_CAPS)
15 #define LU1FAC   lu1fac
16 #define LU6SOL   lu6sol
17 #define M1PAGE   m1page
18 #define M5SETX   m5setx
19 #define M6RDEL   m6rdel
20 #endif
21 
22 /*
23     Dummy symbols that the MINOS files mi25bfac.f and mi15blas.f may require
24 */
25 PETSC_EXTERN void PETSC_STDCALL M1PAGE()
26 {
27   ;
28 }
29 PETSC_EXTERN void PETSC_STDCALL M5SETX()
30 {
31   ;
32 }
33 
34 PETSC_EXTERN void PETSC_STDCALL M6RDEL()
35 {
36   ;
37 }
38 
39 PETSC_EXTERN void PETSC_STDCALL LU1FAC(int *m, int *n, int *nnz, int *size, int *luparm,
40                                  double *parmlu, double *data, int *indc, int *indr,
41                                  int *rowperm, int *colperm, int *collen, int *rowlen,
42                                  int *colstart, int *rowstart, int *rploc, int *cploc,
43                                  int *rpinv, int *cpinv, double *w, int *inform);
44 
45 PETSC_EXTERN void PETSC_STDCALL LU6SOL(int *mode, int *m, int *n, double *rhs, double *x,
46                                  int *size, int *luparm, double *parmlu, double *data,
47                                  int *indc, int *indr, int *rowperm, int *colperm,
48                                  int *collen, int *rowlen, int *colstart, int *rowstart,
49                                  int *inform);
50 
51 extern PetscErrorCode MatDuplicate_LUSOL(Mat,MatDuplicateOption,Mat*);
52 
53 typedef struct  {
54   double *data;
55   int    *indc;
56   int    *indr;
57 
58   int    *ip;
59   int    *iq;
60   int    *lenc;
61   int    *lenr;
62   int    *locc;
63   int    *locr;
64   int    *iploc;
65   int    *iqloc;
66   int    *ipinv;
67   int    *iqinv;
68   double *mnsw;
69   double *mnsv;
70 
71   double elbowroom;
72   double luroom;                /* Extra space allocated when factor fails   */
73   double parmlu[30];            /* Input/output to LUSOL                     */
74 
75   int n;                        /* Number of rows/columns in matrix          */
76   int nz;                       /* Number of nonzeros                        */
77   int nnz;                      /* Number of nonzeros allocated for factors  */
78   int luparm[30];               /* Input/output to LUSOL                     */
79 
80   PetscBool CleanUpLUSOL;
81 
82 } Mat_LUSOL;
83 
84 /*  LUSOL input/Output Parameters (Description uses C-style indexes
85  *
86  *  Input parameters                                        Typical value
87  *
88  *  luparm(0) = nout     File number for printed messages.         6
89  *  luparm(1) = lprint   Print level.                              0
90  *                    < 0 suppresses output.
91  *                    = 0 gives error messages.
92  *                    = 1 gives debug output from some of the
93  *                        other routines in LUSOL.
94  *                   >= 2 gives the pivot row and column and the
95  *                        no. of rows and columns involved at
96  *                        each elimination step in lu1fac.
97  *  luparm(2) = maxcol   lu1fac: maximum number of columns         5
98  *                        searched allowed in a Markowitz-type
99  *                        search for the next pivot element.
100  *                        For some of the factorization, the
101  *                        number of rows searched is
102  *                        maxrow = maxcol - 1.
103  *
104  *
105  *  Output parameters
106  *
107  *  luparm(9) = inform   Return code from last call to any LU routine.
108  *  luparm(10) = nsing    No. of singularities marked in the
109  *                        output array w(*).
110  *  luparm(11) = jsing    Column index of last singularity.
111  *  luparm(12) = minlen   Minimum recommended value for  lena.
112  *  luparm(13) = maxlen   ?
113  *  luparm(14) = nupdat   No. of updates performed by the lu8 routines.
114  *  luparm(15) = nrank    No. of nonempty rows of U.
115  *  luparm(16) = ndens1   No. of columns remaining when the density of
116  *                        the matrix being factorized reached dens1.
117  *  luparm(17) = ndens2   No. of columns remaining when the density of
118  *                        the matrix being factorized reached dens2.
119  *  luparm(18) = jumin    The column index associated with dumin.
120  *  luparm(19) = numl0    No. of columns in initial  L.
121  *  luparm(20) = lenl0    Size of initial  L  (no. of nonzeros).
122  *  luparm(21) = lenu0    Size of initial  U.
123  *  luparm(22) = lenl     Size of current  L.
124  *  luparm(23) = lenu     Size of current  U.
125  *  luparm(24) = lrow     Length of row file.
126  *  luparm(25) = ncp      No. of compressions of LU data structures.
127  *  luparm(26) = mersum   lu1fac: sum of Markowitz merit counts.
128  *  luparm(27) = nutri    lu1fac: triangular rows in U.
129  *  luparm(28) = nltri    lu1fac: triangular rows in L.
130  *  luparm(29) =
131  *
132  *
133  *  Input parameters                                        Typical value
134  *
135  *  parmlu(0) = elmax1   Max multiplier allowed in  L           10.0
136  *                        during factor.
137  *  parmlu(1) = elmax2   Max multiplier allowed in  L           10.0
138  *                        during updates.
139  *  parmlu(2) = small    Absolute tolerance for             eps**0.8
140  *                        treating reals as zero.     IBM double: 3.0d-13
141  *  parmlu(3) = utol1    Absolute tol for flagging          eps**0.66667
142  *                        small diagonals of U.       IBM double: 3.7d-11
143  *  parmlu(4) = utol2    Relative tol for flagging          eps**0.66667
144  *                        small diagonals of U.       IBM double: 3.7d-11
145  *  parmlu(5) = uspace   Factor limiting waste space in  U.      3.0
146  *                        In lu1fac, the row or column lists
147  *                        are compressed if their length
148  *                        exceeds uspace times the length of
149  *                        either file after the last compression.
150  *  parmlu(6) = dens1    The density at which the Markowitz      0.3
151  *                        strategy should search maxcol columns
152  *                        and no rows.
153  *  parmlu(7) = dens2    the density at which the Markowitz      0.6
154  *                        strategy should search only 1 column
155  *                        or (preferably) use a dense LU for
156  *                        all the remaining rows and columns.
157  *
158  *
159  *  Output parameters
160  *
161  *  parmlu(9) = amax     Maximum element in  A.
162  *  parmlu(10) = elmax    Maximum multiplier in current  L.
163  *  parmlu(11) = umax     Maximum element in current  U.
164  *  parmlu(12) = dumax    Maximum diagonal in  U.
165  *  parmlu(13) = dumin    Minimum diagonal in  U.
166  *  parmlu(14) =
167  *  parmlu(15) =
168  *  parmlu(16) =
169  *  parmlu(17) =
170  *  parmlu(18) =
171  *  parmlu(19) = resid    lu6sol: residual after solve with U or U'.
172  *  ...
173  *  parmlu(29) =
174  */
175 
176 #define Factorization_Tolerance       1e-1
177 #define Factorization_Pivot_Tolerance pow(2.2204460492503131E-16, 2.0 / 3.0)
178 #define Factorization_Small_Tolerance 1e-15 /* pow(DBL_EPSILON, 0.8) */
179 
180 #undef __FUNCT__
181 #define __FUNCT__ "MatDestroy_LUSOL"
182 PetscErrorCode MatDestroy_LUSOL(Mat A)
183 {
184   PetscErrorCode ierr;
185   Mat_LUSOL      *lusol=(Mat_LUSOL*)A->spptr;
186 
187   PetscFunctionBegin;
188   if (lusol && lusol->CleanUpLUSOL) {
189     ierr = PetscFree(lusol->ip);CHKERRQ(ierr);
190     ierr = PetscFree(lusol->iq);CHKERRQ(ierr);
191     ierr = PetscFree(lusol->lenc);CHKERRQ(ierr);
192     ierr = PetscFree(lusol->lenr);CHKERRQ(ierr);
193     ierr = PetscFree(lusol->locc);CHKERRQ(ierr);
194     ierr = PetscFree(lusol->locr);CHKERRQ(ierr);
195     ierr = PetscFree(lusol->iploc);CHKERRQ(ierr);
196     ierr = PetscFree(lusol->iqloc);CHKERRQ(ierr);
197     ierr = PetscFree(lusol->ipinv);CHKERRQ(ierr);
198     ierr = PetscFree(lusol->iqinv);CHKERRQ(ierr);
199     ierr = PetscFree(lusol->mnsw);CHKERRQ(ierr);
200     ierr = PetscFree(lusol->mnsv);CHKERRQ(ierr);
201     ierr = PetscFree3(lusol->data,lusol->indc,lusol->indr);CHKERRQ(ierr);
202   }
203   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
204   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
205   PetscFunctionReturn(0);
206 }
207 
208 #undef __FUNCT__
209 #define __FUNCT__  "MatSolve_LUSOL"
210 PetscErrorCode MatSolve_LUSOL(Mat A,Vec b,Vec x)
211 {
212   Mat_LUSOL      *lusol=(Mat_LUSOL*)A->spptr;
213   double         *bb,*xx;
214   int            mode=5;
215   PetscErrorCode ierr;
216   int            i,m,n,nnz,status;
217 
218   PetscFunctionBegin;
219   ierr = VecGetArray(x, &xx);CHKERRQ(ierr);
220   ierr = VecGetArray(b, &bb);CHKERRQ(ierr);
221 
222   m   = n = lusol->n;
223   nnz = lusol->nnz;
224 
225   for (i = 0; i < m; i++) lusol->mnsv[i] = bb[i];
226 
227   LU6SOL(&mode, &m, &n, lusol->mnsv, xx, &nnz,
228          lusol->luparm, lusol->parmlu, lusol->data,
229          lusol->indc, lusol->indr, lusol->ip, lusol->iq,
230          lusol->lenc, lusol->lenr, lusol->locc, lusol->locr, &status);
231 
232   if (status) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"solve failed, error code %d",status);
233 
234   ierr = VecRestoreArray(x, &xx);CHKERRQ(ierr);
235   ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr);
236   PetscFunctionReturn(0);
237 }
238 
239 #undef __FUNCT__
240 #define __FUNCT__ "MatLUFactorNumeric_LUSOL"
241 PetscErrorCode MatLUFactorNumeric_LUSOL(Mat F,Mat A,const MatFactorInfo *info)
242 {
243   Mat_SeqAIJ     *a;
244   Mat_LUSOL      *lusol = (Mat_LUSOL*)F->spptr;
245   PetscErrorCode ierr;
246   int            m, n, nz, nnz, status;
247   int            i, rs, re;
248   int            factorizations;
249 
250   PetscFunctionBegin;
251   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);CHKERRQ(ierr);
252   a    = (Mat_SeqAIJ*)A->data;
253 
254   if (m != lusol->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"factorization struct inconsistent");
255 
256   factorizations = 0;
257   do {
258     /*******************************************************************/
259     /* Check the workspace allocation.                                 */
260     /*******************************************************************/
261 
262     nz  = a->nz;
263     nnz = PetscMax(lusol->nnz, (int)(lusol->elbowroom*nz));
264     nnz = PetscMax(nnz, 5*n);
265 
266     if (nnz < lusol->luparm[12]) {
267       nnz = (int)(lusol->luroom * lusol->luparm[12]);
268     } else if ((factorizations > 0) && (lusol->luroom < 6)) {
269       lusol->luroom += 0.1;
270     }
271 
272     nnz = PetscMax(nnz, (int)(lusol->luroom*(lusol->luparm[22] + lusol->luparm[23])));
273 
274     if (nnz > lusol->nnz) {
275       ierr       = PetscFree3(lusol->data,lusol->indc,lusol->indr);CHKERRQ(ierr);
276       ierr       = PetscMalloc3(nnz,&lusol->data,nnz,&lusol->indc,nnz,&lusol->indr);CHKERRQ(ierr);
277       lusol->nnz = nnz;
278     }
279 
280     /*******************************************************************/
281     /* Fill in the data for the problem.      (1-based Fortran style)  */
282     /*******************************************************************/
283 
284     nz = 0;
285     for (i = 0; i < n; i++) {
286       rs = a->i[i];
287       re = a->i[i+1];
288 
289       while (rs < re) {
290         if (a->a[rs] != 0.0) {
291           lusol->indc[nz] = i + 1;
292           lusol->indr[nz] = a->j[rs] + 1;
293           lusol->data[nz] = a->a[rs];
294           nz++;
295         }
296         rs++;
297       }
298     }
299 
300     /*******************************************************************/
301     /* Do the factorization.                                           */
302     /*******************************************************************/
303 
304     LU1FAC(&m, &n, &nz, &nnz,
305            lusol->luparm, lusol->parmlu, lusol->data,
306            lusol->indc, lusol->indr, lusol->ip, lusol->iq,
307            lusol->lenc, lusol->lenr, lusol->locc, lusol->locr,
308            lusol->iploc, lusol->iqloc, lusol->ipinv,
309            lusol->iqinv, lusol->mnsw, &status);
310 
311     switch (status) {
312     case 0:         /* factored */
313       break;
314 
315     case 7:         /* insufficient memory */
316       break;
317 
318     case 1:
319     case -1:        /* singular */
320       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Singular matrix");
321 
322     case 3:
323     case 4:         /* error conditions */
324       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"matrix error");
325 
326     default:        /* unknown condition */
327       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"matrix unknown return code");
328     }
329 
330     factorizations++;
331   } while (status == 7);
332   F->ops->solve   = MatSolve_LUSOL;
333   F->assembled    = PETSC_TRUE;
334   F->preallocated = PETSC_TRUE;
335   PetscFunctionReturn(0);
336 }
337 
338 #undef __FUNCT__
339 #define __FUNCT__ "MatLUFactorSymbolic_LUSOL"
340 PetscErrorCode MatLUFactorSymbolic_LUSOL(Mat F,Mat A, IS r, IS c,const MatFactorInfo *info)
341 {
342   /************************************************************************/
343   /* Input                                                                */
344   /*     A  - matrix to factor                                            */
345   /*     r  - row permutation (ignored)                                   */
346   /*     c  - column permutation (ignored)                                */
347   /*                                                                      */
348   /* Output                                                               */
349   /*     F  - matrix storing the factorization;                           */
350   /************************************************************************/
351   Mat_LUSOL      *lusol;
352   PetscErrorCode ierr;
353   int            i, m, n, nz, nnz;
354 
355   PetscFunctionBegin;
356   /************************************************************************/
357   /* Check the arguments.                                                 */
358   /************************************************************************/
359 
360   ierr = MatGetSize(A, &m, &n);CHKERRQ(ierr);
361   nz   = ((Mat_SeqAIJ*)A->data)->nz;
362 
363   /************************************************************************/
364   /* Create the factorization.                                            */
365   /************************************************************************/
366 
367   F->ops->lufactornumeric = MatLUFactorNumeric_LUSOL;
368   lusol                   = (Mat_LUSOL*)(F->spptr);
369 
370   /************************************************************************/
371   /* Initialize parameters                                                */
372   /************************************************************************/
373 
374   for (i = 0; i < 30; i++) {
375     lusol->luparm[i] = 0;
376     lusol->parmlu[i] = 0;
377   }
378 
379   lusol->luparm[1] = -1;
380   lusol->luparm[2] = 5;
381   lusol->luparm[7] = 1;
382 
383   lusol->parmlu[0] = 1 / Factorization_Tolerance;
384   lusol->parmlu[1] = 1 / Factorization_Tolerance;
385   lusol->parmlu[2] = Factorization_Small_Tolerance;
386   lusol->parmlu[3] = Factorization_Pivot_Tolerance;
387   lusol->parmlu[4] = Factorization_Pivot_Tolerance;
388   lusol->parmlu[5] = 3.0;
389   lusol->parmlu[6] = 0.3;
390   lusol->parmlu[7] = 0.6;
391 
392   /************************************************************************/
393   /* Allocate the workspace needed by LUSOL.                              */
394   /************************************************************************/
395 
396   lusol->elbowroom = PetscMax(lusol->elbowroom, info->fill);
397   nnz              = PetscMax((int)(lusol->elbowroom*nz), 5*n);
398 
399   lusol->n      = n;
400   lusol->nz     = nz;
401   lusol->nnz    = nnz;
402   lusol->luroom = 1.75;
403 
404   ierr = PetscMalloc(sizeof(int)*n,&lusol->ip);
405   ierr = PetscMalloc(sizeof(int)*n,&lusol->iq);
406   ierr = PetscMalloc(sizeof(int)*n,&lusol->lenc);
407   ierr = PetscMalloc(sizeof(int)*n,&lusol->lenr);
408   ierr = PetscMalloc(sizeof(int)*n,&lusol->locc);
409   ierr = PetscMalloc(sizeof(int)*n,&lusol->locr);
410   ierr = PetscMalloc(sizeof(int)*n,&lusol->iploc);
411   ierr = PetscMalloc(sizeof(int)*n,&lusol->iqloc);
412   ierr = PetscMalloc(sizeof(int)*n,&lusol->ipinv);
413   ierr = PetscMalloc(sizeof(int)*n,&lusol->iqinv);
414   ierr = PetscMalloc(sizeof(double)*n,&lusol->mnsw);
415   ierr = PetscMalloc(sizeof(double)*n,&lusol->mnsv);
416 
417   ierr = PetscMalloc3(nnz,&lusol->data,nnz,&lusol->indc,nnz,&lusol->indr);CHKERRQ(ierr);
418 
419   lusol->CleanUpLUSOL     = PETSC_TRUE;
420   F->ops->lufactornumeric = MatLUFactorNumeric_LUSOL;
421   PetscFunctionReturn(0);
422 }
423 
424 #undef __FUNCT__
425 #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_lusol"
426 PetscErrorCode MatFactorGetSolverPackage_seqaij_lusol(Mat A,const MatSolverPackage *type)
427 {
428   PetscFunctionBegin;
429   *type = MATSOLVERLUSOL;
430   PetscFunctionReturn(0);
431 }
432 
433 #undef __FUNCT__
434 #define __FUNCT__ "MatGetFactor_seqaij_lusol"
435 PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_lusol(Mat A,MatFactorType ftype,Mat *F)
436 {
437   Mat            B;
438   Mat_LUSOL      *lusol;
439   PetscErrorCode ierr;
440   int            m, n;
441 
442   PetscFunctionBegin;
443   ierr = MatGetSize(A, &m, &n);CHKERRQ(ierr);
444   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
445   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
446   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
447   ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
448 
449   ierr     = PetscNewLog(B,&lusol);CHKERRQ(ierr);
450   B->spptr = lusol;
451 
452   B->ops->lufactorsymbolic = MatLUFactorSymbolic_LUSOL;
453   B->ops->destroy          = MatDestroy_LUSOL;
454 
455   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_lusol);CHKERRQ(ierr);
456 
457   B->factortype = MAT_FACTOR_LU;
458   PetscFunctionReturn(0);
459 }
460 
461 #undef __FUNCT__
462 #define __FUNCT__ "MatSolverPackageRegister_Lusol"
463 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_Lusol(void)
464 {
465   PetscErrorCode ierr;
466 
467   PetscFunctionBegin;
468   ierr = MatSolverPackageRegister(MATSOLVERLUSOL,MATSEQAIJ,        MAT_FACTOR_LU,MatGetFactor_seqaij_lusol);CHKERRQ(ierr);
469   PetscFunctionReturn(0);
470 }
471 
472 /*MC
473   MATSOLVERLUSOL - "lusol" - Provides direct solvers (LU) for sequential matrices
474                          via the external package LUSOL.
475 
476   If LUSOL is installed (see the manual for
477   instructions on how to declare the existence of external packages),
478 
479   Works with MATSEQAIJ matrices
480 
481    Level: beginner
482 
483 .seealso: PCLU, PCFactorSetMatSolverPackage(), MatSolverPackage
484 
485 M*/
486