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