xref: /petsc/src/mat/impls/baij/seq/dgefa.c (revision 5f8bbcca2696d61632b174eeae8a5433a128d797)
1 
2 /*
3        This routine was converted by f2c from Linpack source
4              linpack. this version dated 08/14/78
5       cleve moler, university of new mexico, argonne national lab.
6 
7         Does an LU factorization with partial pivoting of a dense
8      n by n matrix.
9 
10        Used by the sparse factorization routines in
11      src/mat/impls/baij/seq
12 
13 */
14 #include <petscsys.h>
15 
16 #undef __FUNCT__
17 #define __FUNCT__ "PetscLINPACKgefa"
18 PetscErrorCode PetscLINPACKgefa(MatScalar *a,PetscInt n,PetscInt *ipvt,PetscBool allowzeropivot,PetscBool *zeropivotdetected)
19 {
20   PetscInt  i__2,i__3,kp1,nm1,j,k,l,ll,kn,knp1,jn1;
21   MatScalar t,*ax,*ay,*aa;
22   MatReal   tmp,max;
23 
24   /* gaussian elimination with partial pivoting */
25 
26   PetscFunctionBegin;
27   if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
28 
29   /* Parameter adjustments */
30   --ipvt;
31   a -= n + 1;
32 
33   /* Function Body */
34   nm1 = n - 1;
35   for (k = 1; k <= nm1; ++k) {
36     kp1  = k + 1;
37     kn   = k*n;
38     knp1 = k*n + k;
39 
40     /* find l = pivot index */
41 
42     i__2 = n - k + 1;
43     aa   = &a[knp1];
44     max  = PetscAbsScalar(aa[0]);
45     l    = 1;
46     for (ll=1; ll<i__2; ll++) {
47       tmp = PetscAbsScalar(aa[ll]);
48       if (tmp > max) { max = tmp; l = ll+1;}
49     }
50     l      += k - 1;
51     ipvt[k] = l;
52 
53     if (a[l + kn] == 0.0) {
54       if (allowzeropivot) {
55         PetscErrorCode ierr;
56         ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",k-1);CHKERRQ(ierr);
57         *zeropivotdetected = PETSC_TRUE;
58       } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
59     }
60 
61     /* interchange if necessary */
62     if (l != k) {
63       t         = a[l + kn];
64       a[l + kn] = a[knp1];
65       a[knp1]   = t;
66     }
67 
68     /* compute multipliers */
69     t    = -1. / a[knp1];
70     i__2 = n - k;
71     aa   = &a[1 + knp1];
72     for (ll=0; ll<i__2; ll++) aa[ll] *= t;
73 
74     /* row elimination with column indexing */
75     ax = aa;
76     for (j = kp1; j <= n; ++j) {
77       jn1 = j*n;
78       t   = a[l + jn1];
79       if (l != k) {
80         a[l + jn1] = a[k + jn1];
81         a[k + jn1] = t;
82       }
83 
84       i__3 = n - k;
85       ay   = &a[1+k+jn1];
86       for (ll=0; ll<i__3; ll++) ay[ll] += t*ax[ll];
87     }
88   }
89   ipvt[n] = n;
90   if (a[n + n * n] == 0.0) {
91     if (allowzeropivot) {
92       PetscErrorCode ierr;
93       ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",n-1);CHKERRQ(ierr);
94       *zeropivotdetected = PETSC_TRUE;
95     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",n-1);
96   }
97   PetscFunctionReturn(0);
98 }
99 
100