Actual source code: dgefa.c
petsc-3.8.4 2018-03-24
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.
7: Does an LU factorization with partial pivoting of a dense
8: n by n matrix.
10: Used by the sparse factorization routines in
11: src/mat/impls/baij/seq
13: */
14: #include <petscsys.h>
16: PETSC_INTERN PetscErrorCode PetscLINPACKgefa(MatScalar *a,PetscInt n,PetscInt *ipvt,PetscBool allowzeropivot,PetscBool *zeropivotdetected)
17: {
18: PetscInt i__2,i__3,kp1,nm1,j,k,l,ll,kn,knp1,jn1;
19: MatScalar t,*ax,*ay,*aa;
20: MatReal tmp,max;
23: if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
25: /* Parameter adjustments */
26: --ipvt;
27: a -= n + 1;
29: /* Function Body */
30: nm1 = n - 1;
31: for (k = 1; k <= nm1; ++k) {
32: kp1 = k + 1;
33: kn = k*n;
34: knp1 = k*n + k;
36: /* find l = pivot index */
38: i__2 = n - k + 1;
39: aa = &a[knp1];
40: max = PetscAbsScalar(aa[0]);
41: l = 1;
42: for (ll=1; ll<i__2; ll++) {
43: tmp = PetscAbsScalar(aa[ll]);
44: if (tmp > max) { max = tmp; l = ll+1;}
45: }
46: l += k - 1;
47: ipvt[k] = l;
49: if (a[l + kn] == 0.0) {
50: if (allowzeropivot) {
52: PetscInfo1(NULL,"Zero pivot, row %D\n",k-1);
53: if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
54: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
55: }
57: /* interchange if necessary */
58: if (l != k) {
59: t = a[l + kn];
60: a[l + kn] = a[knp1];
61: a[knp1] = t;
62: }
64: /* compute multipliers */
65: t = -1. / a[knp1];
66: i__2 = n - k;
67: aa = &a[1 + knp1];
68: for (ll=0; ll<i__2; ll++) aa[ll] *= t;
70: /* row elimination with column indexing */
71: ax = aa;
72: for (j = kp1; j <= n; ++j) {
73: jn1 = j*n;
74: t = a[l + jn1];
75: if (l != k) {
76: a[l + jn1] = a[k + jn1];
77: a[k + jn1] = t;
78: }
80: i__3 = n - k;
81: ay = &a[1+k+jn1];
82: for (ll=0; ll<i__3; ll++) ay[ll] += t*ax[ll];
83: }
84: }
85: ipvt[n] = n;
86: if (a[n + n * n] == 0.0) {
87: if (allowzeropivot) {
89: PetscInfo1(NULL,"Zero pivot, row %D\n",n-1);
90: if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE;
91: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",n-1);
92: }
93: return(0);
94: }