xref: /petsc/src/mat/impls/baij/seq/dgefa.c (revision 59539b86a521cd2148a0c66cdb96b74e0742dafb)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: dgefa.c,v 1.11 1998/12/17 22:10:39 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 (used by the sparse factorization routines in
11      src/mat/impls/baij/seq and src/mat/impls/bdiag/seq
12 
13 */
14 #include "petsc.h"
15 
16 #undef __FUNC__
17 #define __FUNC__ "Linpack_DGEFA"
18 int Linpack_DGEFA(MatScalar *a, int n, int *ipvt)
19 {
20     int        i__2, i__3, kp1, nm1, j, k, l,ll,kn,knp1,jn1;
21     MatScalar  t,*ax,*ay,*aa;
22     MatFloat   tmp,max;
23 
24 /*     gaussian elimination with partial pivoting */
25 
26     PetscFunctionBegin;
27     /* Parameter adjustments */
28     --ipvt;
29     a       -= n + 1;
30 
31     /* Function Body */
32     nm1 = n - 1;
33     for (k = 1; k <= nm1; ++k) {
34 	kp1  = k + 1;
35         kn   = k*n;
36         knp1 = k*n + k;
37 
38 /*        find l = pivot index */
39 
40 	i__2 = n - k + 1;
41         aa = &a[knp1];
42         max = PetscAbsScalar(aa[0]);
43         l = 1;
44         for ( ll=1; ll<i__2; ll++ ) {
45           tmp = PetscAbsScalar(aa[ll]);
46           if (tmp > max) { max = tmp; l = ll+1;}
47         }
48         l += k - 1;
49 	ipvt[k] = l;
50 
51 	if (a[l + kn] == 0.) {
52 	  SETERRQ(k,0,"Zero pivot");
53 	}
54 
55 /*           interchange if necessary */
56 
57 	if (l != k) {
58 	  t = a[l + kn];
59 	  a[l + kn] = a[knp1];
60 	  a[knp1] = t;
61         }
62 
63 /*           compute multipliers */
64 
65 	t = -1. / a[knp1];
66 	i__2 = n - k;
67         aa = &a[1 + knp1];
68         for ( ll=0; ll<i__2; ll++ ) {
69           aa[ll] *= t;
70         }
71 
72 /*           row elimination with column indexing */
73 
74 	ax = aa;
75         for (j = kp1; j <= n; ++j) {
76             jn1 = j*n;
77 	    t = a[l + jn1];
78 	    if (l != k) {
79 	      a[l + jn1] = a[k + jn1];
80 	      a[k + jn1] = t;
81             }
82 
83 	    i__3 = n - k;
84             ay = &a[1+k+jn1];
85             for ( ll=0; ll<i__3; ll++ ) {
86               ay[ll] += t*ax[ll];
87             }
88 	}
89     }
90     ipvt[n] = n;
91     if (a[n + n * n] == 0.) {
92 	SETERRQ(n,0,"Zero pivot,final row");
93     }
94     PetscFunctionReturn(0);
95 }
96 
97