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