1 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. 6 7 Does an LU factorization with partial pivoting of a dense 8 n by n matrix. 9 10 Used by the sparse factorization routines in 11 src/mat/impls/baij/seq 12 13 */ 14 #include <petscsys.h> 15 #include <petsc/private/kernels/blockinvert.h> 16 17 PetscErrorCode PetscLINPACKgefa(MatScalar *a,PetscInt n,PetscInt *ipvt,PetscBool allowzeropivot,PetscBool *zeropivotdetected) 18 { 19 PetscInt i__2,i__3,kp1,nm1,j,k,l,ll,kn,knp1,jn1; 20 MatScalar t,*ax,*ay,*aa; 21 MatReal tmp,max; 22 23 PetscFunctionBegin; 24 if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE; 25 26 /* Parameter adjustments */ 27 --ipvt; 28 a -= n + 1; 29 30 /* Function Body */ 31 nm1 = n - 1; 32 for (k = 1; k <= nm1; ++k) { 33 kp1 = k + 1; 34 kn = k*n; 35 knp1 = k*n + k; 36 37 /* find l = pivot index */ 38 39 i__2 = n - k + 1; 40 aa = &a[knp1]; 41 max = PetscAbsScalar(aa[0]); 42 l = 1; 43 for (ll=1; ll<i__2; ll++) { 44 tmp = PetscAbsScalar(aa[ll]); 45 if (tmp > max) { max = tmp; l = ll+1;} 46 } 47 l += k - 1; 48 ipvt[k] = l; 49 50 if (a[l + kn] == 0.0) { 51 if (allowzeropivot) { 52 PetscErrorCode ierr; 53 ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",k-1);CHKERRQ(ierr); 54 if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE; 55 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1); 56 } 57 58 /* interchange if necessary */ 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 t = -1. / a[knp1]; 67 i__2 = n - k; 68 aa = &a[1 + knp1]; 69 for (ll=0; ll<i__2; ll++) aa[ll] *= t; 70 71 /* row elimination with column indexing */ 72 ax = aa; 73 for (j = kp1; j <= n; ++j) { 74 jn1 = j*n; 75 t = a[l + jn1]; 76 if (l != k) { 77 a[l + jn1] = a[k + jn1]; 78 a[k + jn1] = t; 79 } 80 81 i__3 = n - k; 82 ay = &a[1+k+jn1]; 83 for (ll=0; ll<i__3; ll++) ay[ll] += t*ax[ll]; 84 } 85 } 86 ipvt[n] = n; 87 if (a[n + n * n] == 0.0) { 88 if (allowzeropivot) { 89 PetscErrorCode ierr; 90 ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",n-1);CHKERRQ(ierr); 91 if (zeropivotdetected) *zeropivotdetected = PETSC_TRUE; 92 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",n-1); 93 } 94 PetscFunctionReturn(0); 95 } 96 97