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