Actual source code: dgefa.c
1: /*$Id: dgefa.c,v 1.22 2001/03/23 23:22:07 balay Exp $*/
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 and src/mat/impls/bdiag/seq
13: See also src/inline/ilu.h
14: */
15: #include petsc.h
19: int LINPACKdgefa(MatScalar *a,int n,int *ipvt)
20: {
21: int i__2,i__3,kp1,nm1,j,k,l,ll,kn,knp1,jn1;
22: MatScalar t,*ax,*ay,*aa;
23: MatReal tmp,max;
25: /* gaussian elimination with partial pivoting */
28: /* Parameter adjustments */
29: --ipvt;
30: a -= n + 1;
32: /* Function Body */
33: nm1 = n - 1;
34: for (k = 1; k <= nm1; ++k) {
35: kp1 = k + 1;
36: kn = k*n;
37: knp1 = k*n + k;
39: /* find l = pivot index */
41: i__2 = n - k + 1;
42: aa = &a[knp1];
43: max = PetscAbsScalar(aa[0]);
44: l = 1;
45: for (ll=1; ll<i__2; ll++) {
46: tmp = PetscAbsScalar(aa[ll]);
47: if (tmp > max) { max = tmp; l = ll+1;}
48: }
49: l += k - 1;
50: ipvt[k] = l;
52: if (a[l + kn] == 0.) {
53: SETERRQ(k,"Zero pivot");
54: }
56: /* 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 */
66: t = -1. / a[knp1];
67: i__2 = n - k;
68: aa = &a[1 + knp1];
69: for (ll=0; ll<i__2; ll++) {
70: aa[ll] *= t;
71: }
73: /* 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: }
84: i__3 = n - k;
85: ay = &a[1+k+jn1];
86: for (ll=0; ll<i__3; ll++) {
87: ay[ll] += t*ax[ll];
88: }
89: }
90: }
91: ipvt[n] = n;
92: if (a[n + n * n] == 0.) {
93: SETERRQ(n,"Zero pivot,final row");
94: }
95: return(0);
96: }