xref: /petsc/src/mat/impls/baij/seq/dgefa.c (revision 5f8bbcca2696d61632b174eeae8a5433a128d797)
1be1d678aSKris Buschelman 
2a97b7df5SSatish Balay /*
39fd41dd0SBarry Smith        This routine was converted by f2c from Linpack source
49fd41dd0SBarry Smith              linpack. this version dated 08/14/78
59fd41dd0SBarry Smith       cleve moler, university of new mexico, argonne national lab.
659539b86SBarry Smith 
759539b86SBarry Smith         Does an LU factorization with partial pivoting of a dense
89d8b60e5SBarry Smith      n by n matrix.
99d8b60e5SBarry Smith 
109d8b60e5SBarry Smith        Used by the sparse factorization routines in
11dd882469SBarry Smith      src/mat/impls/baij/seq
1259539b86SBarry Smith 
139fd41dd0SBarry Smith */
14c6db04a5SJed Brown #include <petscsys.h>
159fd41dd0SBarry Smith 
164a2ae208SSatish Balay #undef __FUNCT__
1796b95a6bSBarry Smith #define __FUNCT__ "PetscLINPACKgefa"
18*5f8bbccaSHong Zhang PetscErrorCode PetscLINPACKgefa(MatScalar *a,PetscInt n,PetscInt *ipvt,PetscBool allowzeropivot,PetscBool *zeropivotdetected)
199fd41dd0SBarry Smith {
20690b6cddSBarry Smith   PetscInt  i__2,i__3,kp1,nm1,j,k,l,ll,kn,knp1,jn1;
213f1db9ecSBarry Smith   MatScalar t,*ax,*ay,*aa;
22329f5518SBarry Smith   MatReal   tmp,max;
239fd41dd0SBarry Smith 
249fd41dd0SBarry Smith   /* gaussian elimination with partial pivoting */
259fd41dd0SBarry Smith 
263a40ed3dSBarry Smith   PetscFunctionBegin;
27*5f8bbccaSHong Zhang   if (zeropivotdetected) *zeropivotdetected = PETSC_FALSE;
28*5f8bbccaSHong Zhang 
299fd41dd0SBarry Smith   /* Parameter adjustments */
309fd41dd0SBarry Smith   --ipvt;
3139d66777SBarry Smith   a -= n + 1;
329fd41dd0SBarry Smith 
339fd41dd0SBarry Smith   /* Function Body */
349fd41dd0SBarry Smith   nm1 = n - 1;
3539d66777SBarry Smith   for (k = 1; k <= nm1; ++k) {
369fd41dd0SBarry Smith     kp1  = k + 1;
3739d66777SBarry Smith     kn   = k*n;
3839d66777SBarry Smith     knp1 = k*n + k;
399fd41dd0SBarry Smith 
409fd41dd0SBarry Smith     /* find l = pivot index */
419fd41dd0SBarry Smith 
429fd41dd0SBarry Smith     i__2 = n - k + 1;
4339d66777SBarry Smith     aa   = &a[knp1];
449fd41dd0SBarry Smith     max  = PetscAbsScalar(aa[0]);
459fd41dd0SBarry Smith     l    = 1;
469fd41dd0SBarry Smith     for (ll=1; ll<i__2; ll++) {
479fd41dd0SBarry Smith       tmp = PetscAbsScalar(aa[ll]);
489fd41dd0SBarry Smith       if (tmp > max) { max = tmp; l = ll+1;}
499fd41dd0SBarry Smith     }
509fd41dd0SBarry Smith     l      += k - 1;
519fd41dd0SBarry Smith     ipvt[k] = l;
529fd41dd0SBarry Smith 
53*5f8bbccaSHong Zhang     if (a[l + kn] == 0.0) {
54*5f8bbccaSHong Zhang       if (allowzeropivot) {
55*5f8bbccaSHong Zhang         PetscErrorCode ierr;
56*5f8bbccaSHong Zhang         ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",k-1);CHKERRQ(ierr);
57*5f8bbccaSHong Zhang         *zeropivotdetected = PETSC_TRUE;
58*5f8bbccaSHong Zhang       } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
59*5f8bbccaSHong Zhang     }
609fd41dd0SBarry Smith 
619fd41dd0SBarry Smith     /* interchange if necessary */
629fd41dd0SBarry Smith     if (l != k) {
6339d66777SBarry Smith       t         = a[l + kn];
6439d66777SBarry Smith       a[l + kn] = a[knp1];
6539d66777SBarry Smith       a[knp1]   = t;
669fd41dd0SBarry Smith     }
679fd41dd0SBarry Smith 
689fd41dd0SBarry Smith     /* compute multipliers */
6939d66777SBarry Smith     t    = -1. / a[knp1];
709fd41dd0SBarry Smith     i__2 = n - k;
7139d66777SBarry Smith     aa   = &a[1 + knp1];
7226fbe8dcSKarl Rupp     for (ll=0; ll<i__2; ll++) aa[ll] *= t;
739fd41dd0SBarry Smith 
749fd41dd0SBarry Smith     /* row elimination with column indexing */
7539d66777SBarry Smith     ax = aa;
769fd41dd0SBarry Smith     for (j = kp1; j <= n; ++j) {
778d3e6ddaSBarry Smith       jn1 = j*n;
788d3e6ddaSBarry Smith       t   = a[l + jn1];
799fd41dd0SBarry Smith       if (l != k) {
808d3e6ddaSBarry Smith         a[l + jn1] = a[k + jn1];
818d3e6ddaSBarry Smith         a[k + jn1] = t;
829fd41dd0SBarry Smith       }
839fd41dd0SBarry Smith 
849fd41dd0SBarry Smith       i__3 = n - k;
858d3e6ddaSBarry Smith       ay   = &a[1+k+jn1];
8626fbe8dcSKarl Rupp       for (ll=0; ll<i__3; ll++) ay[ll] += t*ax[ll];
879fd41dd0SBarry Smith     }
889fd41dd0SBarry Smith   }
899fd41dd0SBarry Smith   ipvt[n] = n;
90*5f8bbccaSHong Zhang   if (a[n + n * n] == 0.0) {
91*5f8bbccaSHong Zhang     if (allowzeropivot) {
92*5f8bbccaSHong Zhang       PetscErrorCode ierr;
93*5f8bbccaSHong Zhang       ierr = PetscInfo1(NULL,"Zero pivot, row %D\n",n-1);CHKERRQ(ierr);
94*5f8bbccaSHong Zhang       *zeropivotdetected = PETSC_TRUE;
95*5f8bbccaSHong Zhang     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",n-1);
96*5f8bbccaSHong Zhang   }
973a40ed3dSBarry Smith   PetscFunctionReturn(0);
989fd41dd0SBarry Smith }
999fd41dd0SBarry Smith 
100