xref: /petsc/src/mat/impls/baij/seq/dgefa.c (revision 3a40ed3dce77c081171d005ae1a6ff4bb9d13b6f)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*3a40ed3dSBarry Smith static char vcid[] = "$Id: dgefa.c,v 1.8 1997/07/09 20:55:07 balay Exp bsmith $";
34431cf12SSatish Balay #endif
49fd41dd0SBarry Smith /*
59fd41dd0SBarry Smith        This routine was converted by f2c from Linpack source
69fd41dd0SBarry Smith              linpack. this version dated 08/14/78
79fd41dd0SBarry Smith       cleve moler, university of new mexico, argonne national lab.
89fd41dd0SBarry Smith */
99fd41dd0SBarry Smith #include "petsc.h"
109fd41dd0SBarry Smith 
115615d1e5SSatish Balay #undef __FUNC__
125615d1e5SSatish Balay #define __FUNC__ "Linpack_DGEFA"
139fd41dd0SBarry Smith int Linpack_DGEFA(Scalar *a, int n, int *ipvt)
149fd41dd0SBarry Smith {
1539d66777SBarry Smith     int     i__2, i__3, kp1, nm1, j, k, l,ll,kn,knp1,jn;
166d84be18SBarry Smith     Scalar  t,*aa,*ax,*ay;
176d84be18SBarry Smith     double  tmp,max;
189fd41dd0SBarry Smith 
199fd41dd0SBarry Smith /*     gaussian elimination with partial pivoting */
209fd41dd0SBarry Smith 
21*3a40ed3dSBarry Smith     PetscFunctionBegin;
229fd41dd0SBarry Smith     /* Parameter adjustments */
239fd41dd0SBarry Smith     --ipvt;
2439d66777SBarry Smith     a       -= n + 1;
259fd41dd0SBarry Smith 
269fd41dd0SBarry Smith     /* Function Body */
279fd41dd0SBarry Smith     nm1 = n - 1;
2839d66777SBarry Smith     for (k = 1; k <= nm1; ++k) {
299fd41dd0SBarry Smith 	kp1  = k + 1;
3039d66777SBarry Smith         kn   = k*n;
3139d66777SBarry Smith         knp1 = k*n + k;
329fd41dd0SBarry Smith 
339fd41dd0SBarry Smith /*        find l = pivot index */
349fd41dd0SBarry Smith 
359fd41dd0SBarry Smith 	i__2 = n - k + 1;
3639d66777SBarry Smith         aa = &a[knp1];
379fd41dd0SBarry Smith         max = PetscAbsScalar(aa[0]);
389fd41dd0SBarry Smith         l = 1;
399fd41dd0SBarry Smith         for ( ll=1; ll<i__2; ll++ ) {
409fd41dd0SBarry Smith           tmp = PetscAbsScalar(aa[ll]);
419fd41dd0SBarry Smith           if (tmp > max) { max = tmp; l = ll+1;}
429fd41dd0SBarry Smith         }
439fd41dd0SBarry Smith         l += k - 1;
449fd41dd0SBarry Smith 	ipvt[k] = l;
459fd41dd0SBarry Smith 
4639d66777SBarry Smith 	if (a[l + kn] == 0.) {
47e3372554SBarry Smith 	  SETERRQ(k,0,"Zero pivot");
489fd41dd0SBarry Smith 	}
499fd41dd0SBarry Smith 
509fd41dd0SBarry Smith /*           interchange if necessary */
519fd41dd0SBarry Smith 
529fd41dd0SBarry Smith 	if (l != k) {
5339d66777SBarry Smith 	  t = a[l + kn];
5439d66777SBarry Smith 	  a[l + kn] = a[knp1];
5539d66777SBarry Smith 	  a[knp1] = t;
569fd41dd0SBarry Smith         }
579fd41dd0SBarry Smith 
589fd41dd0SBarry Smith /*           compute multipliers */
599fd41dd0SBarry Smith 
6039d66777SBarry Smith 	t = -1. / a[knp1];
619fd41dd0SBarry Smith 	i__2 = n - k;
6239d66777SBarry Smith         aa = &a[1 + knp1];
639fd41dd0SBarry Smith         for ( ll=0; ll<i__2; ll++ ) {
649fd41dd0SBarry Smith           aa[ll] *= t;
659fd41dd0SBarry Smith         }
669fd41dd0SBarry Smith 
679fd41dd0SBarry Smith /*           row elimination with column indexing */
689fd41dd0SBarry Smith 
6939d66777SBarry Smith 	ax = aa;
709fd41dd0SBarry Smith         for (j = kp1; j <= n; ++j) {
7139d66777SBarry Smith             jn = j*n;
7239d66777SBarry Smith 	    t = a[l + jn];
739fd41dd0SBarry Smith 	    if (l != k) {
7439d66777SBarry Smith 	      a[l + jn] = a[k + jn];
7539d66777SBarry Smith 	      a[k + jn] = t;
769fd41dd0SBarry Smith             }
779fd41dd0SBarry Smith 
789fd41dd0SBarry Smith 	    i__3 = n - k;
7939d66777SBarry Smith             ay = &a[1+k+jn];
809fd41dd0SBarry Smith             for ( ll=0; ll<i__3; ll++ ) {
819fd41dd0SBarry Smith               ay[ll] += t*ax[ll];
829fd41dd0SBarry Smith             }
839fd41dd0SBarry Smith 	}
849fd41dd0SBarry Smith     }
859fd41dd0SBarry Smith     ipvt[n] = n;
869fd41dd0SBarry Smith     if (a[n + n * n] == 0.) {
87e3372554SBarry Smith 	SETERRQ(n,0,"Zero pivot,final row");
889fd41dd0SBarry Smith     }
89*3a40ed3dSBarry Smith     PetscFunctionReturn(0);
909fd41dd0SBarry Smith }
919fd41dd0SBarry Smith 
92