xref: /petsc/src/mat/impls/baij/seq/dgedi.c (revision b2002411cd4bf8262b8aeca5b64ca630feca1e56)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: dgedi.c,v 1.7 1997/10/19 03:26:08 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(Scalar *a,int n,int *ipvt,Scalar *work)
19 {
20     int     i__2,kb, kp1, nm1,i, j, k, l, ll,kn,knp1,jn1;
21     Scalar  t, *aa,*ax,*ay,tmp;
22 
23     PetscFunctionBegin;
24     --work;
25     --ipvt;
26     a       -= n + 1;
27 
28    /*     compute inverse(u) */
29 
30     for (k = 1; k <= n; ++k) {
31         kn           = k*n;
32         knp1         = kn + k;
33 	a[knp1]      = 1.0 / a[knp1];
34 	t            = -a[knp1];
35 	i__2         = k - 1;
36         aa           = &a[1 + kn];
37         for ( ll=0; ll<i__2; ll++ ) aa[ll] *= t;
38 	kp1 = k + 1;
39 	if (n < kp1) continue;
40         ax = aa;
41         for (j = kp1; j <= n; ++j) {
42             jn1 = j*n;
43 	    t = a[k + jn1];
44 	    a[k + jn1] = 0.;
45             ay = &a[1 + jn1];
46             for ( ll=0; ll<k; ll++ ) {
47               ay[ll] += t*ax[ll];
48             }
49 	}
50     }
51 
52    /*    form inverse(u)*inverse(l) */
53 
54     nm1 = n - 1;
55     if (nm1 < 1) {
56 	PetscFunctionReturn(0);
57     }
58     for (kb = 1; kb <= nm1; ++kb) {
59 	k   = n - kb;
60         kn  = k*n;
61 	kp1 = k + 1;
62         aa  = a + kn;
63 	for (i = kp1; i <= n; ++i) {
64 	    work[i] = aa[i];
65 	    aa[i]   = 0.;
66 	}
67 	for (j = kp1; j <= n; ++j) {
68 	    t = work[j];
69             ax = &a[j * n + 1];
70             ay = &a[kn + 1];
71             for ( ll=0; ll<n; ll++ ) {
72               ay[ll] += t*ax[ll];
73             }
74 	}
75 	l = ipvt[k];
76 	if (l != k) {
77             ax = &a[kn + 1];
78             ay = &a[l * n + 1];
79             for ( ll=0; ll<n; ll++ ) {
80               tmp    = ax[ll];
81               ax[ll] = ay[ll];
82               ay[ll] = tmp;
83             }
84 	}
85     }
86     PetscFunctionReturn(0);
87 }
88 
89