xref: /petsc/src/mat/impls/baij/seq/dgedi.c (revision 3f1db9ec2fd39765c6c3a00831044586630c4cca)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: dgedi.c,v 1.8 1998/04/21 17:47:46 bsmith Exp bsmith $";
3 #endif
4 
5 /*
6               This file creating by running f2c
7             linpack. this version dated 08/14/78
8       cleve moler, university of new mexico, argonne national lab.
9 
10       Computes the inverse of a matrix given its factors and pivots
11     calculated by Linpack_DGEFA().
12 */
13 
14 #include "petsc.h"
15 
16 #undef __FUNC__
17 #define __FUNC__ "Linpack_DGEDI"
18 int Linpack_DGEDI(MatScalar *a,int n,int *ipvt,MatScalar *work)
19 {
20     int        i__2,kb, kp1, nm1,i, j, k, l, ll,kn,knp1,jn1;
21     MatScalar  *aa,*ax,*ay,tmp;
22     MatScalar  t;
23 
24     PetscFunctionBegin;
25     --work;
26     --ipvt;
27     a       -= n + 1;
28 
29    /*     compute inverse(u) */
30 
31     for (k = 1; k <= n; ++k) {
32         kn           = k*n;
33         knp1         = kn + k;
34 	a[knp1]      = 1.0 / a[knp1];
35 	t            = -a[knp1];
36 	i__2         = k - 1;
37         aa           = &a[1 + kn];
38         for ( ll=0; ll<i__2; ll++ ) aa[ll] *= t;
39 	kp1 = k + 1;
40 	if (n < kp1) continue;
41         ax = aa;
42         for (j = kp1; j <= n; ++j) {
43             jn1 = j*n;
44 	    t = a[k + jn1];
45 	    a[k + jn1] = 0.;
46             ay = &a[1 + jn1];
47             for ( ll=0; ll<k; ll++ ) {
48               ay[ll] += t*ax[ll];
49             }
50 	}
51     }
52 
53    /*    form inverse(u)*inverse(l) */
54 
55     nm1 = n - 1;
56     if (nm1 < 1) {
57 	PetscFunctionReturn(0);
58     }
59     for (kb = 1; kb <= nm1; ++kb) {
60 	k   = n - kb;
61         kn  = k*n;
62 	kp1 = k + 1;
63         aa  = a + kn;
64 	for (i = kp1; i <= n; ++i) {
65 	    work[i] = aa[i];
66 	    aa[i]   = 0.;
67 	}
68 	for (j = kp1; j <= n; ++j) {
69 	    t = work[j];
70             ax = &a[j * n + 1];
71             ay = &a[kn + 1];
72             for ( ll=0; ll<n; ll++ ) {
73               ay[ll] += t*ax[ll];
74             }
75 	}
76 	l = ipvt[k];
77 	if (l != k) {
78             ax = &a[kn + 1];
79             ay = &a[l * n + 1];
80             for ( ll=0; ll<n; ll++ ) {
81               tmp    = ax[ll];
82               ax[ll] = ay[ll];
83               ay[ll] = tmp;
84             }
85 	}
86     }
87     PetscFunctionReturn(0);
88 }
89 
90