xref: /petsc/src/mat/impls/baij/seq/dgefa.c (revision 9d8b60e591cfb2044a029c66b8ad588a5ba4b21f)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: dgefa.c,v 1.12 1998/12/24 02:54:37 bsmith Exp bsmith $";
3 #endif
4 /*
5        This routine was converted by f2c from Linpack source
6              linpack. this version dated 08/14/78
7       cleve moler, university of new mexico, argonne national lab.
8 
9         Does an LU factorization with partial pivoting of a dense
10      n by n matrix.
11 
12        Used by the sparse factorization routines in
13      src/mat/impls/baij/seq and src/mat/impls/bdiag/seq
14 
15 */
16 #include "petsc.h"
17 
18 #undef __FUNC__
19 #define __FUNC__ "Linpack_DGEFA"
20 int Linpack_DGEFA(MatScalar *a, int n, int *ipvt)
21 {
22     int        i__2, i__3, kp1, nm1, j, k, l,ll,kn,knp1,jn1;
23     MatScalar  t,*ax,*ay,*aa;
24     MatFloat   tmp,max;
25 
26 /*     gaussian elimination with partial pivoting */
27 
28     PetscFunctionBegin;
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.) {
54 	  SETERRQ(k,0,"Zero pivot");
55 	}
56 
57 /*           interchange if necessary */
58 
59 	if (l != k) {
60 	  t = a[l + kn];
61 	  a[l + kn] = a[knp1];
62 	  a[knp1] = t;
63         }
64 
65 /*           compute multipliers */
66 
67 	t = -1. / a[knp1];
68 	i__2 = n - k;
69         aa = &a[1 + knp1];
70         for ( ll=0; ll<i__2; ll++ ) {
71           aa[ll] *= t;
72         }
73 
74 /*           row elimination with column indexing */
75 
76 	ax = aa;
77         for (j = kp1; j <= n; ++j) {
78             jn1 = j*n;
79 	    t = a[l + jn1];
80 	    if (l != k) {
81 	      a[l + jn1] = a[k + jn1];
82 	      a[k + jn1] = t;
83             }
84 
85 	    i__3 = n - k;
86             ay = &a[1+k+jn1];
87             for ( ll=0; ll<i__3; ll++ ) {
88               ay[ll] += t*ax[ll];
89             }
90 	}
91     }
92     ipvt[n] = n;
93     if (a[n + n * n] == 0.) {
94 	SETERRQ(n,0,"Zero pivot,final row");
95     }
96     PetscFunctionReturn(0);
97 }
98 
99