xref: /petsc/src/mat/impls/baij/seq/dgefa.c (revision be1d678a52e6eff2808b2fa31ae986cdbf03c9fe)
1*be1d678aSKris Buschelman #define PETSCMAT_DLL
2*be1d678aSKris Buschelman 
3a97b7df5SSatish Balay /*
49fd41dd0SBarry Smith        This routine was converted by f2c from Linpack source
59fd41dd0SBarry Smith              linpack. this version dated 08/14/78
69fd41dd0SBarry Smith       cleve moler, university of new mexico, argonne national lab.
759539b86SBarry Smith 
859539b86SBarry Smith         Does an LU factorization with partial pivoting of a dense
99d8b60e5SBarry Smith      n by n matrix.
109d8b60e5SBarry Smith 
119d8b60e5SBarry Smith        Used by the sparse factorization routines in
1259539b86SBarry Smith      src/mat/impls/baij/seq and src/mat/impls/bdiag/seq
1359539b86SBarry Smith 
1471c5468dSBarry Smith        See also src/inline/ilu.h
159fd41dd0SBarry Smith */
169fd41dd0SBarry Smith #include "petsc.h"
179fd41dd0SBarry Smith 
184a2ae208SSatish Balay #undef __FUNCT__
194a2ae208SSatish Balay #define __FUNCT__ "LINPACKdgefa"
20690b6cddSBarry Smith PetscErrorCode LINPACKdgefa(MatScalar *a,PetscInt n,PetscInt *ipvt)
219fd41dd0SBarry Smith {
22690b6cddSBarry Smith     PetscInt   i__2,i__3,kp1,nm1,j,k,l,ll,kn,knp1,jn1;
233f1db9ecSBarry Smith     MatScalar  t,*ax,*ay,*aa;
24329f5518SBarry Smith     MatReal    tmp,max;
259fd41dd0SBarry Smith 
269fd41dd0SBarry Smith /*     gaussian elimination with partial pivoting */
279fd41dd0SBarry Smith 
283a40ed3dSBarry Smith     PetscFunctionBegin;
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 
535b8514ebSBarry Smith 	if (a[l + kn] == 0.0) {
5477431f27SBarry Smith 	  SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",k-1);
559fd41dd0SBarry Smith 	}
569fd41dd0SBarry Smith 
579fd41dd0SBarry Smith /*           interchange if necessary */
589fd41dd0SBarry Smith 
599fd41dd0SBarry Smith 	if (l != k) {
6039d66777SBarry Smith 	  t = a[l + kn];
6139d66777SBarry Smith 	  a[l + kn] = a[knp1];
6239d66777SBarry Smith 	  a[knp1] = t;
639fd41dd0SBarry Smith         }
649fd41dd0SBarry Smith 
659fd41dd0SBarry Smith /*           compute multipliers */
669fd41dd0SBarry Smith 
6739d66777SBarry Smith 	t = -1. / a[knp1];
689fd41dd0SBarry Smith 	i__2 = n - k;
6939d66777SBarry Smith         aa = &a[1 + knp1];
709fd41dd0SBarry Smith         for (ll=0; ll<i__2; ll++) {
719fd41dd0SBarry Smith           aa[ll] *= t;
729fd41dd0SBarry Smith         }
739fd41dd0SBarry Smith 
749fd41dd0SBarry Smith /*           row elimination with column indexing */
759fd41dd0SBarry Smith 
7639d66777SBarry Smith 	ax = aa;
779fd41dd0SBarry Smith         for (j = kp1; j <= n; ++j) {
788d3e6ddaSBarry Smith             jn1 = j*n;
798d3e6ddaSBarry Smith 	    t = a[l + jn1];
809fd41dd0SBarry Smith 	    if (l != k) {
818d3e6ddaSBarry Smith 	      a[l + jn1] = a[k + jn1];
828d3e6ddaSBarry Smith 	      a[k + jn1] = t;
839fd41dd0SBarry Smith             }
849fd41dd0SBarry Smith 
859fd41dd0SBarry Smith 	    i__3 = n - k;
868d3e6ddaSBarry Smith             ay = &a[1+k+jn1];
879fd41dd0SBarry Smith             for (ll=0; ll<i__3; ll++) {
889fd41dd0SBarry Smith               ay[ll] += t*ax[ll];
899fd41dd0SBarry Smith             }
909fd41dd0SBarry Smith 	}
919fd41dd0SBarry Smith     }
929fd41dd0SBarry Smith     ipvt[n] = n;
935b8514ebSBarry Smith     if (a[n + n * n] == 0.0) {
9477431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D",n-1);
959fd41dd0SBarry Smith     }
963a40ed3dSBarry Smith     PetscFunctionReturn(0);
979fd41dd0SBarry Smith }
989fd41dd0SBarry Smith 
99